CompositeRMSDForce¶
- class CompositeRMSDForce(referencePositions)¶
Bases:
Force
This is a force whose energy equals the root mean squared deviation (RMSD) between the current coordinates of multiple particle groups and reference structures defined for each group. The composite RMSD is computed after aligning every particle group to its reference while keeping the relative orientations of all groups fixed. Therefore, it is not a mere composition of independent RMSD calculations done in the usual way.
Consider m particle groups, with each group \({\bf g}_j\) containing \(n_j\) particles. The centroid of each group is given by
\[{\bf c}_j = \frac{1}{n_j} \sum_{k \in {\bf g}_j} {\bf r}_k\]where \({\bf r}_k\) is the position of particle \(k\) and the sum is over all particles in the group. Analogously, \({\bf c}_j^{\rm ref}\) is the centroid of the reference structure for group \(j\), whose each particle \(k\) is located at \({\bf r}_k^{\rm ref}\).
The composite RMSD is then defined as
\[d_{\rm crms}({\bf r}) = \sqrt{ \frac{1}{n} \min_{ \bf q \in \mathbb{R}^4 \atop \|{\bf q}\| = 1 } \sum_{j=1}^m \sum_{k \in {\bf g}_j} \left\| {\bf A}({\bf q})\left({\bf r}_k - {\bf c}_j\right) - {\bf r}_k^{\rm ref} + {\bf c}_j^{\rm ref} \right\|^2 }\]where \(n = \sum_{j=1}^m n_j\) is the total number of particles in all groups, \({\bf A}(\bf q)\) is the rotation matrix corresponding to a unit quaternion \({\bf q}\). Note that a single rotation matrix is used to align all groups, thus keeping their relative orientations fixed.
This force is intended for use with openmm.CustomCVForce.
- Parameters:
referencePositions – the reference positions to compute the deviation from. The length of this vector must equal the number of particles in the system, even if not all particles are used in computing the Composite RMSD.
Examples
>>> import openmmcppforces as mmcpp >>> import openmm as mm >>> from openmmtools import testsystems >>> model = testsystems.HostGuestVacuum() >>> host_atoms, guest_atoms = ( ... [a.index for a in r.atoms()] ... for r in model.topology.residues() ... ) >>> force = mmcpp.CompositeRMSDForce(model.positions) >>> _ = force.addGroup(host_atoms) >>> _ = force.addGroup(guest_atoms) >>> force.setForceGroup(1) >>> _ = model.system.addForce(force) >>> context = mm.Context( ... model.system, ... mm.VerletIntegrator(1.0 * unit.femtoseconds), ... mm.Platform.getPlatformByName('Reference'), ... ) >>> context.setPositions(model.positions) >>> state = context.getState(getEnergy=True, groups={1}) >>> state.getPotentialEnergy().value_in_unit(unit.kilojoules_per_mole) 0.0 >>> model.positions[guest_atoms] += 1.0 * unit.nanometers >>> context.setPositions(model.positions) >>> state = context.getState(getEnergy=True, groups={1}) >>> state.getPotentialEnergy().value_in_unit(unit.kilojoules_per_mole) 0.0
- __init__(referencePositions)¶
Methods
addGroup
(particles)Add a group of particles to make part of the composite RMSD calculation.
getForceGroup
(self)Get the force group this Force belongs to.
getGroup
(index)Get the particles of a group included in the composite RMSD calculation.
getName
(self)Get the name of this Force.
Get the number of particle groups included in the composite RMSD calculation.
Get the reference positions to compute the deviation from.
setForceGroup
(self, group)Set the force group this Force belongs to.
setGroup
(index, particles)Set the particles of a group included in the composite RMSD calculation.
setName
(self, name)Set the name of this Force.
setReferencePositions
(positions)Set the reference positions to compute the deviation from.
updateParametersInContext
(context)Update the reference positions and particle groups in a Context to match those stored in this OpenMM::Force object.
Returns whether or not this force makes use of periodic boundary conditions.
Attributes
thisown
The membership flag
- addGroup(particles)¶
Add a group of particles to make part of the composite RMSD calculation.
- Parameters:
particles – the indices of the particles to add
- Returns:
int – the index of the group that was added
- static cast(force)¶
Returns whether or not this force makes use of periodic boundary conditions.
- Returns:
bool – Whether force uses PBC
- getForceGroup(self) int ¶
Get the force group this Force belongs to.
- getGroup(index)¶
Get the particles of a group included in the composite RMSD calculation.
- Parameters:
index – the index of the group whose particles are to be retrieved
- Returns:
List[int] – the indices of the particles in the group
- getName(self) std::string const & ¶
Get the name of this Force. This is an arbitrary, user modifiable identifier. By default it equals the class name, but you can change it to anything useful.
- getNumGroups()¶
Get the number of particle groups included in the composite RMSD calculation.
- getReferencePositions()¶
Get the reference positions to compute the deviation from.
- setForceGroup(self, group)¶
Set the force group this Force belongs to.
- Parameters:
group (int) – the group index. Legal values are between 0 and 31 (inclusive).
- setGroup(index, particles)¶
Set the particles of a group included in the composite RMSD calculation.
- Parameters:
index – the index of the group whose particles are to be set
particles – the indices of the particles to include
- setName(self, name)¶
Set the name of this Force. This is an arbitrary, user modifiable identifier. By default it equals the class name, but you can change it to anything useful.
- setReferencePositions(positions)¶
Set the reference positions to compute the deviation from.
- Parameters:
positions – the reference positions to compute the deviation from. The length of this vector must equal the number of particles in the system, even if not all particles are used in computing the composite RMSD.
- updateParametersInContext(context)¶
Update the reference positions and particle groups in a Context to match those stored in this OpenMM::Force object. This method provides an efficient way to update these parameters in an existing openmm.Context without needing to reinitialize it. Simply call
setReferencePositions()
andsetGroup()
to modify this object’s parameters, then callupdateParametersInContext()
to copy them over to the openmm.Context.
- usesPeriodicBoundaryConditions()¶
Returns whether or not this force makes use of periodic boundary conditions.
- Returns:
bool – Whether force uses PBC