OpenMM
|
A chemical structure instantiated from CHARMM files. More...
Public Member Functions | |
def | __init__ |
Opens and parses a PSF file, then instantiates a CharmmPsfFile instance from the data. More... | |
def | writePsf |
Writes a PSF file from the stored molecule. More... | |
def | loadParameters |
Loads parameters from a parameter set that was loaded via CHARMM RTF, PAR, and STR files. More... | |
def | setBox |
Sets the periodic box boundary conditions. More... | |
def | topology |
Create an OpenMM Topology object from the stored bonded network. More... | |
def | createSystem |
Construct an OpenMM System representing the topology described by the prmtop file. More... | |
def | system |
Return the cached system class – it needs to be initialized via "createSystem" first! More... | |
def | boxLengths |
Return tuple of 3 units. More... | |
def | boxLengths |
def | boxVectors |
Return the box vectors. More... | |
def | boxVectors |
Sets the box vectors. More... | |
def | deleteCmap |
Deletes the CMAP terms from the CHARMM PSF. More... | |
Static Public Attributes | |
int | BOND_FORCE_GROUP = 0 |
int | ANGLE_FORCE_GROUP = 1 |
int | DIHEDRAL_FORCE_GROUP = 2 |
int | UREY_BRADLEY_FORCE_GROUP = 3 |
int | IMPROPER_FORCE_GROUP = 4 |
int | CMAP_FORCE_GROUP = 5 |
int | NONBONDED_FORCE_GROUP = 6 |
int | GB_FORCE_GROUP = 6 |
A chemical structure instantiated from CHARMM files.
Example:
cs = CharmmPsfFile("testfiles/test.psf")
This structure has numerous attributes that are lists of the elements of this structure, including atoms, bonds, torsions, etc. The attributes are
Additional attribute is available if a CharmmParameterSet is loaded into this structure.
The lengths of each of these lists gives the pointers (e.g., natom, nres, etc.)
Example:
cs = CharmmPsfFile("testfiles/test.psf") len(cs.atom_list)
33
len(cs.bond_list)
32
def __init__ | ( | self, | |
psf_name | |||
) |
Opens and parses a PSF file, then instantiates a CharmmPsfFile instance from the data.
psf_name (str) : Name of the PSF file (it must exist)
Exceptions Raised: IOError : If file "psf_name" does not exist CharmmPSFError: If any parsing errors are encountered
References simtk.openmm.app.charmmpsffile.set_molecules().
def boxLengths | ( | self | ) |
Return tuple of 3 units.
References CharmmPsfFile._boxLengths, and CharmmPsfFile.box_vectors.
Referenced by CharmmPsfFile.boxLengths(), CharmmPsfFile.createSystem(), and CharmmPsfFile.topology().
def boxLengths | ( | self, | |
stuff | |||
) |
References CharmmPsfFile.boxLengths().
def boxVectors | ( | self | ) |
Return the box vectors.
References CharmmPsfFile.box_vectors.
Referenced by CharmmPsfFile.boxVectors().
def boxVectors | ( | self, | |
stuff | |||
) |
Sets the box vectors.
References CharmmPsfFile._boxLengths, CharmmPsfFile.box_vectors, and CharmmPsfFile.boxVectors().
def createSystem | ( | self, | |
params, | |||
nonbondedMethod = ff.NoCutoff , |
|||
nonbondedCutoff = 1.0*u.nanometer , |
|||
switchDistance = 0.0*u.nanometer , |
|||
constraints = None , |
|||
rigidWater = True , |
|||
implicitSolvent = None , |
|||
implicitSolventKappa = None , |
|||
implicitSolventSaltConc = 0.0*u.moles/u.liter , |
|||
temperature = 298.15*u.kelvin , |
|||
soluteDielectric = 1.0 , |
|||
solventDielectric = 78.5 , |
|||
removeCMMotion = True , |
|||
hydrogenMass = None , |
|||
ewaldErrorTolerance = 0.0005 , |
|||
flexibleConstraints = True , |
|||
verbose = False |
|||
) |
Construct an OpenMM System representing the topology described by the prmtop file.
You MUST have loaded a parameter set into this PSF before calling createSystem. If not, AttributeError will be raised. ValueError is raised for illegal input.
params | (CharmmParameterSet) The parameter set to use to parametrize this molecule |
nonbondedMethod | (object=NoCutoff) The method to use for nonbonded interactions. Allowed values are NoCutoff, CutoffNonPeriodic, CutoffPeriodic, Ewald, or PME. |
nonbondedCutoff | (distance=1*nanometer) The cutoff distance to use for nonbonded interactions. |
switchDistance | (distance=0*nanometer) The distance at which the switching function is active for nonbonded interactions. If the switchDistance evaluates to boolean False (if it is 0), no switching function will be used. Illegal values will raise a ValueError |
constraints | (object=None) Specifies which bonds or angles should be implemented with constraints. Allowed values are None, HBonds, AllBonds, or HAngles. |
rigidWater | (boolean=True) If true, water molecules will be fully rigid regardless of the value passed for the constraints argument |
implicitSolvent | (object=None) If not None, the implicit solvent model to use. Allowed values are HCT, OBC1, OBC2, or GBn |
implicitSolventKappa | (float=None): Debye screening parameter to model salt concentrations in GB solvent. |
implicitSolventSaltConc | (float=0.0*u.moles/u.liter): Salt concentration for GB simulations. Converted to Debye length `kappa' |
temperature | (float=298.15*u.kelvin): Temperature used in the salt concentration-to-kappa conversion for GB salt concentration term |
soluteDielectric | (float=1.0) The solute dielectric constant to use in the implicit solvent model. |
solventDielectric | (float=78.5) The solvent dielectric constant to use in the implicit solvent model. |
removeCMMotion | (boolean=True) If true, a CMMotionRemover will be added to the System. |
hydrogenMass | (mass=None) The mass to use for hydrogen atoms bound to heavy atoms. Any mass added to a hydrogen is subtracted from the heavy atom to keep their total mass the same. |
ewaldErrorTolerance | (float=0.0005) The error tolerance to use if the nonbonded method is Ewald or PME. |
flexibleConstraints | (bool=True) Are our constraints flexible or not? |
verbose | (bool=False) Optionally prints out a running progress report |
References CharmmPsfFile._get_gb_params(), CharmmPsfFile._system, Context._system, CharmmPsfFile.ANGLE_FORCE_GROUP, CharmmPsfFile.angle_list, CharmmPsfFile.atom_list, CharmmPsfFile.BOND_FORCE_GROUP, CharmmPsfFile.bond_list, CharmmPsfFile.box_vectors, CharmmPsfFile.boxLengths(), CharmmPsfFile.CMAP_FORCE_GROUP, CharmmPsfFile.cmap_list, CharmmPsfFile.DIHEDRAL_FORCE_GROUP, CharmmPsfFile.dihedral_parameter_list, CharmmPsfFile.GB_FORCE_GROUP, CharmmPsfFile.IMPROPER_FORCE_GROUP, CharmmPsfFile.improper_list, CharmmPsfFile.loadParameters(), CharmmPsfFile.NONBONDED_FORCE_GROUP, CharmmPsfFile.UREY_BRADLEY_FORCE_GROUP, and CharmmPsfFile.urey_bradley_list.
Referenced by CharmmPsfFile.topology().
def deleteCmap | ( | self | ) |
Deletes the CMAP terms from the CHARMM PSF.
References CharmmPsfFile.cmap_list.
def loadParameters | ( | self, | |
parmset | |||
) |
Loads parameters from a parameter set that was loaded via CHARMM RTF, PAR, and STR files.
parmset | (CharmmParameterSet) : List of all parameters |
Notes:
If | any parameters that are necessary cannot be found, a MissingParameter exception is raised. |
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 dihedral_parameter_list attribute by adding a separate Dihedral object for each term for types that have a multi-term expansion |
References CharmmPsfFile.atom_list, and CharmmPsfFile.bond_list.
Referenced by CharmmPsfFile.createSystem().
def setBox | ( | self, | |
a, | |||
b, | |||
c, | |||
alpha = 90.0*u.degrees , |
|||
beta = 90.0*u.degrees , |
|||
gamma = 90.0*u.degrees |
|||
) |
Sets the periodic box boundary conditions.
a,b,c | (floats) : Lengths of the periodic cell |
alpha,beta,gamma | (floats, optional) : Angles between the periodic cells. |
References CharmmPsfFile._boxLengths, CharmmPsfFile._topology, and CharmmPsfFile.box_vectors.
def system | ( | self | ) |
Return the cached system class – it needs to be initialized via "createSystem" first!
References CharmmPsfFile._system, and Context._system.
def topology | ( | self | ) |
Create an OpenMM Topology object from the stored bonded network.
References CharmmPsfFile._topology, CharmmPsfFile.atom_list, CharmmPsfFile.bond_list, CharmmPsfFile.box_vectors, CharmmPsfFile.boxLengths(), and CharmmPsfFile.createSystem().
Referenced by DesmondDMSFile.getTopology().
def writePsf | ( | self, | |
dest, | |||
vmd = False |
|||
) |
Writes a PSF file from the stored molecule.
dest | (str or file-like) : The place to write the output PSF file. If it has a "write" attribute, it will be used to print the PSF file. Otherwise, it will be treated like a string and a file will be opened, printed, then closed |
vmd | (bool) : If True, it will write out a PSF in the format that VMD prints it in (i.e., no NUMLP/NUMLPH or MOLNT sections) Example:
|
References CharmmPsfFile.acceptor_list, CharmmPsfFile.angle_list, CharmmPsfFile.atom_list, CharmmPsfFile.bond_list, CharmmPsfFile.cmap_list, CharmmPsfFile.dihedral_list, CharmmPsfFile.donor_list, CharmmPsfFile.flags, CharmmPsfFile.group_list, CharmmPsfFile.improper_list, simtk.openmm.app.charmmpsffile.set_molecules(), CharmmCrdFile.title, CharmmRstFile.title, and CharmmPsfFile.title.
acceptor_list |
Referenced by CharmmPsfFile.writePsf().
|
static |
Referenced by CharmmPsfFile.createSystem().
angle_list |
Referenced by CharmmPsfFile.createSystem(), and CharmmPsfFile.writePsf().
atom_list |
|
static |
Referenced by CharmmPsfFile.createSystem().
bond_list |
box_vectors |
|
static |
Referenced by CharmmPsfFile.createSystem().
cmap_list |
Referenced by CharmmPsfFile.createSystem(), CharmmPsfFile.deleteCmap(), and CharmmPsfFile.writePsf().
|
static |
Referenced by CharmmPsfFile.createSystem().
dihedral_list |
Referenced by CharmmPsfFile.writePsf().
dihedral_parameter_list |
Referenced by CharmmPsfFile.createSystem().
donor_list |
Referenced by CharmmPsfFile.writePsf().
flags |
Referenced by CharmmPsfFile.writePsf().
|
static |
Referenced by CharmmPsfFile.createSystem().
group_list |
Referenced by CharmmPsfFile.writePsf().
|
static |
Referenced by CharmmPsfFile.createSystem().
improper_list |
Referenced by CharmmPsfFile.createSystem(), and CharmmPsfFile.writePsf().
|
static |
Referenced by CharmmPsfFile.createSystem().
residue_list |
title |
Referenced by CharmmPsfFile.writePsf().
|
static |
Referenced by CharmmPsfFile.createSystem().
urey_bradley_list |
Referenced by CharmmPsfFile.createSystem().