ResidueCoordination

class cvpack.ResidueCoordination(residueGroup1, residueGroup2, pbc=True, stepFunction='1/(1+x^6)', thresholdDistance=1.0 nm, normalize=False, weighByMass=True, includeHydrogens=True, name='residue_coordination')[source]

The number of contacts between two disjoint groups of residues:

\[N({\bf r}) = \sum_{i \in {\bf G}_1} \sum_{j \in {\bf G}_2} S\left( \frac{\|{\bf R}_j({\bf r}) - {\bf R}_i({\bf r})\|}{r_0} \right)\]

where \({\bf G}_1\) and \({\bf G}_2\) are the two groups of residues, \({\bf R}_i\) is the centroid of the residue \(i\), \(r_0\) is the threshold distance for defining a contact and \(S(x)\) is a step function equal to \(1\) if a contact is made or equal to \(0\) otherwise. For trajectory analysis, it is fine to make \(S(x) = H(1-x)\), where H is the Heaviside step function. For molecular dynamics, however, \(S(x)\) should be a continuous approximation of \(H(1-x)\) for \(x \geq 0\). By default [5], the following function is used:

\[S(x) = \frac{1-x^6}{1-x^{12}} = \frac{1}{1+x^6}\]
Parameters:
  • residueGroup1 (Iterable[Residue]) – The residues in the first group.

  • residueGroup2 (Iterable[Residue]) – The residues in the second group.

  • pbc (bool) – Whether the system has periodic boundary conditions.

  • stepFunction (str) – The function “step(1-x)” (for analysis only) or a continuous approximation thereof.

  • thresholdDistance (Quantity | Real) – The threshold distance (\(r_0\)) for considering two residues as being in contact.

  • normalize (bool) – Whether the number of contacts should be normalized by the total number of possible contacts.

  • weighByMass (bool) – Whether the centroid of each residue should be weighted by the mass of the atoms in the residue.

  • includeHydrogens (bool) – Whether hydrogen atoms should be included in the centroid calculations.

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

Raises:

ValueError – If the two groups of residues are not disjoint

Example

>>> import cvpack
>>> import itertools as it
>>> import openmm
>>> from openmm import app, unit
>>> from openmmtools import testsystems
>>> model = testsystems.LysozymeImplicit()
>>> group1 = list(it.islice(model.topology.residues(), 125, 142))
>>> print(*[r.name for r in group1])  
TRP ASP GLU ... ASN GLN THR
>>> group2 = list(it.islice(model.topology.residues(), 142, 156))
>>> print(*[r.name for r in group2])  
PRO ASN ARG ... ARG THR GLY
>>> residue_coordination = cvpack.ResidueCoordination(
...     group1,
...     group2,
...     stepFunction="step(1-x)",
...     thresholdDistance=0.6*unit.nanometers,
... )
>>> residue_coordination.addToSystem(model.system)
>>> platform = openmm.Platform.getPlatformByName('Reference')
>>> context = openmm.Context(
...     model.system, openmm.CustomIntegrator(0), platform
... )
>>> context.setPositions(model.positions)
>>> residue_coordination.getValue(context)
26.0 dimensionless
>>> residue_coordination.setReferenceValue(26 * unit.dimensionless)
>>> context.reinitialize(preserveState=True)
>>> residue_coordination.getValue(context)
0.99999... dimensionless

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]

getReferenceValue()[source]

Get the reference value used for normalizing the residue coordination.

Returns:

The reference value.

Return type:

mmunit.Quantity

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
setReferenceValue(value)[source]

Set the reference value used for normalizing the residue coordination.

Parameters:

value (Quantity | Real) – The reference value.