parmed.gromacs package¶
Module contents¶
Contains classes for parsing GROMACS topology and parameter files
-
class
parmed.gromacs.
GromacsGroFile
[source]¶ Bases:
object
Parses and writes Gromacs GRO files
Methods
id_format
(filename)Identifies the file as a GROMACS GRO file
parse
(filename[, skip_bonds])Parses a Gromacs GRO file
write
(struct, dest[, precision, nobox, combine])Write a Gromacs Topology File from a Structure
-
static
id_format
(filename)[source]¶ Identifies the file as a GROMACS GRO file
- Parameters
- filenamestr
Name of the file to check if it is a Gromacs GRO file
- Returns
- is_fmtbool
If it is identified as a Gromacs GRO file, return True. False otherwise
-
static
parse
(filename, skip_bonds=False)[source]¶ Parses a Gromacs GRO file
- Parameters
- filenamestr or file-like
Name of the file or the GRO file object
- skip_bondsbool, optional
If True, skip trying to assign bonds. This can save substantial time when parsing large files with non-standard residue names. However, no bonds are assigned. This is OK if, for instance, the GRO file is being parsed simply for its coordinates. This will also reduce the accuracy of assigned atomic numbers for typical ions. Default is False.
- Returns
- struct
Structure
The Structure instance instantiated with just residues and atoms populated (with coordinates)
- struct
-
static
write
(struct, dest, precision=3, nobox=False, combine=False)[source]¶ Write a Gromacs Topology File from a Structure
- Parameters
- struct
Structure
The structure to write to a Gromacs GRO file (must have coordinates)
- deststr or file-like
The name of a file or a file object to write the Gromacs topology to
- precisionint, optional
The number of decimal places to print in the coordinates. Default 3
- noboxbool, optional
If the system does not have a periodic box defined, and this option is True, no box will be written. If False, the periodic box will be defined to enclose the solute with 0.5 nm clearance on all sides. If periodic box dimensions are defined, this variable has no effect.
- combine‘all’, None, or list of iterables, optional
Equivalent to the combine argument of the GromacsTopologyFile.write method. If None, system atom order may be changed to meet the need for contiguously bonded groups of atoms to be part of a single moleculetype. All other values leave the atom order unchanged. Default is None.
- struct
-
static
-
class
parmed.gromacs.
GromacsTopologyFile
(fname=None, defines=None, parametrize=True, xyz=None, box=None)[source]¶ Bases:
parmed.structure.Structure
Class providing a parser and writer for a GROMACS topology file
- Parameters
- fnamestr
The name of the file to read
- defineslist of str=None
If specified, this is the set of defines to use when parsing the topology file
- parametrizedbool, optional
If True, parameters are assigned after parsing is done from the parametertypes sections. If False, only parameter types defined in the parameter sections themselves are loaded (i.e., on the same line as the parameter was defined). Default is True
- xyzstr or array, optional
The source of atomic coordinates. It can be a string containing the name of a coordinate file from which to fill the coordinates (and optionally the unit cell information), or it can be an array with the coordinates. Default is None
- boxarray, optional
If provided, the unit cell information will be set from this variable. If provided, it must be a collection of 6 floats representing the unit cell dimensions a, b, c, alpha, beta, and gamma, respectively. Default is None.
Notes
If the
xyz
argument is a file name that contains the unit cell information, this unit cell information is set. However, thebox
argument takes precedence and will override values given in the coordinate file unless it has its default value ofNone
.- Attributes
- box
box_vectors
3, 3-element tuple of unit cell vectors that are Quantity objects of
- combining_rule
- coordinates
positions
A list of 3-element Quantity tuples of dimension length representing the atomic positions for every atom in the system.
topology
The OpenMM Topology object.
velocities
A (natom, 3)-shape numpy array with atomic velocities for every atom in
view
Returns an indexable object that can be indexed like a standard
Methods
add_atom
(atom, resname, resnum[, chain, …])Adds a new atom to the Structure, adding a new residue to residues if it has a different name or number as the last residue added and adding it to the atoms list.
add_atom_to_residue
(atom, residue)Adds a new atom to the Structure at the end if the given residue
assign_bonds
(*reslibs)Assigns bonds to all atoms based on the provided residue template libraries.
copy
(cls[, split_dihedrals])Makes a copy of the current structure as an instance of a specified subclass
createSystem
([nonbondedMethod, …])Construct an OpenMM System representing the topology described by the prmtop file.
from_structure
(struct[, copy])Instantiates a GromacsTopologyFile instance from a Structure
get_box
([frame])In some cases, multiple conformations may be stored in the Structure.
get_coordinates
([frame])In some cases, multiple conformations may be stored in the Structure.
has_NBFIX
()Returns whether or not any pairs of atom types have their LJ interactions modified by an NBFIX definition
id_format
(filename)Identifies the file as a GROMACS topology file
is_changed
()Determines if any of the topology has changed for this structure
join_dihedrals
()Joins multi-term torsions into a single term and makes all of the parameters DihedralTypeList instances.
load_dataframe
(df)Loads atomic properties from an input DataFrame
omm_add_constraints
(system, constraints, …)Adds constraints to a given system
omm_angle_force
([constraints, …])Creates an OpenMM HarmonicAngleForce object (or AmoebaAngleForce if the angles are for an Amoeba-parametrized system)
omm_bond_force
([constraints, rigidWater, …])Creates an OpenMM Bond Force object (or AmoebaBondForce if the bonds are for an Amoeba-parametrized system)
omm_cmap_force
()Creates the OpenMM CMAP torsion force
omm_dihedral_force
([split])Creates the OpenMM PeriodicTorsionForce modeling dihedrals
omm_gbsa_force
(implicitSolvent[, …])Creates a Generalized Born force for running implicit solvent calculations
omm_improper_force
()Creates the OpenMM improper torsion force (quadratic bias)
omm_nonbonded_force
([nonbondedMethod, …])Creates the OpenMM NonbondedForce instance
omm_out_of_plane_bend_force
()Creates the Amoeba out-of-plane bend force
omm_pi_torsion_force
()Creates the Amoeba pi-torsion force
omm_rb_torsion_force
()Creates the OpenMM RBTorsionForce for Ryckaert-Bellemans torsions
omm_set_virtual_sites
(system)Sets the virtual sites in a given OpenMM System object from the extra points defined in this system
omm_stretch_bend_force
()Create the OpenMM Amoeba stretch-bend force for this system
omm_torsion_torsion_force
()Create the OpenMM Amoeba coupled-torsion (CMAP) force
omm_trigonal_angle_force
()Creates the Amoeba trigonal-angle force
omm_urey_bradley_force
()Creates the OpenMM Urey-Bradley force
Assign parameters to the current structure.
prune_empty_terms
()Looks through all of the topological lists and gets rid of terms in which at least one of the atoms is None or has an idx attribute set to -1 (indicating that it has been removed from the atoms atom list)
read
(fname[, defines, parametrize])Reads the topology file into the current instance
save
(fname[, format, overwrite])Saves the current Structure in the requested file format.
split
()Split the current Structure into separate Structure instances for each unique molecule.
strip
(selection)Deletes a subset of the atoms corresponding to an atom-based selection.
to_dataframe
()Generates a DataFrame from the current Structure’s atomic properties
unchange
()Toggles all lists so that they do not indicate any changes
update_dihedral_exclusions
()Nonbonded exclusions and exceptions have the following priority:
visualize
(*args, **kwargs)Use nglview for visualization.
write
(dest[, combine, parameters, molfile, itp])Write a Gromacs Topology File from a Structure
write_cif
(dest[, renumber, coordinates, …])Write a PDB file from the current Structure instance
write_pdb
(dest[, renumber, coordinates, …])Write a PDB file from a Structure instance
write_psf
(dest[, vmd])Writes a PSF file from the stored molecule
-
copy
(cls, split_dihedrals=False)[source]¶ Makes a copy of the current structure as an instance of a specified subclass
- Parameters
- clsStructure subclass
The returned object is a copy of this structure as a cls instance
- split_dihedrals
bool
If True, then the Dihedral entries will be split up so that each one is paired with a single DihedralType (rather than a DihedralTypeList)
- Returns
- cls instance
The instance of the Structure subclass cls with a copy of the current Structure’s topology information
-
classmethod
from_structure
(struct, copy=False)[source]¶ Instantiates a GromacsTopologyFile instance from a Structure
- Parameters
- struct
parmed.Structure
The input structure to generate from
- copybool, optional
If True, assign from a copy of
struct
(this is a lot slower). Default is False
- struct
- Returns
- gmxtop
GromacsTopologyFile
The topology file defined by the given struct
- gmxtop
-
static
id_format
(filename)[source]¶ Identifies the file as a GROMACS topology file
- Parameters
- filenamestr
Name of the file to check if it is a gromacs topology file
- Returns
- is_fmtbool
If it is identified as a gromacs topology, return True. False otherwise
-
read
(fname, defines=None, parametrize=True)[source]¶ Reads the topology file into the current instance
-
write
(dest, combine=None, parameters='inline', molfile=None, itp=False)[source]¶ Write a Gromacs Topology File from a Structure
- Parameters
- deststr or file-like
The name of a file or a file object to write the Gromacs topology to
- combine‘all’, None, or list of iterables, optional
If None, no molecules are combined into a single moleculetype. If ‘all’, all molecules are combined into a single moleculetype. Otherwise, the list of molecule indices (start from 0) will control which atoms are combined into single moleculetype’s. Default is None
- parameters‘inline’ or str or file-like object, optional
This specifies where parameters should be printed. If ‘inline’ (default), the parameters are written on the same lines as the valence terms are defined on. Any other string is interpreted as a filename for an ITP that will be written to and then included at the top of dest. If it is a file-like object, parameters will be written there. If parameters is the same as
dest
, then the parameter types will be written to the same topologyfile.- molfileNone or str of file-like object, optional
If specified as other than None, the molecules will be written to a separate file that is included in the main topology file. The name of this file will be the provided srting. If None or the same as the ``dest’’, the molecules will be written into the body of the topology file. If it is a file-like object, the molecules will be written there. Using this option can make it easier to combine multiple molecules into the same topology. This will change where the following topology sections are written: moleculetype, atoms, bonds, pairs, angles, dihedrals, cmap, settles, virtual_sites2, virtual_sites3 and exclusions.
- itpbool, optional
If True the following topology sections are not written: defaults, atomtypes, nonbond_params, bondtypes, pairtypes, angletypes, dihedraltypes, cmaptypes, system and molecules Thus only the individual molecules will be written in a stand-alone fashion, i.e. an itp-file. If True the molfile parameter will be set to None
- Raises
- ValueError if the same molecule number appears in multiple combine lists
- TypeError if the dest input cannot be parsed
- ValueError if the combine, parameters, or molfile input cannot be parsed