“nbsphinx-toctree”: { “maxdepth”: 2 }

Introduction to BLUES

In this Jupyter Notebook, we will cover the following topics:

  • YAML configuration

  • Setting up a system for BLUES

  • Advanced options (HMR, restraints, freezing)

  • Configuring reporters

  • Running a BLUES simulation

Background

Coupling MD simulations with random NCMC moves for enhanced sampling of ligand binding modes via BLUES

BLUES is an approach that combines molecular dynamics (MD) simulations and the Non-equilibrium Candidate Monte Carlo (NCMC) framework to enhance ligand binding mode sampling (Github)

During a MD simulation, BLUES will perform a random rotation of the bound ligand and then allow the system to relax through alchemically scaling off/on the ligand-receptor interactions. BLUES enables us to sample alternative ligand binding modes, that would normally take very long simulations to capture in a traditional MD simulation because of the large gap in timescales between atomistic motions and biological motions.

YAML Configuration

BLUES can be configured in either pure python through dictionaries of the appropriate parameters or using a YAML file (which is converted to a dict under the hood). Below we will walk through the keywords for configuring BLUES. An example YAML configuration file can be found in rotmove_cuda.yaml.

Note: Code blocks in this notebook denoted with —— indicate a section from a YAML configuration file.

Input/Output

------
output_dir: .
outfname: t4-toluene
logger_level: info #critical, error, warning, info, debug
-----

Output files

Specify the directory you want all the simulation output files to be saved to with ``output_dir``. By default, BLUES will save them in the current directory that you’re running BLUES in. The parameter ``outfname`` will be used for the filename prefix for all output files (e.g. t4-toluene.nc, t4-toleune.log). The level of verbosity can be controlled by ``logger_level``. The default ``logger_level`` is set to info and valid choices are critical, error, warning, info or debug.

Input files to generate the structure of your system

  • Input a Parameter/topology file and a Coordinate file, which will be used to generate the ParmEd Structure.

  • The ParmEd Structure is a chemical structure composed of atoms, bonds, angles, torsions, and other topological features.

To see a full list of supported file formats: https://parmed.github.io/ParmEd/html/readwrite.html

----------------
structure:
  filename: tests/data/eqToluene.prmtop
  xyz: tests/data/eqToluene.inpcrd
  restart: t4-toluene_2.rst7
----------------

BLUES simulations all begin from a parmed.Structure. Keywords nested under ``structure`` are for generating the parmed.Structure by calling parmed.load_file() under the hood. The ``filename`` keyword should point to the file containing the parameters for your system. For example, if coming from AMBER you would specify the .prmtop file. The ``xyz`` keyword is intended for specifying the coordinates of your system (e.g. .inpcrd or .pdb).

Restart files

BLUES supports “soft” restarts of simulations from AMBER restart files ``.rst7`` via parmed.amber.Rst7. “Soft” restart implies the simulation will begin from the saved positions, velocities, and box vectors but does not store any internal data such as the states of random number generators. It should be noted that velocities are re-initialized at every BLUES iteration, so storing these is not so important.

System configuration

From the YAML file, there is a section dedicated for generating an openmm.System from a parmed.Structure. For definitions of ``system`` keywords and valid options see parmed.Structure.createSystem()

Below we provide an example for generating a system in a cubic box with explicit solvent.

--------
system:
  nonbondedMethod: PME
  nonbondedCutoff: 10 * angstroms
  constraints: HBonds
  rigidWater: True
  removeCMMotion: True
  ewaldErrorTolerance: 0.005
  flexibleConstraints: True
  splitDihedrals: False
---------

(Optional) Hydrogen mass repartitioning

BLUES has the option to use the hydrogen mass repartitioning scheme HMR to allow use of longer time steps in the simulation. Simply provide the keyword ``hydrogenMass`` in the YAML file like below:

--------
system:
  nonbondedMethod: PME
  nonbondedCutoff: 10 * angstroms
  constraints: HBonds
  rigidWater: True
  removeCMMotion: True
  ewaldErrorTolerance: 0.005
  flexibleConstraints: True
  splitDihedrals: False
  hydrogenMass: 3.024 * daltons
---------

