openxps

class openxps.openxps.AuxiliaryVariable(id, unit, boundary, minval, maxval, mass, **parameters)[source]

Bases: object

An extended phase-space variable whose dynamics is coupled to that of one of more collective variables through a potential energy term.

Parameters
  • id (str) – a valid identifier string.

  • unit (UnitOrStr) – the unity of measurement of this auxiliary variable.

  • boundary (str) – boundary conditions to be applied to this auxiliary variable. Valid options are “periodic” and “reflective”.

  • minval (QuantityOrFloat) – the minimum allowable value for this auxiliary variable.

  • maxval (QuantityOrFloat) – the maximum allowable value for this auxiliary variable.

  • mass (QuantityOrFloat) – the mass assigned to this auxiliary variable, whose unit of measurement must be compatible with dalton*(nanometer/unit)**2, where unit is the auxiliary variable’s own unit (see above).

Keyword Arguments

**parameters – Names and values of additional parameters for this auxiliary variable.

Example

>>> import math
>>> import openxps as xps
>>> import openmm as mm
>>> from openmm import unit
>>> cv = xps.CollectiveVariable('psi', mm.CustomTorsionForce('theta'), 2*math.pi, 'radians')
>>> index = cv.force.addTorsion(0, 1, 2, 3, [])
>>> mass = 50*unit.dalton*(unit.nanometer/unit.radians)**2
>>> limit = 180*unit.degrees
>>> xps.AuxiliaryVariable('s_psi', 'radian', 'periodic', -limit, limit, mass)
<s_psi in [-3.141592653589793 rad, 3.141592653589793 rad], periodic, m=50 nm**2 Da/(rad**2)>
class openxps.openxps.CollectiveVariable(id, force, unit, period)[source]

Bases: object

A scalar-valued function of the particle coordinates, evaluated by means of an OpenMM Force object.

Quoting OpenMM’s CustomCVForce manual entry:

“Each collective variable is defined by a Force object. The Force’s potential energy is computed, and that becomes the value of the variable. This provides enormous flexibility in defining collective variables, especially by using custom forces. Anything that can be computed as a potential function can also be used as a collective variable.”

Parameters
  • id (str) – a valid identifier string.

  • force (mm.Force) – an OpenMM Force object whose energy function is used to evaluate this collective variable.

  • unit (UnitOrStr) – the unity of measurement of this collective variable.

  • period (Union[QuantityOrFloat, None]) – the period of this collective variable, if it is periodic. If period=None, it will be considered as non-periodic.

Raises
  • TypeError – if period has a unit of measurement that is incompatible with unit.

  • ValueError – if id is not a valid identifier string (like a Python variable name)

Example

>>> import openmm as mm
>>> import openxps as xps
>>> dihedral_angle = mm.CustomTorsionForce('theta')
>>> index = dihedral_angle.addTorsion(0, 1, 2, 3, [])
>>> xps.CollectiveVariable('psi', dihedral_angle, 'radians', 360*mm.unit.degrees)
psi: CustomTorsionForce (period=6.283185307179586 rad)
>>> distance = mm.CustomBondForce('10*r')
>>> index = distance.addBond(1, 18, [])
>>> xps.CollectiveVariable('distance', distance, 'angstroms', None)
distance: CustomBondForce (non-periodic, unit=A)
effective_mass(system, positions)[source]

Computes the effective mass of the collective variable for a given system and a given set of particle coordinates.

The effective mass of a collective variable \(q(\mathbf{r})\) is defined as [1]:

\[m_\mathrm{eff} = \left( \sum_{j=1}^N \frac{1}{m_j} \left\|\frac{dq}{d\mathbf{r}_j}\right\|^2 \right)^{-1}\]
Parameters
  • system (openmm.openmm.System) – the system for which the collective variable will be evaluated.

  • positions (List[openmm.vec3.Vec3]) – a list whose size equals the number of particles in the system and which contains the coordinates of these particles.

Return type

Union[openmm.unit.quantity.Quantity, float]

Example

>>> import openxps as xps
>>> from openmm import app
>>> model = xps.AlanineDipeptideModel()
>>> force_field = app.ForceField('amber03.xml')
>>> system = force_field.createSystem(model.topology)
>>> print(model.phi.effective_mass(system, model.positions))
0.0479588726559707 nm**2 Da/(rad**2)
>>> print(model.psi.effective_mass(system, model.positions))
0.05115582071188152 nm**2 Da/(rad**2)
evaluate(system, positions)[source]

Computes the value of the collective variable for a given system and a given set of particle coordinates.

Parameters
  • system (openmm.openmm.System) – the system for which the collective variable will be evaluated.

  • positions (List[openmm.vec3.Vec3]) – a list whose size equals the number of particles in the system and which contains the coordinates of these particles.

Return type

Union[openmm.unit.quantity.Quantity, float]

Example

>>> import openxps as xps
>>> from openmm import app
>>> model = xps.AlanineDipeptideModel()
>>> force_field = app.ForceField('amber03.xml')
>>> system = force_field.createSystem(model.topology)
>>> print(model.phi.evaluate(system, model.positions))
3.141592653589793 rad
>>> print(model.psi.evaluate(system, model.positions))
3.141592653589793 rad