grand.samplers module

Description

This module is written to execute GCMC moves with water molecules in OpenMM, via a series of Sampler objects.

Marley Samways Ollie Melling

class grand.samplers.BaseGrandCanonicalMonteCarloSampler(system, topology, temperature, ghostFile='gcmc-ghost-wats.txt', log='gcmc.log', dcd=None, rst=None, overwrite=False)

Bases: object

Base class for carrying out GCMC moves in OpenMM. All other Sampler objects are derived from this

__init__(system, topology, temperature, ghostFile='gcmc-ghost-wats.txt', log='gcmc.log', dcd=None, rst=None, overwrite=False)

Initialise the object to be used for sampling water insertion/deletion moves

Parameters
  • system (simtk.openmm.System) – System object to be used for the simulation

  • topology (simtk.openmm.app.Topology) – Topology object for the system to be simulated

  • temperature (simtk.unit.Quantity) – Temperature of the simulation, must be in appropriate units

  • ghostFile (str) – Name of a file to write out the residue IDs of ghost water molecules. This is useful if you want to visualise the sampling, as you can then remove these waters from view, as they are non-interacting. Default is ‘gcmc-ghost-wats.txt’

  • log (str) – Log file to write out

  • dcd (str) – Name of the DCD file to write the system out to

  • rst (str) – Name of the restart file to write out (.pdb or .rst7)

  • overwrite (bool) – Overwrite any data already present

adjustSpecificWater(atoms, new_lambda)

Adjust the coupling of a specific water molecule, by adjusting the lambda value

Parameters
  • atoms (list) – List of the atom indices of the water to be adjusted

  • new_lambda (float) – Value to set lambda to for this particle

customiseForces()

Create a CustomNonbondedForce to handle water-water interactions and modify the original NonbondedForce to ignore water interactions

deleteGhostWaters(ghostResids=None, ghostFile=None)

Switch off nonbonded interactions involving the ghost molecules initially added This function should be executed before beginning the simulation, to prevent any explosions.

Parameters
  • context (simtk.openmm.Context) – Current context of the simulation

  • ghostResids (list) – List of residue IDs corresponding to the ghost waters added

  • ghostFile (str) – File containing residue IDs of ghost waters. Will switch off those on the last line. This will be useful in restarting simulations

Returns

context – Updated context, with ghost waters switched off

Return type

simtk.openmm.Context

getWaterParameters(water_resname='HOH')

Get the non-bonded parameters for each of the atoms in the water model used

Parameters

water_resname (str) – Name of the water residues

Returns

wat_params – List of dictionaries containing the charge, sigma and epsilon for each water atom

Return type

list

getWaterResids(water_resname='HOH')

Get the residue IDs of all water molecules in the system

Parameters

water_resname (str) – Name of the water residues

Returns

resid_list – List of residue ID numbers

Return type

list

getWaterStatusResids(value)

Get a list of resids which have a particular status value

Parameters

value (int) – Value of the water status. 0: ghost, 1: GCMC water, 2: Non-tracked water

Returns

resids – List of residues which match that status

Return type

numpy.array

getWaterStatusValue(resid)

Get the status value of a particular resid

Parameters

resid (int) – Residue to update the status for

Returns

value – Value of the water status. 0: ghost, 1: GCMC water, 2: Non-tracked water

Return type

int

move(context, n=1)

Returns an error if someone attempts to execute a move with the parent object Parameters are designed to match the signature of the inheriting classes

Parameters
  • context (simtk.openmm.Context) – Current context of the simulation

  • n (int) – Number of moves to execute

raiseError(error_msg)

Make it nice and easy to report an error in a consisent way - also easier to manage error handling in future

Parameters

error_msg (str) – Message describing the error

report(simulation)

Function to report any useful data

Parameters

simulation (simtk.openmm.app.Simulation) – Simulation object being used

reset()

Reset counted values (such as number of total or accepted moves) to zero

setWaterStatus(resid, new_value)

Set the status of a perticular water to a particular value

Parameters
  • resid (int) – Residue to update the status for

  • new_value (int) – New value of the water status. 0: ghost, 1: GCMC water, 2: Non-tracked water

writeGhostWaterResids()

Write out a comma-separated list of the residue IDs of waters which are non-interacting, so that they can be removed from visualisations. It is important to execute this function when writing to trajectory files, so that each line in the ghost water file corresponds to a frame in the trajectory