If using HMR, you can set the timestep (``dt``) for the simulation to 4fs and ``constraints`` should be set to either ``HBonds`` or ``AllBonds``.

(Optional) Alchemical system configuration

Nested under the system parameters, you can modify parameters for the alchemical system. Below are the default settings and are not required to be specified in the YAML configuration. Modifications to these parameters are for advanced users.

-------
system:
  nonbondedMethod: PME
  nonbondedCutoff: 10 * angstroms
  constraints: HBonds
  rigidWater: True
  removeCMMotion: True
  ewaldErrorTolerance: 0.005
  flexibleConstraints: True
  splitDihedrals: False
  alchemical:
    # Sterics
    softcore_alpha: 0.5
    softcore_a : 1
    softcore_b : 1
    softcore_c : 6

    # Electrostatics
    softcore_beta : 0.0
    softcore_d : 1
    softcore_e : 1
    softcore_f : 2

    annihilate_electrostatics : True
    annihilate_sterics : False
-------

For further details on alchemical parameters see: http://getyank.org/0.16.2/yamlpages/options.html

Simulation Configuration

The keywords for configuring the Simulations for BLUES are explained below:

  • ``dt``: timestep

  • ``nIter``: number of iterations or proposed moves

  • ``nstepsMD``: number of MD steps

  • ``nstepsNC``: number NCMC steps

The configuration below will run BLUES in NVT with 2fs timesteps for 10 iterations. The MD and NCMC simulation will run 10,000 steps per iteration.

----------
simulation:
  platform: CUDA
  dt: 0.002 * picoseconds
  friction: 1 * 1/picoseconds
  temperature: 300 * kelvin
  nIter: 10
  nstepsMD: 10000
  nstepsNC: 10000
----------

NPT Simulation

To run BLUES in NPT, simply specify a ``pressure``:

----------
simulation:
  platform: CUDA
  dt: 0.002 * picoseconds
  friction: 1 * 1/picoseconds
  temperature: 300 * kelvin
  nIter: 10
  nstepsMD: 10000
  nstepsNC: 10000
  pressure: 1 * atmospheres
----------

(Optional) Additional relaxation steps in the NCMC simulation

Keywords ``nprop`` and ``propLambda`` allow you to add additional relxation steps between a set range in the lambda schedule for the alchemical process in the NCMC simulation. Setting ``propLambda`` to 0.3 will select a lambda range of -/+ 0.3 from the midpoint (0.5), giving [0.2, 0.8]. During the alchemical process, when lambda is between 0.2 to 0.8, ``nprop`` controls the number of additional relaxation steps to add at each lambda step (change in lambda). Additional relaxation steps has been show to increase acceptance proposed NCMC moves.

----------
simulation:
  platform: CUDA
  dt: 0.002 * picoseconds
  friction: 1 * 1/picoseconds
  temperature: 300 * kelvin
  nIter: 10
  nstepsMD: 10000
  nstepsNC: 10000
  pressure: 1 * atmospheres
  nprop: 3
  propLambda: 0.3
----------

(Optional) Platform Properties

If you need to modify platform properties for the simulation, you can set the keyword ``properties`` like below:

Example: OpenCL in single precision on GPU device 2

Note: works for running on the GPU on MacBook Pro 2017

----------
simulation:
  platform: OpenCL
  properties:
    OpenCLPrecision: single
    OpenCLDeviceIndex: 2
  dt: 0.002 * picoseconds
  friction: 1 * 1/picoseconds
  temperature: 300 * kelvin
  nIter: 10
  nstepsMD: 10000
  nstepsNC: 10000
  pressure: 1 * atmospheres
----------

Example: CUDA in double precision on GPU device 0

----------
simulation:
  platform: CUDA
  properties:
    CudaPrecision: double
    CudaDeviceIndex: 0
  dt: 0.002 * picoseconds
  friction: 1 * 1/picoseconds
  temperature: 300 * kelvin
  nIter: 10
  nstepsMD: 10000
  nstepsNC: 10000
  pressure: 1 * atmospheres
----------

Reporter Configuration

