parmed.openmm package¶
Submodules¶
Module contents¶
This is a collection of all of the OpenMM functionality supported in ParmEd
- 
class parmed.openmm.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.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
- simulationapp.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.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
- simulationapp.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.OpenMMParameterSet(*filenames)[source]¶
- Bases: - parmed.parameters.ParameterSet,- parmed.charmm.parameters.CharmmImproperMatchingMixin- Class storing parameters from an OpenMM parameter set - Parameters
- filenamesstr, list of str, file-like, or list of file-like; optional
- Either the name of a file or a list of filenames from which parameters should be parsed. 
 
 - See also - Notes - Order is important in the list of files provided. The parameters are loaded in the order they are provided, and any parameters that are specified in multiple places are overwritten (that is, the last occurrence is the parameter type that is used) - Attributes
- combining_rule
 
 - Methods - condense([do_dihedrals])- This function goes through each of the parameter type dicts and eliminates duplicate types. - from_parameterset(params[, copy, …])- Instantiates an OpenMMParameterSet from another ParameterSet (or subclass). - from_structure(struct[, …])- Extracts known parameters from a Structure instance - id_format(filename)- Identifies the file type as either an Amber-style frcmod or parm.dat file. - match_improper_type(a1, a2, a3, a4)- Matches an improper type based on atom type names - typeify_templates()- Assign atom types to atom names in templates - write(dest[, provenance, write_unused, …])- Write the parameter set to an XML file for use with OpenMM - 
classmethod from_parameterset(params, copy=False, remediate_residues=True, unique_atom_types=False)[source]¶
- Instantiates an OpenMMParameterSet from another ParameterSet (or subclass). The main thing this feature is responsible for is converting lower-case atom type names into all upper-case and decorating the name to ensure each atom type name is unique. - Parameters
- paramsparmed.parameters.ParameterSet
- ParameterSet containing the list of parameters to be converted to as OpenMM-compatible parameter set 
- copybool, optional, default=False
- If True, the returned parameter set is a deep copy of - params. If False, the returned parameter set is a shallow copy. Default False.
- remediate_residuesbool, optional, default=True
- If True, will remove non-chemical bonds and drop Residue definitions that are missing parameters 
- unique_atom_typesbool
- If True, a unique OpenMM atom type will be created for every atom of every residue. In this case, the - AtomTypeobjects correspond to atom classes rather than atom types (in the OpenMM terminology).
 
- params
- Returns
- new_paramsOpenMMParameterSet
- OpenMMParameterSet with the same parameters as that defined in the input parameter set 
 
 
 - 
static id_format(filename)[source]¶
- Identifies the file type as either an Amber-style frcmod or parm.dat file. - Parameters
- filenamestr
- Name of the file to check format for 
 
- Returns
- is_fmtbool
- True if it is an Amber-style parameter file. False otherwise. 
 
 
 - 
write(dest, provenance=None, write_unused=True, separate_ljforce=False, improper_dihedrals_ordering='default', charmm_imp=False, skip_duplicates=True)[source]¶
- Write the parameter set to an XML file for use with OpenMM - Parameters
- deststr or file-like
- The name of the file or the file-like object (with a - writeattribute) to which the XML file will be written
- provenancedict, optional
- If present, the XML file will be tagged with the available fields. Keys of the dictionary become XML etree.Element tags, the values of the dictionary must be instances of any of: - str / unicode (Py2) or str (Py3) - one XML element with this content is written - list - one XML element per each item of the list is written, all these XML elements use the same tag (key in provenance dict) - dict - one of the keys of this dict must be the same as the key of of the provenance dict under which this dict is nested. The value under this key becomes the content of the XML element. Remaining keys and their values are used to construct attributes of the XML element. Note that OrderedDict’s should be used to ensure appropriate order of the XML elements and their attributes. Default is no provenance. Example (unordered): provenance = {‘Reference’ : [‘Nature’, ‘Cell’], - ‘Source’ : {‘Source’: ‘leaprc.ff14SB’, sourcePackage : ‘AmberTools’, sourcePackageVersion : ‘15’}, ‘User’ : ‘Mark’} 
- write_unusedbool
- If False: a) residue templates using unavailable atom types will not be written, b) atom types that are not used in any of the residue templates remaining and parameters including those atom types will not be written. A ParameterWarning is issued if any such residues are found in a). 
- separate_ljforcebool
- If True will use a separate LennardJonesForce to create a CostumNonbondedForce to compute L-J interactions. It will set sigma to 1 and epsilon to 0 in the NonbondedForce so that the NonbondedForce only calculates the electrostatic contribution. It should be set to True when converting a CHARMM force field file that doesn’t have pair-specific L-J modifications (NBFIX in CHARMM) so that the ffxml conversion is compatible with the main charmm36.xml file. Note: —- When pair-specific L-J modifications are present (NBFIX in CHARMM), this behavior is always present and this flag is ignored. 
- improper_dihedrals_orderingstr
- The ordering to use when assigning improper torsions in OpenMM. Default is ‘default’, other option is ‘amber’ 
- charmm_imp: bool
- If True, will check for existence of IMPR in each residue and patch template, and write out the explicit improper definition without wildcards in the ffxml file. 
- skip_duplicatesbool
- If True: residues which appear identical to an existing residue will not be written. This is usually the best choice. The most common reason for setting skip_duplicates to false is if you have different parametrizations for stereoisomers. Note that since OpenMM’s residue hashing is not aware of chirality, if you wish to use the results in simulations you will need to explicitly provide the template names for affected residues. 
 
 - Notes - The generated XML file will have the XML tag - DateGeneratedadded to the provenance information set to the current date. Therefore, you should not provide this information in- provenance(it will be removed if it is provided).
 
