ShortestDistance

class cvpack.ShortestDistance(group1, group2, numAtoms, sigma=0.01 nm, magnitude=0.2 nm, cutoffDistance=0.5 nm, switchDistance=0.4 nm, pbc=True, name='shortest_distance')[source]

A smooth approximation of the shortest distance between two atom groups:

\[r_{\rm min}({\bf r}) = \frac{ \sum_{i \in {\bf g}_1} \sum_{j \in {\bf g}_2} r_{ij} e^{-\frac{r_{ij}^2}{2 \sigma^2}} }{ \sum_{i \in {\bf g}_1} \sum_{j \in {\bf g}_2} e^{-\frac{r_{ij}^2}{2 \sigma^2}} }\]

where \(r_{ij} = \|{\bf r}_j - {\bf r}_i\|\) is the distance between atoms \(i\) and \(j\) and \(\sigma\) is a parameter that controls the degree of approximation. The smaller the value of \(\sigma\), the closer the approximation to the true shortest distance.

In practice, a cutoff distance \(r_c\) is applied to the atom pairs so that the collective variable is computed only for pairs of atoms separated by a distance less than \(r_c\). A switching function is also applied to smoothly turn off the collective variable starting from a distance \(r_s < r_c\).

Note

Atoms are allowed to be in both groups. In this case, terms for which \(i = j\) are ignored.

Parameters:
  • group1 (Sequence[int]) – The indices of the atoms in the first group.

  • group2 (Sequence[int]) – The indices of the atoms in the second group.

  • numAtoms (int) – The total number of atoms in the system, including those that are not in any of the groups.

  • sigma (Quantity) – The parameter that controls the degree of approximation.

  • magnitude (Quantity) – The expected order of magnitude of the shortest distance. This parameter does not affect the value of the collective variable, but a good estimate can improve numerical stability.

  • cutoffDistance (Quantity) – The cutoff distance for evaluating atom pairs.

  • switchDistance (Quantity) – The distance at which the switching function starts to be applied.

  • pbc (bool) – Whether to consider periodic boundary conditions in distance calculations.

  • name (str) – The name of the collective variable.

Example

>>> import cvpack
>>> import numpy as np
>>> from openmm import unit
>>> from openmmtools import testsystems
>>> model = testsystems.HostGuestVacuum()
>>> group1, group2 = [], []
>>> for residue in model.topology.residues():
...     group = group1 if residue.name == "B2" else group2
...     group.extend(atom.index for atom in residue.atoms())
>>> coords = model.positions.value_in_unit(mmunit.nanometers)
>>> np.linalg.norm(
...     coords[group1, None, :] - coords[None, group2, :],
...     axis=-1,
... ).min()
0.17573...
>>> num_atoms = model.system.getNumParticles()
>>> min_dist = cvpack.ShortestDistance(group1, group2, num_atoms)
>>> min_dist.addToSystem(model.system)
>>> platform = openmm.Platform.getPlatformByName("Reference")
>>> integrator = openmm.VerletIntegrator(1.0 * mmunit.femtoseconds)
>>> context = openmm.Context(model.system, integrator, platform)
>>> context.setPositions(model.positions)
>>> min_dist.getValue(context)
0.17573... nm

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
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

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