We provide functionality to configure a recommended set of reporters from the YAML file. These are used to record information for either the MD or NCMC simulation. Below are the keywords for each reporter. Each reporter will require the ``reportInterval`` keyword to specify the frequency to store the simulation data: - ``state`` : State data reporter for OpenMM simulations, but it is a little more generalized. Writes to a .ene file. - For full list of parameters see parmed.openmm.reporters.StateDataReporter - ``traj_netcdf`` : Customized AMBER NetCDF (.nc) format reporter - ``restart`` : Restart AMBER NetCDF (.rst7) format reporter - ``progress`` : Write to a file (.prog), the progress report of how many steps has been done, how fast the simulation is running, and how much time is left (similar to the mdinfo file in Amber). File is overwritten at each reportInterval. - For full list of parameters see parmed.openmm.reporters.ProgressReporter - ``stream`` : Customized version of openmm.app.StateDataReporter. This will instead stream/print the information to the terminal as opposed to writing to a file. - takes the same parameters as the openmm.app.StateDataReporter

To attach them to the MD simulation. You nest the reporter keywords under the keyword ``md_reporters`` like below. To attach the reporters to NCMC simulation, use the ``ncmc_reporters`` keyword instead.

------
md_reporters:
  state:
    reportInterval: 250
  traj_netcdf:
    reportInterval: 250
  restart:
    reportInterval: 1000
  progress:
    totalSteps: 10000
    reportInterval: 10
  stream:
    title: md
    reportInterval: 250
    totalSteps: 10000
    step: True
    speed: True
    progress: True
    remainingTime: True
----

In the above example, we are using the ``stream`` reporter to print the speeds on the intergrator at regular intervals. This may be a bit redudant with the ``progress`` reporter if are running the job remotely and don’t need the information streamed to terminal.

(Optional) Advanced options to the ``traj_netcdf`` reporter

The ``traj_netcdf`` reporter can store additional information that may be useful for the NCMC simulation or record at specific frames. In the example below, we will store the first, midpoint (when the move is applied), and last frame of each NCMC iteration, along with the ``alchemicalLambda`` step and the ``protocolWork``.

----
ncmc_reporters:
  traj_netcdf:
    frame_indices: [1, 0.5, -1]
    alchemicalLambda: True
    protocolWork: True
-----

To access the numerical data stored in the NetCDF file:

from netCDF4 import Dataset

f = Dataset("t4-toluene-ncmc.nc")
print(f.variables['alchemicalLambda'][:])
print(f.variables['protocolWork'][:])

>>> [ 0.001  0.5    1.   ]
>>> [  0.03706791  30.72696877  25.708498  ]

Running a BLUES simulation

Below we will provide an example for running an NPT BLUES simulation which applies random rotational moves to the toluene ligand in T4-lysozyme from a YAML configuration.

[1]:
yaml_cfg = """
output_dir: .
outfname: t4-toluene
logger_level: info

structure:
  filename: ../blues/tests/data/eqToluene.prmtop
  xyz: ../blues/tests/data/eqToluene.inpcrd

system:
  nonbondedMethod: PME
  nonbondedCutoff: 10 * angstroms
  constraints: HBonds
  rigidWater: True
  removeCMMotion: True
  hydrogenMass: 3.024 * daltons
  ewaldErrorTolerance: 0.005
  flexibleConstraints: True
  splitDihedrals: False

freeze:
  freeze_center: ':LIG'
  freeze_solvent: ':WAT,Cl-'
  freeze_distance: 5 * angstroms

simulation:
  platform: CUDA
  properties:
    CudaPrecision: single
    CudaDeviceIndex: 0
  dt: 0.004 * picoseconds
  friction: 1 * 1/picoseconds
  pressure: 1 * atmospheres
  temperature: 300 * kelvin
  nIter: 5
  nstepsMD: 1000
  nstepsNC: 1000

md_reporters:
  state:
    reportInterval: 250
  traj_netcdf:
    reportInterval: 250
  restart:
    reportInterval: 1000
  stream:
    title: md
    reportInterval: 250
    totalSteps: 5000 # nIter * nstepsMD
    step: True
    speed: True
    progress: True
    remainingTime: True

ncmc_reporters:
  stream:
    title: ncmc
    reportInterval: 250
    totalSteps: 1000 # Use nstepsNC
    step: True
    speed: True
    progress: True
    remainingTime: True
"""

