parmed.openmm.reporters module¶
This module contains a handful of extra reporter classes for OpenMM simulations
-
class
parmed.openmm.reporters.
EnergyMinimizerReporter
(f, volume=False, **kwargs)[source]¶ Bases:
parmed.openmm.reporters.StateDataReporter
This class acts as a simple energy reporter class for minimizations. This is not meant to be used as a reporter for OpenMM’s molecular dynamics routines, but instead passed a simulation object for printing out single-point energies.
- Parameters
- fstr or file-like
File name or object to write energies to
- volumebool, optional
If True, write the system volume (default is False)
Methods
describeNextReport
(*args, **kwargs)Disable this reporter inside MD
finalize
()Closes any open file
report
(simulation[, frame])Print out the current energy
-
class
parmed.openmm.reporters.
MdcrdReporter
(**kwargs)[source]¶ Bases:
object
MdcrdReporter prints a trajectory in ASCII Amber format. This reporter will be significantly slower than binary file reporters (like DCDReporter or NetCDFReporter).
- Parameters
- filestr
Name of the file to write the trajectory to
- reportIntervalint
Number of steps between writing trajectory frames
- crdsbool=True
Write coordinates to this trajectory file?
- velsbool=False
Write velocities to this trajectory file?
- frcsbool=False
Write forces to this trajectory file?
Notes
You can only write one of coordinates, forces, or velocities to a mdcrd file.
Methods
describeNextReport
(simulation)Get information about the next report this object will generate.
finalize
()Closes any open file
report
(simulation, state)Generate a report.
-
describeNextReport
(simulation)[source]¶ Get information about the next report this object will generate.
- Parameters
- simulation
app.Simulation
The Simulation to generate a report for
- simulation
- Returns
- nsteps, pos, vel, frc, eneint, bool, bool, bool, bool
nsteps is the number of steps until the next report pos, vel, frc, and ene are flags indicating whether positions, velocities, forces, and/or energies are needed from the Context
-
class
parmed.openmm.reporters.
NetCDFReporter
(**kwargs)[source]¶ Bases:
object
NetCDFReporter prints a trajectory in NetCDF format
Methods
describeNextReport
(simulation)Get information about the next report this object will generate.
finalize
()Closes any open file
report
(simulation, state)Generate a report.
-
describeNextReport
(simulation)[source]¶ Get information about the next report this object will generate.
- Parameters
- simulation
app.Simulation
The Simulation to generate a report for
- simulation
- Returns
- nsteps, pos, vel, frc, eneint, bool, bool, bool, bool
nsteps is the number of steps until the next report pos, vel, frc, and ene are flags indicating whether positions, velocities, forces, and/or energies are needed from the Context
-
-
class
parmed.openmm.reporters.
ProgressReporter
(**kwargs)[source]¶ Bases:
parmed.openmm.reporters.StateDataReporter
A class that prints out a progress report of how much MD (or minimization) has been done, how fast the simulation is running, and how much time is left (similar to the mdinfo file in Amber)
- Parameters
- fstr
The file name of the progress report file (overwritten each time)
- reportIntervalint
The step interval between which to write frames
- totalStepsint
The total number of steps that will be run in the simulation (used to estimate time remaining)
- potentialEnergybool, optional
Whether to print the potential energy (default True)
- kineticEnergybool, optional
Whether to print the kinetic energy (default True)
- totalEnergybool, optional
Whether to print the total energy (default True)
- temperaturebool, optional
Whether to print the system temperature (default True)
- volumebool, optional
Whether to print the system volume (default False)
- densitybool, optional
Whether to print the system density (default False)
- systemMassfloat or
unit.Quantity
The mass of the system used when reporting density (useful in instances where masses are set to 0 to constrain their positions)
Notes
In addition to the above,
ProgressReporter
also accepts arguments forStateDataReporter
Methods
describeNextReport
(simulation)Get information about the next report this object will generate.
finalize
()Closes any open file
report
(simulation, state)Generate a report and predict the time to completion (and current/overall MD performance)
-
describeNextReport
(simulation)[source]¶ Get information about the next report this object will generate.
- Parameters
- simulation
app.Simulation
The Simulation to generate a report for
- simulation
- Returns
- nsteps, pos, vel, frc, eneint, bool, bool, bool, bool
nsteps is the number of steps until the next report pos, vel, frc, and ene are flags indicating whether positions, velocities, forces, and/or energies are needed from the Context
-
class
parmed.openmm.reporters.
RestartReporter
(**kwargs)[source]¶ Bases:
object
Use a reporter to handle writing restarts at specified intervals.
- Parameters
- filestr
Name of the file to write the restart to.
- reportIntervalint
Number of steps between writing restart files
- write_multiplebool=False
Either write a separate restart each time (appending the step number in the format .# to the file name given above) if True, or overwrite the same file each time if False
- netcdfbool=False
Use the Amber NetCDF restart file format
- write_velocitiesbool=True
Write velocities to the restart file. You can turn this off for passing in, for instance, a minimized structure.
Methods
describeNextReport
(simulation)Get information about the next report this object will generate.
finalize
()No-op here
report
(sim, state)Generate a report.
-
describeNextReport
(simulation)[source]¶ Get information about the next report this object will generate.
- Parameters
- simulation
app.Simulation
The Simulation to generate a report for
- simulation
- Returns
- nsteps, pos, vel, frc, eneint, bool, bool, bool, bool
nsteps is the number of steps until the next report pos, vel, frc, and ene are flags indicating whether positions, velocities, forces, and/or energies are needed from the Context
-
class
parmed.openmm.reporters.
StateDataReporter
(**kwargs)[source]¶ Bases:
object
This class acts as a state data reporter for OpenMM simulations, but it is a little more generalized. Notable differences are:
It allows the units of the output to be specified, with defaults being those used in Amber (e.g., kcal/mol, angstroms^3, etc.)
It will write to any object containing a ‘write’ method; not just files. This allows, for instance, writing to a GUI window that implements the desired ‘write’ attribute.
Most of this code is copied from the OpenMM StateDataReporter class, with the above-mentioned changes made.
- Parameters
- fstr or file-like
Destination to write the state data (file name or file object)
- reportIntervalint
Number of steps between state data reports
- stepbool, optional
Print out the step number (Default True)
- timebool, optional
Print out the simulation time (Defaults True)
- potentialEnergybool, optional
Print out the potential energy of the structure (Default True)
- kineticEnergybool, optional
Print out the kinetic energy of the structure (Default True)
- totalEnergybool, optional
Print out the total energy of the system (Default True)
- temperaturebool, optional
Print out the temperature of the system (Default True)
- volumebool, optional
Print out the volume of the unit cell. If the system is not periodic, the value is meaningless (Default False)
- densitybool, optional
Print out the density of the unit cell. If the system is not periodic, the value is meaningless (Default False)
- separatorstr, optional
The string to separate data fields (Default ‘,’)
- systemMassfloat, optional
If not None, the density will be computed from this mass, since setting a mass to 0 is used to constrain the position of that particle. (Default None)
- energyUnitunit, optional
The units to print energies in (default unit.kilocalories_per_mole)
- timeUnitunit, optional
The units to print time in (default unit.picoseconds)
- volumeUnitunit, optional
The units print volume in (default unit.angstroms**3)
- densityUnitunit, optional
The units to print density in (default unit.grams/unit.item/unit.milliliter)
Methods
describeNextReport
(simulation)Get information about the next report this object will generate.
finalize
()Closes any open file
report
(simulation, state)Generate a report.
-
describeNextReport
(simulation)[source]¶ Get information about the next report this object will generate.
- Parameters
- simulation
app.Simulation
The simulation to generate a report for
- simulation
- Returns
- nsteps, pos, vel, frc, eneint, bool, bool, bool, bool
nsteps is the number of steps until the next report pos, vel, frc, and ene are flags indicating whether positions, velocities, forces, and/or energies are needed from the Context