hps.core package
Submodules
hps.core.geometry module
- class hps.core.geometry.geometry[source]
Bases:
object
A class to hold methods for calculating geometrical values given sets of atom coordinates.
- static bond(coord1: Quantity, coord2: Quantity)[source]
Calculate the distance length between two (x,y,z) quantity coordinates.
- Parameters:
coord1 (openmm.unit.quantity.Quantity array) – Vector for the first coordinate.
coord2 (openmm.unit.quantity.Quantity array) – Vector for the second coordinate.
- Returns:
Quantity (value and unit) of the distance length in nanometers.
- Return type:
simtk.unit.quantity.Quantity
- static position2Array(position: Quantity, output_unit)[source]
Converts an OpenMM position object quantity into a numpy array.
- Parameters:
position (openmm.unit.quantity.Quantity) – Array containing quantity objects [e.g. (x,y,z) array returned from positions].
output_unit (openmm.unit.nanometer) – Unit in which to return the items of the array.
- Returns:
A numpy array containing the quantity values in unit of nm, converted to float.
- Return type:
numpy.ndarray
hps.core.models module
- class hps.core.models.models[source]
Bases:
object
A class to hold functions for the automated generation of default hps models.
- static buildHPSModel(structure_file: str, minimize: bool = False, model: str = 'hps_urry', box_dimension: Optional[Any] = None)[source]
Creates an alpha-carbon only
hpsOpenMM
system class object with default initialized parameters.- Parameters:
structure_file (string [required]) – Path to the input structure file.
minimize (boolean (False)) – If True, the initial structure will undergo energy minimization.
model (string [Optional, hps_urry]) – HPS scale. Available options are ‘hps_urry’, ‘hps_ss’, ‘hps_kr’, and ‘mpipi’.
box_dimension (float or array (None)) – If box_dimension is supplied, a PBC will be used. If a float is given, a cubic box will be used. If an array of (3,1) is given, a rectangular box with the given dimension will be used. If not specified, PBC will not be used.
- Returns:
hps – Initialized hpsOpenMM.system class with default options for defining a coarse-grained CA force field.
- Return type:
hpsOpenMM.system
hps.core.system module
- class hps.core.system.system(structure_path: str, model: str = 'hps_urry')[source]
Bases:
object
A class for generating coarse-grained (CG) systems that can be simulated using the OpenMM interface. It allows for the creation of both default and custom CG systems and easy modification of their parameters.
- Parameters:
structure_path (str) – The path to the input PDB or CIF file. Required.
model (str, optional, default='hps_urry') – The hydropathy scale to be used. Currently, two models are supported: ‘hps_kr’ and ‘hps_urry’.
- structure
An object that holds the information of the parsed PDB or CIF file.
- Type:
openmm.app.pdbfile.PDBFile
oropenmm.app.pdbxfile.PDBxFile
- topology
The topology of the model, as generated by OpenMM.
- positions
unit.quantity.Quantity
The atomic positions of the model.
- particles_massfloat or list
The mass of each particle. If a float is provided, all particles will have the same mass. If a list is provided, per-particle masses will be assigned.
- particles_chargelist
The charge of each particle.
- Type:
openmm.app.topology.Topology
- positions
- rf_sigma
The sigma parameter used in the pairwise force object. This is the vdw radius of beads.
- Type:
float
- atoms
A list of the current atoms in the model. The items are
openmm.app.topology.atoms
initialised classes.- Type:
list
- n_atoms
The total number of atoms in the model.
- Type:
int
- bonds
A dictionary that uses bonds (2-tuple of
openmm.app.topology.bonds
objects) present in the model as keys and their forcefield properties as values.- Type:
collections.OrderedDict
- bonds_indexes
A list containing the zero-based indexes of the atoms defining the bonds in the model.
- Type:
list
- n_bonds
The total number of bonds in the model.
- Type:
int
- bonded_exclusions_index
The exclusion rule for nonbonded force. =1 for hps_kr and hps_urry, =3 for hps_ss
- Type:
int
- harmonicBondForce
The
openmm.HarmonicBondForce
object that implements a harmonic bond potential between pairs of particles, that depends quadratically on their distance.- Type:
openmm.HarmonicBondForce
- n_angles
The total number of angles in the model.
- Type:
int
- gaussianAngleForce
The
openmm.CustomAngleForce
object that implements a Gaussian angle bond potential between pairs of three particles.- Type:
openmm.CustomAngleForce
- n_torsions
Total number of torsion angles in the model.
- Type:
int
- gaussianTorsionForce
Stores the OpenMM
CustomTorsionForce
initialised-class. Implements a Gaussian torsion angle bond potential between pairs of four particles.- Type:
openmm.CustomTorsionForce
- yukawaForce
Stores the OpenMM
CustomNonbondedForce
initialized-class. Implements the Debye-Huckle potential.- Type:
openmm.CustomNonbondedForce
- ashbaugh_HatchForce
Stores the OpenMM
CustomNonbondedForce
initialized-class. Implements the pairwise short-range potential.- Type:
openmm.CustomNonbondedForce
- forceGroups
A dict that uses force names as keys and their corresponding force as values.
- Type:
collections.OrderedDict
- system
Stores the OpenMM System initialised class. It stores all the forcefield information for the hps model.
- Type:
openmm.System
- loadForcefieldFromFile()
Loads forcefield parameters from a force field file written with the
dumpForceFieldData()
method.
- static _setParameters(term, parameters)[source]
General function to set up or change force field parameters. protected method, can be called only inside class system.
- Parameters:
term (dict) – Dictionary object containing the set of degrees of freedom (DOF) to set up attributes to (e.g.
bonds
attribute)parameters (integer or float or list) – Value(s) for the specific forcefield parameters. If integer or float, sets up the same value for all the DOF in terms. If a list is given, sets a unique parameter for each DOF.
- Return type:
None
- addAshbaughHatchForces(use_pbc: bool) None [source]
HPS family used Ashbaugh-Hatch Potential instead of Wang-Frenkel Creates a nonbonded force term for pairwise interaction (customize LJ 12-6 potential).
Creates an
openmm.CustomNonbondedForce()
object with the parameters sigma and epsilon given to this method. The custom non-bonded force is initialized with the formula: (note: hps here is \(\lambda_{ij}^{0}\) in the paper)Unlike
BondForce
class, where we specify index for atoms pair to add bond, it means that number of bondForces may differ from number of particle.NonBondedForce
is added to all particles, hence we don’t need to pass theatom index
.\[ \begin{align}\begin{aligned}\Phi_{i,j}^{vdw}(r) = step(2^{1/6}\sigma_{ij}-r) \times \left( 4\epsilon\left[\left(\frac{\sigma_{ij}}{r}\right)^{12}- \left(\frac{\sigma_{ij}}{r}\right)^{6}\right]+(1-\lambda_{ij})\epsilon\right)\\+ \left[1-step(2^{1/6}\sigma_{ij}-r)\right]\times\left[(\lambda_{ij})\times 4\epsilon \left[\left(\frac{\sigma_{ij}}{r}\right)^{12}-\left(\frac{\sigma_{ij}}{r}\right)^6\right]\right]\end{aligned}\end{align} \]Here, \(\sigma= \frac{(\sigma_1+\sigma_2)}{2}; \lambda_{ij}^{0}=\frac{(\lambda_i+\lambda_j)}{2}; \epsilon = 0.8368 kj/mol\)
The force object is stored at the
ashbaugh_HatchForce
attribute.- epsilonfloat
Value of the epsilon constant in the energy function.
- sigmafloat or list
Value of the sigma constant (in nm) in the energy function. If float the same sigma value is used for every particle. If list a unique parameter is given for each particle.
- cutofffloat
The cutoff distance (in nm) being used for the non-bonded interactions.
- Parameters:
use_pbc (bool. Whether use PBC, cutoff periodic boundary condition) –
- Return type:
None
- addGaussianAngleForces() None [source]
Add Gaussian functional form of angle. Note that in openMM log is neutral logarithm.
Angle potential take Gaussian functional form in hps-ss model.
\[U_{angle}(\theta) = \frac{-1}{\gamma} \ln{[e^{-\gamma[k_\alpha( \theta-\theta_\alpha)^2+\epsilon_\alpha]} +e^{-\gamma k_\beta(\theta-\theta_\beta)^2}]}\]Angle potential is taken from reference:
- addGaussianTorsionForces() None [source]
Torsion potential in hps-ss model takes the form:
\[U_{torsion}(\theta) = -\ln\left[ U_{torsion, \alpha}(\theta, \epsilon_d) + U_{torsion, \beta}(\theta, \epsilon_d)\right]\]where,
\[ \begin{align}\begin{aligned}U_{torsion, \alpha}(\theta, \epsilon_d) &= e^{-k_{\alpha, 1}(\theta-\theta_{\alpha,1})^2-\epsilon_d} + e^{-k_{\alpha, 2}(\theta-\theta_{\alpha,2})^4 + e_0} + e^{-k_{\alpha, 2}(\theta-\theta_{\alpha,2}+2\pi)^4 + e_0} \\\ U_{torsion, \beta}(\theta, \epsilon_d) &= e^{-k_{\beta,1}(\theta-\theta_{\beta,1})^2+e_1+\epsilon_d} + e^{-k_{\beta,1}(\theta-\theta_{\beta,1}-2\pi)^2+e_1+\epsilon_d} \\\ &+ e^{-k_{\beta,2}(\theta-\theta_{\beta,2})^4+e_2} + e^{-k_{\beta,2}(\theta-\theta_{\beta,2}-2\pi)^4+e_2}\end{aligned}\end{align} \]
- addHarmonicBondForces() None [source]
Creates a harmonic bonded force term for each bond in the main class using their defined forcefield parameters.
Creates an
openmm.HarmonicBondForce()
object with the bonds and parameters set up in the “bonds” dictionary attribute. The force object is stored at theharmonicBondForce
attribute.openMM uses harmonic bond that has energy term of form:
\[E= \frac{1}{2}k(r-r_0)^2\]The force parameters must be contained in self.bonds as follows:
- self.bonds is a dictionary:
The keys are 2-tuples for two atom items in
self.topology.atoms
attribute.- The values are a 2-tuple of parameters in the following order:
first -> bond0 (quantity)
second -> k (float) (measured in unit of kj/mol/nm^2)
- Return type:
None
- addParticles() None [source]
Add particles to the system OpenMM class instance.
Add a particle to the system for each atom in it. The mass of each particle is set up with the values in the
particles_mass
attribute.
- addSystemForces() None [source]
Add forces to the system OpenMM class instance. It also save names for the added forces to include them in the reporter class.
Adds generated forces to the system, also adding a force group to the
forceGroups
attribute dictionary.
- addYukawaForces(use_pbc: bool) None [source]
Creates a nonbonded force term for electrostatic interaction DH potential.
Creates an
openmm.CustomNonbondedForce()
object with the parameters sigma and epsilon given to this method. The custom non-bonded force is initialized with the formula:\[energy = f \times \frac{q_1q_2}{\epsilon_r \times r}\times e^{(-r/lD)}\]where \(f=\frac{1}{4\pi\epsilon_0}=138.935458\) is the factor for short to convert dimensionless in calculation to \(kj.nm/(mol\times e^2)\) unit.
\(\epsilon_r=80\): Dielectric constant of water at 100mM mono-valent ion
The force object is stored at the
yukawaForce
attribute.- Parameters:
use_pbc ((bool) whether use PBC, cutoff periodic boundary condition) –
- Return type:
None
- add_Wang_Frenkel_Forces(use_pbc: bool)[source]
MPIPI model. using TabulatedFunction for pair interaction. More information about TabulatedFUnction can be found here: http://docs.openmm.org/7.2.0/api-c++/generated/OpenMM.Discrete2DFunction.html
- checkBondDistances(threshold: float = 0.5) None [source]
Searches for large bond distances for the atom pairs defined in the ‘bonds’ attribute. It raises an error when large bonds are found.
- Parameters:
threshold ((float, default=0.5 nm)) – Threshold to check for large bond distances.
- Return type:
None
- checkLargeForces(minimize: bool = False, threshold: float = 10) None [source]
Prints the hps system energies of the input configuration of the system. It optionally checks for large forces acting upon all particles in the hps system and iteratively minimizes the system configuration until no forces larger than a threshold are found.
- Parameters:
threshold ((float, default=10)) – Threshold to check for large forces.
minimize ((bool, default= False)) – Whether to iteratively minimize the system until all forces are lower or equal to the threshold value.
- Return type:
None
- createSystemObject(check_bond_distances: bool = True, minimize: bool = False, check_large_forces: bool = True, force_threshold: float = 10.0, bond_threshold: float = 0.5) None [source]
Creates OpenMM system object adding particles, masses and forces. It also groups the added forces into Force-Groups for the hpsReporter class.
Creates an
openmm.System()
object using the force field parameters given to the ‘system’ class. It adds particles, forces and creates a force group for each force object. Optionally the method can check for large bond distances (default) and minimize the atomic positions if large forces are found in any atom (default False).- Parameters:
minimize (boolean (False)) – Whether to minimize the system if large forces are found.
check_bond_distances (boolean (True)) – Whether to check for large bond distances.
check_large_forces (boolean (False)) – Whether to print force summary of force groups
force_threshold (float (10.0)) – Threshold to check for large forces.
bond_threshold (float (0.5)) – Threshold to check for large bond distances.
- Return type:
None
- dumpForceFieldData(output_file: str) None [source]
Writes to a file the parameters of the forcefield.
Writes a file containing the current forcefield parameters in the CG system.
- Parameters:
output_file (string [requires]) – name of the output file.
- Return type:
None
- dumpStructure(output_file: str) None [source]
Writes a structure file of the system in its current state.
Writes a PDB file containing the currently defined CG system atoms and its positions.
- Parameters:
output_file (string) – name of the PDB output file.
- Return type:
None
- dumpTopology(output_file: str) None [source]
Writes a topology file of the system in PSF format, this is used for visualization and post-analysis.
Writes a file containing the current topology in the hpsOpenMM system. This file contains topology of system, used in visualization and analysis.
Here, we used
parmed
to loadopenMM topology
andopenMM system
to createStructure
object inparmed
. Because parmed doesn’t automatically recognizecharge
,mass
of atoms by their name. We need to setcharge
,mass
back to residues properties.- Parameters:
output_file (string [requires]) – name of the output PSF file.
- Return type:
None
- getAngles()[source]
Reads bonds from topology, adds them to the main class and sorts them into a dictionary to store their forcefield properties.
Angles are built from bond (bondedTo instance). This function modify
self.angles
,self.angles_indexes
andself.n_angles
- getAtoms()[source]
Reads atoms from topology, adds them to the main class and sorts them into a dictionary to store their forcefield properties.
After getCAlphaOnly, C-alpha atoms are stored on
self.topology only
. We need to add them to atoms attribute and system also. Addsatoms
in theOpenMM topology
instance to thehpsOpenMM system
class.- Return type:
None
- getBonds(except_chains=None)[source]
Reads bonds from topology, adds them to the main class and sorts them into a dictionary to store their forcefield properties.
Adds
bonds
in theOpenMM topology
instance to thehpsOpenMM system
class.- Parameters:
except_chains (String [optional]) –
- Return type:
None
- getCAlphaOnly() None [source]
Filter in only alpha carbon atoms from the input structure and updates the topology object to add new bonds between them. Used specially for creating alpha-carbon (CA) coarse-grained models.
Keeps in the
hps system
only the alpha carbon atoms from theOpenMM topology
.- Return type:
None
- getTorsions()[source]
Reads bonds from topology, adds them to the main class and sorts them into a dictionary to store their forcefield properties.
Torsion Angles are built from angles and bondedTo instance. This function modify
self.torsions
,self.torsions_indexes
andself.n_torsions
- setBondForceConstants()[source]
Change the forcefield parameters for bonded terms.
Set the harmonic bond constant force parameters. The input can be a float, to set the same parameter for all force interactions, or a list, to define a unique parameter for each force interaction.
- Return type:
None
- setCAChargePerResidueType()[source]
Sets the charge of the alpha carbon atoms to characteristic charge of their corresponding amino acid residue.
- Return type:
None
- setCAHPSPerResidueType()[source]
Sets alpha carbon atoms to their residue hydropathy scale. Used specially for modifying alpha-carbon (CA) coarse-grained models.
Sets the HPS model of the alpha carbon atoms using corresponding scale.
- Return type:
None
- setCAMassPerResidueType()[source]
Sets alpha carbon atoms to their average residue mass. Used specially for modifying alpha-carbon (CA) coarse-grained models.
Sets the masses of the alpha carbon atoms to the average mass of its amino acid residue.
- Return type:
None
- setCARadiusPerResidueType()[source]
Sets alpha carbon atoms to their average residue mass. Used specially for modifying alpha-carbon (CA) coarse-grained models.
Sets the excluded volume radii of the alpha carbon atoms to characteristic radii of their corresponding amino acid residue.
- Return type:
None
- setParticlesCharge(particles_charge)[source]
Set the charge of the particles in the system. The input can be a float, to set the same charge for all particles, or a list, to define a unique charge for each particle.
- Parameters:
particles_charge (float or list) – Charge values to add for the particles in the hpsOpenMM system class.
- Return type:
None
- setParticlesHPS(particles_hps)[source]
Set the hydropathy scale of the particles in the system. The input can be a float, to set the same hydropathy for all particles, or a list, to define a unique hydropathy for each particle.
- Parameters:
particles_hps (float or list) – HPS scale values to add for the particles in the hpsOpenMM system class.
- Return type:
None
- setParticlesMass(particles_mass)[source]
Change the mass parameter for each atom in the system.
Set the masses of the particles in the system. The input can be a float, to set the same mass for all particles, or a list, to define a unique mass for each particle.
- Parameters:
particles_mass (float or list) – Mass(es) values to add for the particles in the hpsOpenMM system class.
- Return type:
None
- setParticlesRadii(particles_radii)[source]
Change the excluded volume radius parameter for each atom in the system.
Set the radii of the particles in the system. The input can be a float, to set the same radius for all particles, or a list, to define a unique radius for each particle.
- Parameters:
particles_radii (float or list) – Radii values to add for the particles in the hpsOpenMM system class.
- Return type:
None
Module contents
core package of the hpsOpenMM package that contains the main hpsOpenMM classes.
The hpsOpenMM.core package contains the three hpsOpenMM main classes:
geometry
models
system
The first class, geometry, contains methods to calculate the geometrical parameters from the input structures.
The second class, models, allows to easily set up predefined CG models. The final class, system, is the main class that holds all the methods to define, modify and create CG to be simulated with OpenMM.