Import the following BLUES modules required for the following steps.

  • Make sure to specify the type of move that you want to import from blues.moves

  • Available moves can be veiwed in moves.py, supported/tested moves include:

    • RandomLigandRotationMove

    • SideChainMove

[2]:
from blues.moves import RandomLigandRotationMove
from blues.engine import MoveEngine
from blues.simulation import *
from blues.settings import *
[3]:
#Read in the YAML file
cfg = Settings(yaml_cfg).asDict()
#Shortcut to access `parmed.Structure` from dict
structure = cfg['Structure']
./t4-toluene

Below is what the resulting configuration dictionary looks like (formatted into JSON for readability).

[4]:
import json
print(json.dumps(cfg, sort_keys=True, indent=2, skipkeys=True, default=str))
{
  "Logger": "<logging.RootLogger object at 0x7f96b1bff940>",
  "Structure": "../blues/tests/data/eqToluene.prmtop",
  "freeze": {
    "freeze_center": ":LIG",
    "freeze_distance": "5.0 A",
    "freeze_solvent": ":WAT,Cl-"
  },
  "logger_level": "info",
  "md_reporters": [
    "<parmed.openmm.reporters.StateDataReporter object at 0x7f966832c128>",
    "<blues.reporters.NetCDF4Reporter object at 0x7f966574e898>",
    "<parmed.openmm.reporters.RestartReporter object at 0x7f966574e860>",
    "<blues.reporters.BLUESStateDataReporter object at 0x7f966574e828>"
  ],
  "ncmc_reporters": [
    "<blues.reporters.BLUESStateDataReporter object at 0x7f966574e908>"
  ],
  "outfname": "./t4-toluene",
  "output_dir": ".",
  "simulation": {
    "dt": "0.004 ps",
    "friction": "1.0 /ps",
    "md_trajectory_interval": 250,
    "moveStep": 500,
    "nIter": 5,
    "nprop": 1,
    "nstepsMD": 1000,
    "nstepsNC": 1000,
    "outfname": "./t4-toluene",
    "platform": "CUDA",
    "pressure": "1.0 atm",
    "propLambda": 0.3,
    "propSteps": 1000,
    "properties": {
      "CudaDeviceIndex": 0,
      "CudaPrecision": "single"
    },
    "temperature": "300.0 K",
    "verbose": false
  },
  "structure": {
    "filename": "../blues/tests/data/eqToluene.prmtop",
    "xyz": "../blues/tests/data/eqToluene.inpcrd"
  },
  "system": {
    "constraints": "HBonds",
    "ewaldErrorTolerance": 0.005,
    "flexibleConstraints": true,
    "hydrogenMass": "3.024 Da",
    "nonbondedCutoff": "10.0 A",
    "nonbondedMethod": "PME",
    "removeCMMotion": true,
    "rigidWater": true,
    "splitDihedrals": false,
    "verbose": false
  },
  "verbose": false
}

Selecting a move and initialize the move engine

Here, we initialize the ``RandomLigandRotationMove`` from the ``blues.moves`` module which proposes random rotations on the toluene ligand. We select the toluene ligand by providing the residue name LIG and the parmed.Structure to select the atoms from. If we begin BLUES from our YAML configuration, the parmed.Structure for our system is generated from the call to startup(). We can access it at the top level with cfg['Structure']

After initialization of the selected move, we pass the move object to the ``MoveEngine`` from the ``blues.engine`` module. The ``MoveEngine`` controls what types of moves will be performed during the NCMC protocol and with a given probability. This will be more useful when we use multiple move types.

[5]:
#Initialize the move class and pass it to the engine
ligand = RandomLigandRotationMove(structure, 'LIG')
ligand_mover = MoveEngine(ligand)

Generating the Systems for openMM: ``SystemFactory``

Next, we must generate the openmm.System from the parmed.Structure by calling the ``SystemFactory`` class from the blues.simulation module. The class must be initialized by providing 3 required arguments:

structure : parmed.Structure
    A chemical structure composed of atoms, bonds, angles, torsions, and
    other topological features.
atom_indices : list of int
    Atom indicies of the move or designated for which the nonbonded forces
    (both sterics and electrostatics components) have to be alchemically
    modified.
