STORMM Source Documentation
|
A struct to hold information relating to an Amber topology. This struct's member functions are limited to getters for its private data. Because the primary means of constructing a topology will be complex, i.e. meticulous parsing of a file or expansion of an existing topology based on some new information, the constructors will need to be buried within wrapper functions that perform such procedures. This struct contains all information contained within an Amber prmtop-format file, nothng more, as described on: More...
#include <atomgraph.h>
Public Member Functions | |
AtomGraph () | |
The blank constructor makes a blank AtomGraph, which is used by the general- purpose file-based constructor to delegate initialization. | |
~AtomGraph ()=default | |
The default destructor is adequate. | |
AtomGraph (const AtomGraph &original) | |
The copy constructor takes a lot of tedium upon itself to prevent more complexity in working with AtomGraph objects downstream. | |
AtomGraph & | operator= (const AtomGraph &other) |
Copy constructor. | |
AtomGraph (AtomGraph &&original) | |
The move constructor makes prodigious use of std::move for each string and Hybrid member variable. | |
AtomGraph & | operator= (AtomGraph &&other) |
Move constructor. | |
void | buildFromPrmtop (const std::string &file_name, const ExceptionResponse policy=ExceptionResponse::WARN, double coulomb_constant_in=amber_ancient_bioq, double default_elec14_screening=amber_default_elec14_screen, double default_vdw14_screening=amber_default_vdw14_screen, double charge_rounding_tol=default_charge_rounding_tol, double charge_discretization=default_charge_precision_inc) |
Build an AtomGraph form a file. This is called by the general-purpose constructor or also by the developer after instantiating an empty object. | |
std::string | getTitle () const |
Get the title of the topology (may be blank). | |
std::string | getFileName () const |
Get the file name where a topology originated (may be blank, indicating that the topology was produced by some other means). | |
int | getAtomCount () const |
Get the number of atoms, using the topology's dedicated private variable rather than a list of input dimensions that was mostly stored for developers most familiar with Amber. | |
int | getResidueCount () const |
Get the number of residues. | |
int | getMoleculeCount () const |
Get the number of separate molecules in the system. | |
MoleculeKind | getMoleculeKind (int mol_index) const |
Get the kind of a specific molecule wihtin the topology. | |
int | getOrganicMoleculeCount () const |
Get the number of organic molecules in the system, based on the criteria of greater than or equal to eight total real atoms and at least one carbon atom. | |
int | getOrganicAtomCount () const |
Get the number of organic atoms in the system, based on the criteria of greater than or equal to eight total real atoms and at least one carbon atom. | |
int | getWaterResidueSize () const |
Get the number of particles (which may include virtual sites) in a single water molecule within the system. | |
int | getWaterResidueCount () const |
Get the number of water residues in the system. | |
int | getLargestResidueSize () const |
Get the number of separate molecules in the system. | |
int | getLastSoluteResidue () const |
Get the final solute residue, indexed according to C/C++ array indexing conventions. | |
int | getLastSoluteAtom () const |
Get the final solute atom, indexed according to C/C++ array indexing conventions. | |
int | getFirstSolventMolecule () const |
Get the first solvent molecule, indexed according to C/C++ array indexing conventions. | |
int | getFirstSolventAtom () const |
Get the first solvent atom, indexed according to C/C++ array indexing conventions. | |
int | getLargestMoleculeSize () const |
Get the largest molecule's size. | |
double | getTotalMass () const |
Get the total mass of all atoms in the topology (this is a computation, not a stored value). | |
int | getDegreesOfFreedom () const |
Get the number of degrees of freedom, without consideration to geometric constraints. | |
int | getConstrainedDegreesOfFreedom () const |
Get the number of degrees of freedom after geometric constraints have been applied. | |
std::vector< int > | getParticlesPerMolecule () const |
Get the sizes of all molecules in the system. Rather redundant with the molecule limits being accessible, but convenient. | |
int | getParticlesPerMolecule (int index) const |
std::string | getFullAtomName (int index) const |
Get the full name (including the residue and structural residue number) of an atom, as a string. | |
std::vector< std::vector< char4 > > | getAtomTypeNameTable () const |
Get a table with the names of all unique atom types, arranged according to their Lennard-Jones indices as they appear in the topology (in ascending order). | |
int | getUreyBradleyTermCount () const |
Get the number of Urey-Bradley terms in the system. | |
int | getCharmmImprTermCount () const |
Get the number of CHARMM improper terms. | |
int | getCmapTermCount () const |
Get the number of CMAP terms (the number of CMAP parameter surfaces is different) | |
int | getUreyBradleyParameterCount () const |
Return the number of unique Urey-Bradley parameter sets. | |
int | getCharmmImprParameterCount () const |
Return the number of unique CHARMM improper parameter sets. | |
int | getCmapSurfaceCount () const |
Return the number of unique CMAP surfaces. | |
int | getCmapDimension (int index) const |
Return the dimension of a particular CMAP surface (all CMAPs are assumed to be square) | |
template<typename T> | |
UreyBradleyTerm< T > | getUreyBradleyTerm (int index) const |
Collect the atom indices and actual parameters to implement one of the system's Urey-Bradley terms. The special-purpose, unguarded struct that holds the result also provides the original parameter index, if of interest. This is for computing the potential energy and forces in the most straightforward manner possible, albeit slow. | |
template<typename T> | |
CharmmImprTerm< T > | getCharmmImprTerm (int index) const |
Collect the atom indices and actual parameters to implement one of the system's CHARMM improper terms. The special-pupose, unguarded struct that holds the result also provides the original parameter index, if of interest. This is for computing the potential energy and forces in the most straightforward manner possible, albeit slow. | |
template<typename T> | |
CmapTerm< T > | getCmapTerm (int index) const |
Collect the atom indices and a pointer to the surface to implement one of the system's CMAP terms. The special-purpose, unguarded struct that holds the result also provides the original parameter index, if of interest. | |
int | getBondTermCount () const |
Get the number of bond stretching terms in the system. | |
int | getAngleTermCount () const |
Get the number of angle bending terms in the system. | |
int | getDihedralTermCount () const |
Get the number of dihedral torsion terms in the system. | |
int | getBondParameterCount () const |
Return the number of unique bond stretching parameter sets. | |
int | getAngleParameterCount () const |
Return the number of unique bond angle bending parameter sets. | |
int | getDihedralParameterCount () const |
Return the number of unique dihedral parameter sets. | |
template<typename T> | |
BondTerm< T > | getBondTerm (int index) const |
Collect the atom indices and actual parameters for implement one of the system's harmonic bond stretching terms. | |
template<typename T> | |
AngleTerm< T > | getAngleTerm (int index) const |
Collect the atom indices and actual parameters for implement one of the system's harmonic bond angle bending terms. | |
template<typename T> | |
DihedralTerm< T > | getDihedralTerm (int index) const |
Collect the atom indices and actual parameters for implement one of the system's cosine-based dihedral torsion terms. | |
int | getVirtualSiteCount () const |
Get the number of virtual sites in the system. | |
int | getVirtualSiteParameterSetCount () const |
Get the number of unique virtual site frames in the system. For example, a system of 1024 TIP-4P waters would have one unique frame, applied to 1024 different groups of atoms. | |
int | getVirtualSiteIndex (int atom_index, ExceptionResponse policy=ExceptionResponse::DIE) const |
Get the virtual site index of an atom index. If the atom index is not a virtual site, this function will throw an error. | |
VirtualSiteKind | getVirtualSiteFrameType (int index) const |
Get the frame type of a particular virtual site. | |
int | getVirtualSiteFrameAtom (int index, int nfrm) const |
Get the general topology index of one of a virtual site's frame atoms. | |
template<typename T> | |
T | getVirtualSiteFrameDimension (int index, int ndim) |
Get the dimensions (could be length in A, could be scaling factor for a distance or cross product) of one virtual site's frame. | |
int | getChargeTypeCount () const |
Get the number of charge types in the system. | |
int | getLJTypeCount () const |
Get the number of atom types in the system. | |
int | getTotalExclusions () const |
Get the number of exclusions in the entire topology. | |
UnitCellType | getUnitCellType () const |
Get the simulation cell type (isolated boundary conditions, rectilinear, triclinic) | |
ImplicitSolventModel | getImplicitSolventModel () const |
Get the implicit solvent model type, i.e. some flavor of GB or NONE. | |
double | getDielectricConstant () const |
Get the dielectric constant (part of the implicit solvent setup) | |
double | getSaltConcentration () const |
Get the salt concentration (part of the implicit solvent setup) | |
double | getCoulombConstant () const |
Get the fundamental coulomb constant that this system uses in electrostatic calculations. | |
double | getElectrostatic14Screening () const |
Get the general 1:4 screening factor on electrostatic interactions. | |
double | getVanDerWaals14Screening () const |
Get the general 1:4 screening factor on van-der Waals interactions. | |
std::string | getPBRadiiSet () const |
Get the PB Radii set. | |
bool | usesShake () const |
Get the status of SHAKE and RATTLE constraints in the system. | |
bool | usesSettle () const |
Get the status of SETTLE rigid water constraints in the system. | |
std::string | getWaterResidueName () const |
Get the water residue name. | |
int | getRigidWaterCount () const |
Get the number of rigid waters in the system. This will check through the SETTLE groups and confirm the number that are, by their chemistry, water. | |
int | getRigidWaterAtomCount () const |
Get the number of atoms in a rigid water molecule of the system, if the system has rigid waters. This function will return zero if the system has no rigid waters. | |
int | getSettleGroupCount () const |
Get the number of SETTLE groups in the system. | |
int | getBondConstraintCount () const |
Get the number of bond constraints. | |
int | getConstraintGroupCount () const |
Get the total number of constrained groups. | |
std::vector< int > | getConstraintGroupAtoms (int index) const |
Get the atoms of a single constraint group. | |
int | getConstraintGroupTotalSize () const |
Get the total size of the constrained group atoms list. | |
int | getNonrigidParticleCount () const |
Get the number of non-rigid particles in the system. | |
std::vector< int > | getAtomExclusions (int index) const |
Get the exclusion list for a particular atom. This will concatenate 1:1, 1:2, 1:3, and 1:4 exclusions, without regard to which is which. To get a specific type of exclusion, use one of the getNonbondedXXExclusions() functions. All functions of this sort will return a complete list such that if exclusions for all atoms are considered in sequence, all exclusions will be double-counted. | |
std::vector< int > | getNonbonded11Exclusions (int index) const |
Get the 1:1 exclusions for a particular atom. This list will only be populated if the atom is a virtual site or has virtual sites that claim it as their parent. | |
std::vector< int > | getNonbonded12Exclusions (int index) const |
Get all 1:2 exclusions for a particular atom. | |
std::vector< int > | getNonbonded13Exclusions (int index) const |
Get all 1:3 exclusions for a particular atom. | |
std::vector< int > | getNonbonded14Exclusions (int index) const |
Get all 1:4 exclusions for a particular atom. | |
template<typename T> | |
T | getChargeParameter (int index) const |
Get the charge parameter for a specific atom. The quantity is derived from tables at the specified level of precision (float or double), and equals the value obtained from getPartialCharge<T>(index), although that function reaches into a separate array with a great deal of redundant information. | |
std::vector< char4 > | matchOverflowAtomName (const std::string &query) const |
Match a long atom name and return its key in the standard list of char4 atom names. | |
std::vector< char4 > | matchOverflowAtomType (const std::string &query) const |
Match a long atom type and return its key in the standard list of char4 atom types. | |
std::vector< char4 > | matchOverflowResidueName (const std::string &query) const |
Match a long residue name and return its key among the standard char4 residue names. | |
ValenceKit< double > | getDoublePrecisionValenceKit (HybridTargetLevel tier=HybridTargetLevel::HOST) const |
Get an abstract for calculating valence terms in double precision. | |
ValenceKit< float > | getSinglePrecisionValenceKit (HybridTargetLevel tier=HybridTargetLevel::HOST) const |
Get a single-precision abstract for calculating valence terms. | |
NonbondedKit< double > | getDoublePrecisionNonbondedKit (HybridTargetLevel tier=HybridTargetLevel::HOST) const |
Get a double-precision abstract for calculating non-bonded energy and forces. | |
NonbondedKit< float > | getSinglePrecisionNonbondedKit (HybridTargetLevel tier=HybridTargetLevel::HOST) const |
Get a single-precision abstract for calculating non-bonded energy and forces. | |
ImplicitSolventKit< double > | getDoublePrecisionImplicitSolventKit (HybridTargetLevel tier=HybridTargetLevel::HOST) const |
Get a double-precision abstract for calculating implicit solvent energy and forces. | |
ImplicitSolventKit< float > | getSinglePrecisionImplicitSolventKit (HybridTargetLevel tier=HybridTargetLevel::HOST) const |
Get a single-precision abstract for calculating implicit solvent energy and forces. | |
ChemicalDetailsKit | getChemicalDetailsKit (HybridTargetLevel tier=HybridTargetLevel::HOST) const |
Get an abstract of the chemical details of the system, such as atom name and type or residue names and limits. | |
VirtualSiteKit< double > | getDoublePrecisionVirtualSiteKit (HybridTargetLevel tier=HybridTargetLevel::HOST) const |
Get a double-precision abstract for placing virtual sites and transmitting their forces to frame atoms with mass. | |
VirtualSiteKit< float > | getSinglePrecisionVirtualSiteKit (HybridTargetLevel tier=HybridTargetLevel::HOST) const |
Get a single-precision abstract for placing virtual sites and transmitting their forces to frame atoms with mass. | |
ConstraintKit< double > | getDoublePrecisionConstraintKit (HybridTargetLevel tier=HybridTargetLevel::HOST) const |
Get an abstract for managing constraints based on double-precision parameters. | |
ConstraintKit< float > | getSinglePrecisionConstraintKit (HybridTargetLevel tier=HybridTargetLevel::HOST) const |
Get an abstract for managing constraints based on single-precision parameters. While the overall implementation for constraints cannot be implemented entirely in single precision, having the underlying parameters in this form is not, in principle, a drvier of energy drift. It will simply set the system into a slightly different geometry, which is probably still closer to the double-precision result than many cycles of SHAKE or RATTLE would get it anyway. | |
void | setSource (const std::string &new_source) |
Set the source file for the topology. Primary use is to impart a destination file for a topology to be written to disk, but can also be used to "disguise" one topology from another even if both are from the same original file. | |
void | setImplicitSolventModel (ImplicitSolventModel igb_in, double dielectric_in=80.0, double saltcon_in=0.0, AtomicRadiusSet radii_set=AtomicRadiusSet::NONE, ExceptionResponse policy=ExceptionResponse::WARN) |
Set the implicit solvent model. This will leverage a lot of hard-coded constants to fill out some allocated but otherwise blank arrays and impart one particular Generalized Born or other implicit solvent method. | |
void | printToFile (const std::string &output_file={}, TopologyKind output_style=TopologyKind::AMBER, PrintSituation expectation=PrintSituation::OPEN_NEW, ExceptionResponse pr_policy=ExceptionResponse::DIE) const |
Print out a topology in the specified format. | |
void | printAmberFormat (const std::string &output_file={}, PrintSituation expectation=PrintSituation::OPEN_NEW, ExceptionResponse pr_policy=ExceptionResponse::DIE) const |
Print the contents of the topology in the AMBER file format. Descriptions of input parameters follow from printToFile(), above. | |
AtomGraph (const std::string &file_name, ExceptionResponse policy=ExceptionResponse::WARN, TopologyKind engine_format=TopologyKind::AMBER) | |
The general-purpose constructor for file-based topology creation. | |
AtomGraph (const std::string &file_name, ExceptionResponse policy, TopologyKind engine_format, double coulomb_constant_in, double default_elec14_screening, double default_vdw14_screening, double charge_rounding_tol, double charge_discretization, ApplyConstraints use_bond_constraints_in=ApplyConstraints::NO, ApplyConstraints use_settle_in=ApplyConstraints::NO) | |
AtomGraph (const std::vector< AtomGraph * > &agv, const std::vector< int > &counts, MoleculeOrdering arrangement=MoleculeOrdering::REORDER_ALL, ExceptionResponse policy=ExceptionResponse::DIE) | |
AtomGraph (const AtomGraph &ag_a, int n_a, const AtomGraph &ag_b, int n_b, MoleculeOrdering arrangement=MoleculeOrdering::REORDER_ALL, ExceptionResponse policy=ExceptionResponse::DIE) | |
AtomGraph (const AtomGraph &ag_a, const AtomGraph &ag_b, int n_b, MoleculeOrdering arrangement=MoleculeOrdering::REORDER_ALL, ExceptionResponse policy=ExceptionResponse::DIE) | |
AtomGraph (const AtomGraph &ag_a, const AtomGraph &ag_b, MoleculeOrdering arrangement=MoleculeOrdering::REORDER_ALL, ExceptionResponse policy=ExceptionResponse::DIE) | |
AtomGraph (const AtomGraph &original, const std::vector< int > &atom_subset, ExceptionResponse policy=ExceptionResponse::DIE, double charge_rounding_tol=default_charge_rounding_tol, double charge_discretization=default_charge_precision_inc, ApplyConstraints use_bond_constraints_in=ApplyConstraints::NO, ApplyConstraints use_settle_in=ApplyConstraints::NO) | |
AtomGraph (const AtomGraph &original, const std::vector< bool > &mask, ExceptionResponse policy=ExceptionResponse::DIE, double charge_rounding_tol=default_charge_rounding_tol, double charge_discretization=default_charge_precision_inc, ApplyConstraints use_bond_constraints_in=ApplyConstraints::NO, ApplyConstraints use_settle_in=ApplyConstraints::NO) | |
int | getDescriptor (TopologyDescriptor choice) const |
Get a descriptor from within the array of topology descriptors. If this topology were read from an Amber-format file, all descriptors are taken from the preamble in its first ~15 lines. | |
int | getDescriptor (SanderDescriptor choice) const |
const Hybrid< int > & | getResidueLimits () const |
Get residue limits for the beginning and ends of up to N residues. | |
int2 | getResidueLimits (int index) const |
std::vector< int > | getResidueIndex () const |
Get the index (topology index, not natural / structure-informed residue number) of a particular residue. | |
int | getResidueIndex (int atom_index) const |
std::vector< int > | getStructuralAtomNumber () const |
Get the (arbitrarily defined) atom numbers for atoms in the system. For example, a PDB file could ahve its own atom numbering scheme. This provides a way to transcribe that into the AtomGraph object. | |
std::vector< int > | getStructuralAtomNumber (int low_index, int high_index) const |
int | getStructuralAtomNumber (int index) const |
std::vector< int > | getResidueNumber () const |
Get the residue numbers for atoms in the system. | |
std::vector< int > | getResidueNumber (int low_index, int high_index) const |
int | getResidueNumber (int index) const |
std::vector< int > | getMoleculeLimits () const |
Get the molecule limits, indices into molecule_contents, for every molecule in the system. This allows that molecules not be contiguous in the topology atom indexing, but will be recoverable by molecule_contents. | |
int2 | getMoleculeLimits (int index) const |
std::vector< int > | getAtomicNumber () const |
Get the atomic numbers for atoms in the system. | |
std::vector< int > | getAtomicNumber (int low_index, int high_index) const |
int | getAtomicNumber (int index) const |
std::vector< bool > | getAtomMobility () const |
Get the mobile atom mask for the system, as a boolean vector translated from the bitmask stored internally. | |
std::vector< bool > | getAtomMobility (int low_index, int high_index) const |
bool | getAtomMobility (int index) const |
std::vector< uint > | getAtomMobilityMask () const |
Get the raw mobility bitmask over a particular range. | |
std::vector< uint > | getAtomMobilityMask (int low_index, int high_index) const |
void | modifyAtomMobility (MobilitySetting movement) |
Change an atom's mobility in the topology. | |
void | modifyAtomMobility (int low_index, int high_index, MobilitySetting movement) |
void | modifyAtomMobility (int index, MobilitySetting movement) |
void | modifyAtomMobility (const std::vector< int > &mask, MobilitySetting movement) |
std::vector< int > | getMoleculeMembership () const |
Get molecule membership for atoms within the system. | |
std::vector< int > | getMoleculeMembership (int low_index, int high_index) const |
int | getMoleculeMembership (int index) const |
std::vector< int > | getMoleculeContents () const |
Get the contents of a one or more molecules in the system. | |
std::vector< int > | getMoleculeContents (int index) const |
template<typename T> | |
std::vector< T > | getPartialCharge () const |
Get the atomic partial charges of atoms in the system. | |
template<typename T> | |
std::vector< T > | getPartialCharge (int low_index, int high_index) const |
template<typename T> | |
T | getPartialCharge (int index) const |
template<typename T> | |
std::vector< T > | getAtomicMass (MassForm rep=MassForm::ORDINARY) const |
Get the masses of atoms in the system. The developer may also choose to take the pre-computed inverse masses (this can be useful when dealing with extra points, as inverses of zero-mass particles must also be set to zero by a conditional statement). | |
template<typename T> | |
std::vector< T > | getAtomicMass (int low_index, int high_index, MassForm rep=MassForm::ORDINARY) const |
template<typename T> | |
T | getAtomicMass (int index, MassForm rep=MassForm::ORDINARY) const |
const Hybrid< char4 > & | getAtomName () const |
Get the names of atoms in the system. | |
std::vector< char4 > | getAtomName (int low_index, int high_index) const |
char4 | getAtomName (int index) const |
const Hybrid< char4 > & | getAtomType () const |
Get the types of atoms in the system. | |
std::vector< char4 > | getAtomType (int low_index, int high_index) const |
char4 | getAtomType (int index) const |
const Hybrid< char4 > & | getResidueName () const |
Get the names of residues in the system. | |
char4 | getResidueName (int index) const |
std::vector< int > | findVirtualSites () const |
Get a list of the general topological indices of one or more virtual sites. | |
std::vector< int2 > | findVirtualSites (int low_index, int high_index) const |
int | findVirtualSites (int index) const |
std::vector< int > | getChargeIndex () const |
Get the charge type index of atoms in the system. | |
std::vector< int > | getChargeIndex (int low_index, int high_index) const |
int | getChargeIndex (int index) const |
std::vector< int > | getLennardJonesIndex () const |
Get the Lennard-Jones type index of atoms in the system. | |
std::vector< int > | getLennardJonesIndex (int low_index, int high_index) const |
int | getLennardJonesIndex (int index) const |
template<typename T> | |
T | getLennardJonesSigma (int index_a, int index_b) const |
Get the Lennard-Jones interaction sigma. The quantity is derived from tables at the specified level of precision (float or double), but computed in double precision based on that information before it is returned in the template data type. | |
template<typename T> | |
T | getLennardJonesSigma (int index_a) const |
template<typename T> | |
std::vector< T > | getLennardJonesSigma () const |
template<typename T> | |
T | getLennardJonesEpsilon (int index_a, int index_b) const |
Get the Lennard-Jones interaction epsilon value. The quantity is derived from tables at the specified level of precision (float or double), but computed in double precision based on that information before it is returned in the template data type. | |
template<typename T> | |
T | getLennardJonesEpsilon (int index_a) const |
template<typename T> | |
std::vector< T > | getLennardJonesEpsilon () const |
template<typename T> | |
std::vector< T > | getAtomPBRadius () const |
Get the Poisson-Boltzmann radius of a particular atom. | |
template<typename T> | |
std::vector< T > | getAtomPBRadius (int low_index, int high_index) const |
template<typename T> | |
T | getAtomPBRadius (int index) const |
template<typename T> | |
std::vector< T > | getGBScreeningFactor () const |
Get GB screening factors for one or more atoms. | |
template<typename T> | |
std::vector< T > | getGBScreeningFactor (int low_index, int high_index) const |
template<typename T> | |
T | getGBScreeningFactor (int index) const |
const AtomGraph * | getSelfPointer () const |
Get a pointer to the object itself. | |
AtomGraph * | getSelfPointer () |
void | setBondParameters (double new_keq, double new_leq, bool set_keq, bool set_leq, int parm_idx, ExceptionResponse policy=ExceptionResponse::SILENT) |
Alter the parameters for a bond parameter set. This will alter the parameters in the topology and all bond terms that use the parameters will then be altered. This function does not add or remove bond terms or parameter sets, nor does it change the mapping of bond terms to parameter sets. | |
void | setBondParameters (double new_keq, double new_leq, bool set_keq, bool set_leq, char4 type_a, char4 type_b, ExceptionResponse policy=ExceptionResponse::SILENT) |
void | setAngleParameters (double new_keq, double new_teq, bool set_keq, bool set_teq, int parm_idx, ExceptionResponse policy=ExceptionResponse::SILENT) |
Alter the parameters for an angle parameter set. This will alter the parameters in the topology and all angle terms that use the parameters will then be use the new parameters. This function does not add or remove angle terms or parameter sets, nor does it change the mapping of angle terms to parameter sets. | |
void | setAngleParameters (double new_keq, double new_teq, bool set_keq, bool set_teq, char4 type_a, char4 type_b, char4 type_c, ExceptionResponse policy=ExceptionResponse::SILENT) |
void | setDihedralParameters (double new_amplitude, double new_phase_angle, bool set_amplitude, bool set_phase_angle, int parm_idx, ExceptionResponse policy=ExceptionResponse::SILENT) |
Alter the parameters for a cosine-based dihedral angle parameter set. This will alter the parameters in the topology and all dihedral terms that use the parameters will then be use the new parameters. This function does not add or remove dihedral terms or parameter sets, nor does it change the mapping of dihedral terms to parameter sets. | |
void | setDihedralParameters (double new_amplitude, double new_phase_angle, bool set_amplitude, bool set_phase_angle, char4 type_a, char4 type_b, char4 type_c, char4 type_d, double periodicity, ExceptionResponse policy=ExceptionResponse::SILENT) |
void | setUreyBradleyParameters (double new_keq, double new_leq, bool set_keq, bool set_leq, int parm_idx, ExceptionResponse policy=ExceptionResponse::SILENT) |
Alter the parameters for a Urey-Bradley parameter set. | |
void | setUreyBradleyParameters (double new_keq, double new_leq, bool set_keq, bool set_leq, char4 type_a, char4 type_b, char4 type_c, ExceptionResponse policy=ExceptionResponse::SILENT) |
void | setCharmmImprParameters (double new_stiffness, double new_phase_angle, bool set_stiffness, bool set_phase_angle, int parm_idx, ExceptionResponse policy=ExceptionResponse::SILENT) |
Alter the characteristics of a CHARMM improper dihedral parameter set. | |
void | setCharmmImprParameters (double new_stiffness, double new_phase_angle, bool set_stiffness, bool set_phase_angle, char4 type_a, char4 type_b, char4 type_c, char4 type_d, ExceptionResponse policy=ExceptionResponse::SILENT) |
void | setWaterResidueName (const char4 new_name) |
Set the name of the water residue. | |
void | setWaterResidueName (const std::string &new_name) |
Friends | |
class | AtomGraphStage |
A struct to hold information relating to an Amber topology. This struct's member functions are limited to getters for its private data. Because the primary means of constructing a topology will be complex, i.e. meticulous parsing of a file or expansion of an existing topology based on some new information, the constructors will need to be buried within wrapper functions that perform such procedures. This struct contains all information contained within an Amber prmtop-format file, nothng more, as described on:
http://ambermd.org/FileFormats.php
The design is intended to be both performant as well as accessible to developers.
stormm::topology::AtomGraph::AtomGraph | ( | const std::string & | file_name, |
ExceptionResponse | policy = ExceptionResponse::WARN, | ||
TopologyKind | engine_format = TopologyKind::AMBER ) |
The general-purpose constructor for file-based topology creation.
Overloaded:
file_name | Name of the source file |
policy | The alert level to raise if a problem is encountered |
engine_format | Format of the topology file to read |
coulomb_constant_in | Value of Coulomb's constant, in kcal-A/mol-e^2 |
default_elec14_screening | The 1:4 electrostatic scaling to apply in absense of any other indications |
default_vdw14_screening | The 1:4 van-der Waals scaling to apply in absense of any other indications |
charge_rounding_tol | The maximum tolerance at which to initiate charge rounding |
charge_discretization | Increment with which to discretize charges |
ag_a | The first of two topologies to combine |
ag_b | The second of two topologies to combine |
stormm::topology::AtomGraph::AtomGraph | ( | const AtomGraph & | original | ) |
stormm::topology::AtomGraph::AtomGraph | ( | AtomGraph && | original | ) |
void stormm::topology::AtomGraph::buildFromPrmtop | ( | const std::string & | file_name, |
const ExceptionResponse | policy = ExceptionResponse::WARN, | ||
double | coulomb_constant_in = amber_ancient_bioq, | ||
double | default_elec14_screening = amber_default_elec14_screen, | ||
double | default_vdw14_screening = amber_default_vdw14_screen, | ||
double | charge_rounding_tol = default_charge_rounding_tol, | ||
double | charge_discretization = default_charge_precision_inc ) |
Build an AtomGraph form a file. This is called by the general-purpose constructor or also by the developer after instantiating an empty object.
file_name | Name of the source file |
policy | Indicates the alert level to raise if a problem is encountered |
std::vector< int > stormm::topology::AtomGraph::findVirtualSites | ( | ) | const |
Get a list of the general topological indices of one or more virtual sites.
Overloaded:
low_index | The start of a stretch of atoms in the general topology |
high_index | The end of a stretch of atoms in the general topology |
index | A specific virtual site index |
AngleTerm< T > stormm::topology::AtomGraph::getAngleTerm | ( | int | index | ) | const |
Collect the atom indices and actual parameters for implement one of the system's harmonic bond angle bending terms.
index | The bond angle term to query, indexed according to the list of bond angle terms |
std::vector< int > stormm::topology::AtomGraph::getAtomExclusions | ( | int | index | ) | const |
Get the exclusion list for a particular atom. This will concatenate 1:1, 1:2, 1:3, and 1:4 exclusions, without regard to which is which. To get a specific type of exclusion, use one of the getNonbondedXXExclusions() functions. All functions of this sort will return a complete list such that if exclusions for all atoms are considered in sequence, all exclusions will be double-counted.
index | The atom of interest |
std::vector< T > stormm::topology::AtomGraph::getAtomicMass | ( | MassForm | rep = MassForm::ORDINARY | ) | const |
Get the masses of atoms in the system. The developer may also choose to take the pre-computed inverse masses (this can be useful when dealing with extra points, as inverses of zero-mass particles must also be set to zero by a conditional statement).
Overloaded:
rep | The representation of the mass to take: ORDINARY or INVERSE |
low_index | The start of a stretch of atoms |
high_index | The end of a stretch of atoms |
index | A specific atom index |
std::vector< int > stormm::topology::AtomGraph::getAtomicNumber | ( | ) | const |
Get the atomic numbers for atoms in the system.
Overloaded:
low_index | Index of the first atom for which to read a Z number |
high_index | Index of the last atom for which to read a Z number |
index | Index of a specific atom |
std::vector< bool > stormm::topology::AtomGraph::getAtomMobility | ( | ) | const |
Get the mobile atom mask for the system, as a boolean vector translated from the bitmask stored internally.
Overloaded:
low_index | Index of the first atom for which to read the mobility |
high_index | Index of the last atom for which to read the mobility |
index | Index of a specific atom |
std::vector< uint > stormm::topology::AtomGraph::getAtomMobilityMask | ( | ) | const |
Get the raw mobility bitmask over a particular range.
Overloaded:
low_index | Index of the first atom for which to read the mobility |
high_index | Index of the last atom for which to read the mobility |
Get the names of atoms in the system.
Overloaded:
low_index | The start of a stretch of atoms |
high_index | The end of a stretch of atoms |
index | A specific atom index |
std::vector< T > stormm::topology::AtomGraph::getAtomPBRadius | ( | ) | const |
Get the Poisson-Boltzmann radius of a particular atom.
Overloaded:
low_index | The start of a stretch of atoms |
high_index | The end of a stretch of atoms |
index | A specific atom index |
Get the types of atoms in the system.
Overloaded:
low_index | The start of a stretch of atoms |
high_index | The end of a stretch of atoms |
index | A specific atom index |
BondTerm< T > stormm::topology::AtomGraph::getBondTerm | ( | int | index | ) | const |
Collect the atom indices and actual parameters for implement one of the system's harmonic bond stretching terms.
index | The bond term to query, indexed according to the list of bond terms |
std::vector< int > stormm::topology::AtomGraph::getChargeIndex | ( | ) | const |
Get the charge type index of atoms in the system.
Overloaded:
low_index | The start of a stretch of atoms |
high_index | The end of a stretch of atoms |
index | A specific atom index |
CharmmImprTerm< T > stormm::topology::AtomGraph::getCharmmImprTerm | ( | int | index | ) | const |
Collect the atom indices and actual parameters to implement one of the system's CHARMM improper terms. The special-pupose, unguarded struct that holds the result also provides the original parameter index, if of interest. This is for computing the potential energy and forces in the most straightforward manner possible, albeit slow.
index | The CHARMM improper term to query, indexed according to the list of CHARMM improper terms (not the general atom list) |
ChemicalDetailsKit stormm::topology::AtomGraph::getChemicalDetailsKit | ( | HybridTargetLevel | tier = HybridTargetLevel::HOST | ) | const |
Get an abstract of the chemical details of the system, such as atom name and type or residue names and limits.
tier | Indicator of whether pointers should target the CPU or the GPU layer |
int stormm::topology::AtomGraph::getCmapDimension | ( | int | index | ) | const |
Return the dimension of a particular CMAP surface (all CMAPs are assumed to be square)
index | The CMAP surface of interest (along the list of available CMAP surfaces, not the list of CMAP terms) |
CmapTerm< T > stormm::topology::AtomGraph::getCmapTerm | ( | int | index | ) | const |
Collect the atom indices and a pointer to the surface to implement one of the system's CMAP terms. The special-purpose, unguarded struct that holds the result also provides the original parameter index, if of interest.
index | The CMAP term to query, indexed according to the list of CMAP terms (not the general atom list) |
std::vector< int > stormm::topology::AtomGraph::getConstraintGroupAtoms | ( | int | index | ) | const |
Get the atoms of a single constraint group.
index | The index of the constraint group within the topology |
int stormm::topology::AtomGraph::getDescriptor | ( | TopologyDescriptor | choice | ) | const |
Get a descriptor from within the array of topology descriptors. If this topology were read from an Amber-format file, all descriptors are taken from the preamble in its first ~15 lines.
Overloaded:
choice | Index of the descriptor in the list. Call this function using the enum classes TopologyDescriptors or SanderDescriptors (see above) for easy, annotated access. |
DihedralTerm< T > stormm::topology::AtomGraph::getDihedralTerm | ( | int | index | ) | const |
Collect the atom indices and actual parameters for implement one of the system's cosine-based dihedral torsion terms.
index | The dihedral term to query, indexed according to the list of dihedral terms |
ConstraintKit< double > stormm::topology::AtomGraph::getDoublePrecisionConstraintKit | ( | HybridTargetLevel | tier = HybridTargetLevel::HOST | ) | const |
Get an abstract for managing constraints based on double-precision parameters.
tier | Indicator of whether pointers should target the CPU or the GPU layer |
ImplicitSolventKit< double > stormm::topology::AtomGraph::getDoublePrecisionImplicitSolventKit | ( | HybridTargetLevel | tier = HybridTargetLevel::HOST | ) | const |
Get a double-precision abstract for calculating implicit solvent energy and forces.
tier | Indicator of whether pointers should target the CPU or the GPU layer |
NonbondedKit< double > stormm::topology::AtomGraph::getDoublePrecisionNonbondedKit | ( | HybridTargetLevel | tier = HybridTargetLevel::HOST | ) | const |
Get a double-precision abstract for calculating non-bonded energy and forces.
tier | Indicator of whether pointers should target the CPU or the GPU layer |
ValenceKit< double > stormm::topology::AtomGraph::getDoublePrecisionValenceKit | ( | HybridTargetLevel | tier = HybridTargetLevel::HOST | ) | const |
Get an abstract for calculating valence terms in double precision.
tier | Indicator of whether pointers should target the CPU or the GPU layer |
VirtualSiteKit< double > stormm::topology::AtomGraph::getDoublePrecisionVirtualSiteKit | ( | HybridTargetLevel | tier = HybridTargetLevel::HOST | ) | const |
Get a double-precision abstract for placing virtual sites and transmitting their forces to frame atoms with mass.
tier | Indicator of whether pointers should target the CPU or the GPU layer |
std::string stormm::topology::AtomGraph::getFullAtomName | ( | int | index | ) | const |
Get the full name (including the residue and structural residue number) of an atom, as a string.
index | The topological index of the atom of interest |
std::vector< T > stormm::topology::AtomGraph::getGBScreeningFactor | ( | ) | const |
Get GB screening factors for one or more atoms.
Overloaded:
low_index | The start of a stretch of atoms |
high_index | The end of a stretch of atoms |
index | A specific atom index |
T stormm::topology::AtomGraph::getLennardJonesEpsilon | ( | int | index_a, |
int | index_b ) const |
Get the Lennard-Jones interaction epsilon value. The quantity is derived from tables at the specified level of precision (float or double), but computed in double precision based on that information before it is returned in the template data type.
Overloaded:
index_a | Index for the first atom |
index_b | Index for the second atom |
std::vector< int > stormm::topology::AtomGraph::getLennardJonesIndex | ( | ) | const |
Get the Lennard-Jones type index of atoms in the system.
Overloaded:
low_index | The start of a stretch of atoms |
high_index | The end of a stretch of atoms |
index | A specific atom index |
T stormm::topology::AtomGraph::getLennardJonesSigma | ( | int | index_a, |
int | index_b ) const |
Get the Lennard-Jones interaction sigma. The quantity is derived from tables at the specified level of precision (float or double), but computed in double precision based on that information before it is returned in the template data type.
Overloaded:
index_a | Index for the first atom |
index_b | Index for the second atom |
std::vector< int > stormm::topology::AtomGraph::getMoleculeContents | ( | ) | const |
Get the contents of a one or more molecules in the system.
Overloaded:
index | Index of a specific molecule |
MoleculeKind stormm::topology::AtomGraph::getMoleculeKind | ( | int | mol_index | ) | const |
Get the kind of a specific molecule wihtin the topology.
mol_index | The index of the molecule of interesta |
std::vector< int > stormm::topology::AtomGraph::getMoleculeLimits | ( | ) | const |
Get the molecule limits, indices into molecule_contents, for every molecule in the system. This allows that molecules not be contiguous in the topology atom indexing, but will be recoverable by molecule_contents.
Overloaded:
index | Index of a specific molecule |
std::vector< int > stormm::topology::AtomGraph::getMoleculeMembership | ( | ) | const |
Get molecule membership for atoms within the system.
Overloaded:
low_index | Index of the first atom for which to read the mobility |
high_index | Index of the last atom for which to read the mobility |
index | Index of a specific atom |
std::vector< int > stormm::topology::AtomGraph::getNonbonded11Exclusions | ( | int | index | ) | const |
Get the 1:1 exclusions for a particular atom. This list will only be populated if the atom is a virtual site or has virtual sites that claim it as their parent.
index | The atom of interest |
std::vector< int > stormm::topology::AtomGraph::getNonbonded12Exclusions | ( | int | index | ) | const |
Get all 1:2 exclusions for a particular atom.
index | The atom of interest |
std::vector< int > stormm::topology::AtomGraph::getNonbonded13Exclusions | ( | int | index | ) | const |
Get all 1:3 exclusions for a particular atom.
index | The atom of interest |
std::vector< int > stormm::topology::AtomGraph::getNonbonded14Exclusions | ( | int | index | ) | const |
Get all 1:4 exclusions for a particular atom.
index | The atom of interest |
std::vector< T > stormm::topology::AtomGraph::getPartialCharge | ( | ) | const |
Get the atomic partial charges of atoms in the system.
Overloaded:
low_index | The start of a stretch of atoms |
high_index | The end of a stretch of atoms |
index | A specific atom index |
std::vector< int > stormm::topology::AtomGraph::getParticlesPerMolecule | ( | ) | const |
Get the sizes of all molecules in the system. Rather redundant with the molecule limits being accessible, but convenient.
Overloaded:
index | Index of a specific residue |
std::vector< int > stormm::topology::AtomGraph::getResidueIndex | ( | ) | const |
Get the index (topology index, not natural / structure-informed residue number) of a particular residue.
Overloaded:
atom_index | Index of the atom in question |
const Hybrid< int > & stormm::topology::AtomGraph::getResidueLimits | ( | ) | const |
Get residue limits for the beginning and ends of up to N residues.
Overloaded:
index | Index of a specific residue |
Get the names of residues in the system.
Overloaded:
index | A specific residue index |
std::vector< int > stormm::topology::AtomGraph::getResidueNumber | ( | ) | const |
Get the residue numbers for atoms in the system.
Overloaded:
low_index | Index of the first atom for which to read a residue number |
high_index | Index of the last atom for which to read a residue number |
index | Index of a specific atom |
const AtomGraph * stormm::topology::AtomGraph::getSelfPointer | ( | ) | const |
ConstraintKit< float > stormm::topology::AtomGraph::getSinglePrecisionConstraintKit | ( | HybridTargetLevel | tier = HybridTargetLevel::HOST | ) | const |
Get an abstract for managing constraints based on single-precision parameters. While the overall implementation for constraints cannot be implemented entirely in single precision, having the underlying parameters in this form is not, in principle, a drvier of energy drift. It will simply set the system into a slightly different geometry, which is probably still closer to the double-precision result than many cycles of SHAKE or RATTLE would get it anyway.
tier | Indicator of whether pointers should target the CPU or the GPU layer |
ImplicitSolventKit< float > stormm::topology::AtomGraph::getSinglePrecisionImplicitSolventKit | ( | HybridTargetLevel | tier = HybridTargetLevel::HOST | ) | const |
Get a single-precision abstract for calculating implicit solvent energy and forces.
tier | Indicator of whether pointers should target the CPU or the GPU layer |
NonbondedKit< float > stormm::topology::AtomGraph::getSinglePrecisionNonbondedKit | ( | HybridTargetLevel | tier = HybridTargetLevel::HOST | ) | const |
Get a single-precision abstract for calculating non-bonded energy and forces.
tier | Indicator of whether pointers should target the CPU or the GPU layer |
ValenceKit< float > stormm::topology::AtomGraph::getSinglePrecisionValenceKit | ( | HybridTargetLevel | tier = HybridTargetLevel::HOST | ) | const |
Get a single-precision abstract for calculating valence terms.
tier | Indicator of whether pointers should target the CPU or the GPU layer |
VirtualSiteKit< float > stormm::topology::AtomGraph::getSinglePrecisionVirtualSiteKit | ( | HybridTargetLevel | tier = HybridTargetLevel::HOST | ) | const |
Get a single-precision abstract for placing virtual sites and transmitting their forces to frame atoms with mass.
tier | Indicator of whether pointers should target the CPU or the GPU layer |
std::vector< int > stormm::topology::AtomGraph::getStructuralAtomNumber | ( | ) | const |
Get the (arbitrarily defined) atom numbers for atoms in the system. For example, a PDB file could ahve its own atom numbering scheme. This provides a way to transcribe that into the AtomGraph object.
Overloaded:
low_index | Index of the first atom for which to read a residue number |
high_index | Index of the last atom for which to read a residue number |
index | Index of a specific atom |
UreyBradleyTerm< T > stormm::topology::AtomGraph::getUreyBradleyTerm | ( | int | index | ) | const |
Collect the atom indices and actual parameters to implement one of the system's Urey-Bradley terms. The special-purpose, unguarded struct that holds the result also provides the original parameter index, if of interest. This is for computing the potential energy and forces in the most straightforward manner possible, albeit slow.
index | The Urey-Bradley term to query, indexed according to the list of Urey-Bradley terms (not the general atom list) |
int stormm::topology::AtomGraph::getVirtualSiteFrameAtom | ( | int | index, |
int | nfrm ) const |
Get the general topology index of one of a virtual site's frame atoms.
index | The virtual site of interest, indexed according to its place in the list of virtual sites (not the general topology atom list) |
nfrm | The frame atom to pull out (values of 1:4 for one, two, three, or four-atom frames are acceptable) |
T stormm::topology::AtomGraph::getVirtualSiteFrameDimension | ( | int | index, |
int | ndim ) |
Get the dimensions (could be length in A, could be scaling factor for a distance or cross product) of one virtual site's frame.
index | The virtual site of interest, indexed according to its place in the list of virtual sites (not the general topology atom list) |
ndim | The frame dimension to pull out (values of 1-3 for frames with one, two, or three details are acceptable) |
VirtualSiteKind stormm::topology::AtomGraph::getVirtualSiteFrameType | ( | int | index | ) | const |
Get the frame type of a particular virtual site.
index | The virtual site of interest, indexed according to its place in the list of virtual sites (not the general topology atom list) |
int stormm::topology::AtomGraph::getVirtualSiteIndex | ( | int | atom_index, |
ExceptionResponse | policy = ExceptionResponse::DIE ) const |
Get the virtual site index of an atom index. If the atom index is not a virtual site, this function will throw an error.
atom_index | Index of the atom in question |
policy | Action to take if the atom turns out not to be a virtual site at all |
std::vector< char4 > stormm::topology::AtomGraph::matchOverflowAtomName | ( | const std::string & | query | ) | const |
Match a long atom name and return its key in the standard list of char4 atom names.
query | The extended atom name to match |
std::vector< char4 > stormm::topology::AtomGraph::matchOverflowAtomType | ( | const std::string & | query | ) | const |
Match a long atom type and return its key in the standard list of char4 atom types.
query | The extended atom type name to match |
std::vector< char4 > stormm::topology::AtomGraph::matchOverflowResidueName | ( | const std::string & | query | ) | const |
Match a long residue name and return its key among the standard char4 residue names.
query | The extended residue name to match |
void stormm::topology::AtomGraph::modifyAtomMobility | ( | MobilitySetting | movement | ) |
Change an atom's mobility in the topology.
Overloaded:
movement | The desired atom mobility (OFF, ON, or TOGGLE) |
low_index | The start of a stretch of atoms |
high_index | The end of a stretch of atoms |
index | A specific atom index |
Move constructor.
other | Another AtomGraph to form the basis for re-assigning members of this one |
Copy constructor.
other | Another AtomGraph to form the basis for re-assigning members of this one |
void stormm::topology::AtomGraph::printToFile | ( | const std::string & | output_file = {}, |
TopologyKind | output_style = TopologyKind::AMBER, | ||
PrintSituation | expectation = PrintSituation::OPEN_NEW, | ||
ExceptionResponse | pr_policy = ExceptionResponse::DIE ) const |
Print out a topology in the specified format.
output_file | The name of the file to print. If an empty string is given, the source file will be scheduled for overwriting, pending a user-specifiable directive as to whether overwriting is enabled. |
output_style | The format of topology to write, each format being readable by specific third-party molecular mechancis packages |
expectation | Indicates what state the output file is expected to be found in, if such a file may be found to already exist at all |
pr_policy | Specifies the course of action to take if errors are encountered |
void stormm::topology::AtomGraph::setAngleParameters | ( | double | new_keq, |
double | new_teq, | ||
bool | set_keq, | ||
bool | set_teq, | ||
int | parm_idx, | ||
ExceptionResponse | policy = ExceptionResponse::SILENT ) |
Alter the parameters for an angle parameter set. This will alter the parameters in the topology and all angle terms that use the parameters will then be use the new parameters. This function does not add or remove angle terms or parameter sets, nor does it change the mapping of angle terms to parameter sets.
Overloaded:
new_keq | The new stiffness constant to apply |
new_teq | The new angular equilibrium constant to apply |
set_keq | Whether to change the stiffness constant |
set_teq | Whether to change the equilibrium |
parm_idx | Index of the angle parameter set from within the topology's table (one cannot depend on this index to be consistent across topologies of different systems making use of the same force field) |
type_a | Atom type of the first atom in the angle (this will be consistent across multiple topologies using the same force field) |
type_b | Atom type of the second atom in the angle |
type_c | Atom type of the third atom in the angle |
policy | Behavior if an angle parameter set matching the index or atom types is not found in the topology |
void stormm::topology::AtomGraph::setBondParameters | ( | double | new_keq, |
double | new_leq, | ||
bool | set_keq, | ||
bool | set_leq, | ||
int | parm_idx, | ||
ExceptionResponse | policy = ExceptionResponse::SILENT ) |
Alter the parameters for a bond parameter set. This will alter the parameters in the topology and all bond terms that use the parameters will then be altered. This function does not add or remove bond terms or parameter sets, nor does it change the mapping of bond terms to parameter sets.
Overloaded:
new_keq | The new stiffness constant to apply |
new_leq | The new equilibrium constant to apply |
set_keq | Whether to change the stiffness constant |
set_leq | Whether to change the equilibrium |
parm_idx | Index of the bond parameter set from within the topology's table (one cannot depend on this index to be consistent across topologies of different systems making use of the same force field) |
type_a | Atom type of the first atom in the bond (this will be consistent across multiple topologies using the same force field) |
type_b | Atom type of the second atom in the bond |
policy | Behavior if a bond parameter set matching the index or atom types is not found in the topology |
void stormm::topology::AtomGraph::setCharmmImprParameters | ( | double | new_stiffness, |
double | new_phase_angle, | ||
bool | set_stiffness, | ||
bool | set_phase_angle, | ||
int | parm_idx, | ||
ExceptionResponse | policy = ExceptionResponse::SILENT ) |
Alter the characteristics of a CHARMM improper dihedral parameter set.
Overloaded:
new_stiffness | The new stiffness constant to apply |
new_phase_angle | The new phase angle constant to apply |
set_stiffness | Whether to change the stiffness constant |
set_phase_angle | Whether to change the equilibrium |
parm_idx | Index of the CHARMM improper dihedral parameter set from within the topology's table (one cannot depend on this index to be consistent across topologies of different systems making use of the same force field) |
type_a | Atom type of the first atom in the improper (this will be consistent across multiple topologies using the same force field) |
type_b | Atom type of the second atom in the improper |
type_c | Atom type of the third atom in the improper |
type_d | Atom type of the fourth atom in the improper |
policy | Behavior if a CHARMM improper dihedral parameter set matching the index or atom types and periodicity is not found in the topology |
void stormm::topology::AtomGraph::setDihedralParameters | ( | double | new_amplitude, |
double | new_phase_angle, | ||
bool | set_amplitude, | ||
bool | set_phase_angle, | ||
int | parm_idx, | ||
ExceptionResponse | policy = ExceptionResponse::SILENT ) |
Alter the parameters for a cosine-based dihedral angle parameter set. This will alter the parameters in the topology and all dihedral terms that use the parameters will then be use the new parameters. This function does not add or remove dihedral terms or parameter sets, nor does it change the mapping of dihedral terms to parameter sets.
Overloaded:
new_amplitude | The new cosine series term amplitude constant to apply |
new_phase_angle | The new phase angle to apply |
set_amplitude | Whether to change the amplitude constant |
set_phase_angle | Whether to change the equilibrium |
parm_idx | Index of the dihedral angle parameter set from within the topology's table (one cannot depend on this index to be consistent across topologies of different systems making use of the same force field) |
type_a | Atom type of the first atom in the dihedral (this will be consistent across multiple topologies using the same force field) |
type_b | Atom type of the second atom in the dihedral |
type_c | Atom type of the third atom in the dihedral |
type_d | Atom type of the fourth atom in the dihedral |
periodicity | Periodicity of the dihedral parameter set (this disambiguates dihedral parameters that apply on top of one another to the same four atom types) |
policy | Behavior if a dihedral parameter set matching the index or atom types and periodicity is not found in the topology |
void stormm::topology::AtomGraph::setImplicitSolventModel | ( | ImplicitSolventModel | igb_in, |
double | dielectric_in = 80.0, | ||
double | saltcon_in = 0.0, | ||
AtomicRadiusSet | radii_set = AtomicRadiusSet::NONE, | ||
ExceptionResponse | policy = ExceptionResponse::WARN ) |
Set the implicit solvent model. This will leverage a lot of hard-coded constants to fill out some allocated but otherwise blank arrays and impart one particular Generalized Born or other implicit solvent method.
igb_in | The implicit solvent model to impart onto the system |
dielectric_in | The desired dielectric constant |
saltcon_in | The intended salt concentration (affects the GB decay parameter Kappa) |
radii_set | Radii to impart to the topology (this is often coupled to the choice of implicit solvent model, but for purposes of experimentation or new model development might be flexible) |
policy | Indicator of what to do if the topology's PB radii to not meet the implicit solvent model requirements, or there is some other problem |
void stormm::topology::AtomGraph::setUreyBradleyParameters | ( | double | new_keq, |
double | new_leq, | ||
bool | set_keq, | ||
bool | set_leq, | ||
int | parm_idx, | ||
ExceptionResponse | policy = ExceptionResponse::SILENT ) |
Alter the parameters for a Urey-Bradley parameter set.
Overloaded:
new_keq | The new stiffness constant to apply |
new_leq | The new equilibrium constant to apply |
set_keq | Whether to change the stiffness constant |
set_leq | Whether to change the equilibrium |
parm_idx | Index of the bond parameter set from within the topology's table (one cannot depend on this index to be consistent across topologies of different systems making use of the same force field) |
type_a | Atom type of the first atom in the Urey-Bradley term, the I atom (this will be consistent across multiple topologies using the same force field) |
type_b | Atom type of the second atom in the Urey-Bradley term (the K atom) |
policy | Behavior if a Urey-Bradley parameter set matching the index or atom types is not found in the topology |
void stormm::topology::AtomGraph::setWaterResidueName | ( | const char4 | new_name | ) |
Set the name of the water residue.
Overloaded:
new_name | The new water residue name to use |
|
friend |
Allow the AtomGraphStage to access private members directly to modify and guide construction of topologies used in STORMM's production calculations.