- 
class parmed.openmm.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, - ProgressReporteralso accepts arguments for- StateDataReporter- 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
- simulationapp.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.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
- simulationapp.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.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
- simulationapp.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.XmlFile[source]¶
- Bases: - object- Wrapper for parsing OpenMM-serialized objects. Supports serialized State, System, Integrator, and ForceField objects. - Methods - id_format(filename)- Identifies the file type as an XML file - parse(*args, **kwargs)- Parses XML file and returns deserialized object. - 
static id_format(filename)[source]¶
- Identifies the file type as an XML file - Parameters
- filenamestr
- Name of the file to check format for 
 
- Returns
- is_fmtbool
- True if it is an XML format, False otherwise 
 
 
 - 
static parse(*args, **kwargs)[source]¶
- Parses XML file and returns deserialized object. The return value depends on the serialized object, summarized below - System : returns simtk.openmm.System 
- State : returns simtk.openmm.State 
- Integrator : returns simtk.openmm.Integrator subclass 
- ForceField : returns simtk.openmm.app.ForceField 
 - Parameters
- filenamestr or file-like
- The file name or file object containing the XML-serialized object 
 
- Returns
- objSystem, State, Integrator, or ForceField
- The deserialized object 
 
 - Notes - OpenMM requires the entire contents of this file read into memory. As a result, this function may require a significant amount of memory. 
 
- 
static 
- 
parmed.openmm.energy_decomposition(structure, context, nrg=Unit({BaseUnit(base_dim=BaseDimension('amount'), name='mole', symbol='mol'): - 1.0, ScaledUnit(factor=4.184, master=kilojoule, name='kilocalorie', symbol='kcal'): 1.0}))[source]¶
- This computes the energy of every force group in the given structure and computes the energy for each force group for the given Context. Note, the context must have positions already assigned. - Parameters
- structureStructure
- This should be the Structure object from which the System object in the Context was created 
- contextmm.Context
- The OpenMM context set up for computing forces and energies 
- nrgenergy unit, optional
- The unit to convert all energies into. Default is kcal/mol 
 
- Returns
- dict {str:float}
- A dictionary mapping the name of the force group (taken from the attribute names of the format XXX_FORCE_GROUP in the structure object) with the energy of that group in 
 
 
- 
parmed.openmm.energy_decomposition_system(structure, system, platform=None, nrg=Unit({BaseUnit(base_dim=BaseDimension('amount'), name='mole', symbol='mol'): - 1.0, ScaledUnit(factor=4.184, master=kilojoule, name='kilocalorie', symbol='kcal'): 1.0}))[source]¶
- This function computes the energy contribution for all of the different force groups. - Parameters
- structureStructure
- The Structure with the coordinates for this system 
- systemmm.System
- The OpenMM System object to get decomposed energies from 
- platformstr
- The platform to use. Options are None (default), ‘CPU’, ‘Reference’, ‘CUDA’, and ‘OpenCL’. None will pick default OpenMM platform for this installation and computer 
- nrgenergy unit, optional
- The unit to convert all energies into. Default is kcal/mol 
 
- Returns
- energieslist of tuple
- Each entry is a tuple with the name of the force followed by its contribution 
 
 
- 
parmed.openmm.load_topology(topology, system=None, xyz=None, box=None, condense_atom_types=True)[source]¶
- Creates a - parmed.structure.Structureinstance from an OpenMM Topology, optionally filling in parameters from a System- Parameters
- topologysimtk.openmm.app.Topology
- The Topology instance with the list of atoms and bonds for this system 
- systemsimtk.openmm.Systemor str, optional
- If provided, parameters from this System will be applied to the Structure. If a string is given, it will be interpreted as the file name of an XML-serialized System, and it will be deserialized into a System before used to supply parameters 
- xyzstr or array of float
- Name of a file containing coordinate information or an array of coordinates. If file has unit cell information, it also uses that information unless - box(below) is also specified
- boxarray of 6 floats
- Unit cell dimensions 
- condense_atom_typesbool, default=True
- If True, create unique atom types based on de-duplicating properties. If False, create one atom type for each atom in the system, even if its properties match an existing atom type. 
 
- topology
- Returns
- structStructure
- The structure from the provided topology 
 
- struct
- Raises
- OpenMMWarning if parameters are found that cannot be interpreted or
- processed by ParmEd
- TypeError if there are any mismatches between the provided topology and
- system (e.g., they have different numbers of atoms)
- IOError if system is a string and it is not an existing file
 
 - Notes - Due to its flexibility with CustomForces, it is entirely possible that the functional form of the potential will be unknown to ParmEd. This function will try to use the energy expression to identify supported potential types that are implemented as CustomForce objects. In particular, quadratic improper torsions, when recognized, will be extracted. - Other CustomForces, including the CustomNonbondedForce used to implement NBFIX (off-diagonal L-J modifications) and the 12-6-4 potential, will not be processed and will result in an unknown functional form. - If an OpenMM Atom.id attribute is populated by a non-integer, it will be used to name the corresponding ParmEd AtomType object.