HelixRMSDContent

class cvpack.HelixRMSDContent(residues, numAtoms, thresholdRMSD=0.08 nm, stepFunction='(1+x^4)/(1+x^4+x^8)', normalize=False, name='helix_rmsd_content')[source]

The alpha-helix RMSD content of a sequence of n residues [3]:

\[\alpha_{\rm rmsd}({\bf r}) = \sum_{i=1}^{n-5} S\left( \frac{r_{\rm rmsd}({\bf g}_i, {\bf g}_{\rm ref})}{r_0} \right)\]

where \({\bf g}_i\) represents a group of atoms selected from residues \(i\) to \(i+5\) of the sequence, \({\bf g}_{\rm ref}\) represents the same atoms in an ideal alpha-helix configuration, \(r_{\rm rmsd}({\bf g}, {\bf g}_{\rm ref})\) is the root-mean-square distance (RMSD) between groups \({\bf g}\) and \({\bf g}_{\rm ref}\), \(r_0\) is a threshold RMSD value, and \(S(x)\) is a smooth step function whose default form is

\[S(x) = \frac{1 + x^4}{1 + x^4 + x^8}\]

The residues must be a contiguous sequence from a single chain, ordered from the N-terminus to the C-terminus.

Every group \({\bf g}_{i,j}\) is formed by the N, \({\rm C}_\alpha\), \({\rm C}_\beta\), C, and O atoms of the six residues involvend, thus comprising 30 atoms in total. In the case of glycine, the missing \({\rm C}_\beta\) atom is replaced by the corresponding H atom. The RMSD is then defined as

\[r_{\rm rmsd}({\bf g}, {\bf g}_{\rm ref}) = \sqrt{\frac{1}{30} \sum_{j=1}^{30} \left\| \hat{\bf r}_j({\bf g}) - {\bf A}({\bf g})\hat{\bf r}_j({\bf g}_{\rm ref}) \right\|^2}\]

where \(\hat{\bf r}_k({\bf g})\) is the position of the \(k\)-th atom in a group \({\bf g}\) relative to the group’s center of geometry and \({\bf A}({\bf g})\) is the rotation matrix that minimizes the RMSD between \({\bf g}\) and \({\bf g}_{\rm ref}\).

Optionally, the alpha-helix RMSD content can be normalized to the range \([0, 1]\). This is done by dividing its value by \(N_{\rm groups} = n - 5\).

This collective variable was introduced in Ref. [3]. The default step function shown above is identical to the one in the original paper, but written in a numerically safe form. The ideal alpha-helix configuration is the same used for the collective variable ALPHARMSD in PLUMED v2.8.1 .

Note

The present implementation is limited to \(1 \leq N_{\rm groups} \leq 1024\).

Parameters:
  • residues (Sequence[Residue]) – The residue sequence to be used in the calculation.

  • numAtoms (int) – The total number of atoms in the system (required by OpenMM).

  • thresholdRMSD (Quantity | Real) – The threshold RMSD value for considering a group of residues as a close match to an ideal alpha-helix.

  • stepFunction (str) – The form of the step function \(S(x)\).

  • normalize (bool) – Whether to normalize the collective variable to the range \([0, 1]\).

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

Example

>>> import itertools as it
>>> import cvpack
>>> import openmm
>>> from openmm import app, unit
>>> from openmmtools import testsystems
>>> model = testsystems.LysozymeImplicit()
>>> residues = list(it.islice(model.topology.residues(), 59, 80))
>>> print(*[r.name for r in residues])
LYS ASP GLU ... ILE LEU ARG
>>> helix_content = cvpack.HelixRMSDContent(
...     residues, model.system.getNumParticles()
... )
>>> helix_content.getNumResidueBlocks()
16
>>> helix_content.addToSystem(model.system)
>>> normalized_helix_content = cvpack.HelixRMSDContent(
...     residues, model.system.getNumParticles(), normalize=True
... )
>>> normalized_helix_content.addToSystem(model.system)
>>> platform = openmm.Platform.getPlatformByName('Reference')
>>> integrator = openmm.VerletIntegrator(0)
>>> context = openmm.Context(model.system, integrator, platform)
>>> context.setPositions(model.positions)
>>> helix_content.getValue(context)
15.98... dimensionless
>>> normalized_helix_content.getValue(context)
0.998... 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]

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