HelixTorsionContent¶
- class cvpack.HelixTorsionContent(residues, pbc=False, phiReference=-63.8 deg, psiReference=-41.1 deg, tolerance=25 deg, halfExponent=3, normalize=False, name='helix_torsion_content')[source]¶
The alpha-helix Ramachandran content of a sequence of n residues:
\[\alpha_{\phi,\psi}({\bf r}) = \frac{1}{2} \sum_{i=2}^{n-1} \left[ B_m\left( \frac{\phi_i({\bf r}) - \phi_{\rm ref}}{\theta_{\rm tol}} \right) + B_m\left( \frac{\psi_i({\bf r}) - \psi_{\rm ref}}{\theta_{\rm tol}} \right) \right]\]where \(\phi_i({\bf r})\) and \(\psi_i({\bf r})\) are the Ramachandran dihedral angles of residue \(i\), \(\phi_{\rm ref}\) and \(\psi_{\rm ref}\) are their reference values in an alpha helix [4], and \(\theta_{\rm tol}\) is the threshold tolerance around these refenrences. \(B_m(x)\) is a smooth boxcar function given by
\[B_m(x) = \frac{1}{1 + x^{2m}}\]where \(m\) is an integer parameter that controls its steepness. Note that \(x\) needs to be elevated to an even power for \(B_m(x)\) to be an even function.
Optionally, this collective variable can be normalized to the range \([0, 1]\).
Note
The \(\phi\) and \(\psi\) angles of the first and last residues are not considered. They are used to compute the dihedral angles of the second and penultimate residues, respectively.
The residues must be a contiguous sequence from a single chain, ordered from the N- to the C-terminus. Due to an OpenMM limitation, the maximum supported number of residues is 37.
- Parameters:
residues (Sequence[Residue]) – The residues in the sequence.
pbc (bool) – Whether to use periodic boundary conditions.
phiReference (Quantity | Real) – The reference value of the phi dihedral angle in an alpha helix.
psiReference (Quantity | Real) – The reference value of the psi dihedral angle in an alpha helix.
tolerance (Quantity | Real) – The threshold tolerance around the reference values.
halfExponent (int) – The parameter \(m\) of the boxcar function.
normalize (bool) – Whether to normalize the collective variable to the range \([0, 1]\).
name (str) – The name of the collective variable.
- Raises:
ValueError – If some residue does not contain a \(\phi\) or \(\psi\) angle.
Example
>>> import cvpack >>> import openmm >>> from openmm import app, unit >>> from openmmtools import testsystems >>> model = testsystems.LysozymeImplicit() >>> residues = list(model.topology.residues())[59:80] >>> print(*[r.name for r in residues]) LYS ASP GLU ... ILE LEU ARG >>> helix_content = cvpack.HelixTorsionContent(residues) >>> 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) 17.452... 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