ExtraDOF¶
- class openxps.ExtraDOF(name, unit, mass, bounds)[source]¶
Extra degree of freedom for extended phase-space simulations with OpenMM.
- Parameters:
name (str) – The name of the context parameter to be turned into an extra degree of freedom.
unit (openmm.unit.unit.Unit) – The unity of measurement of this extra degree of freedom. It must be compatible with OpenMM’s MD unit system (mass in
dalton, distance innanometer, angle inradian, time inpicosecond, temperature inkelvin, energy inkilojoules_per_mol). If the extra degree of freedom does not have a unit, usedimensionless.mass (openmm.unit.quantity.Quantity) – The mass assigned to this extra degree of freedom, whose unit of measurement must be compatible with
dalton*(nanometer/unit)**2, whereunitis the extra degree of freedom’s own unit (see above).bounds (openxps.bounds.Bounds | None) – The boundary condition to be applied to this extra degree of freedom. It must be an instance of
openxps.bounds.Periodic,openxps.bounds.Reflective, orNone(for unbounded variables). If it is notNone, its unit of measurement must be compatible with the extra degree of freedom’s own unit.
Example
>>> import openxps as xps >>> import yaml >>> from openmm import unit >>> dv = xps.ExtraDOF( ... "psi", ... unit.radian, ... 3 * unit.dalton*(unit.nanometer/unit.radian)**2, ... xps.bounds.Periodic(-180, 180, unit.degree) ... ) >>> dv ExtraDOF(name='psi', unit=rad, mass=3 nm**2 Da/(rad**2), bounds=...) >>> dv.bounds Periodic(lower=-3.14159..., upper=3.14159..., unit=rad) >>> assert yaml.safe_load(yaml.safe_dump(dv)) == dv
Methods
- createCollectiveVariable(particle, name=None)[source]¶
Returns a cvpack.OpenMMForceWrapper object associating this extra degree of freedom with the x coordinate of a particle in an OpenMM system.
- Parameters:
particle (int) – The index of the particle whose x coordinate will be associated with this extra degree of freedom.
name (str | None) – The name of the context parameter to be used in the OpenMM system. If
None, the name of the extra degree of freedom will be used.
- Returns:
The collective variable object representing this extra degree of freedom.
- Return type:
cvpack.OpenMMForceWrapper
Example
>>> import openxps as xps >>> from openmm import unit >>> xdof = xps.ExtraDOF( ... "psi", ... unit.radian, ... 3 * unit.dalton*(unit.nanometer/unit.radian)**2, ... xps.bounds.Periodic(-180, 180, unit.degree) ... ) >>> cv = xdof.createCollectiveVariable(0)
- distanceTo(other)[source]¶
Returns a Lepton expression representing the distance between this extra degree of freedom and another variable.
- Parameters:
other (CollectiveVariable) – The other variable to which the distance will be calculated.
- Returns:
A string representing the distance between the two variables.
- Return type:
str
Example
>>> import openxps as xps >>> import pytest >>> from openmm import unit >>> xdof = xps.ExtraDOF( ... "psi0", ... unit.radian, ... 3 * unit.dalton*(unit.nanometer/unit.radian)**2, ... xps.bounds.Periodic(-180, 180, unit.degree) ... ) >>> psi = cvpack.Torsion(6, 8, 14, 16, name="psi") >>> xdof.distanceTo(psi) '(psi-psi0-6.28318...*floor(0.5+(psi-psi0)/6.28318...))' >>>
- isPeriodic()[source]¶
Returns whether this extra degree of freedom is periodic.
- Returns:
Trueif this extra degree of freedom is periodic,Falseotherwise.- Return type:
bool
Example
>>> import openxps as xps >>> from openmm import unit >>> dv = xps.ExtraDOF( ... "psi", ... unit.radian, ... 3 * unit.dalton*(unit.nanometer/unit.radian)**2, ... xps.bounds.Periodic(-180, 180, unit.degree) ... ) >>> dv.isPeriodic() True