utils¶
-
atomsmm.utils.
countDegreesOfFreedom
(system)[source]¶ Counts the number of degrees of freedom in a system, given by:
\[N_\mathrm{DOF} = 3N_\mathrm{moving particles} - 3 - N_\mathrm{constraints}\]- Parameters
system (openmm.System) – The system whose degrees of freedom will be summed up.
-
atomsmm.utils.
findNonbondedForce
(system, position=0)[source]¶ Searches for a NonbondedForce object in an OpenMM system.
- Parameters
system (openmm.System) – The system to which the wanted NonbondedForce object is attached.
position (int, optional, default=0) – The position index of the wanted force among the NonbondedForce objects attached to the system.
- Returns
int – The index of the wanted NonbondedForce object.
-
atomsmm.utils.
hijackForce
(system, index)[source]¶ Extracts a Force object from an OpenMM system.
Warning
Side-effect: the passed system object will no longer have the hijacked Force object in its force list.
-
atomsmm.utils.
splitPotentialEnergy
(system, topology, positions, **globals)[source]¶ Computes the potential energy of a system, with possible splitting into contributions of all Force objects attached to the system.
- Parameters
system (openmm.System) – The system whose energy is to be computed.
topology (openmm.app.topology.Topology) – The topological information about a system.
positions (list(tuple)) – A list of 3D vectors containing the positions of all atoms.
- Keyword Arguments
for global context variables. (Values) –
- Returns
unit.Quantity or dict(str, unit.Quantity) – The total potential energy or a dict containing all potential energy terms.
-
atomsmm.utils.
evaluateForce
(force, positions, boxVectors=None)[source]¶ Computes the value of a Force object for a given set of particle coordinates and box vectors. Whether periodic boundary conditions will be used or not depends on the corresponding attribute of the Force object specified as the collective variable.
- Parameters
positions (list(openmm.Vec3)) – A list whose length equals the number of particles in the system and which contains the coordinates of these particles.
boxVectors (list(openmm.Vec3), optional, default=None) – A list with three vectors which describe the edges of the simulation box.
- Returns
unit.Quantity
Example
>>> import afed >>> from simtk import unit >>> model = afed.AlanineDipeptideModel() >>> psi_angle, _ = model.getDihedralAngles() >>> psi = afed.DrivenCollectiveVariable('psi', psi_angle, unit.radians, period=360*unit.degrees) >>> psi.evaluate(model.getPositions()) Quantity(value=3.141592653589793, unit=radian)