MetaCVWriter

class cvpack.reporting.MetaCVWriter(metaCV, values=(), emasses=(), parameters=(), derivatives=())[source]

A custom writer for reporting meta-collective variable data.

Parameters:
  • metaCV (MetaCollectiveVariable) – The meta-collective variable whose associated values will be reported.

  • innerValues – The names of the inner variables whose values will be reported.

  • innerMasses – The names of the inner variables whose effective masses will be reported.

  • parameters (Sequence[str]) – The names of the parameters whose values will be reported.

  • derivatives (Sequence[str]) – The names of the parameters with respect to which the derivatives of the meta-collective variable will be reported.

Examples

>>> import cvpack
>>> import openmm
>>> from math import pi
>>> from cvpack.reporting import StateDataReporter, MetaCVWriter
>>> 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 = StateDataReporter(
...     stdout,
...     100,
...     writers=[
...         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","phi (rad)",...,"d[umbrella]/d[psi0] (kJ/(mol rad))"
100,2.36849...,40.3718...
200,2.88515...,27.9109...
300,2.43112...,-12.743...
400,2.96786...,3.97688...
500,2.58383...,41.8782...
600,2.72482...,25.2626...
700,2.55836...,25.3424...
800,2.71046...,11.3498...
900,2.43913...,37.3804...
1000,2.7584...,31.1599...

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]