class grand.samplers.GCMCSphereSampler(system, topology, temperature, adams=None, excessChemicalPotential=Quantity(value=- 6.09, unit=kilocalorie / mole), standardVolume=Quantity(value=30.345, unit=angstrom ** 3), adamsShift=0.0, ghostFile='gcmc-ghost-wats.txt', referenceAtoms=None, sphereRadius=None, sphereCentre=None, log='gcmc.log', dcd=None, rst=None, overwrite=False)

Bases: grand.samplers.BaseGrandCanonicalMonteCarloSampler

Base class for carrying out GCMC moves in OpenMM, using a GCMC sphere to sample the system

__init__(system, topology, temperature, adams=None, excessChemicalPotential=Quantity(value=- 6.09, unit=kilocalorie / mole), standardVolume=Quantity(value=30.345, unit=angstrom ** 3), adamsShift=0.0, ghostFile='gcmc-ghost-wats.txt', referenceAtoms=None, sphereRadius=None, sphereCentre=None, log='gcmc.log', dcd=None, rst=None, overwrite=False)

Initialise the object to be used for sampling water insertion/deletion moves

Parameters
  • system (simtk.openmm.System) – System object to be used for the simulation

  • topology (simtk.openmm.app.Topology) – Topology object for the system to be simulated

  • temperature (simtk.unit.Quantity) – Temperature of the simulation, must be in appropriate units

  • adams (float) – Adams B value for the simulation (dimensionless). Default is None, if None, the B value is calculated from the box volume and chemical potential

  • excessChemicalPotential (simtk.unit.Quantity) – Excess chemical potential of the system that the simulation should be in equilibrium with, default is -6.09 kcal/mol. This should be the hydration free energy of water, and may need to be changed for specific simulation parameters.

  • standardVolume (simtk.unit.Quantity) – Standard volume of water - corresponds to the volume per water molecule in bulk. The default value is 30.345 A^3

  • adamsShift (float) – Shift the B value from Bequil, if B isn’t explicitly set. Default is 0.0

  • ghostFile (str) – Name of a file to write out the residue IDs of ghost water molecules. This is useful if you want to visualise the sampling, as you can then remove these waters from view, as they are non-interacting. Default is ‘gcmc-ghost-wats.txt’

  • referenceAtoms (list) – List containing dictionaries describing the atoms to use as the centre of the GCMC region Must contain ‘name’ and ‘resname’ as keys, and optionally ‘resid’ (recommended) and ‘chain’ e.g. [{‘name’: ‘C1’, ‘resname’: ‘LIG’, ‘resid’: ‘123’}]

  • sphereRadius (simtk.unit.Quantity) – Radius of the spherical GCMC region

  • sphereCentre (simtk.unit.Quantity) – Coordinates around which the GCMC sphere is based

  • log (str) – Log file to write out

  • dcd (str) – Name of the DCD file to write the system out to

  • rst (str) – Name of the restart file to write out (.pdb or .rst7)

  • overwrite (bool) – Overwrite any data already present

deleteRandomWater()

Choose a random water to be deleted

Returns

  • delete_water (int) – Resid of the water to delete

  • atom_indices (list) – List of the atom IDs for this molecule

deleteWatersInGCMCSphere()

Function to delete all of the waters currently present in the GCMC region This may be useful the plan is to generate a water distribution for this region from scratch. If so, it would be recommended to interleave the GCMC sampling with coordinate propagation, as this will converge faster.

Parameters

context (simtk.openmm.Context) – Current context of the system. Only needs to be supplied if the context has changed since the last update

Returns

context – Updated context after deleting the relevant waters

Return type

simtk.openmm.Context

getReferenceAtomIndices(ref_atoms)

Get the index of the atom used to define the centre of the GCMC box

Parameters

ref_atoms (list) – List of dictionaries containing the atom name, residue name and (optionally) residue ID and chain, as marked by keys ‘name’, ‘resname’, ‘resid’ and ‘chain’

Returns

atom_indices – Indices of the atoms chosen

Return type

list

getSphereCentre()

Update the coordinates of the sphere centre Need to make sure it isn’t affected by the reference atoms being split across PBCs

initialise(context, ghostResids=[])

Prepare the GCMC sphere for simulation by loading the coordinates from a Context object.

Parameters
  • context (simtk.openmm.Context) – Current context of the simulation

  • ghostResids (list) – List of residue IDs corresponding to the ghost waters added

insertRandomWater()

Translate a random ghost to a random point in the GCMC sphere to allow subsequent insertion

