ExtensionWriter

class openxps.ExtensionWriter(context, potential=False, kinetic=False, total=False, temperature=False)[source]

A custom writer for reporting state data from an extension context.

Parameters:
  • context (ExtendedSpaceContext) – The extended space context whose extension-context state data will be reported.

  • potential (bool) – If True, the potential energy of the extension context will be reported.

  • kinetic (bool) – If True, the kinetic energy of the extension context will be reported.

  • total (bool) – If True, the total energy of the extension context will be reported.

  • temperature (bool) – If True, the temperature of the extension context will be reported.

Example

>>> import openxps as xps
>>> from math import pi
>>> from sys import stdout
>>> import cvpack
>>> import openmm
>>> from openmm import app, 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,
... )
>>> integrator = openmm.LangevinMiddleIntegrator(
...     300 * unit.kelvin, 1 / unit.picosecond, 4 * unit.femtosecond
... )
>>> integrator.setRandomNumberSeed(1234)
>>> platform = openmm.Platform.getPlatformByName("Reference")
>>> simulation = app.Simulation(
...     model.topology, model.system, integrator, platform
... )
>>> mass = 3 * unit.dalton*(unit.nanometer/unit.radian)**2
>>> phi0 = xps.ExtraDOF("phi0", unit.radian, mass, xps.bounds.CIRCULAR)
>>> context = xps.ExtendedSpaceContext(
...     simulation.context, [phi0], umbrella_potential
... )
>>> context.setPositions(model.positions)
>>> context.setVelocitiesToTemperature(300 * unit.kelvin, 1234)
>>> context.setExtraValues([180 * unit.degree])
>>> context.setExtraVelocitiesToTemperature(300 * unit.kelvin, 1234)
>>> reporter = cvpack.reporting.StateDataReporter(
...     stdout,
...     10,
...     step=True,
...     kineticEnergy=True,
...     writers=[xps.ExtensionWriter(context, kinetic=True)],
... )
>>> simulation.reporters.append(reporter)
>>> simulation.step(100)  
#"Step","Kinetic Energy (kJ/mole)","Extension Kinetic Energy (kJ/mole)"
10,60.512...,1.7013...
20,75.765...,2.5089...
30,61.116...,1.3375...
40,52.359...,0.4791...
50,61.382...,0.7065...
60,48.674...,0.6520...
70,60.771...,1.3525...
80,46.518...,2.0280...
90,66.111...,0.9597...
100,60.94...,0.9695...

Methods

getHeaders()[source]

Gets a list of strigs containing the headers to be added to the report.

Return type:

List[str]

getValues(simulation)[source]

Gets a list of floats containing the values to be added to the report.

Parameters:

simulation (Simulation) – The simulation object.

Return type:

List[float]

initialize(simulation)[source]

Initializes the writer. This method is called before the first report and can be used to perform any necessary setup.

Parameters:

simulation (Simulation) – The simulation object.