parmed.charmm package¶
Submodules¶
Module contents¶
This package contains code for reading CHARMM structure files for setting up a simulation with CHARMM; specifically PSF, PAR, RTF, and STR files
- PARParameter file (PRM) – this contains all of the force field
parameters (e.g., bond lengths and strengths) for all of the atom types
- RTFResidue Topology File – this file contains the residue
connectivity tables as well as a definition of all of the atom types. Also contains an internal coordinate representation of the residues
- PSFProtein Structure File – this is the main file type in CHARMM
simulations that defines all of the residues in a system as well as the atom types and connectivity between the atoms
- STRStream file – Source of additional information and CHARMM commands
that can contain RTF and PAR information. Allows users to define additional parameters without ‘contaminating’ the original force field parameter files
-
class
parmed.charmm.
CharmmCrdFile
(fname)[source]¶ Bases:
object
Reads and parses a CHARMM coordinate file (.crd) into its components, namely the coordinates, CHARMM atom types, resid, resname, etc.
- Parameters
- fnamestr
Name of the restart file to parse
- Attributes
- natomint
Number of atoms in the system
- resnamelist of str
List of all residue names in the system
- coordinatesnp.ndarray with shape (1, natom, 3)
Atomic coordinates in a numpy array
positions
natom x 3 distance QuantityAtomic coordinates with units attached to them with the shape (natom, 3)
Methods
id_format
(filename)Identifies the file type as a CHARMM coordinate file
write
(struct, dest)Writes a CHARMM coordinate file from a structure
-
property
box
¶
-
property
coordinates
¶
-
static
id_format
(filename)[source]¶ Identifies the file type as a CHARMM coordinate file
- Parameters
- filenamestr
Name of the file to check format for
- Returns
- is_fmtbool
True if it is a CHARMM coordinate file
-
property
positions
¶ Atomic coordinates with units attached to them with the shape (natom, 3)
-
static
write
(struct, dest)[source]¶ Writes a CHARMM coordinate file from a structure
- Parameters
- struct
parmed.structure.Structure
The input structure to write the CHARMM coordinate file from
- deststr or file-like object
The file name or file object to write the coordinate file to
- struct
-
class
parmed.charmm.
CharmmParameterSet
(*args)[source]¶ Bases:
parmed.parameters.ParameterSet
,parmed.charmm.parameters.CharmmImproperMatchingMixin
Stores a parameter set defined by CHARMM files. It stores the equivalent of the information found in the MASS section of the CHARMM topology file (TOP/RTF) and all of the information in the parameter files (PAR)
- Parameters
- *filenamesvariable length arguments of str
The list of topology, parameter, and stream files to load into the parameter set. The following file type suffixes are recognized:
.rtf, .top – Residue topology file .par, .prm – Parameter file .str – Stream file .inp – If “par” is in the file name, it is a parameter file. If
“top” is in the file name, it is a topology file. Otherwise, ValueError is raised.
See also
- 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 a CharmmParameterSet from another ParameterSet (or subclass).
from_structure
(struct)Extracts known parameters from a Structure instance
load_set
([tfile, pfile, sfiles])Instantiates a CharmmParameterSet from a Topology file and a Parameter file (or just a Parameter file if it has all information)
match_improper_type
(a1, a2, a3, a4)Matches an improper type based on atom type names
read_parameter_file
(pfile[, comments])Reads all of the parameters from a parameter file.
read_stream_file
(sfile)Reads RTF and PAR sections from a stream file and dispatches the sections to read_topology_file or read_parameter_file
read_topology_file
(tfile)Reads _only_ the atom type definitions from a topology file.
typeify_templates
()Assign atom types to atom names in templates
write
([top, par, str])Write a CHARMM parameter set to a file
-
classmethod
from_parameterset
(params, copy=False)[source]¶ Instantiates a CharmmParameterSet 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
- params
parmed.parameters.ParameterSet
ParameterSet containing the list of parameters to be converted to a CHARMM-compatible set
- copybool, optional
If True, the returned parameter set is a deep copy of
params
. If False, the returned parameter set is a shallow copy, and the original set may be modified if any lower-case atom type names are present. Default is False.
- params
- Returns
- new_paramsCharmmParameterSet
The parameter set whose atom type names are converted to all upper-case
-
classmethod
from_structure
(struct)[source]¶ Extracts known parameters from a Structure instance
- Parameters
- struct
parmed.structure.Structure
The parametrized
Structure
instance from which to extract parameters into a ParameterSet
- struct
- Returns
- params
ParameterSet
The parameter set with all parameters defined in the Structure
- params
Notes
The parameters here are copies of the ones in the Structure, so modifying the generated ParameterSet will have no effect on
struct
. Furthermore, the first occurrence of each parameter will be used. If future ones differ, they will be silently ignored, since this is expected behavior in some instances (like with Gromacs topologies in the ff99sb-ildn force field).Dihedrals are a little trickier. They can be multi-term, which can be represented either as a single entry in dihedrals with a type of DihedralTypeList or multiple entries in dihedrals with a DihedralType parameter type. In this case, the parameter is constructed from either the first DihedralTypeList found or the first DihedralType of each periodicity found if no matching DihedralTypeList is found.
-
classmethod
load_set
(tfile=None, pfile=None, sfiles=None)[source]¶ Instantiates a CharmmParameterSet from a Topology file and a Parameter file (or just a Parameter file if it has all information)
- Parameters
- tfilestr
The name of the Topology (RTF/TOP) file to parse
- pfilestr
The name of the Parameter (PAR/PRM) file to parse
- sfileslist(str)
Iterable of stream (STR) file names
- Returns
- New CharmmParameterSet populated with parameters found in the provided
- files
Notes
The RTF file is read first (if provided), followed by the PAR file, followed by the list of stream files (in the order they are provided). Parameters in each stream file will overwrite those that came before (or simply append to the existing set if they are different)
-
read_parameter_file
(pfile, comments=None)[source]¶ Reads all of the parameters from a parameter file. Versions 36 and later of the CHARMM force field files have an ATOMS section defining all of the atom types. Older versions need to load this information from the RTF/TOP files.
- Parameters
- pfilestr or list of lines
Name of the CHARMM parameter file to read or list of lines to parse as a file
- commentslist of str, optional
List of comments on each of the pfile lines (if pfile is a list of lines)
Notes
The atom types must all be loaded by the end of this routine. Either supply a PAR file with atom definitions in them or read in a RTF/TOP file first. Failure to do so will result in a raised RuntimeError.
-
read_stream_file
(sfile)[source]¶ Reads RTF and PAR sections from a stream file and dispatches the sections to read_topology_file or read_parameter_file
- Parameters
- sfilestr or CharmmStreamFile
Stream file to parse
-
read_topology_file
(tfile)[source]¶ Reads _only_ the atom type definitions from a topology file. This is unnecessary for versions 36 and later of the CHARMM force field.
- Parameters
- tfilestr
Name of the CHARMM topology file to read
-
write
(top=None, par=None, str=None)[source]¶ Write a CHARMM parameter set to a file
- Parameters
- topstr or file-like object, optional
If provided, the atom types will be written to this file in RTF format.
- parstr or file-like object, optional
If provided, the parameters will be written to this file in PAR format. Either this or the
str
argument must be provided- strstr or file-like object, optional
If provided, the atom types and parameters will be written to this file as separate RTF and PAR cards that can be read as a CHARMM stream file. Either this or the
par
argument must be provided
- Raises
- ValueError if both par and str are None
-
class
parmed.charmm.
CharmmPsfFile
(**kwargs)[source]¶ Bases:
parmed.structure.Structure
A chemical
Structure
instantiated from CHARMM files.- Parameters
- psf_namestr, optional
Name of the PSF file (it must exist)
- Raises
- IOErrorIf file
psf_name
does not exist - CharmmPsfErrorIf any parsing errors are encountered
- IOErrorIf file
- 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.
Clear the cmap list to prevent any CMAP parameters from being used
copy
(cls[, split_dihedrals])Makes a copy of the current structure as an instance of a specified subclass
createSystem
([params])Creates an OpenMM System object from the CHARMM PSF file.
from_structure
(struct[, copy])Instantiates a CharmmPsfFile from an input Structure instance.
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
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
load_parameters
(parmset[, copy_parameters])Loads parameters from a parameter set that was loaded via CHARMM RTF, PAR, and STR files.
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
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)
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_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
-
createSystem
(params=None, *args, **kwargs)[source]¶ Creates an OpenMM System object from the CHARMM PSF file. This is a shortcut for calling load_parameters followed by Structure.createSystem. If params is not None, load_parameters will be called on that parameter set, and Structure.createSystem will be called with the remaining args and kwargs
- Parameters
- paramsCharmmParameterSet=None
If not None, this parameter set will be loaded
See also
parmed.structure.Structure.createSystem()
In addition to params, this method also takes all arguments for
parmed.structure.Structure.createSystem()
-
classmethod
from_structure
(struct, copy=False)[source]¶ Instantiates a CharmmPsfFile from an input Structure instance. This method makes sure all atom types have uppercase-only names
- Parameters
- struct
parmed.structure.Structure
The input structure to convert to a CharmmPsfFile instance
- copybool, optional
If True, a copy of all items are made. Otherwise, the resulting CharmmPsfFile is a shallow copy
- struct
- Returns
- psf
CharmmPsfFile
CHARMM PSF file
- psf
- Raises
- ValueError if the functional form is not recognized or cannot be
- implemented through the PSF and parameter/stream files
Notes
If copy is False, the original object may have its atom type names changed if any of them have lower-case letters
-
load_parameters
(parmset, copy_parameters=True)[source]¶ Loads parameters from a parameter set that was loaded via CHARMM RTF, PAR, and STR files.
- Parameters
- parmset
CharmmParameterSet
List of all parameters
- copy_parametersbool, optional, default=True
If False, parmset will not be copied.
Not copying parmset will cause ParameterSet and Structure to share references to types. If you modify the original parameter set, the references in Structure list_types will be silently modified. However, if you change any reference in the parameter set, then that reference will no longer be shared with structure.
Example where the reference in ParameterSet is changed. The following will NOT modify the parameters in the psf:
psf.load_parameters(parmset, copy_parameters=False) parmset.angle_types[('a1', 'a2', a3')] = AngleType(1, 2)
The following WILL change the parameter in the psf because the reference has not been changed in
ParameterSet
:psf.load_parameters(parmset, copy_parameters=False) a = parmset.angle_types[('a1', 'a2', 'a3')] a.k = 10 a.theteq = 100
Extra care should be taken when trying this with dihedral_types. Since dihedral_type is a Fourier sequence, ParameterSet stores DihedralType for every term in DihedralTypeList. Therefore, the example below will STILL modify the type in the
Structure
list_types:parmset.dihedral_types[('a', 'b', 'c', 'd')][0] = DihedralType(1, 2, 3)
This assigns a new instance of DihedralType to an existing DihedralTypeList that ParameterSet and Structure are tracking and the shared reference is NOT changed.
Use with caution!
- parmset
- Raises
- ParameterError if any parameters cannot be found
Notes
If any dihedral or improper parameters cannot be found, I will try inserting wildcards (at either end for dihedrals and as the two central atoms in impropers) and see if that matches. Wild-cards will apply ONLY if specific parameters cannot be found.
This method will expand the dihedrals attribute by adding a separate Dihedral object for each term for types that have a multi-term expansion
-
class
parmed.charmm.
CharmmRstFile
(fname)[source]¶ Bases:
object
Reads and parses data, velocities and coordinates from a CHARMM restart file (.rst) of file name ‘fname’ into class attributes
- Parameters
- fnamestr
Name of the restart file to parse
- Attributes
- natomint
Number of atoms in the system
- resnamelist of str
Names of all residues in the system
- coordinatesnp.ndarray shape(1, natom, 3)
List of all coordinates in the format [x1, y1, z1, x2, y2, z2, …]
- coordinatesoldnp.ndarray shape(1, natom, 3)
List of all old coordinates in the format [x1, y1, z1, x2, y2, z2, …]
velocities
np.ndarray shape(1, natom, 3)Atomic velocities in Angstroms/picoseconds
positions
natom x 3 distance QuantityAtomic positions with units
positionsold
natom x 3 distance QuantityOld atomic positions with units
Methods
id_format
(filename)Identifies the file type as a CHARMM restart file
scan
-
property
box
¶
-
property
coordinates
¶
-
property
coordinatesold
¶
-
static
id_format
(filename)[source]¶ Identifies the file type as a CHARMM restart file
- Parameters
- filenamestr
Name of the file to check format for
- Returns
- is_fmtbool
True if it is a CHARMM restart file
-
property
positions
¶ Atomic positions with units
-
property
positionsold
¶ Old atomic positions with units
-
property
velocities
¶ Atomic velocities in Angstroms/picoseconds