StateDataReporter¶
- class cvpack.reporting.StateDataReporter(file, reportInterval, writers=(), **kwargs)[source]¶
An extended version of OpenMM’s openmm.app.StateDataReporter class that includes custom writers for reporting additional simulation data.
A custom writer is an object that includes the following methods:
getHeaders: returns a list of strings containing the headers to be added to the report. It must have the following signature:
def getHeaders(self) -> List[str]: pass
getValues: returns a list of floats containing the values to be added to the report at a given time step. It must have the following signature:
def getValues(self, simulation: openmm.app.Simulation) -> List[float]: pass
initialize (optional): performs any necessary setup before the first report. If present, it must have the following signature:
def initialize(self, simulation: openmm.app.Simulation) -> None: pass
- Parameters:
file (str | TextIOBase) – The file to write to. This can be a file name or a file object.
reportInterval (int) – The interval (in time steps) at which to report data.
writers (Sequence[CustomWriter]) – A sequence of custom writers.
**kwargs – Additional keyword arguments to be passed to the StateDataReporter constructor.
Examples
>>> import cvpack >>> import openmm >>> from math import pi >>> from cvpack import reporting >>> from openmm import app, unit >>> from sys import stdout >>> from openmmtools import testsystems >>> model = testsystems.AlanineDipeptideVacuum() >>> phi = cvpack.Torsion(6, 8, 14, 16, name="phi") >>> psi = cvpack.Torsion(8, 14, 16, 18, name="psi") >>> umbrella = cvpack.MetaCollectiveVariable( ... f"0.5*kappa*(min(dphi,{2*pi}-dphi)^2+min(dpsi,{2*pi}-dpsi)^2)" ... "; dphi=abs(phi-phi0); dpsi=abs(psi-psi0)", ... [phi, psi], ... unit.kilojoules_per_mole, ... name="umbrella", ... kappa=100 * unit.kilojoules_per_mole/unit.radian**2, ... phi0=5*pi/6 * unit.radian, ... psi0=-5*pi/6 * unit.radian, ... ) >>> reporter = reporting.StateDataReporter( ... stdout, ... 100, ... writers=[ ... reporting.CVWriter(umbrella, value=True, emass=True), ... reporting.MetaCVWriter( ... umbrella, ... values=["phi", "psi"], ... emasses=["phi", "psi"], ... parameters=["phi0", "psi0"], ... derivatives=["phi0", "psi0"], ... ), ... ], ... step=True, ... ) >>> integrator = openmm.LangevinIntegrator( ... 300 * unit.kelvin, ... 1 / unit.picosecond, ... 2 * unit.femtosecond, ... ) >>> integrator.setRandomNumberSeed(1234) >>> umbrella.addToSystem(model.system) >>> simulation = app.Simulation(model.topology, model.system, integrator) >>> simulation.context.setPositions(model.positions) >>> simulation.context.setVelocitiesToTemperature(300 * unit.kelvin, 5678) >>> simulation.reporters.append(reporter) >>> simulation.step(1000) #"Step","umbrella (kJ/mol)",...,"d[umbrella]/d[psi0] (kJ/(mol rad))" 100,11.26...,40.371... 200,7.463...,27.910... 300,2.558...,-12.74... 400,6.199...,3.9768... 500,8.827...,41.878... 600,3.761...,25.262... 700,3.388...,25.342... 800,1.071...,11.349... 900,8.586...,37.380... 1000,5.84...,31.159...