PathInRMSDSpace¶
- class cvpack.PathInRMSDSpace(metric, milestones, numAtoms, sigma, name=None)[source]¶
A metric of the system’s progress (\(s\)) or deviation (\(z\)) with respect to a path defined by a sequence of \(n\) milestones defined as reference structures [6]:
\[s({\bf r}) = \frac{ \dfrac{\sum_{i=1}^n i w_i({\bf r})}{\sum_{i=1}^n w_i({\bf r})} - 1 }{n-1} \quad \text{or} \quad z({\bf r}) = - 2 \sigma ^2 \ln \sum_{i=1}^n w_i({\bf r})\]with \(w_i({\bf r})\) being a Gaussian kernel centered at the \(i\)-th milestone, i.e.,
\[w_i({\bf r}) = \exp\left(\ -\frac{d^2_{\rm rms}({\bf r},{\bf r}^{\rm ref}_i)}{2 \sigma^2} \right)\]where \(d_{\rm rms}({\bf r},{\bf r}^{\rm ref}_i)\) is the root-mean-square distance between the current system state and the \(i\)-th reference structure and \(\sigma\) sets the width of the kernels.
Note
The kernel width \(\sigma\) is related to the parameter \(\lambda\) of Ref. [6] by \(\sigma = \frac{1}{\sqrt{2\lambda}}\).
- Parameters:
metric (Metric) – The path-related metric to compute. Use
cvpack.path.progress
for computing \(s({\bf r})\) orcvpack.path.deviation
for computing \(z({\bf r})\).milestones (Sequence[Dict[int, Quantity | ndarray | Vec3]]) – A sequence of reference structures, each represented as a dictionary mapping atom indices to coordinate vectors.
numAtoms (int) – The total number of atoms in the system, including those that are not in any of the reference structures.
sigma (Quantity) – The width of the Gaussian kernels in nanometers
name (str | None) – The name of the collective variable. If not provided, it is set to “path_progress_in_rmsd_space” or “path_deviation_in_rmsd_space” depending on the metric.
- Raises:
ValueError – The number of milestones is less than 2
ValueError – If the metric is not cvpack.path.progress or cvpack.path.deviation
Examples
>>> import cvpack >>> import networkx as nx >>> import numpy as np >>> import openmm >>> from openmm import unit >>> from openmmtools import testsystems >>> from scipy.spatial.transform import Rotation >>> model = testsystems.AlanineDipeptideVacuum() >>> atom1, atom2 = 8, 14 >>> graph = model.mdtraj_topology.to_bondgraph() >>> nodes = list(graph.nodes) >>> graph.remove_edge(nodes[atom1], nodes[atom2]) >>> movable = list(nx.connected_components(graph))[1] >>> x = model.positions / model.positions.unit >>> x0 = x[atom1, :] >>> vector = x[atom2, :] - x0 >>> vector /= np.linalg.norm(vector) >>> rotation = Rotation.from_rotvec((np.pi / 6) * vector) >>> atoms = [nodes.index(atom) for atom in movable] >>> frames = [x.copy()] >>> for _ in range(6): ... x[atoms, :] = x0 + rotation.apply(x[atoms, :] - x0) ... frames.append(x.copy()) >>> milestones = [ ... {i: row for i, row in enumerate(frame)} ... for frame in frames ... ] >>> s, z = [ ... cvpack.PathInRMSDSpace( ... metric, milestones, len(x), 0.5 * unit.angstrom ... ) ... for metric in (cvpack.path.progress, cvpack.path.deviation) ... ] >>> s.addToSystem(model.system) >>> z.addToSystem(model.system) >>> context = openmm.Context(model.system, openmm.VerletIntegrator(0.001)) >>> context.setPositions(model.positions) >>> s.getValue(context) 0.172... dimensionless >>> z.getValue(context) -0.004... nm**2
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