ExtendedSpaceContext¶
- class openxps.ExtendedSpaceContext(context, extra_dofs, coupling_potential, integrator_template=None, biasing_potential=None)[source]¶
Wraps an openmm.Context object to include extra degrees of freedom (DOFs) and allow for extended phase-space (XPS) simulations.
Note: The system and integrator attached to the context are modified in-place.
A provided cvpack.MetaCollectiveVariable is added to the system to couple the physical DOFs and the extra ones. The integrator’s
stepmethod is replaced with a custom function that advances the physical and extension systems in tandem.- Parameters:
context (Context) – The original OpenMM context containing the physical system.
extra_dofs (Sequence[ExtraDOF]) – A group of extra degrees of freedom to be included in the XPS simulation.
coupling_potential (MetaCollectiveVariable) – A meta-collective variable defining the potential energy term that couples the physical system to the extra DOFs. It must have units of
kilojoules_per_moleor equivalent.integrator_template (Integrator | None) – An openmm.Integrator object to be used as a template for the algorithm that advances the extra DOFs. If not provided, the physical system’s integrator is used as a template.
biasing_potential (BiasingPotential | None) – A bias potential applied to the extra DOFs. If not provided, no bias is applied.
Example
>>> import openxps as xps >>> from math import pi >>> import cvpack >>> import openmm >>> from openmm import unit >>> from openmmtools import testsystems >>> model = testsystems.AlanineDipeptideVacuum() >>> umbrella_potential = cvpack.MetaCollectiveVariable( ... f"0.5*kappa*min(delta,{2*pi}-delta)^2; delta=abs(phi-phi0)", ... [cvpack.Torsion(6, 8, 14, 16, name="phi")], ... unit.kilojoule_per_mole, ... kappa=1000 * unit.kilojoule_per_mole / unit.radian**2, ... phi0=pi*unit.radian, ... ) >>> temp = 300 * unit.kelvin >>> integrator = openmm.LangevinMiddleIntegrator( ... temp, 1 / unit.picosecond, 4 * unit.femtosecond ... ) >>> platform = openmm.Platform.getPlatformByName("Reference") >>> mass = 3 * unit.dalton*(unit.nanometer/unit.radian)**2 >>> phi0 = xps.ExtraDOF("phi0", unit.radian, mass, xps.bounds.CIRCULAR) >>> height = 2 * unit.kilojoule_per_mole >>> sigma = 18 * unit.degree >>> context = xps.ExtendedSpaceContext( ... openmm.Context(model.system, integrator, platform), ... [phi0], ... umbrella_potential, ... biasing_potential=xps.MetadynamicsBias( ... [phi0], [sigma], height, temp, 10, [100] ... ), ... ) >>> context.setPositions(model.positions) >>> context.setVelocitiesToTemperature(temp) >>> context.setExtraValues([180 * unit.degree]) >>> context.setExtraVelocitiesToTemperature(temp) >>> context.getIntegrator().step(100) >>> context.getExtraValues() (Quantity(value=..., unit=radian),) >>> context.addBiasKernel() >>> state = context.getExtensionContext().getState(getEnergy=True) >>> state.getPotentialEnergy() Quantity(value=..., unit=kilojoule/mole)
Methods
- getExtensionContext()[source]¶
Get a reference to the OpenMM context containing the extension system.
- Returns:
The context containing the extension system.
- Return type:
mm.Context
- getExtraDOFs()[source]¶
Get the extra degrees of freedom included in the extended phase-space system.
- Returns:
A tuple containing the extra degrees of freedom.
- Return type:
t.Tuple[ExtraDOF]
- getExtraValues()[source]¶
Get the values of the extra degrees of freedom.
- Returns:
A tuple containing the values of the extra degrees of freedom.
- Return type:
t.Tuple[mmunit.Quantity]
- getExtraVelocities()[source]¶
Get the velocities of the extra degrees of freedom.
- Returns:
A tuple containing the velocities of the extra degrees of freedom.
- Return type:
t.Tuple[mmunit.Quantity]
- setExtraValues(values)[source]¶
Set the values of the extra degrees of freedom.
- Parameters:
values (Iterable[Quantity]) – A sequence of quantities containing the values and units of all extra degrees of freedom.
- setExtraVelocities(velocities)[source]¶
Set the velocities of the extra degrees of freedom.
- Parameters:
velocities (Iterable[Quantity]) – A dictionary containing the velocities of the extra degrees of freedom.