Returns

  • new_positions (simtk.unit.Quantity) – Positions following the ‘insertion’ of the ghost water

  • insert_water (int) – Residue ID of the water to insert

  • atom_indices (list) – List of the atom IDs for this molecule

report(simulation)

Function to report any useful data

Parameters

simulation (simtk.openmm.app.Simulation) – Simulation object being used

updateGCMCSphere(state)

Update the relevant GCMC-sphere related parameters. This also involves monitoring which water molecules are in/out of the region

Parameters

state (simtk.openmm.State) – Current State

class grand.samplers.GCMCSystemSampler(system, topology, temperature, adams=None, excessChemicalPotential=Quantity(value=- 6.09, unit=kilocalorie / mole), standardVolume=Quantity(value=30.345, unit=angstrom ** 3), adamsShift=0.0, boxVectors=None, ghostFile='gcmc-ghost-wats.txt', log='gcmc.log', dcd=None, rst=None, overwrite=False)

Bases: grand.samplers.BaseGrandCanonicalMonteCarloSampler

Base class for carrying out GCMC moves in OpenMM, sampling the whole system with GCMC

__init__(system, topology, temperature, adams=None, excessChemicalPotential=Quantity(value=- 6.09, unit=kilocalorie / mole), standardVolume=Quantity(value=30.345, unit=angstrom ** 3), adamsShift=0.0, boxVectors=None, ghostFile='gcmc-ghost-wats.txt', log='gcmc.log', dcd=None, rst=None, overwrite=False)

Initialise the object to be used for sampling water insertion/deletion moves

Parameters
  • system (simtk.openmm.System) – System object to be used for the simulation

  • topology (simtk.openmm.app.Topology) – Topology object for the system to be simulated

  • temperature (simtk.unit.Quantity) – Temperature of the simulation, must be in appropriate units

  • adams (float) – Adams B value for the simulation (dimensionless). Default is None, if None, the B value is calculated from the box volume and chemical potential

  • excessChemicalPotential (simtk.unit.Quantity) – Excess chemical potential of the system that the simulation should be in equilibrium with, default is -6.09 kcal/mol. This should be the hydration free energy of water, and may need to be changed for specific simulation parameters.

  • standardVolume (simtk.unit.Quantity) – Standard volume of water - corresponds to the volume per water molecule in bulk. The default value is 30.345 A^3

  • adamsShift (float) – Shift the B value from Bequil, if B isn’t explicitly set. Default is 0.0

  • boxVectors (simtk.unit.Quantity) – Box vectors for the simulation cell

  • ghostFile (str) – Name of a file to write out the residue IDs of ghost water molecules. This is useful if you want to visualise the sampling, as you can then remove these waters from view, as they are non-interacting. Default is ‘gcmc-ghost-wats.txt’

  • log (str) – Log file to write out

  • dcd (str) – Name of the DCD file to write the system out to

  • rst (str) – Name of the restart file to write out (.pdb or .rst7)

  • overwrite (bool) – Overwrite any data already present

deleteRandomWater()

Choose a random water to be deleted

Returns

  • delete_water (int) – Resid of the water to delete

  • atom_indices (list) – List of the atom IDs for this molecule

initialise(context, ghostResids)

Prepare the GCMC sphere for simulation by loading the coordinates from a Context object.

Parameters
  • context (simtk.openmm.Context) – Current context of the simulation

  • ghostResids (list) – List of residue IDs corresponding to the ghost waters added

insertRandomWater()

Translate a random ghost to a random point in the simulation box to allow subsequent insertion

Returns

  • new_positions (simtk.unit.Quantity) – Positions following the ‘insertion’ of the ghost water

  • insert_water (int) – Residue ID of the water to insert

  • atom_indices (list) – List of the atom IDs for this molecule

class grand.samplers.NonequilibriumGCMCSphereSampler(system, topology, temperature, integrator, adams=None, excessChemicalPotential=Quantity(value=- 6.09, unit=kilocalorie / mole), standardVolume=Quantity(value=30.345, unit=angstrom ** 3), adamsShift=0.0, nPertSteps=1, nPropStepsPerPert=1, timeStep=Quantity(value=2, unit=femtosecond), lambdas=None, ghostFile='gcmc-ghost-wats.txt', referenceAtoms=None, sphereRadius=None, sphereCentre=None, log='gcmc.log', dcd=None, rst=None, overwrite=False)

Bases: grand.samplers.GCMCSphereSampler

Class to carry out GCMC moves in OpenMM, using nonequilibrium candidate Monte Carlo (NCMC) to boost acceptance rates