config : dict, parameters for generating the `openmm.System for the MD
    and NCMC simulation. For complete parameters, see docs for `generateSystem`
    and `generateAlchSystem`

Upon initialization, this class will create the system for the MD simulation and the NCMC simulation. They can be accessed through the attributes ``systems.md`` or ``systems.alch``. Any modifications to either of these Systems should be done within the context of this object. Once the systems are passed into an openmm.Simulation, you will not be able to modify the system easily.

[6]:
systems = SystemFactory(structure, ligand.atom_indices, cfg['system'])
INFO: Adding bonds...
INFO: Adding angles...
INFO: Adding dihedrals...
INFO: Adding Ryckaert-Bellemans torsions...
INFO: Adding Urey-Bradleys...
INFO: Adding improper torsions...
INFO: Adding CMAP torsions...
INFO: Adding trigonal angle terms...
INFO: Adding out-of-plane bends...
INFO: Adding pi-torsions...
INFO: Adding stretch-bends...
INFO: Adding torsion-torsions...
INFO: Adding Nonbonded force...

(Optional) Applying restraints or freezing atoms

The ``SystemFactory`` class also provides functionality for restraining or freezing the atoms. Use extreme caution when freezing/restraining atoms. You should consider if freezing/restraining should be applied to BOTH the MD and alchemical system.

Selections for either restraining or freezing atoms in your system use the Amber mask syntax.

Positional restraints: ``SystemFactory.restrain_positions()``

To apply positional restraints, you can call ``SystemFactory.restrain_positions()``. You can specify the parameters/selection for applying positional restraints in the YAML file.

------
restraints:
    selection: '@CA,C,N'
    weight: 5
------

restrain keywords: - ``selection``: Specify what to apply positional restraints to using Amber mask syntax. Default = ‘@CA,C,N’ - ``weight``: Restraint weight for xyz atom restraints in kcal/(mol A^2). Default = 5

From the YAML example above, we would be applying positional restraints to the backbone atoms of the protein. If applying restraints, you most likely will want to apply it to BOTH the MD and alchemical systems like below:

systems.md = SystemFactory.restrain_positions(structure, systems.md, **cfg['restraints'])
systems.alch = SystemFactory.restrain_positions(structure, systems.alch, **cfg['restraints'])

Freezing selected atoms: ``SystemFactory.freeze_atoms()``

To freeze a selection, call ``SystemFactory.freeze_atoms()``. Atoms that have a mass of zero will be ignored by the integrator and will not change positions during the simulation, effectively they are frozen.

To freeze atoms using a given selection string. Use the keyword ``freeze_selection``.

-----
freeze:
  freeze_selection: ':LIG'
------

From the YAML example above, we would be freezing only the atoms belonging to the residue LIG. Although freezing the ligand in this example wouldn’t be very useful. It would be applied like

systems.md = SystemFactory.freeze_atoms(structure, systems.md, **cfg['freeze'])

Freezing atoms around a selection: ``SystemFactory.freeze_radius()``

Alternatively, you can choose to freeze atoms around a given selection. To do so, call ``SystemFactory.freeze_radius()``. For example, you may want to freeze atoms that are 5 angstroms away from the ligand and include the solvent.

-----
freeze:
  freeze_center: ':LIG'
  freeze_solvent: ':WAT,Cl-'
  freeze_distance: 5 * angstroms
------
  • ``freeze_center``: Specifies the center of the object for freezing, masses will be zeroed. Default = ‘:LIG’

  • ``freeze_solvent``: select which solvent atoms should have their masses zeroed. Default = ‘:HOH,NA,CL’

  • ``freeze_distance``: Distance ( in angstroms) to select atoms for retaining their masses. Atoms outside the set distance will have their masses set to 0.0. Default = 5.0

We often utilize this type of freezing to speed up the alchemical process during the NCMC simulation while leaving them completely free in the MD simulation for proper relaxation.

systems.alch = SystemFactory.freeze_radius(structure, systems.alch, **cfg['freeze'])

In this notebook example, our YAML config indicates we will be freezing around the ligand (keyword: ``freeze_center``). So we will call the ``freeze_radius`` function.

[7]:
systems.alch = SystemFactory.freeze_radius(structure, systems.alch, **cfg['freeze'])
INFO: Freezing 22065 atoms 5.0 Angstroms from ':LIG' on <simtk.openmm.openmm.System; proxy of <Swig Object of type 'OpenMM::System *' at 0x7f96aa74fab0> >

Generating the OpenMM Simulations: ``SimulationFactory``

Now that we have generated our openmm.System for the MDand frozen the solvent around the ligand. We are now ready to create the set of simulations for running BLUES. We do this by calling the ``SimulationFactory`` class . The expected parameters are:

systems : blues.simulation.SystemFactory object
    The object containing the MD and alchemical openmm.Systems
move_engine : blues.engine.MoveEngine object
    MoveProposal object which contains the dict of moves performed
    in the NCMC simulation.
config : dict of parameters for the simulation (i.e timestep, temperature, etc.)
md_reporters : list of Reporter objects for the MD openmm.Simulation
ncmc_reporters : list of Reporter objects for the NCMC openmm.Simulation

If you wish to use your own openmm reporters, simply pass them into the arguments as a list of Reporter objects. Since we have configured our reporters from the YAML file, we can pass them into the arguments ``md_reporters`` and ``ncmc_reporters``.

[8]:
# List of MD reporters
cfg['md_reporters']
[8]:
[<parmed.openmm.reporters.StateDataReporter at 0x7f966832c128>,
 <blues.reporters.NetCDF4Reporter at 0x7f966574e898>,
 <parmed.openmm.reporters.RestartReporter at 0x7f966574e860>,
 <blues.reporters.BLUESStateDataReporter at 0x7f966574e828>]
[9]:
# List of NCMC reporters
cfg['ncmc_reporters']
[9]:
[<blues.reporters.BLUESStateDataReporter at 0x7f966574e908>]
[10]:
simulations = SimulationFactory(systems, ligand_mover, cfg['simulation'],
                                cfg['md_reporters'], cfg['ncmc_reporters'])
INFO: Adding MonteCarloBarostat with 1.0 atm. MD simulation will be 300.0 K NPT.
WARNING: NCMC simulation will NOT have pressure control. NCMC will use pressure from last MD state.
INFO: OpenMM(7.1.1.dev-c1a64aa) simulation generated for CUDA platform
system = Linux
node = titanpascal
release = 4.13.0-41-generic
version = #46~16.04.1-Ubuntu SMP Thu May 3 10:06:43 UTC 2018
machine = x86_64
processor = x86_64
DeviceIndex = 0
DeviceName = TITAN Xp
UseBlockingSync = true
Precision = single
UseCpuPme = false
CudaCompiler = /usr/local/cuda-8.0/bin/nvcc
TempDirectory = /tmp
CudaHostCompiler =
DisablePmeStream = false
DeterministicForces = false

If you would like to access the MD or NCMC simulation. You can access them as attributes to the ``SimulationFactory`` class with ``simulations.md`` or ``simulations.ncmc``. This will allow you to do things like energy minimize the system or run a few steps of regular dynamics before running the hybrid (MD+NCMC) BLUES approach. The NCMC simulation will automatically be synced to the state of the MD simulation when running the BLUES simulation.

[11]:
# Energy minimization
state = simulations.md.context.getState(getPositions=True, getEnergy=True)
print('Pre-Minimized energy = {}'.format(state.getPotentialEnergy().in_units_of(unit.kilocalorie_per_mole)))

simulations.md.minimizeEnergy(maxIterations=0)
state = simulations.md.context.getState(getPositions=True, getEnergy=True)
print('Minimized energy = {}'.format(state.getPotentialEnergy().in_units_of(unit.kilocalorie_per_mole)))
Pre-Minimized energy = -69057.34671058532 kcal/mol
Minimized energy = -87007.38938198877 kcal/mol
[12]:
# Running only the MD simulation
simulations.md.step(500)
#"Progress (%)" "Step"  "Speed (ns/day)"        "Time Remaining"
md: 5.0%        250     0       --
md: 10.0%       500     271     0:05

Run the BLUES Simulation

To run the full BLUES simulation, where we apply NCMC moves and follow-up with the MD simulation, we simply pass the ``SimulationFactory`` object to the ``Simulation`` class and call the ``run()`` function which takes ``nIter``, ``nstepsNC``, ``nstepsMD`` as arguments (or we can pass it the simulation configuration from the YAML on initialization of the class).

[13]:
blues = BLUESSimulation(simulations, cfg['simulation'])
blues.run()
INFO: Total BLUES Simulation Time = 40.0 ps (8.0 ps/Iter)
Total Force Evaluations = 10000
Total NCMC time = 20.0 ps (4.0 ps/iter)
Total MD time = 20.0 ps (4.0 ps/iter)
Trajectory Interval = 2.0 ps/frame (4.0 frames/iter)
INFO: Running 5 BLUES iterations...
INFO: BLUES Iteration: 0
INFO: Advancing 1000 NCMC switching steps...
#"Progress (%)" "Step"  "Speed (ns/day)"        "Time Remaining"
ncmc: 25.0%     250     0       --
ncmc: 50.0%     500     110     0:01
Performing RandomLigandRotationMove...
ncmc: 75.0%     750     93.3    0:00
ncmc: 100.0%    1000    93      0:00
NCMC MOVE REJECTED: work_ncmc -17.61007128922719 < -3.476065340494845
Advancing 1000 MD steps...
md: 15.0%       750     31.9    0:46
md: 20.0%       1000    43.1    0:32
md: 25.0%       1250    54.8    0:23
md: 30.0%       1500    65.5    0:18
BLUES Iteration: 1
Advancing 1000 NCMC switching steps...
ncmc: 25.0%     250     57.2    --
ncmc: 50.0%     500     62.6    0:13
Performing RandomLigandRotationMove...
ncmc: 75.0%     750     65.4    0:03
ncmc: 100.0%    1000    69.1    0:00
NCMC MOVE REJECTED: work_ncmc -28.930984747001787 < -2.602259467723195
Advancing 1000 MD steps...
md: 35.0%       1750    45.5    0:24
md: 40.0%       2000    50.5    0:20
md: 45.0%       2250    56.3    0:16
md: 50.0%       2500    61.9    0:13
BLUES Iteration: 2
Advancing 1000 NCMC switching steps...
ncmc: 25.0%     250     57.9    --
ncmc: 50.0%     500     60.1    0:25
Performing RandomLigandRotationMove...
ncmc: 75.0%     750     62.4    0:06
ncmc: 100.0%    1000    64.6    0:00
NCMC MOVE REJECTED: work_ncmc -9.178847171172448 < -0.3081492900307722
Advancing 1000 MD steps...
md: 55.0%       2750    49.7    0:15
md: 60.0%       3000    53.1    0:13
md: 65.0%       3250    56.7    0:10
md: 70.0%       3500    60.4    0:08
BLUES Iteration: 3
Advancing 1000 NCMC switching steps...
ncmc: 25.0%     250     58.4    --
ncmc: 50.0%     500     60.4    0:37
Performing RandomLigandRotationMove...
ncmc: 75.0%     750     61.6    0:09
ncmc: 100.0%    1000    63.6    0:00
NCMC MOVE REJECTED: work_ncmc -6.022089748510081 < -0.28447339331708726
Advancing 1000 MD steps...
md: 75.0%       3750    52.6    0:08
md: 80.0%       4000    55      0:06
md: 85.0%       4250    58      0:04
md: 90.0%       4500    60.8    0:02
BLUES Iteration: 4
Advancing 1000 NCMC switching steps...
ncmc: 25.0%     250     58.6    --
ncmc: 50.0%     500     59.8    0:49
Performing RandomLigandRotationMove...
ncmc: 75.0%     750     60.7    0:12
ncmc: 100.0%    1000    61.9    0:00
NCMC MOVE REJECTED: work_ncmc -63.43688730563919 < -0.5417610940298393
Advancing 1000 MD steps...
md: 95.0%       4750    53.3    0:01
md: 100.0%      5000    55.2    0:00
md: 105.0%      5250    57.5    23:59:59
md: 110.0%      5500    59.8    23:59:58
Acceptance Ratio: 0.0
nIter: 5
[ ]: