MetaCollectiveVariable

class cvpack.MetaCollectiveVariable(function, variables, unit, periodicBounds=None, name='meta_collective_variable', **parameters)[source]

A collective variable that is a function of other collective variables:

\[f({\bf r}) = F\left(c_1({\bf r}), c_2({\bf r}), \ldots, c_n({\bf r})\right)\]

where \(F(c_1,c_2,\ldots,c_n)\) is a user-defined function and \(c_i({\bf r})\) is the value of the \(i\)-th collective variable at a given configuration \({\bf r}\).

The function \(F\) is defined as a string and can be any expression supported by the openmm.CustomCVForce class. If the expression contains named parameters, the value of each parameter can be passed in one of two ways:

  1. By a semicolon-separated definition in the function string, such as described in the openmm.CustomCompoundBondForce documentation. In this case, the parameter value will be the same for all groups of atoms.

  2. By a scalar passed as a keyword argument to the AtomicFunction constructor. In this case, the parameter will apply to all atom groups and will become available for on-the-fly modification during a simulation via the setParameter method of an openmm.Context object. Warning: other collective variables or openmm.Force objects in the same system will share the same values of equal-named parameters.

Parameters:
  • function (str) – The function to be evaluated. It must be a valid openmm.CustomCVForce expression.

  • variables (t.Iterable[str]) – A sequence of CollectiveVariable instances that represent the collective variables on which the meta-collective variable depends. The name of each collective variable must be unique and match a symbol used in the function.

  • unit (mmunit.Unit) – The unit of measurement of the collective variable. It must be compatible with the MD unit system (mass in daltons, distance in nanometers, time in picoseconds, temperature in kelvin, energy in kilojoules_per_mol, angle in radians). If the collective variables does not have a unit, use unit.dimensionless

  • periodicBounds (t.Optional[VectorQuantity]) – The periodic bounds of the collective variable if it is periodic, or None if it is not periodic.

Keyword Arguments:

**parameters – The named parameters of the function, if any. They will become settable context parameters when this meta-collective variable is added to an openmm.System. The passed objects must be scalar quantities. Their values will be converted to OpenMM’s MD unit system to serve as default values for the context parameters.

Example

>>> import cvpack
>>> import math
>>> import openmm
>>> from openmm import unit
>>> from openmmtools import testsystems
>>> model = testsystems.AlanineDipeptideVacuum()
>>> phi = cvpack.Torsion(6, 8, 14, 16, name="phi")
>>> driving_force = cvpack.MetaCollectiveVariable(
...     f"0.5*kappa*min(delta,{2*math.pi}-delta)^2"
...     "; delta=abs(phi-phi0)",
...     [phi],
...     unit.kilojoules_per_mole,
...     kappa = 1e3 * unit.kilojoules_per_mole/unit.radian**2,
...     phi0 = 120 * unit.degrees
... )
>>> driving_force.addToSystem(model.system)
>>> integrator = openmm.VerletIntegrator(0)
>>> platform = openmm.Platform.getPlatformByName('Reference')
>>> context = openmm.Context(model.system, integrator, platform)
>>> context.setPositions(model.positions)
>>> driving_force.getParameterDefaultValues()
{'kappa': 1000.0 kJ/(mol rad**2), 'phi0': 2.094... rad}
>>> driving_force.getParameterValues(context)
{'kappa': 1000.0 kJ/(mol rad**2), 'phi0': 2.094... rad}
>>> driving_force.getValue(context)
548.3... kJ/mol
>>> driving_force.getInnerValues(context)
{'phi': 3.14... rad}
>>> driving_force.getInnerEffectiveMasses(context)
{'phi': 0.05119... nm**2 Da/(rad**2)}
>>> driving_force.getParameterDerivatives(context)
{'kappa': 0.548... rad**2, 'phi0': -1047.19... kJ/(mol rad)}

Methods

addToSystem(system, setUnusedForceGroup=True)

Add this collective variable to an openmm.System.

Parameters:
  • system (System) – The system to which this collective variable should be added

  • setUnusedForceGroup (bool) – If True, the force group of this collective variable will be set to the first available force group in the system

getEffectiveMass(context, allowReinitialization=False)

Compute the effective mass of this collective variable at a given openmm.Context.

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

\[m_\mathrm{eff}({\bf r}) = \left( \sum_{i=1}^N \frac{1}{m_i} \left\| \frac{dq}{d{\bf r}_i} \right\|^2 \right)^{-1}\]

If this collective variable share the force group with other forces, then evaluating its effective mass requires reinitializing the openmm.Context twice at each call. This is inefficient and should be avoided. To allow this behavior, the allowReinitialization parameter must be set to True.

Note

By adding this collective variable to the system using the addToSystem() method, the force group of this collective variable is set to an available force group in the system by default.

Parameters:
  • context (Context) – The context at which this collective variable’s effective mass should be evaluated.

  • allowReinitialization (bool) – If True, the context will be reinitialized if necessary.

Returns:

The effective mass of this collective variable at the given context.

Return type:

Quantity

Example

In this example, we compute the effective masses of the backbone dihedral angles and the radius of gyration of an alanine dipeptide molecule in water:

>>> import cvpack
>>> import openmm
>>> from openmmtools import testsystems
>>> model = testsystems.AlanineDipeptideExplicit()
>>> top = model.mdtraj_topology
>>> backbone_atoms = top.select("name N C CA and resid 1 2")
>>> phi = cvpack.Torsion(*backbone_atoms[0:4])
>>> psi = cvpack.Torsion(*backbone_atoms[1:5])
>>> radius_of_gyration = cvpack.RadiusOfGyration(
...     top.select("not water")
... )
>>> for cv in [phi, psi, radius_of_gyration]:
...     cv.addToSystem(model.system)
>>> context = openmm.Context(
...     model.system, openmm.VerletIntegrator(0)
... )
>>> context.setPositions(model.positions)
>>> phi.getEffectiveMass(context)
0.05119... nm**2 Da/(rad**2)
>>> psi.getEffectiveMass(context)
0.05186... nm**2 Da/(rad**2)
>>> radius_of_gyration.getEffectiveMass(context)
30.946... Da
getInnerEffectiveMasses(context)[source]

Get the effective masses of the collective variables on which the meta-collective variable depends. The effective masses are calculated from the forces acting on the particles that represent the collective variables.

Parameters:

context (Context) – The context in which the collective variables will be evaluated.

Returns:

A dictionary with the names of the collective variables as keys and their effective masses as values.

Return type:

Dict[str, Quantity]

getInnerValues(context)[source]

Get the values of the collective variables on which the meta-collective variable depends. The values are returned as a dictionary with the names of the collective variables as keys.

Parameters:

context (Context) – The context in which the collective variables will be evaluated.

Returns:

A dictionary with the names of the collective variables as keys and their values as values.

Return type:

Dict[str, Quantity]

getInnerVariables()[source]

Get the collective variables on which the meta-collective variable depends.

Returns:

A tuple with the collective variables.

Return type:

Tuple[CollectiveVariable]

getMassUnit()

Get the unit of measurement of the effective mass of this collective variable.

Return type:

Unit

getName()

Get the name of this collective variable.

Return type:

str

getParameterDefaultValues()[source]

Get the default values of the named parameters of this meta-collective variable.

Returns:

A dictionary with the names of the named parameters as keys and their default values as values.

Return type:

Dict[str, Quantity]

getParameterDerivatives(context, allowReinitialization=False)[source]

Get the derivatives of the named parameters of this meta-collective variable.

Returns:

A dictionary with the names of the named parameters as keys and their derivatives as values.

Return type:

Dict[str, Quantity]

getParameterValues(context)[source]

Get the values of the named parameters of this meta-collective variable. The values are returned as a dictionary with the names of the parameters as keys.

Parameters:

context (Context) – The context in which the named parameters will be evaluated.

Returns:

A dictionary with the names of the named parameters as keys and their values as values.

Return type:

Dict[str, Quantity]

getPeriodicBounds()

Get the periodic bounds of this collective variable.

Returns:

The periodic bounds of this collective variable or None if it is not periodic

Return type:

Union[Quantity, None]

getUnit()

Get the unit of measurement of this collective variable.

Return type:

Unit

getValue(context, allowReinitialization=False)

Evaluate this collective variable at a given openmm.Context.

If this collective variable share the force group with other forces, then its evaluation requires reinitializing the openmm.Context twice at each call. This is inefficient and should be avoided. To allow this behavior, the allowReinitialization parameter must be set to True.

Note

By adding this collective variable to the system using the addToSystem() method, the force group of this collective variable is set to an available force group in the system by default.

Parameters:
  • context (Context) – The context at which this collective variable should be evaluated.

  • allowReinitialization (bool) – If True, the context will be reinitialized if necessary.

Returns:

The value of this collective variable at the given context.

Return type:

Quantity

Example

In this example, we compute the values of the backbone dihedral angles and the radius of gyration of an alanine dipeptide molecule in water:

>>> import cvpack
>>> import openmm
>>> from openmmtools import testsystems
>>> model = testsystems.AlanineDipeptideExplicit()
>>> top = model.mdtraj_topology
>>> backbone_atoms = top.select("name N C CA and resid 1 2")
>>> phi = cvpack.Torsion(*backbone_atoms[0:4])
>>> psi = cvpack.Torsion(*backbone_atoms[1:5])
>>> radius_of_gyration = cvpack.RadiusOfGyration(
...     top.select("not water")
... )
>>> for cv in [phi, psi, radius_of_gyration]:
...     cv.addToSystem(model.system)
>>> context = openmm.Context(
...     model.system, openmm.VerletIntegrator(0)
... )
>>> context.setPositions(model.positions)
>>> phi.getValue(context)
3.1415... rad
>>> psi.getValue(context)
3.1415... rad
>>> radius_of_gyration.getValue(context)
0.29514... nm