__init__(system, topology, temperature, integrator, adams=None, excessChemicalPotential=Quantity(value=- 6.09, unit=kilocalorie / mole), standardVolume=Quantity(value=30.345, unit=angstrom ** 3), adamsShift=0.0, nPertSteps=1, nPropStepsPerPert=1, timeStep=Quantity(value=2, unit=femtosecond), lambdas=None, ghostFile='gcmc-ghost-wats.txt', referenceAtoms=None, sphereRadius=None, sphereCentre=None, log='gcmc.log', dcd=None, rst=None, overwrite=False)

Initialise the object to be used for sampling NCMC-enhanced water insertion/deletion moves

Parameters
  • system (simtk.openmm.System) – System object to be used for the simulation

  • topology (simtk.openmm.app.Topology) – Topology object for the system to be simulated

  • temperature (simtk.unit.Quantity) – Temperature of the simulation, must be in appropriate units

  • integrator (simtk.openmm.CustomIntegrator) – Integrator to use to propagate the dynamics of the system. Currently want to make sure that this is the customised Langevin integrator found in openmmtools which uses BAOAB (VRORV) splitting.

  • adams (float) – Adams B value for the simulation (dimensionless). Default is None, if None, the B value is calculated from the box volume and chemical potential

  • excessChemicalPotential (simtk.unit.Quantity) – Excess chemical potential of the system that the simulation should be in equilibrium with, default is -6.09 kcal/mol. This should be the hydration free energy of water, and may need to be changed for specific simulation parameters.

  • standardVolume (simtk.unit.Quantity) – Standard volume of water - corresponds to the volume per water molecule in bulk. The default value is 30.345 A^3

  • adamsShift (float) – Shift the B value from Bequil, if B isn’t explicitly set. Default is 0.0

  • nPertSteps (int) – Number of pertubation steps over which to shift lambda between 0 and 1 (or vice versa).

  • nPropStepsPerPert (int) – Number of propagation steps to carry out for

  • timeStep (simtk.unit.Quantity) – Time step to use for non-equilibrium integration during the propagation steps

  • lambdas (list) – Series of lambda values corresponding to the pathway over which the molecules are perturbed

  • ghostFile (str) – Name of a file to write out the residue IDs of ghost water molecules. This is useful if you want to visualise the sampling, as you can then remove these waters from view, as they are non-interacting. Default is ‘gcmc-ghost-wats.txt’

  • referenceAtoms (list) – List containing dictionaries describing the atoms to use as the centre of the GCMC region Must contain ‘name’ and ‘resname’ as keys, and optionally ‘resid’ (recommended) and ‘chain’ e.g. [{‘name’: ‘C1’, ‘resname’: ‘LIG’, ‘resid’: ‘123’}]

  • sphereRadius (simtk.unit.Quantity) – Radius of the spherical GCMC region

  • sphereCentre (simtk.unit.Quantity) – Coordinates around which the GCMC sphere is based

  • log (str) – Name of the log file to write out

  • dcd (str) – Name of the DCD file to write the system out to

  • rst (str) – Name of the restart file to write out (.pdb or .rst7)

  • overwrite (bool) – Indicates whether to overwrite already existing data

deletionMove()

Carry out a nonequilibrium deletion move for a random water molecule

insertionMove()

Carry out a nonequilibrium insertion move for a random water molecule

move(context, n=1)

Carry out a nonequilibrium GCMC move

Parameters
  • context (simtk.openmm.Context) – Current context of the simulation

  • n (int) – Number of moves to execute

reset()

Reset counted values (such as number of total or accepted moves) to zero

class grand.samplers.NonequilibriumGCMCSystemSampler(system, topology, temperature, integrator, adams=None, excessChemicalPotential=Quantity(value=- 6.09, unit=kilocalorie / mole), standardVolume=Quantity(value=30.345, unit=angstrom ** 3), adamsShift=0.0, nPertSteps=1, nPropStepsPerPert=1, timeStep=Quantity(value=2, unit=femtosecond), boxVectors=None, ghostFile='gcmc-ghost-wats.txt', log='gcmc.log', dcd=None, rst=None, overwrite=False, lambdas=None)

Bases: grand.samplers.GCMCSystemSampler

Class to carry out GCMC moves in OpenMM, using nonequilibrium candidate Monte Carlo (NCMC) to boost acceptance rates

__init__(system, topology, temperature, integrator, adams=None, excessChemicalPotential=Quantity(value=- 6.09, unit=kilocalorie / mole), standardVolume=Quantity(value=30.345, unit=angstrom ** 3), adamsShift=0.0, nPertSteps=1, nPropStepsPerPert=1, timeStep=Quantity(value=2, unit=femtosecond), boxVectors=None, ghostFile='gcmc-ghost-wats.txt', log='gcmc.log', dcd=None, rst=None, overwrite=False, lambdas=None)

Initialise the object to be used for sampling NCMC-enhanced water insertion/deletion moves

Parameters
  • system (simtk.openmm.System) – System object to be used for the simulation

  • topology (simtk.openmm.app.Topology) – Topology object for the system to be simulated

  • temperature (simtk.unit.Quantity) – Temperature of the simulation, must be in appropriate units

  • integrator (simtk.openmm.CustomIntegrator) – Integrator to use to propagate the dynamics of the system. Currently want to make sure that this is the customised Langevin integrator found in openmmtools which uses BAOAB (VRORV) splitting.

  • adams (float) – Adams B value for the simulation (dimensionless). Default is None, if None, the B value is calculated from the box volume and chemical potential

  • excessChemicalPotential (simtk.unit.Quantity) – Excess chemical potential of the system that the simulation should be in equilibrium with, default is -6.09 kcal/mol. This should be the hydration free energy of water, and may need to be changed for specific simulation parameters.

  • standardVolume (simtk.unit.Quantity) – Standard volume of water - corresponds to the volume per water molecule in bulk. The default value is 30.345 A^3

  • adamsShift (float) – Shift the B value from Bequil, if B isn’t explicitly set. Default is 0.0

  • nPertSteps (int) – Number of pertubation steps over which to shift lambda between 0 and 1 (or vice versa).

  • nPropStepsPerPert (int) – Number of propagation steps to carry out for

  • timeStep (simtk.unit.Quantity) – Time step to use for non-equilibrium integration during the propagation steps

  • lambdas (list) – Series of lambda values corresponding to the pathway over which the molecules are perturbed

  • boxVectors (simtk.unit.Quantity) – Box vectors for the simulation cell

  • ghostFile (str) – Name of a file to write out the residue IDs of ghost water molecules. This is useful if you want to visualise the sampling, as you can then remove these waters from view, as they are non-interacting. Default is ‘gcmc-ghost-wats.txt’

  • log (str) – Name of the log file to write out

  • dcd (str) – Name of the DCD file to write the system out to

  • rst (str) – Name of the restart file to write out (.pdb or .rst7)

  • overwrite (bool) – Indicates whether to overwrite already existing data

deletionMove()

Carry out a nonequilibrium deletion move for a random water molecule

insertionMove()

Carry out a nonequilibrium insertion move for a random water molecule

move(context, n=1)

Carry out a nonequilibrium GCMC move

Parameters
  • context (simtk.openmm.Context) – Current context of the simulation

  • n (int) – Number of moves to execute

reset()

Reset counted values (such as number of total or accepted moves) to zero

class grand.samplers.StandardGCMCSphereSampler(system, topology, temperature, adams=None, excessChemicalPotential=Quantity(value=- 6.09, unit=kilocalorie / mole), standardVolume=Quantity(value=30.345, unit=angstrom ** 3), adamsShift=0.0, ghostFile='gcmc-ghost-wats.txt', referenceAtoms=None, sphereRadius=None, sphereCentre=None, log='gcmc.log', dcd=None, rst=None, overwrite=False)

Bases: grand.samplers.GCMCSphereSampler

Class to carry out instantaneous GCMC moves in OpenMM

__init__(system, topology, temperature, adams=None, excessChemicalPotential=Quantity(value=- 6.09, unit=kilocalorie / mole), standardVolume=Quantity(value=30.345, unit=angstrom ** 3), adamsShift=0.0, ghostFile='gcmc-ghost-wats.txt', referenceAtoms=None, sphereRadius=None, sphereCentre=None, log='gcmc.log', dcd=None, rst=None, overwrite=False)

Initialise the object to be used for sampling instantaneous water insertion/deletion moves

