ufedmm

class ufedmm.ufedmm.CollectiveVariable(id, force)[source]

Bases: object

A function of the particle coordinates, evaluated by means of an openmm.Force object.

Quoting the manual entry for openmm.CustomCVForce:

“Each collective variable is defined by a Force object. The Force’s potential energy is computed, and that becomes the value of the variable. This provides enormous flexibility in defining collective variables, especially by using custom forces. Anything that can be computed as a potential function can also be used as a collective variable.”

Parameters:
  • id (str) – A valid identifier string for this collective variable.

  • force (openmm.Force) – An openmm.Force object whose energy function is used to evaluate this collective variable.

Example

>>> import openmm
>>> import ufedmm
>>> from openmm import unit
>>> cv = ufedmm.CollectiveVariable('psi', openmm.CustomTorsionForce('theta'))
>>> cv.force.addTorsion(0, 1, 2, 3, [])
0
effective_mass(system, positions, cv_unit=None)[source]

Computes the effective mass of the collective variable for a given system and a given set of particle coordinates.

The effective mass of a collective variable \(q(\mathbf{r})\) is defined as [1]:

\[m_\mathrm{eff} = \left( \sum_{j=1}^N \frac{1}{m_j} \left\| \frac{dq}{d\mathbf{r}_j} \right\|^2 \right)^{-1}\]
Parameters:
  • system (openmm.System) – The system for which the collective variable will be evaluated.

  • positions (list of openmm.Vec3) – A list whose size equals the number of particles in the system and which contains the coordinates of these particles.

Keyword Arguments:

cv_unit (unit.Unit, default=None) – The unity of measurement of the collective variable. If this is None, then a numerical value is returned based on the OpenMM default units.

Return type:

float or unit.Quantity

Example

>>> import ufedmm
>>> from openmm import unit
>>> model = ufedmm.AlanineDipeptideModel()
>>> model.phi.effective_mass(model.system, model.positions)
0.04795887...
>>> model.psi.effective_mass(model.system, model.positions)
0.05115582...
evaluate(system, positions, cv_unit=None)[source]

Computes the value of the collective variable for a given system and a given set of particle coordinates.

Parameters:
  • system (openmm.System) – The system for which the collective variable will be evaluated.

  • positions (list of openmm.Vec3) – A list whose size equals the number of particles in the system and which contains the coordinates of these particles.

Keyword Arguments:

cv_unit (unit.Unit, default=None) – The unity of measurement of the collective variable. If this is None, then a numerical value is returned based on the OpenMM default units.

Return type:

float or unit.Quantity

Example

>>> import ufedmm
>>> from openmm import unit
>>> model = ufedmm.AlanineDipeptideModel()
>>> model.phi.evaluate(model.system, model.positions)
3.141592653589793
>>> model.psi.evaluate(model.system, model.positions)
3.141592653589793
class ufedmm.ufedmm.DynamicalVariable(id, min_value, max_value, mass, temperature, colvars, potential, periodic=True, sigma=None, grid_size=None, **parameters)[source]

Bases: object

An extended phase-space variable, whose dynamics is coupled to that of one of more collective variables of a system.

The coupling occurs in the form of a potential energy term involving this dynamical variable and its associated collective variables.

The default potential is a harmonic driving of the type:

\[V(s, \mathbf r) = \frac{\kappa}{2} [s - q(\mathbf r)]^2\]

where \(s\) is the new dynamical variable, \(q(\mathbf r)\) is its associated collective variable, and \(kappa\) is a force constant.

Parameters:
  • id (str) – A valid identifier string for this dynamical variable.

  • min_value (float or unit.Quantity) – The minimum allowable value for this dynamical variable.

  • max_value (float or unit.Quantity) – The maximum allowable value for this dynamical variable.

  • mass (float or unit.Quantity) – The mass assigned to this dynamical variable, whose unit of measurement must be compatible with unit.dalton*(unit.nanometers/X)**2, where X is the unit of measurement of the dynamical variable itself.

  • temperature (float or unit.Quantity) – The temperature of the heat bath attached to this variable.

  • colvars (CollectiveVariable or list thereof) – Either a single colective variable or a list.

  • potential (float or unit.Quantity or str) – Either the value of the force constant of a harmonic driving potential or an algebraic expression giving the energy of the system as a function of this dynamical variable and its associated collective variable. Such expression can also contain a set of global parameters, whose values must be passed as keyword arguments (see below).

Keyword Arguments:
  • periodic (bool, default=True) – Whether the collective variable is periodic with period L=max_value-min_value.

  • sigma (float or unit.Quantity, default=None) – The standard deviation. If this is None, then no bias will be considered.

  • grid_size (int, default=None) – The grid size. If this is None and sigma is finite, then a convenient value will be automatically chosen.

  • **parameters – Names and values of global parameters present in the algebraic expression defined as potential (see above).

Example

>>> import openmm
>>> import ufedmm
>>> from openmm import unit
>>> cv = ufedmm.CollectiveVariable('psi', openmm.CustomTorsionForce('theta'))
>>> cv.force.addTorsion(0, 1, 2, 3, [])
0
>>> mass = 50*unit.dalton*(unit.nanometer/unit.radians)**2
>>> K = 1000*unit.kilojoules_per_mole/unit.radians**2
>>> Ts = 1500*unit.kelvin
>>> ufedmm.DynamicalVariable(
...     's_psi', -180*unit.degrees, 180*unit.degrees, mass, Ts, cv, K
... )
<s_psi in [-3.141592653589793, 3.141592653589793], periodic, m=50, T=1500>
evaluate(x, Lx)[source]

Computes the value of this dynamical variable for a given x coordinate and a given length of the simulation box in the x direction.

Parameters:
  • x (float or unit.Quantity) – The x coordinate.

  • Lx (float or unit.Quantity) – The length of the simulation box in the x direction.

Return type:

float

class ufedmm.ufedmm.ExtendedSpaceContext(variables, system, *args, **kwargs)[source]

Bases: Context

An extension of the openmm.Context class.

Parameters:
  • variables (list(DynamicalVariable)) – The dynamical variables to be added to the system’s phase space.

  • system (openmm.System) – The openmm.System which will be simulated

  • integrator (openmm.Integrator) – The openmm.Integrator which will be used to simulate the openmm.System.

  • platform (openmm.Platform) – The openmm.Platform to use for calculations.

  • properties (dict(str: str)) – A set of values for platform-specific properties. Keys are the property names.

Variables:
  • variables (list(DynamicalVariable)) – The dynamical variables added to the system’s phase space.

  • driving_forces (list[openmm.CustomCVForce]) – A list of openmm.CustomCVForce objects responsible for evaluating the potential energy terms which couples the extra dynamical variables to their associated collective variables.

getState(**kwargs)[source]

Returns a ExtendedSpaceState object.

Keyword Arguments:

**kwargs – See getState.

setPeriodicBoxVectors(a, b, c)[source]

Set the vectors defining the axes of the periodic box.

Warning

Only orthorhombic boxes are allowed.

Parameters:
  • a (openmm.Vec3) – The vector defining the first edge of the periodic box.

  • b (openmm.Vec3) – The vector defining the second edge of the periodic box.

  • c (openmm.Vec3) – The vector defining the third edge of the periodic box.

setPositions(positions, extended_positions=None)[source]

Sets the positions of all particles and extended-space variables in this context. If the latter are not provided, then suitable values are automatically determined from the particle positions.

Parameters:

positions (list of openmm.Vec3) – The positions of physical particles.

Keyword Arguments:

extended_positions (list of float or unit.Quantity) – The positions of extended-space particles.

setVelocitiesToTemperature(temperature, randomSeed=None)[source]

Sets the velocities of all particles in the system to random values chosen from a Maxwell-Boltzmann distribution at a given temperature.

Warning

The velocities of the extended-space variables are set to values consistent with the distribution at its own specified temperature.

Parameters:

temperature (float or unit.Quantity) – The temperature of the system.

Keyword Arguments:

randomSeed (int, default=None) – A seed for the random number generator.

class ufedmm.ufedmm.ExtendedSpaceSimulation(variables, topology, system, integrator, platform=None, platformProperties=None)[source]

Bases: Simulation

A simulation involving extended phase-space variables.

Parameters:
  • variables (list of DynamicalVariable) – The dynamical variables to be added to the system’s phase space.

  • topology (openmm.app.Topology) – A Topology describing the system to be simulated.

  • system (openmm.System) – The openmm.System object to be simulated.

  • integrator (openmm.Integrator) – The OpenMM Integrator to use for simulating the system dynamics.

Keyword Arguments:
  • platform (openmm.Platform, default=None) – If not None, the openmm.Platform to use.

  • platformProperties (dict, default=None) – If not None, a set of platform-specific properties to pass to the Context’s constructor.

add_periodic_task(task, force_group=0)[source]

Adds a task to be executed periodically along this simulation.

Parameters:

task (PeriodicTask) – A PeriodicTask object.

loadCheckpoint(file)[source]

Load a checkpoint file that was created with saveCheckpoint().

Parameters:

file (string or file) – a File-like object to load the checkpoint from, or alternatively a filename

saveCheckpoint(file)[source]

Save a checkpoint of the simulation to a file.

The output is a binary file that contains a complete representation of the current state of the Simulation. It includes both publicly visible data such as the particle positions and velocities, and also internal data such as the states of random number generators. Reloading the checkpoint will put the Simulation back into precisely the same state it had before, so it can be exactly continued.

A checkpoint file is highly specific to the Simulation it was created from. It can only be loaded into another Simulation that has an identical System, uses the same Platform and OpenMM version, and is running on identical hardware. If you need a more portable way to resume simulations, consider using saveState() instead.

Parameters:

file (string or file) – a File-like object to write the checkpoint to, or alternatively a filename

step(steps)[source]

Executed a given number of simulation steps.

Parameters:

steps (int) – The number of steps to be executed.

class ufedmm.ufedmm.ExtendedSpaceState(variables, state)[source]

Bases: State

An extension of the openmm.State class.

getDynamicalVariables()[source]

Gets the values of the extended-space dynamical variables.

Return type:

list(float)

Raises:

Exception – If positions were not requested in the context.getState() call.

Example

>>> import ufedmm
>>> import numpy as np
>>> model = ufedmm.AlanineDipeptideModel()
>>> integrator = ufedmm.CustomIntegrator(300, 0.001)
>>> args = [-np.pi, np.pi, 50, 1500, model.phi, 1000]
>>> s_phi = ufedmm.DynamicalVariable('s_phi', *args)
>>> s_psi = ufedmm.DynamicalVariable('s_psi', *args)
>>> context = ufedmm.ExtendedSpaceContext(
...     [s_phi, s_psi], model.system, integrator
... )
>>> context.setPositions(model.positions)
>>> context.getState(getPositions=True).getDynamicalVariables()
[-3.141592653589793, -3.141592653589793]
getPositions(asNumpy=False, extended=False, split=True)[source]

Gets the positions of all physical particles and optionally also gets the positions of the extra particles from which the extended-space dynamical variables are computed.

Keyword Arguments:
  • asNumpy (bool, default=False) – Whether to return Numpy arrays instead of lists of openmm.Vec3.

  • extended (bool, default=False) – Whether to include the positions of the extra particles from which the extended-space dynamical variables are computed.

  • split (bool, default=True) – Whether to return the positions of the physical particles and the positions of the extra particles separately.

Returns:

  • list(openmm.Vec3) or numpy.ndarray – If extended=False or split=False.

  • list(openmm.Vec3) or numpy.ndarray, list(float) or numpy.ndarray – If extended=True and split=True.

Raises:

Exception – If positions were not requested in the context.getState() call.

getVelocities(asNumpy=False, extended=False)[source]

Gets the velocities of all physical particles and optionally also gets the velocities of the extra particles associated to the extended- space dynamical variables.

Keyword Arguments:
  • asNumpy (bool, default=False) – Whether to return Numpy arrays instead of lists of openmm.Vec3.

  • extended (bool, default=False) – Whether to include the velocities of the extra particles.

Returns:

  • list(openmm.Vec3) – If asNumpy=False and extended=False.

  • numpy.ndarray – If asNumpy=True and extended=False.

  • list(openmm.Vec3), list(float) – If asNumpy=False and extended=True.

  • numpy.ndarray, numpy.ndarray – If asNumpy=True and extended=True.

Raises:

Exception – If velocities were not requested in the context.getState() call.

class ufedmm.ufedmm.UnifiedFreeEnergyDynamics(variables, temperature, height=None, frequency=None, bias_factor=None, grid_expansion=20, enforce_gridless=False, buffer_size=100)[source]

Bases: object

A Unified Free-Energy Dynamics (UFED) setup.

Parameters:
  • variables (list of DynamicalVariable) – The variables.

  • temperature (float or unit.Quantity) – The temperature.

Keyword Arguments:
  • height (float or unit.Quantity, default=None) – The height.

  • frequency (int, default=None) – The frequency.

  • bias_factor (float, default=None) – Scales the height of the hills added to the metadynamics bias potential. If it is None, then the hills will have a constant height over time. For a bias factor to be applicable, all bias variables must be at the same temperature T. The extended-space dynamical variables are sampled as if the effective temperature were T*bias_factor.

  • grid_expansion (int, default=20) – The grid expansion.

  • enforce_gridless (bool, default=False) – If this is True, gridless metadynamics is enforced even for 1D to 3D problems.

  • buffer_size (int, default=100) – The buffer size.

Example

>>> import ufedmm
>>> from openmm import unit
>>> model = ufedmm.AlanineDipeptideModel(water='tip3p')
>>> mass = 50*unit.dalton*(unit.nanometer/unit.radians)**2
>>> Ks = 1000*unit.kilojoules_per_mole/unit.radians**2
>>> Ts = 1500*unit.kelvin
>>> limit = 180*unit.degrees
>>> s_phi = ufedmm.DynamicalVariable(
...     's_phi', -limit, limit, mass, Ts, model.phi, Ks
... )
>>> s_psi = ufedmm.DynamicalVariable(
...     's_psi', -limit, limit, mass, Ts, model.psi, Ks
... )
>>> ufedmm.UnifiedFreeEnergyDynamics([s_phi, s_psi], 300*unit.kelvin)
<variables=[s_phi, s_psi], temperature=300, height=None, frequency=None>
simulation(topology, system, integrator, platform=None, platformProperties=None)[source]

Returns a ExtendedSpaceSimulation object.

Warning

If the temperature of any driving parameter is different from the particle-system temperature, then the passed integrator must be a openmm.CustomIntegrator object containing a per-dof variable kT whose content is the Boltzmann constant times the temperature associated to each degree of freedom. This is true for all integrators available in ufedmm.integrators, which are subclasses of ufedmm.integrators.CustomIntegrator.

Parameters:
  • topology (openmm.app.Topology) – The topology.

  • system (openmm.System) – The system.

  • integrator – The integrator.

Keyword Arguments:
  • platform (openmm.Platform, default=None) – The platform.

  • platformProperties (dict, default=None) – The platform properties.

Example

>>> import ufedmm
>>> from openmm import unit
>>> model = ufedmm.AlanineDipeptideModel(water='tip3p')
>>> mass = 50*unit.dalton*(unit.nanometer/unit.radians)**2
>>> Ks = 1000*unit.kilojoules_per_mole/unit.radians**2
>>> Ts = 1500*unit.kelvin
>>> dt = 2*unit.femtoseconds
>>> gamma = 10/unit.picoseconds
>>> limit = 180*unit.degrees
>>> s_phi = ufedmm.DynamicalVariable(
...    's_phi', -limit, limit, mass, Ts, model.phi, Ks
... )
>>> s_psi = ufedmm.DynamicalVariable(
...    's_psi', -limit, limit, mass, Ts, model.psi, Ks
... )
>>> ufed = ufedmm.UnifiedFreeEnergyDynamics([s_phi, s_psi], 300*unit.kelvin)
>>> integrator = ufedmm.GeodesicLangevinIntegrator(
...    300*unit.kelvin, gamma, dt
... )
>>> simulation = ufed.simulation(model.topology, model.system, integrator)