Parameters
  • system (simtk.openmm.System) – System object to be used for the simulation

  • topology (simtk.openmm.app.Topology) – Topology object for the system to be simulated

  • temperature (simtk.unit.Quantity) – Temperature of the simulation, must be in appropriate units

  • adams (float) – Adams B value for the simulation (dimensionless). Default is None, if None, the B value is calculated from the box volume and chemical potential

  • excessChemicalPotential (simtk.unit.Quantity) – Excess chemical potential of the system that the simulation should be in equilibrium with, default is -6.09 kcal/mol. This should be the hydration free energy of water, and may need to be changed for specific simulation parameters.

  • standardVolume (simtk.unit.Quantity) – Standard volume of water - corresponds to the volume per water molecule in bulk. The default value is 30.345 A^3

  • adamsShift (float) – Shift the B value from Bequil, if B isn’t explicitly set. Default is 0.0

  • ghostFile (str) – Name of a file to write out the residue IDs of ghost water molecules. This is useful if you want to visualise the sampling, as you can then remove these waters from view, as they are non-interacting. Default is ‘gcmc-ghost-wats.txt’

  • referenceAtoms (list) – List containing dictionaries describing the atoms to use as the centre of the GCMC region Must contain ‘name’ and ‘resname’ as keys, and optionally ‘resid’ (recommended) and ‘chain’ e.g. [{‘name’: ‘C1’, ‘resname’: ‘LIG’, ‘resid’: ‘123’}]

  • sphereRadius (simtk.unit.Quantity) – Radius of the spherical GCMC region

  • sphereCentre (simtk.unit.Quantity) – Coordinates around which the GCMC sphere is based

  • log (str) – Name of the log file to write out

  • dcd (str) – Name of the DCD file to write the system out to

  • rst (str) – Name of the restart file to write out (.pdb or .rst7)

  • overwrite (bool) – Indicates whether to overwrite already existing data

deletionMove()

Carry out a random water deletion move on the current system

insertionMove()

Carry out a random water insertion move on the current system

move(context, n=1)

Execute a number of GCMC moves on the current system

Parameters
  • context (simtk.openmm.Context) – Current context of the simulation

  • n (int) – Number of moves to execute

class grand.samplers.StandardGCMCSystemSampler(system, topology, temperature, adams=None, excessChemicalPotential=Quantity(value=- 6.09, unit=kilocalorie / mole), standardVolume=Quantity(value=30.345, unit=angstrom ** 3), adamsShift=0.0, boxVectors=None, ghostFile='gcmc-ghost-wats.txt', log='gcmc.log', dcd=None, rst=None, overwrite=False)

Bases: grand.samplers.GCMCSystemSampler

Class to carry out instantaneous GCMC moves in OpenMM

__init__(system, topology, temperature, adams=None, excessChemicalPotential=Quantity(value=- 6.09, unit=kilocalorie / mole), standardVolume=Quantity(value=30.345, unit=angstrom ** 3), adamsShift=0.0, boxVectors=None, ghostFile='gcmc-ghost-wats.txt', log='gcmc.log', dcd=None, rst=None, overwrite=False)

Initialise the object to be used for sampling instantaneous water insertion/deletion moves

Parameters
  • system (simtk.openmm.System) – System object to be used for the simulation

  • topology (simtk.openmm.app.Topology) – Topology object for the system to be simulated

  • temperature (simtk.unit.Quantity) – Temperature of the simulation, must be in appropriate units

  • adams (float) – Adams B value for the simulation (dimensionless). Default is None, if None, the B value is calculated from the box volume and chemical potential

  • excessChemicalPotential (simtk.unit.Quantity) – Excess chemical potential of the system that the simulation should be in equilibrium with, default is -6.09 kcal/mol. This should be the hydration free energy of water, and may need to be changed for specific simulation parameters.

  • standardVolume (simtk.unit.Quantity) – Standard volume of water - corresponds to the volume per water molecule in bulk. The default value is 30.345 A^3

  • adamsShift (float) – Shift the B value from Bequil, if B isn’t explicitly set. Default is 0.0

  • boxVectors (simtk.unit.Quantity) – Box vectors for the simulation cell

  • ghostFile (str) – Name of a file to write out the residue IDs of ghost water molecules. This is useful if you want to visualise the sampling, as you can then remove these waters from view, as they are non-interacting. Default is ‘gcmc-ghost-wats.txt’

  • log (str) – Name of the log file to write out

  • dcd (str) – Name of the DCD file to write the system out to

  • rst (str) – Name of the restart file to write out (.pdb or .rst7)

  • overwrite (bool) – Indicates whether to overwrite already existing data

deletionMove()

Carry out a random water deletion move on the current system

insertionMove()

Carry out a random water insertion move on the current system

move(context, n=1)

Execute a number of GCMC moves on the current system

Parameters
  • context (simtk.openmm.Context) – Current context of the simulation

  • n (int) – Number of moves to execute