Molmodel
Classes | Namespaces | Modules | Enumerator | Functions | Friends

Molecular Mechanics in Molmodel

Once you have built a system of Compounds (molecules) in Molmodel, and modeled them with the degrees of freedom you want, you can apply molecular forces using DuMM, the Molmodel molecular mechanics force field. More...

Classes

class  SimTK::DuMM::CustomBondStretch
 Abstract base class for custom bond stretch terms, that is, functional forms that apply a distance-based force along a 1-2 bond between a pair of atoms. More...
class  SimTK::DuMMForceFieldSubsystem
 This is a concrete subsystem that provides basic molecular mechanics functionality for coarse-grained molecules built in the SimTK framework. More...
class  SimTK::Amber99ForceSubsystem
 This class is just a DuMMForceFieldSubsystem for which the constructor pre-loads the definitions of the Amber99 force field. More...
class  SimTK::DuMM::IncludedAtomIndex
 This is the unique integer index that DuMM assigns to atoms that are to be included in force calculations. These represent a subset of all the atoms and you can map from an IncludedAtomIndex to the corresponding AtomIndex. More...
class  SimTK::DuMM::NonbondAtomIndex
 This is the unique integer index that DuMM assigns to included atoms that are involved in nonbonded force calculations (Coulomb, van der Waals, and/or GBSA). These represent a subset of all the included atoms and you can map from a NonbondAtomIndex to the corresponding IncludedAtomIndex. More...
class  SimTK::DuMM::BondIndex
 This is the unique integer index that DuMM assigns to bonds as they are added via addBond(). More...
class  SimTK::DuMM::ClusterIndex
 This is the unique integer index that DuMM assigns to clusters as they are created via createCluster(). More...
class  SimTK::DuMM::AtomClassIndex
 This is a unique integer associated with each "atom class", assigned by the user when the atom class is first introduced via defineAtomClass(). More...
class  SimTK::DuMM::ChargedAtomTypeIndex
 This is a unique integer associated with each "charged atom type", assigned by the user when the charged atom type is first introduced via defineChargedAtomType(). More...

Namespaces

namespace  SimTK::DuMM
 

This namespace is used for symbols which are useful in conjunction with Molmodel's DuMMForceFieldSubsystem.


Modules

 User defined molecular mechanics force terms
 

These classes are used with DuMMForceFieldSubsystem to permit user-defined functional forms for bonded terms.


 Handy conversion constants
 

These are a set of multiplicative factors for use in converting among the various commonly used molecular mechanics units.


Functions

 SimTK::DuMM::SimTK_DEFINE_UNIQUE_INDEX_TYPE (AtomIndex)
 SimTK::DuMM::SimTK_DEFINE_UNIQUE_INDEX_TYPE (IncludedAtomIndex)
 SimTK::DuMM::SimTK_DEFINE_UNIQUE_INDEX_TYPE (NonbondAtomIndex)
 SimTK::DuMM::SimTK_DEFINE_UNIQUE_INDEX_TYPE (BondIndex)
 SimTK::DuMM::SimTK_DEFINE_UNIQUE_INDEX_TYPE (ClusterIndex)
 SimTK::DuMM::SimTK_DEFINE_UNIQUE_INDEX_TYPE (AtomClassIndex)
 SimTK::DuMM::SimTK_DEFINE_UNIQUE_INDEX_TYPE (ChargedAtomTypeIndex)
virtual Real SimTK::DuMM::CustomBondStretch::calcEnergy (Real distance) const =0
virtual Real SimTK::DuMM::CustomBondStretch::calcForce (Real distance) const =0
virtual Real SimTK::DuMM::CustomBondBend::calcEnergy (Real bendAngle) const =0
virtual Real SimTK::DuMM::CustomBondBend::calcTorque (Real bendAngle) const =0
virtual Real SimTK::DuMM::CustomBondTorsion::calcEnergy (Real dihedralAngle) const =0
virtual Real SimTK::DuMM::CustomBondTorsion::calcTorque (Real dihedralAngle) const =0
 SimTK::DuMMForceFieldSubsystem::DuMMForceFieldSubsystem ()
 SimTK::DuMMForceFieldSubsystem::DuMMForceFieldSubsystem (MolecularMechanicsSystem &)
void SimTK::DuMMForceFieldSubsystem::setTraceOpenMM (bool shouldTrace)
 OBSOLETE NAME -- just use setTracing().
void SimTK::DuMMForceFieldSubsystem::defineIncompleteAtomClass (DuMM::AtomClassIndex classIx, const char *name, int elementNumber, int valence)
void SimTK::DuMMForceFieldSubsystem::defineIncompleteAtomClass_KA (DuMM::AtomClassIndex classIx, const char *name, int elementNumber, int valence)
void SimTK::DuMMForceFieldSubsystem::setAtomClassVdwParameters (DuMM::AtomClassIndex atomClassIx, Real vdwRadiusInNm, Real vdwWellDepthInKJPerMol)
void SimTK::DuMMForceFieldSubsystem::setAtomClassVdwParameters_KA (DuMM::AtomClassIndex atomClassIx, Real radiusInAng, Real wellDepthInKcal)
bool SimTK::DuMMForceFieldSubsystem::isValidAtomClass (DuMM::AtomClassIndex) const
void SimTK::DuMMForceFieldSubsystem::defineIncompleteChargedAtomType (DuMM::ChargedAtomTypeIndex typeIx, const char *name, DuMM::AtomClassIndex classIx)
void SimTK::DuMMForceFieldSubsystem::defineIncompleteChargedAtomType_KA (DuMM::ChargedAtomTypeIndex typeIx, const char *name, DuMM::AtomClassIndex classIx)
void SimTK::DuMMForceFieldSubsystem::setChargedAtomTypeCharge (DuMM::ChargedAtomTypeIndex, Real charge)
void SimTK::DuMMForceFieldSubsystem::setChargedAtomTypeCharge_KA (DuMM::ChargedAtomTypeIndex chargedAtomTypeIx, Real charge)

Friends

class SimTK::DuMMForceFieldSubsystem::MolecularMechanicsSystem

Define particular molecules

Methods in this group are used to define the atoms and bonds in the particular set of molecules being simulated.
DuMM::AtomIndex SimTK::DuMMForceFieldSubsystem::addAtom (DuMM::ChargedAtomTypeIndex chargedAtomTypeIx)
 Add a new atom to the model.
DuMM::BondIndex SimTK::DuMMForceFieldSubsystem::addBond (DuMM::AtomIndex atom1Ix, DuMM::AtomIndex atom2Ix)
 Declare that there is a covalent bond between two atoms.
DuMM::AtomIndex SimTK::DuMMForceFieldSubsystem::getBondAtom (DuMM::BondIndex bond, int which) const
 For a given 1-2 bond, return the atoms which are connected by that bond.
int SimTK::DuMMForceFieldSubsystem::getNumAtoms () const
 How many atoms are currently in the model?
int SimTK::DuMMForceFieldSubsystem::getNumBonds () const
 How many 1-2 bonds are currently in the model?
Real SimTK::DuMMForceFieldSubsystem::getAtomMass (DuMM::AtomIndex atomIx) const
 Obtain the mass in Daltons (g/mol) of the atom indicated by the given AtomIndex.
int SimTK::DuMMForceFieldSubsystem::getAtomElement (DuMM::AtomIndex atomIx) const
 Obtain the element (by atomic number) of the atom indicated by the given AtomIndex.
Real SimTK::DuMMForceFieldSubsystem::getAtomRadius (DuMM::AtomIndex atomIx) const
 Obtain the van der Waals radius of the atom indicated by the given AtomIndex.
MobilizedBodyIndex SimTK::DuMMForceFieldSubsystem::getAtomBody (DuMM::AtomIndex atomIx) const
 Obtain the Simbody MobilizedBodyIndex of the rigid body on which a particular atom has been fixed.
Vec3 SimTK::DuMMForceFieldSubsystem::getAtomStationOnBody (DuMM::AtomIndex atomIx) const
 Obtain the station at which a particular atom is fixed on its body.
Vec3 SimTK::DuMMForceFieldSubsystem::getAtomStationInCluster (DuMM::AtomIndex atomIx, DuMM::ClusterIndex clusterIx) const
 Obtain the station at which a particular atom is fixed within a particular Cluster (an atom can be in more than one Cluster).
Vec3 SimTK::DuMMForceFieldSubsystem::getElementDefaultColor (int atomicNumber) const
 For display purposes, return the RGB value of a suggested color for an element given by atomic number.
Vec3 SimTK::DuMMForceFieldSubsystem::getAtomDefaultColor (DuMM::AtomIndex atomIx) const
 For display purposes, return the RGB value of a suggested color with which to display a particular atom.

Define clusters and bodies

Methods in this group control the grouping of atoms into rigid clusters and the placement of such clusters onto the rigid bodies of the underlying multibody system.Note: we use the term "station" to refer to a fixed location with respect to a cluster or body frame, that is, a point "stationary" on that cluster or body.
DuMM::ClusterIndex SimTK::DuMMForceFieldSubsystem::createCluster (const char *clusterName)
 Create an empty Cluster (rigid group of atoms).
void SimTK::DuMMForceFieldSubsystem::placeAtomInCluster (DuMM::AtomIndex atomIx, DuMM::ClusterIndex clusterIx, const Vec3 &station)
 Place an existing atom at a particular station in the local frame of a Cluster.
void SimTK::DuMMForceFieldSubsystem::placeClusterInCluster (DuMM::ClusterIndex childClusterIndex, DuMM::ClusterIndex parentClusterIndex, const Transform &X_PC)
 Place a Cluster (the child) in another Cluster (the parent).
MassProperties SimTK::DuMMForceFieldSubsystem::calcClusterMassProperties (DuMM::ClusterIndex clusterIx, const Transform &X_BC=Transform()) const
 Calculate the composite mass properties of a Cluster, either in its own reference frame C or in reference frame B with the Cluster placed relative to B using the indicated Transform X_BC.
void SimTK::DuMMForceFieldSubsystem::attachClusterToBody (DuMM::ClusterIndex clusterIx, MobilizedBodyIndex body, const Transform &X_BC=Transform())
 Place a Cluster's local frame C at a particular location and orientation with respect to a MobilizedBody's frame B.
void SimTK::DuMMForceFieldSubsystem::attachAtomToBody (DuMM::AtomIndex atomIx, MobilizedBodyIndex body, const Vec3 &station=Vec3(0))
 Place an individual atom at a particular station on a body without an intervening cluster.
MobilizedBodyIndex SimTK::DuMMForceFieldSubsystem::getClusterBody (DuMM::ClusterIndex clusterIx) const
 Find the MobilizedBody on which a Cluster has been fixed in place.
Transform SimTK::DuMMForceFieldSubsystem::getClusterPlacementOnBody (DuMM::ClusterIndex clusterIx) const
 Find where on its body a Cluster has been placed by returning the transform X_BC giving the orientation and position of Cluster frame C in body frame B.
Transform SimTK::DuMMForceFieldSubsystem::getClusterPlacementInCluster (DuMM::ClusterIndex childClusterIndex, DuMM::ClusterIndex parentClusterIndex) const
 Find where on parent cluster P a child cluster C has been placed, by returning the transform X_PC.

Atom/bond exclusion methods (advanced users only!)

Methods in this group give you fine control over which atoms in the defined molecules are actually used in the calculation of nonbonded forces, and which bonds are included in the calculation of bonded forces.By default, any atom or bond is included if it can contribute to forces that affect motion.
Warning:
You can cause strange effects with these methods if you're not careful -- for example excluding atoms while using GBSA will cause odd surface area calculations to be performed. Don't use these methods unless you know exactly what you're doing!

During realizeTopology() DuMM studies the defined molecules and force field terms to determine (for example) which bonds cross bodies and which atoms can generate forces. Then the list of included atoms and bonds is finalized based on the instructions you give by calling these methods.

void SimTK::DuMMForceFieldSubsystem::clearIncludedNonbondAtomList ()
 Clear the list of atoms to be included in nonbonded force calculations. This is usually called prior to including selected groups of atoms.
void SimTK::DuMMForceFieldSubsystem::clearIncludedBondList ()
 Clear the list of bonds to be included in bonded force calculations. This is usually called prior to including selected bonds.
void SimTK::DuMMForceFieldSubsystem::resetIncludedNonbondAtomListToDefault ()
 Restore the included nonbond atoms list to its DuMM-selected default. This will include all the atoms that can affect interbody nonbonded forces, including Coulomb, Van der Waals, and solvent (GBSA) forces.
void SimTK::DuMMForceFieldSubsystem::resetIncludedBondListToDefault ()
 Request that DuMM set the included bonds list to all bonds that can generate forces that affect motion, meaning bonds whose bonded force terms have a non-zero amplitude and that involve atoms from two or more bodies.
void SimTK::DuMMForceFieldSubsystem::includeNonbondAtom (DuMM::AtomIndex atom)
 Add the given atom to the included nonbond atom list if it isn't already there. This does not include bonded terms that involve this atom; you have to request that explicitly with includeAllInterbodyBondsForOneAtom().
void SimTK::DuMMForceFieldSubsystem::includeAllNonbondAtomsForOneBody (MobilizedBodyIndex mobod)
 Add to the included nonbond atom list all the atoms that are fixed to the given body.
void SimTK::DuMMForceFieldSubsystem::includeAllInterbodyBondsForOneAtom (DuMM::AtomIndex atom)
 Given an atom, include all the bonded force terms that (a) include this atom, and (b) involve a body other than the one to which this atom is fixed. This does not include the atom in nonbonded force calculations; you have to request that explicitly with includeNonbondAtom().
void SimTK::DuMMForceFieldSubsystem::includeAllInterbodyBondsWithBothAtoms (DuMM::AtomIndex atom1, DuMM::AtomIndex atom2)
 Given two atoms that may appear together in some bonded force term, include all the bonded terms that involve both of them.
void SimTK::DuMMForceFieldSubsystem::includeAllInterbodyBondsWithBothAtoms (DuMM::BondIndex bond)
 Given a bond (that is, a 1-2 connection between atoms), extract the two atoms from it and then make sure all interbody bond terms that involve both of those atoms (in any position) are included.
void SimTK::DuMMForceFieldSubsystem::includeAllInterbodyBondsForOneBody (MobilizedBodyIndex mobod)
 Include any bonds that are needed to calculate interbody bonded forces that involve any atoms fixed to the given mobilized body. This will involve a few atoms on neighboring bodies.
void SimTK::DuMMForceFieldSubsystem::includeAllInterbodyBondsBetweenTwoBodies (MobilizedBodyIndex mobod1, MobilizedBodyIndex mobod2)
 Include any bonds that are needed to calculate interbody bonded forces that act between the given pair of mobilized bodies. Any bonded term that involves both an atom from mobod1 and an atom from mobod2 will be included.
int SimTK::DuMMForceFieldSubsystem::getNumIncludedAtoms () const
 Atoms to be included in force calculations are numbered from 0 to getNumIncludedAtoms()-1 after realizeTopology() has been called.
DuMM::AtomIndex SimTK::DuMMForceFieldSubsystem::getAtomIndexOfIncludedAtom (DuMM::IncludedAtomIndex incAtomIx) const
 Given a DuMM::IncludedAtomIndex, return the corresponding DuMM::AtomIndex. You must already have called realizeTopology().
int SimTK::DuMMForceFieldSubsystem::getNumNonbondAtoms () const
 A subset of the included atoms are used in nonbond calculations. These are numbered from 0 to getNumNonbondAtoms()-1 after realizeTopology() has been called.
DuMM::IncludedAtomIndex SimTK::DuMMForceFieldSubsystem::getIncludedAtomIndexOfNonbondAtom (DuMM::NonbondAtomIndex nonbondAtomIx) const
 Given a DuMM::NonbondAtomIndex, return the corresponding DuMM::IncludedAtomIndex. You must already have called realizeTopology().
DuMM::AtomIndex SimTK::DuMMForceFieldSubsystem::getAtomIndexOfNonbondAtom (DuMM::NonbondAtomIndex nonbondAtomIx) const
 Given a DuMM::NonbondAtomIndex, return the corresponding DuMM::AtomIndex. You must already have called realizeTopology().

Define atom categories

An AtomClass is used to collect together a set of properties which are expected to be shared by many individual atoms.The properties are: the element (as atomic number), expected valence, and van der Waals parameters. Charge is not included in AtomClass but in a second more detailed classification level called ChargedAtomType.
void SimTK::DuMMForceFieldSubsystem::defineAtomClass (DuMM::AtomClassIndex atomClassIx, const char *atomClassName, int atomicNumber, int expectedValence, Real vdwRadiusInNm, Real vdwWellDepthInKJ)
 Define a new atom class for this force field, for identifying atoms of a particular element, number of bonds, and van der Waals parameters.
void SimTK::DuMMForceFieldSubsystem::defineAtomClass_KA (DuMM::AtomClassIndex atomClassIx, const char *atomClassName, int element, int valence, Real vdwRadiusInAng, Real vdwWellDepthInKcal)
 Same routine as defineAtomClass() but in Kcal/Angstrom (KA) unit system, that is, radius (still not sigma) is in Angstroms, and well depth in kcal/mol.
void SimTK::DuMMForceFieldSubsystem::defineAtomClass_KA (int atomClassIx, const char *atomClassName, int element, int valence, Real vdwRadiusInAng, Real vdwWellDepthInKcal)
 Obsolete method -- use the other signature.
bool SimTK::DuMMForceFieldSubsystem::hasAtomClass (DuMM::AtomClassIndex) const
 Check whether an atom class has been defined using this index.
bool SimTK::DuMMForceFieldSubsystem::hasAtomClass (const String &atomClassName) const
 Check whether an atom class has been defined using this name.
DuMM::AtomClassIndex SimTK::DuMMForceFieldSubsystem::getAtomClassIndex (const String &atomClassName) const
 Obtain the atom class index corresponding to this atom class name.
DuMM::AtomClassIndex SimTK::DuMMForceFieldSubsystem::getNextUnusedAtomClassIndex () const
 Obtain an atom class index that is numerically larger than the largest currently-defined atom class index.
DuMM::AtomClassIndex SimTK::DuMMForceFieldSubsystem::getAtomClassIndex (DuMM::AtomIndex atomIx) const
 Get the index number of the atom class associated with this atom.
Real SimTK::DuMMForceFieldSubsystem::getVdwRadius (DuMM::AtomClassIndex atomClassIx) const
 Get the van der Waals radius shared by all atoms that belong to the indicated atom class. See comments for this group for a precise definition of what this means; there are ambiguities so don't assume you already know.
Real SimTK::DuMMForceFieldSubsystem::getVdwWellDepth (DuMM::AtomClassIndex atomClassIx) const
 Get the van der Waals energy well depth shared by all atoms that belong to the indicated atom class. See comments for this group for a precise definition of what this means; there are ambiguities so don't assume you already know.
void SimTK::DuMMForceFieldSubsystem::defineChargedAtomType (DuMM::ChargedAtomTypeIndex atomTypeIx, const char *atomTypeName, DuMM::AtomClassIndex atomClassIx, Real partialChargeInE)
 Define a new ChargedAtomType for this force field, for identifying atoms of a particular AtomClass that have a particular partial charge.
void SimTK::DuMMForceFieldSubsystem::defineChargedAtomType_KA (DuMM::ChargedAtomTypeIndex atomTypeIx, const char *atomTypeName, DuMM::AtomClassIndex atomClassIx, Real partialChargeInE)
 This is an alternate name for defineChargedAtomType() but since partial charge uses the same unit (e) in both the MD and KA unit systems, this is identical to defineChargedAtomType().
void SimTK::DuMMForceFieldSubsystem::defineChargedAtomType_KA (int atomTypeIx, const char *atomTypeName, int atomClassIx, Real partialChargeInE)
 Obsolete method -- use the other signature.
bool SimTK::DuMMForceFieldSubsystem::hasChargedAtomType (DuMM::ChargedAtomTypeIndex) const
 Check whether a charged atom type has been defined using this index.
bool SimTK::DuMMForceFieldSubsystem::hasChargedAtomType (const String &chargedTypeName) const
 Check whether a charged atom type has been defined using this name.
DuMM::ChargedAtomTypeIndex SimTK::DuMMForceFieldSubsystem::getChargedAtomTypeIndex (const String &chargedTypeName) const
 Obtain the charged atom type index corresponding to this charged atom type name.
DuMM::ChargedAtomTypeIndex SimTK::DuMMForceFieldSubsystem::getNextUnusedChargedAtomTypeIndex () const
 Obtain a charged atom type index that is numerically larger than the largest currently-defined charged atom type index.

Control force field nonbonded behavior in special circumstances

These methods permit setting overall force field behavior in special circumstances, including van der Waals mixing behavior for dissimilar atom pairs and scaling of non-bonded terms for closely-bonded atoms.
const char * SimTK::DuMMForceFieldSubsystem::getVdwMixingRuleName (VdwMixingRule) const
 Obtain a human-readable name for one of our van der Waals mixing rules.
void SimTK::DuMMForceFieldSubsystem::setVdwMixingRule (VdwMixingRule)
 Set the van der Waals mixing rule -- our default is Waldman-Hagler.
VdwMixingRule SimTK::DuMMForceFieldSubsystem::getVdwMixingRule () const
 Get the van der Waals mixing rule currently in effect.
void SimTK::DuMMForceFieldSubsystem::setVdw12ScaleFactor (Real)
 default 0
void SimTK::DuMMForceFieldSubsystem::setVdw13ScaleFactor (Real)
 default 0
void SimTK::DuMMForceFieldSubsystem::setVdw14ScaleFactor (Real)
 default 1
void SimTK::DuMMForceFieldSubsystem::setVdw15ScaleFactor (Real)
 default 1
void SimTK::DuMMForceFieldSubsystem::setCoulomb12ScaleFactor (Real)
 default 0
void SimTK::DuMMForceFieldSubsystem::setCoulomb13ScaleFactor (Real)
 default 0
void SimTK::DuMMForceFieldSubsystem::setCoulomb14ScaleFactor (Real)
 default 1
void SimTK::DuMMForceFieldSubsystem::setCoulomb15ScaleFactor (Real)
 default 1

Tinker biotypes and pre-defined force field parameter sets

DuMM understands Tinker-format parameter files that can be used to load in a whole force field description.This requires assigning Tinker Biotypes to the atoms in the molecule.
void SimTK::DuMMForceFieldSubsystem::loadAmber99Parameters ()
 Use Amber99 force field parameters.
void SimTK::DuMMForceFieldSubsystem::populateFromTinkerParameterFile (std::istream &)
 Load force field parameters from a TINKER format force field parameter file. (Only the Amber99 force field has been tested.)
void SimTK::DuMMForceFieldSubsystem::setBiotypeChargedAtomType (DuMM::ChargedAtomTypeIndex chargedAtomTypeIndex, BiotypeIndex biotypeIx)
 Associate a Tinker Biotype with a ChargedAtomType in this subsystem.
DuMM::ChargedAtomTypeIndex SimTK::DuMMForceFieldSubsystem::getBiotypeChargedAtomType (BiotypeIndex biotypeIx) const
 Get charged atom type index in this force field associated with a particular Biotype.
std::ostream & SimTK::DuMMForceFieldSubsystem::generateBiotypeChargedAtomTypeSelfCode (std::ostream &os) const
 Generate C++ code from the current contents of this DuMM force field object.

Bond stretch terms

Bond stretch parameters (between 2 atom classes).You can use the standard, built-in functional form (a harmonic) or define your own.
void SimTK::DuMMForceFieldSubsystem::defineBondStretch (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, Real stiffnessInKJperNmSq, Real nominalLengthInNm)
 Define a harmonic bond stretch term between two atom classes using the built-in functional form.
void SimTK::DuMMForceFieldSubsystem::defineBondStretch_KA (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, Real stiffnessInKcalPerAngSq, Real nominalLengthInAng)
 Same as defineBondStretch() except that for convenience this takes stiffness in (kcal/mol)/A^2, and nominal length is in A (angstroms).
void SimTK::DuMMForceFieldSubsystem::defineBondStretch_KA (int class1, int class2, Real stiffnessInKcalPerAngSq, Real nominalLengthInAng)
 Same as defineBondStretch_KA() but takes integer class arguments for backwards compatibility.
void SimTK::DuMMForceFieldSubsystem::defineCustomBondStretch (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::CustomBondStretch *bondStretchTerm)
 Define a custom bond stretch term to be applied to the indicated pair of atom classes whenever they are found in a 1-2 bond.

Bond bending terms

Bond bending parameters (for 3 atom classes bonded 1-2-3).You can use the standard, built-in functional form (a harmonic based on the 1-2-3 angle) or define your own.
void SimTK::DuMMForceFieldSubsystem::defineBondBend (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3, Real stiffnessInKJPerRadSq, Real nominalAngleInDeg)
 Define a harmonic bond bending term applying to a triple of atom classes using the built-in functional form.
void SimTK::DuMMForceFieldSubsystem::defineBondBend_KA (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3, Real stiffnessInKcalPerRadSq, Real nominalAngleInDeg)
 Same as defineBondBend() except that for convenience this takes stiffness in (kcal/mol)/rad^2 (nominal angle is still in degrees).
void SimTK::DuMMForceFieldSubsystem::defineBondBend_KA (int class1, int class2, int class3, Real stiffnessInKcalPerRadSq, Real nominalAngleInDeg)
 Same as defineBondBend_KA() but takes integer class arguments for backwards compatibility.
void SimTK::DuMMForceFieldSubsystem::defineCustomBondBend (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3, DuMM::CustomBondBend *bondBendTerm)
 Define a custom bond bend term to be applied to the indicated triple of atom classes whenever they are found in a 1-2-3 bonded sequence.

Bond torsion terms

Bond torsion (dihedral) parameters (for 4 atom classes bonded 1-2-3-4).You can use the standard, built-in functional form (a combination of sinusoids) or define your own. Bond torsion terms produce energy and a scalar torque that are dependent only on the rotation angle about the 2-3 bond in the specified atom class sequence.Theory Label the four bonded atoms r-x-y-s. Rotation occurs about the axis v=y-x, that is, a vector from x to y. We define a torsion angle theta using the "polymer convention" rather than the IUPAC one which is 180 degrees different. Ours is like this:

            r                         r      s
  theta=0    \             theta=180   \    / 
              x--y                      x--y
                  \
                   s

The sign convention is the same for IUPAC and polymer: A positive angle is defined by considering r-x fixed in space. Then using the right hand rule around v (that is, thumb points from x to y) a positive rotation rotates y->s in the direction of your fingers.We use a periodic energy function like this:

      E(theta) = sum E_n(1 + cos(n*theta - theta0_n))

where n is the periodicity, E_n is the amplitude (kJ/mol) for term n, and theta0_n is the phase offset for term n. The torque term (applied about the v axis) is then

      T(theta) = -[sum -n*E_n*sin(n*theta - theta0_n)]
Note that in the functional forms above the amplitude paramter E_n could be considered a "half amplitude" since the total excursion of the energy and torque functions is 2*E_n. When entering this parameter, be sure you understand whether the data you have represents the (half) amplitude as above (most common) or the full excursion range of the function, in which case you should provide only 1/2 the excursion as the value for DuMM's amplitude parameter.For a custom torsion bond, DuMM will provide the angle theta to the user-written energy and torque routines, expecting a scalar value E(theta) and T(theta) back as above. However, the functional form is specified by the user in that case, and it is the user's responsibility to ensure that the returned torque is the negative gradient of the energy with respect to the angle parameter theta.In either case, DuMM takes care of translating the returned scalar torque into an appropriate mobility torque when possible, or an appropriate set of forces on each of the atoms r-x-y-s such that the desired pure torque is realized as the resultant of the forces.
void SimTK::DuMMForceFieldSubsystem::defineBondTorsion (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3, DuMM::AtomClassIndex class4, int periodicity, Real ampInKJ, Real phaseInDegrees)
 Define bond torsion terms applying to a quadruple of atom classes using the built-in functional form, a combination of sinusoids with different periods (see full discussion in the group header for bond torsion terms).
void SimTK::DuMMForceFieldSubsystem::defineBondTorsion (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3, DuMM::AtomClassIndex class4, int periodicity1, Real amp1InKJ, Real phase1InDegrees, int periodicity2, Real amp2InKJ, Real phase2InDegrees)
 Same as defineBondTorsion() but permits two torsion terms (with different periods) to be specified simultaneously.
void SimTK::DuMMForceFieldSubsystem::defineBondTorsion (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3, DuMM::AtomClassIndex class4, int periodicity1, Real amp1InKJ, Real phase1InDegrees, int periodicity2, Real amp2InKJ, Real phase2InDegrees, int periodicity3, Real amp3InKJ, Real phase3InDegrees)
 Same as defineBondTorsion() but permits three torsion terms (with different periods) to be specified simultaneously.
void SimTK::DuMMForceFieldSubsystem::defineBondTorsion_KA (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3, DuMM::AtomClassIndex class4, int periodicity1, Real amp1InKcal, Real phase1InDegrees)
 Same as defineBondTorsion() but permits takes the amplitude in kcal/mol (but note that this is converted immediately to our MD unit system of kJ/mol).
void SimTK::DuMMForceFieldSubsystem::defineBondTorsion_KA (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3, DuMM::AtomClassIndex class4, int periodicity1, Real amp1InKcal, Real phase1InDegrees, int periodicity2, Real amp2InKcal, Real phase2InDegrees)
 Same as defineBondTorsion_KA() but permits two torsion terms (with different periods) to be specified simultaneously.
void SimTK::DuMMForceFieldSubsystem::defineBondTorsion_KA (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3, DuMM::AtomClassIndex class4, int periodicity1, Real amp1InKcal, Real phase1InDegrees, int periodicity2, Real amp2InKcal, Real phase2InDegrees, int periodicity3, Real amp3InKcal, Real phase3InDegrees)
 Same as defineBondTorsion_KA() but permits three torsion terms (with different periods) to be specified simultaneously.
void SimTK::DuMMForceFieldSubsystem::defineBondTorsion_KA (int class1, int class2, int class3, int class4, int periodicity1, Real amp1InKcal, Real phase1InDegrees)
 Same as defineBondTorsion_KA() but takes integer class arguments for backwards compatibility.
void SimTK::DuMMForceFieldSubsystem::defineBondTorsion_KA (int class1, int class2, int class3, int class4, int periodicity1, Real amp1InKcal, Real phase1InDegrees, int periodicity2, Real amp2InKcal, Real phase2InDegrees)
 Same as defineBondTorsion_KA() but takes integer class arguments for backwards compatibility.
void SimTK::DuMMForceFieldSubsystem::defineBondTorsion_KA (int class1, int class2, int class3, int class4, int periodicity1, Real amp1InKcal, Real phase1InDegrees, int periodicity2, Real amp2InKcal, Real phase2InDegrees, int periodicity3, Real amp3InKcal, Real phase3InDegrees)
 Same as defineBondTorsion_KA() but takes integer class arguments for backwards compatibility.
void SimTK::DuMMForceFieldSubsystem::defineCustomBondTorsion (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3, DuMM::AtomClassIndex class4, DuMM::CustomBondTorsion *bondTorsionTerm)
 Define a custom bond torsion term to be applied to the indicated quadruple of atom classes whenever they are found in a 1-2-3-4 bonded sequence.

Amber-style improper torsions

As with normal torsions, (see defineBondTorsion()), only one term may have a given periodicity.The amplitudes are in kJ/mol. The third atom is the central one to which the other three are bonded; this is not the same in reverse order.
void SimTK::DuMMForceFieldSubsystem::defineAmberImproperTorsion (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3, DuMM::AtomClassIndex class4, int periodicity, Real ampInKJ, Real phaseInDegrees)
 Provide one torsion term in MD units, using kilojoules for amplitude.
void SimTK::DuMMForceFieldSubsystem::defineAmberImproperTorsion (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3, DuMM::AtomClassIndex class4, int periodicity1, Real amp1InKJ, Real phase1InDegrees, int periodicity2, Real amp2InKJ, Real phase2InDegrees)
 Provide two torsion terms in MD units, using kilojoules/mole for amplitude.
void SimTK::DuMMForceFieldSubsystem::defineAmberImproperTorsion (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3, DuMM::AtomClassIndex class4, int periodicity1, Real amp1InKJ, Real phase1InDegrees, int periodicity2, Real amp2InKJ, Real phase2InDegrees, int periodicity3, Real amp3InKJ, Real phase3InDegrees)
 Provide three torsion terms in MD units, using kilojoules/mole for amplitude.
void SimTK::DuMMForceFieldSubsystem::defineAmberImproperTorsion_KA (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3, DuMM::AtomClassIndex class4, int periodicity1, Real amp1InKcal, Real phase1InDegrees)
 Provide one torsion term in KA units, using kilocalories/mole for amplitude.
void SimTK::DuMMForceFieldSubsystem::defineAmberImproperTorsion_KA (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3, DuMM::AtomClassIndex class4, int periodicity1, Real amp1InKcal, Real phase1InDegrees, int periodicity2, Real amp2InKcal, Real phase2InDegrees)
 Provide two torsion terms in KA units, using kilocalories/mole for amplitude.
void SimTK::DuMMForceFieldSubsystem::defineAmberImproperTorsion_KA (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3, DuMM::AtomClassIndex class4, int periodicity1, Real amp1InKcal, Real phase1InDegrees, int periodicity2, Real amp2InKcal, Real phase2InDegrees, int periodicity3, Real amp3InKcal, Real phase3InDegrees)
 Provide three torsion terms in KA units, using kilocalories/mole for amplitude.

GBSA implicit solvation

These methods are used to set GBSA terms.
void SimTK::DuMMForceFieldSubsystem::setSolventDielectric (Real)
void SimTK::DuMMForceFieldSubsystem::setSoluteDielectric (Real)
Real SimTK::DuMMForceFieldSubsystem::getSolventDielectric () const
Real SimTK::DuMMForceFieldSubsystem::getSoluteDielectric () const
void SimTK::DuMMForceFieldSubsystem::setGbsaIncludeAceApproximation (bool)
void SimTK::DuMMForceFieldSubsystem::setGbsaIncludeAceApproximationOn ()
void SimTK::DuMMForceFieldSubsystem::setGbsaIncludeAceApproximationOff ()

Global scale factors

These non-physical parameters can be used to weaken or disable (or magnify) individual force field terms.These are always 1 for correct implementation of any force field; other values are primarily useful for testing the effects of individual terms on results or performance. Set to 0 to disable the corresponding term altogether.
void SimTK::DuMMForceFieldSubsystem::setVdwGlobalScaleFactor (Real)
 scale all van der Waals terms
void SimTK::DuMMForceFieldSubsystem::setCoulombGlobalScaleFactor (Real)
 scale all Coulomb terms
void SimTK::DuMMForceFieldSubsystem::setGbsaGlobalScaleFactor (Real)
 scale all GBSA terms
void SimTK::DuMMForceFieldSubsystem::setBondStretchGlobalScaleFactor (Real)
 scale all built-in bond stretch terms
void SimTK::DuMMForceFieldSubsystem::setBondBendGlobalScaleFactor (Real)
 scale all built-in bond bending terms
void SimTK::DuMMForceFieldSubsystem::setBondTorsionGlobalScaleFactor (Real)
 scale all built-in bond torsion terms
void SimTK::DuMMForceFieldSubsystem::setAmberImproperTorsionGlobalScaleFactor (Real)
 scale all improper torsion terms
void SimTK::DuMMForceFieldSubsystem::setCustomBondStretchGlobalScaleFactor (Real)
 scale all custom bond stretch terms
void SimTK::DuMMForceFieldSubsystem::setCustomBondBendGlobalScaleFactor (Real)
 scale all custom bond bending terms
void SimTK::DuMMForceFieldSubsystem::setCustomBondTorsionGlobalScaleFactor (Real)
 scale all custom bond torsion terms
Real SimTK::DuMMForceFieldSubsystem::getVdwGlobalScaleFactor () const
 get current scale factor for van der Waals terms
Real SimTK::DuMMForceFieldSubsystem::getCoulombGlobalScaleFactor () const
 get current scale factor for Coulomb terms
Real SimTK::DuMMForceFieldSubsystem::getGbsaGlobalScaleFactor () const
 get current scale factor for GBSA terms
Real SimTK::DuMMForceFieldSubsystem::getBondStretchGlobalScaleFactor () const
 get current scale factor for built-in bond stretch terms
Real SimTK::DuMMForceFieldSubsystem::getBondBendGlobalScaleFactor () const
 get current scale factor for built-in bond bending terms
Real SimTK::DuMMForceFieldSubsystem::getBondTorsionGlobalScaleFactor () const
 get current scale factor for built-in bond torsion terms
Real SimTK::DuMMForceFieldSubsystem::getAmberImproperTorsionGlobalScaleFactor () const
 get current scale factor for improper torsion terms
Real SimTK::DuMMForceFieldSubsystem::getCustomBondStretchGlobalScaleFactor () const
 get current scale factor for custom bond stretch terms
Real SimTK::DuMMForceFieldSubsystem::getCustomBondBendGlobalScaleFactor () const
 get current scale factor for custom bond bending terms
Real SimTK::DuMMForceFieldSubsystem::getCustomBondTorsionGlobalScaleFactor () const
 get current scale factor for custom bond torsion terms
void SimTK::DuMMForceFieldSubsystem::setAllGlobalScaleFactors (Real s)
 Set all the global scale factors to the same value.

Computational options

These methods control how DuMM performs its computations.
void SimTK::DuMMForceFieldSubsystem::setTracing (bool)
 For debugging, you can ask DuMM to dump some information to std::clog about its attempts to figure out the best way to use the available hardware.
void SimTK::DuMMForceFieldSubsystem::setUseMultithreadedComputation (bool)
 Enable or disable multithreaded computation (enabled by default).
bool SimTK::DuMMForceFieldSubsystem::getUseMultithreadedComputation () const
 Is multithreaded computation enabled?
void SimTK::DuMMForceFieldSubsystem::setNumThreadsRequested (int)
 Request how many threads we want to use if multithreading is enabled; zero (the default) means let DuMM choose, which will likely be one thread per processor.
int SimTK::DuMMForceFieldSubsystem::getNumThreadsRequested () const
 What was the last value passed to setNumThreadsRequested()? The default is zero meaning DuMM chooses the number of threads.
bool SimTK::DuMMForceFieldSubsystem::isUsingMultithreadedComputation () const
 Is DuMM using the multithreaded code? This could return true even if there is just one thread, if you forced it with setNumThreadsToUse().
int SimTK::DuMMForceFieldSubsystem::getNumThreadsInUse () const
 Find out how many threads DuMM is actually using; will be zero until after realizeTopology().
void SimTK::DuMMForceFieldSubsystem::setUseOpenMMAcceleration (bool)
 This determines whether we use OpenMM GPU acceleration if it is available.
bool SimTK::DuMMForceFieldSubsystem::getUseOpenMMAcceleration () const
 Return the current setting of the flag set by setUseOpenMMAcceleration().
void SimTK::DuMMForceFieldSubsystem::setAllowOpenMMReference (bool)
 This allows us to use OpenMM even if only the Reference platform is available.
bool SimTK::DuMMForceFieldSubsystem::getAllowOpenMMReference () const
 Return the current setting of the flag set by setAllowOpenMMReference().
bool SimTK::DuMMForceFieldSubsystem::isUsingOpenMM () const
 Return true if DuMM is currently using OpenMM for its computations.
std::string SimTK::DuMMForceFieldSubsystem::getOpenMMPlatformInUse () const
 Return the OpenMM Platform currently in use, or the empty string if we're not using OpenMM.

Bookkeeping, debugging, and internal-use-only methods

Hopefully you won't need these.
long long SimTK::DuMMForceFieldSubsystem::getForceEvaluationCount () const
 How many times has the forcefield been evaluated?
void SimTK::DuMMForceFieldSubsystem::dump () const
 Produce an ugly but comprehensive dump of the contents of DuMM's internal data structures, sent to std::cout (stdout).
void SimTK::DuMMForceFieldSubsystem::dumpCForceFieldParameters (std::ostream &os, const String &methodName="loadParameters") const
 Generate C++ code to reproduce forceField parameters presently in memory.
void SimTK::DuMMForceFieldSubsystem::loadTestMoleculeParameters ()
 Load test parameters.

Detailed Description

Once you have built a system of Compounds (molecules) in Molmodel, and modeled them with the degrees of freedom you want, you can apply molecular forces using DuMM, the Molmodel molecular mechanics force field.


Function Documentation

SimTK::DuMM::SimTK_DEFINE_UNIQUE_INDEX_TYPE ( AtomIndex  )
SimTK::DuMM::SimTK_DEFINE_UNIQUE_INDEX_TYPE ( IncludedAtomIndex  )
SimTK::DuMM::SimTK_DEFINE_UNIQUE_INDEX_TYPE ( NonbondAtomIndex  )
SimTK::DuMM::SimTK_DEFINE_UNIQUE_INDEX_TYPE ( BondIndex  )
SimTK::DuMM::SimTK_DEFINE_UNIQUE_INDEX_TYPE ( ClusterIndex  )
SimTK::DuMM::SimTK_DEFINE_UNIQUE_INDEX_TYPE ( AtomClassIndex  )
SimTK::DuMM::SimTK_DEFINE_UNIQUE_INDEX_TYPE ( ChargedAtomTypeIndex  )
virtual Real SimTK::DuMM::CustomBondStretch::calcEnergy ( Real  distance) const [pure virtual, inherited]
virtual Real SimTK::DuMM::CustomBondStretch::calcForce ( Real  distance) const [pure virtual, inherited]
virtual Real SimTK::DuMM::CustomBondBend::calcEnergy ( Real  bendAngle) const [pure virtual, inherited]
virtual Real SimTK::DuMM::CustomBondBend::calcTorque ( Real  bendAngle) const [pure virtual, inherited]
virtual Real SimTK::DuMM::CustomBondTorsion::calcEnergy ( Real  dihedralAngle) const [pure virtual, inherited]
virtual Real SimTK::DuMM::CustomBondTorsion::calcTorque ( Real  dihedralAngle) const [pure virtual, inherited]
SimTK::DuMMForceFieldSubsystem::DuMMForceFieldSubsystem ( ) [inherited]
SimTK::DuMMForceFieldSubsystem::DuMMForceFieldSubsystem ( MolecularMechanicsSystem ) [explicit, inherited]
DuMM::AtomIndex SimTK::DuMMForceFieldSubsystem::addAtom ( DuMM::ChargedAtomTypeIndex  chargedAtomTypeIx) [inherited]

Add a new atom to the model.

The AtomIndex number is returned; you don't get to pick your own. Use the AtomIndex to identify this particular atom subsequently.

DuMM::BondIndex SimTK::DuMMForceFieldSubsystem::addBond ( DuMM::AtomIndex  atom1Ix,
DuMM::AtomIndex  atom2Ix 
) [inherited]

Declare that there is a covalent bond between two atoms.

Note that these are AtomIndex numbers, not AtomClasses or ChargedAtomTypes.

DuMM::AtomIndex SimTK::DuMMForceFieldSubsystem::getBondAtom ( DuMM::BondIndex  bond,
int  which 
) const [inherited]

For a given 1-2 bond, return the atoms which are connected by that bond.

You select one atom at a time by setting parameter which to 0 or 1. 0 will return the lower-numbered AtomIndex, regardless of the order in which the atoms were specified to addBond().

int SimTK::DuMMForceFieldSubsystem::getNumAtoms ( ) const [inherited]

How many atoms are currently in the model?

int SimTK::DuMMForceFieldSubsystem::getNumBonds ( ) const [inherited]

How many 1-2 bonds are currently in the model?

Real SimTK::DuMMForceFieldSubsystem::getAtomMass ( DuMM::AtomIndex  atomIx) const [inherited]

Obtain the mass in Daltons (g/mol) of the atom indicated by the given AtomIndex.

int SimTK::DuMMForceFieldSubsystem::getAtomElement ( DuMM::AtomIndex  atomIx) const [inherited]

Obtain the element (by atomic number) of the atom indicated by the given AtomIndex.

Real SimTK::DuMMForceFieldSubsystem::getAtomRadius ( DuMM::AtomIndex  atomIx) const [inherited]

Obtain the van der Waals radius of the atom indicated by the given AtomIndex.

MobilizedBodyIndex SimTK::DuMMForceFieldSubsystem::getAtomBody ( DuMM::AtomIndex  atomIx) const [inherited]

Obtain the Simbody MobilizedBodyIndex of the rigid body on which a particular atom has been fixed.

An exception will be thrown if this atom is not fixed to any body.

Vec3 SimTK::DuMMForceFieldSubsystem::getAtomStationOnBody ( DuMM::AtomIndex  atomIx) const [inherited]

Obtain the station at which a particular atom is fixed on its body.

An exception will be thrown if this atom is not fixed to any body.

Vec3 SimTK::DuMMForceFieldSubsystem::getAtomStationInCluster ( DuMM::AtomIndex  atomIx,
DuMM::ClusterIndex  clusterIx 
) const [inherited]

Obtain the station at which a particular atom is fixed within a particular Cluster (an atom can be in more than one Cluster).

An exception will be thrown if this atom is not fixed to the cluster.

Vec3 SimTK::DuMMForceFieldSubsystem::getElementDefaultColor ( int  atomicNumber) const [inherited]

For display purposes, return the RGB value of a suggested color for an element given by atomic number.

For example, if the atomicNumber is 8 (Oxygen) the suggested color will be Red (1,0,0).

Vec3 SimTK::DuMMForceFieldSubsystem::getAtomDefaultColor ( DuMM::AtomIndex  atomIx) const [inherited]

For display purposes, return the RGB value of a suggested color with which to display a particular atom.

DuMM::ClusterIndex SimTK::DuMMForceFieldSubsystem::createCluster ( const char *  clusterName) [inherited]

Create an empty Cluster (rigid group of atoms).

The Cluster index number is returned; you don't get to pick your own. The name is just for display; you must use the index to reference the Cluster. Every Cluster has its own reference frame C.

void SimTK::DuMMForceFieldSubsystem::placeAtomInCluster ( DuMM::AtomIndex  atomIx,
DuMM::ClusterIndex  clusterIx,
const Vec3 station 
) [inherited]

Place an existing atom at a particular station in the local frame of a Cluster.

It is fine for an atom to be in more than one Cluster as long as only one of them ends up attached to a body.

void SimTK::DuMMForceFieldSubsystem::placeClusterInCluster ( DuMM::ClusterIndex  childClusterIndex,
DuMM::ClusterIndex  parentClusterIndex,
const Transform X_PC 
) [inherited]

Place a Cluster (the child) in another Cluster (the parent).

The child's local frame C is placed at a given Transform with respect to the parent's frame P. All the atoms in the child Cluster maintain their relative positioning.

MassProperties SimTK::DuMMForceFieldSubsystem::calcClusterMassProperties ( DuMM::ClusterIndex  clusterIx,
const Transform X_BC = Transform() 
) const [inherited]

Calculate the composite mass properties of a Cluster, either in its own reference frame C or in reference frame B with the Cluster placed relative to B using the indicated Transform X_BC.

void SimTK::DuMMForceFieldSubsystem::attachClusterToBody ( DuMM::ClusterIndex  clusterIx,
MobilizedBodyIndex  body,
const Transform X_BC = Transform() 
) [inherited]

Place a Cluster's local frame C at a particular location and orientation with respect to a MobilizedBody's frame B.

All the atoms within the Cluster will become fixed to the body while maintaining the same relative positions as they had in the Cluster.

void SimTK::DuMMForceFieldSubsystem::attachAtomToBody ( DuMM::AtomIndex  atomIx,
MobilizedBodyIndex  body,
const Vec3 station = Vec3(0) 
) [inherited]

Place an individual atom at a particular station on a body without an intervening cluster.

MobilizedBodyIndex SimTK::DuMMForceFieldSubsystem::getClusterBody ( DuMM::ClusterIndex  clusterIx) const [inherited]

Find the MobilizedBody on which a Cluster has been fixed in place.

Transform SimTK::DuMMForceFieldSubsystem::getClusterPlacementOnBody ( DuMM::ClusterIndex  clusterIx) const [inherited]

Find where on its body a Cluster has been placed by returning the transform X_BC giving the orientation and position of Cluster frame C in body frame B.

Transform SimTK::DuMMForceFieldSubsystem::getClusterPlacementInCluster ( DuMM::ClusterIndex  childClusterIndex,
DuMM::ClusterIndex  parentClusterIndex 
) const [inherited]

Find where on parent cluster P a child cluster C has been placed, by returning the transform X_PC.

void SimTK::DuMMForceFieldSubsystem::clearIncludedNonbondAtomList ( ) [inherited]

Clear the list of atoms to be included in nonbonded force calculations. This is usually called prior to including selected groups of atoms.

Any previous calls to include or exclude atoms from the nonbonded list will be forgotten.

void SimTK::DuMMForceFieldSubsystem::clearIncludedBondList ( ) [inherited]

Clear the list of bonds to be included in bonded force calculations. This is usually called prior to including selected bonds.

Any previous calls to include or exclude bonds will be forgotten.

void SimTK::DuMMForceFieldSubsystem::resetIncludedNonbondAtomListToDefault ( ) [inherited]

Restore the included nonbond atoms list to its DuMM-selected default. This will include all the atoms that can affect interbody nonbonded forces, including Coulomb, Van der Waals, and solvent (GBSA) forces.

Any previous calls to include or exclude nonbond atoms will be forgotten.

void SimTK::DuMMForceFieldSubsystem::resetIncludedBondListToDefault ( ) [inherited]

Request that DuMM set the included bonds list to all bonds that can generate forces that affect motion, meaning bonds whose bonded force terms have a non-zero amplitude and that involve atoms from two or more bodies.

Any previous calls to include or exclude bonds will be forgotten.

void SimTK::DuMMForceFieldSubsystem::includeNonbondAtom ( DuMM::AtomIndex  atom) [inherited]

Add the given atom to the included nonbond atom list if it isn't already there. This does not include bonded terms that involve this atom; you have to request that explicitly with includeAllInterbodyBondsForOneAtom().

See also:
includeAllNonbondAtomsForOneBody()
void SimTK::DuMMForceFieldSubsystem::includeAllNonbondAtomsForOneBody ( MobilizedBodyIndex  mobod) [inherited]

Add to the included nonbond atom list all the atoms that are fixed to the given body.

This produces the same result as if includeNonbondAtom() were called on each of this body's atoms individually. Nothing happens if there are no atoms on the given body.

See also:
includeNonbondAtom()
void SimTK::DuMMForceFieldSubsystem::includeAllInterbodyBondsForOneAtom ( DuMM::AtomIndex  atom) [inherited]

Given an atom, include all the bonded force terms that (a) include this atom, and (b) involve a body other than the one to which this atom is fixed. This does not include the atom in nonbonded force calculations; you have to request that explicitly with includeNonbondAtom().

void SimTK::DuMMForceFieldSubsystem::includeAllInterbodyBondsWithBothAtoms ( DuMM::AtomIndex  atom1,
DuMM::AtomIndex  atom2 
) [inherited]

Given two atoms that may appear together in some bonded force term, include all the bonded terms that involve both of them.

Note that even if these atoms are on the same body, some of the bonded terms involving them might still contain atoms which are on other bodies. An alternate signature allows you to specify the two atoms by supplying a Bond instead.

void SimTK::DuMMForceFieldSubsystem::includeAllInterbodyBondsWithBothAtoms ( DuMM::BondIndex  bond) [inherited]

Given a bond (that is, a 1-2 connection between atoms), extract the two atoms from it and then make sure all interbody bond terms that involve both of those atoms (in any position) are included.

See the other signature of this method for details.

void SimTK::DuMMForceFieldSubsystem::includeAllInterbodyBondsForOneBody ( MobilizedBodyIndex  mobod) [inherited]

Include any bonds that are needed to calculate interbody bonded forces that involve any atoms fixed to the given mobilized body. This will involve a few atoms on neighboring bodies.

This does not include the involved atoms in the nonbonded atom list; you have to do that explicitly.

See also:
includeAllInterbodyBondsBetweenTwoBodies()
void SimTK::DuMMForceFieldSubsystem::includeAllInterbodyBondsBetweenTwoBodies ( MobilizedBodyIndex  mobod1,
MobilizedBodyIndex  mobod2 
) [inherited]

Include any bonds that are needed to calculate interbody bonded forces that act between the given pair of mobilized bodies. Any bonded term that involves both an atom from mobod1 and an atom from mobod2 will be included.

Nothing will happen if there is no bond that connects these two bodies.

See also:
includeAllInterbodyBondsForOneBody
int SimTK::DuMMForceFieldSubsystem::getNumIncludedAtoms ( ) const [inherited]

Atoms to be included in force calculations are numbered from 0 to getNumIncludedAtoms()-1 after realizeTopology() has been called.

Use DuMM::IncludedAtomIndex to refer to included atoms by their sequential index numbers.

DuMM::AtomIndex SimTK::DuMMForceFieldSubsystem::getAtomIndexOfIncludedAtom ( DuMM::IncludedAtomIndex  incAtomIx) const [inherited]

Given a DuMM::IncludedAtomIndex, return the corresponding DuMM::AtomIndex. You must already have called realizeTopology().

int SimTK::DuMMForceFieldSubsystem::getNumNonbondAtoms ( ) const [inherited]

A subset of the included atoms are used in nonbond calculations. These are numbered from 0 to getNumNonbondAtoms()-1 after realizeTopology() has been called.

Use the unique integer type DuMM::NonbondAtomIndex to refer to included nonbond atoms by their sequential index numbers. Note that some or all of these nonbond atoms may also be involved in bonded force calculations.

DuMM::IncludedAtomIndex SimTK::DuMMForceFieldSubsystem::getIncludedAtomIndexOfNonbondAtom ( DuMM::NonbondAtomIndex  nonbondAtomIx) const [inherited]

Given a DuMM::NonbondAtomIndex, return the corresponding DuMM::IncludedAtomIndex. You must already have called realizeTopology().

See getAtomIndexOfNonbondAtom() if you want its DuMM::AtomIndex instead.

DuMM::AtomIndex SimTK::DuMMForceFieldSubsystem::getAtomIndexOfNonbondAtom ( DuMM::NonbondAtomIndex  nonbondAtomIx) const [inline, inherited]

Given a DuMM::NonbondAtomIndex, return the corresponding DuMM::AtomIndex. You must already have called realizeTopology().

See getIncludedAtomIndexOfNonbondAtom() if you want its DuMM::IncludedAtomIndex instead.

void SimTK::DuMMForceFieldSubsystem::defineAtomClass ( DuMM::AtomClassIndex  atomClassIx,
const char *  atomClassName,
int  atomicNumber,
int  expectedValence,
Real  vdwRadiusInNm,
Real  vdwWellDepthInKJ 
) [inline, inherited]

Define a new atom class for this force field, for identifying atoms of a particular element, number of bonds, and van der Waals parameters.

You must assign a unique index number and name, either of which can be used to reference this AtomClass subsequently. Typically the number and name are defined by the force field you are using.

Parameters:
[in]atomClassIxA unique integer index identifying this AtomClass.
[in]atomClassNameA unique string identifying this AtomClass.
[in]atomicNumberThe kind of element associated with this AtomClass, that is, the number of protons in the nucleus of atoms that are members of this class.
[in]expectedValenceThe number of bonds expected for atoms in this class; this is a defining feature of the atom class so must be satisfied exactly. For example, this would be 4 for a class representing a tetrahedral carbon atom.
[in]vdwRadiusInNmVan der Waals radius is given as Rmin, the radius (in nm, not A!) at which the energy well minimum is seen (actually it is 1/2 the distance between atom centers for a pair of atoms of this class). This is not Sigma, which we define as the radius (half distance) at which the energy crosses zero, that is, a little closer together than when the energy well is at maximum depth. See warning below.
[in]vdwWellDepthInKJThis is the minimum energy value for the van der Waals force as a positive number in kJ/mol (not kcal/mol!). This is the value used in the definition of the vdwRadiusInNm parameter.
Warning:
It is also common in force fields that van der Waals size terms are given by diameter rather than radius -- be sure you know what convention was used by your source so that you can convert to radius here. To convert for a Lennard-Jones force: Rmin = 2^(1/6)Sigma. The radius is in nm, the well depth in kJ/mol.
See also:
defineAtomClass_KA() to supply values in Angstroms and kcals.
void SimTK::DuMMForceFieldSubsystem::defineAtomClass_KA ( DuMM::AtomClassIndex  atomClassIx,
const char *  atomClassName,
int  element,
int  valence,
Real  vdwRadiusInAng,
Real  vdwWellDepthInKcal 
) [inline, inherited]

Same routine as defineAtomClass() but in Kcal/Angstrom (KA) unit system, that is, radius (still not sigma) is in Angstroms, and well depth in kcal/mol.

void SimTK::DuMMForceFieldSubsystem::defineAtomClass_KA ( int  atomClassIx,
const char *  atomClassName,
int  element,
int  valence,
Real  vdwRadiusInAng,
Real  vdwWellDepthInKcal 
) [inline, inherited]

Obsolete method -- use the other signature.

Same routine as defineAtomClass_KA() but for backwards compatibility and compactness of expression in some contexts, this one accepts an integer for the class index rather than requiring the safer DuMM::AtomClassIndex type.

bool SimTK::DuMMForceFieldSubsystem::hasAtomClass ( DuMM::AtomClassIndex  ) const [inherited]

Check whether an atom class has been defined using this index.

bool SimTK::DuMMForceFieldSubsystem::hasAtomClass ( const String atomClassName) const [inherited]

Check whether an atom class has been defined using this name.

DuMM::AtomClassIndex SimTK::DuMMForceFieldSubsystem::getAtomClassIndex ( const String atomClassName) const [inherited]

Obtain the atom class index corresponding to this atom class name.

DuMM::AtomClassIndex SimTK::DuMMForceFieldSubsystem::getNextUnusedAtomClassIndex ( ) const [inherited]

Obtain an atom class index that is numerically larger than the largest currently-defined atom class index.

DuMM::AtomClassIndex SimTK::DuMMForceFieldSubsystem::getAtomClassIndex ( DuMM::AtomIndex  atomIx) const [inherited]

Get the index number of the atom class associated with this atom.

Real SimTK::DuMMForceFieldSubsystem::getVdwRadius ( DuMM::AtomClassIndex  atomClassIx) const [inherited]

Get the van der Waals radius shared by all atoms that belong to the indicated atom class. See comments for this group for a precise definition of what this means; there are ambiguities so don't assume you already know.

Real SimTK::DuMMForceFieldSubsystem::getVdwWellDepth ( DuMM::AtomClassIndex  atomClassIx) const [inherited]

Get the van der Waals energy well depth shared by all atoms that belong to the indicated atom class. See comments for this group for a precise definition of what this means; there are ambiguities so don't assume you already know.

void SimTK::DuMMForceFieldSubsystem::defineChargedAtomType ( DuMM::ChargedAtomTypeIndex  atomTypeIx,
const char *  atomTypeName,
DuMM::AtomClassIndex  atomClassIx,
Real  partialChargeInE 
) [inline, inherited]

Define a new ChargedAtomType for this force field, for identifying atoms of a particular AtomClass that have a particular partial charge.

You must assign a unique index number and name, either of which can be used to reference this ChargedAtomType subsequently. Typically the number and name are defined by the force field you are using.

Parameters:
[in]atomTypeIxA unique integer index identifying this ChargedAtomType.
[in]atomTypeNameA unique string identifying this ChargedAtomType.
[in]atomClassIxThe index number of a previously defined AtomClass.
[in]partialChargeInEA partial charge in units of e (charge on a proton); this value is the same in both the MD and KA units system.
void SimTK::DuMMForceFieldSubsystem::defineChargedAtomType_KA ( DuMM::ChargedAtomTypeIndex  atomTypeIx,
const char *  atomTypeName,
DuMM::AtomClassIndex  atomClassIx,
Real  partialChargeInE 
) [inline, inherited]

This is an alternate name for defineChargedAtomType() but since partial charge uses the same unit (e) in both the MD and KA unit systems, this is identical to defineChargedAtomType().

void SimTK::DuMMForceFieldSubsystem::defineChargedAtomType_KA ( int  atomTypeIx,
const char *  atomTypeName,
int  atomClassIx,
Real  partialChargeInE 
) [inline, inherited]

Obsolete method -- use the other signature.

Same routine as defineChargedAtomType_KA() but for backwards compatibility and compactness of expression in some contexts, this one accepts an integer for the charged atom type index and the atom class index rather than requiring the safer DuMM::ChargedAtomTypeIndex and DuMM::AtomClassIndex types.

bool SimTK::DuMMForceFieldSubsystem::hasChargedAtomType ( DuMM::ChargedAtomTypeIndex  ) const [inherited]

Check whether a charged atom type has been defined using this index.

bool SimTK::DuMMForceFieldSubsystem::hasChargedAtomType ( const String chargedTypeName) const [inherited]

Check whether a charged atom type has been defined using this name.

DuMM::ChargedAtomTypeIndex SimTK::DuMMForceFieldSubsystem::getChargedAtomTypeIndex ( const String chargedTypeName) const [inherited]

Obtain the charged atom type index corresponding to this charged atom type name.

DuMM::ChargedAtomTypeIndex SimTK::DuMMForceFieldSubsystem::getNextUnusedChargedAtomTypeIndex ( ) const [inherited]

Obtain a charged atom type index that is numerically larger than the largest currently-defined charged atom type index.

const char* SimTK::DuMMForceFieldSubsystem::getVdwMixingRuleName ( VdwMixingRule  ) const [inherited]

Obtain a human-readable name for one of our van der Waals mixing rules.

void SimTK::DuMMForceFieldSubsystem::setVdwMixingRule ( VdwMixingRule  ) [inherited]

Set the van der Waals mixing rule -- our default is Waldman-Hagler.

VdwMixingRule SimTK::DuMMForceFieldSubsystem::getVdwMixingRule ( ) const [inherited]

Get the van der Waals mixing rule currently in effect.

See also:
getVdwMixingRuleName() for a human-readable version.
void SimTK::DuMMForceFieldSubsystem::setVdw12ScaleFactor ( Real  ) [inherited]

default 0

void SimTK::DuMMForceFieldSubsystem::setVdw13ScaleFactor ( Real  ) [inherited]

default 0

void SimTK::DuMMForceFieldSubsystem::setVdw14ScaleFactor ( Real  ) [inherited]

default 1

void SimTK::DuMMForceFieldSubsystem::setVdw15ScaleFactor ( Real  ) [inherited]

default 1

void SimTK::DuMMForceFieldSubsystem::setCoulomb12ScaleFactor ( Real  ) [inherited]

default 0

void SimTK::DuMMForceFieldSubsystem::setCoulomb13ScaleFactor ( Real  ) [inherited]

default 0

void SimTK::DuMMForceFieldSubsystem::setCoulomb14ScaleFactor ( Real  ) [inherited]

default 1

void SimTK::DuMMForceFieldSubsystem::setCoulomb15ScaleFactor ( Real  ) [inherited]

default 1

void SimTK::DuMMForceFieldSubsystem::loadAmber99Parameters ( ) [inherited]

Use Amber99 force field parameters.

This duplicates the Tinker Amber99 parameter set in pre-built code, so you don't need to load the parameters from a file.

void SimTK::DuMMForceFieldSubsystem::populateFromTinkerParameterFile ( std::istream &  ) [inherited]

Load force field parameters from a TINKER format force field parameter file. (Only the Amber99 force field has been tested.)

void SimTK::DuMMForceFieldSubsystem::setBiotypeChargedAtomType ( DuMM::ChargedAtomTypeIndex  chargedAtomTypeIndex,
BiotypeIndex  biotypeIx 
) [inherited]

Associate a Tinker Biotype with a ChargedAtomType in this subsystem.

A Biotype association is required for every Biotype found on atoms in the current System.

Parameters:
chargedAtomTypeIndexPrexisting charged atom type index in this subsystem
biotypeIxPreexisting BiotypeIndex defined in the Biotype class
DuMM::ChargedAtomTypeIndex SimTK::DuMMForceFieldSubsystem::getBiotypeChargedAtomType ( BiotypeIndex  biotypeIx) const [inherited]

Get charged atom type index in this force field associated with a particular Biotype.

std::ostream& SimTK::DuMMForceFieldSubsystem::generateBiotypeChargedAtomTypeSelfCode ( std::ostream &  os) const [inherited]

Generate C++ code from the current contents of this DuMM force field object.

void SimTK::DuMMForceFieldSubsystem::defineBondStretch ( DuMM::AtomClassIndex  class1,
DuMM::AtomClassIndex  class2,
Real  stiffnessInKJperNmSq,
Real  nominalLengthInNm 
) [inherited]

Define a harmonic bond stretch term between two atom classes using the built-in functional form.

You are only allowed to specify a built-in term once per pair of classes, regardless of the order they are listed. It is permitted to redefine the same term as long as the parameters are the same both times, otherwise this method will throw an exception if (class1,class2) or (class2,class1) has already been assigned. Stiffness (energy per length^2) must be given in (kJ/mol)/nm^2, and nominal bond length is in nanometers. CAUTION: Energy is kx^2 using this definition, while force is -2kx; note factor of 2 in force -- conventions vary.

Parameters:
[in]class1,class2The atom class pair to which this term applies; order doesn't matter.
[in]stiffnessInKJperNmSqThe bond stiffness in kilojoules per nm^2.
[in]nominalLengthInNmBond length at which energy and force are zero, in nm.
See also:
defineBondStretch_KA() for kcal-angstrom units
defineCustomBondStretch() to define your own functional form
void SimTK::DuMMForceFieldSubsystem::defineBondStretch_KA ( DuMM::AtomClassIndex  class1,
DuMM::AtomClassIndex  class2,
Real  stiffnessInKcalPerAngSq,
Real  nominalLengthInAng 
) [inline, inherited]

Same as defineBondStretch() except that for convenience this takes stiffness in (kcal/mol)/A^2, and nominal length is in A (angstroms).

Note that these are immediately converted to our standard MD units internally.

See also:
defineBondStretch() for a complete description
Parameters:
[in]class1,class2The atom class pair to which this term applies; order doesn't matter.
[in]stiffnessInKcalPerAngSqThe bond stiffness in kilocalories per Angstrom^2.
[in]nominalLengthInAngBond length at which energy and force are zero, in Angstrom.
void SimTK::DuMMForceFieldSubsystem::defineBondStretch_KA ( int  class1,
int  class2,
Real  stiffnessInKcalPerAngSq,
Real  nominalLengthInAng 
) [inline, inherited]

Same as defineBondStretch_KA() but takes integer class arguments for backwards compatibility.

void SimTK::DuMMForceFieldSubsystem::defineCustomBondStretch ( DuMM::AtomClassIndex  class1,
DuMM::AtomClassIndex  class2,
DuMM::CustomBondStretch bondStretchTerm 
) [inherited]

Define a custom bond stretch term to be applied to the indicated pair of atom classes whenever they are found in a 1-2 bond.

Pass in a pointer to a newly- allocated CustomBondStretch object; the subsystem takes over ownership of the object so it should NOT be deleted by the caller.

Parameters:
[in]class1,class2The atom class pair to which this term applies; order doesn't matter.
[in]bondStretchTermPointer to a heap-allocated object of a concrete type derived from DuMM::CustomBondStretch.
See also:
DuMM::CustomBondStretch
defineBondStretch() to use the built-in harmonic functional form
void SimTK::DuMMForceFieldSubsystem::defineBondBend ( DuMM::AtomClassIndex  class1,
DuMM::AtomClassIndex  class2,
DuMM::AtomClassIndex  class3,
Real  stiffnessInKJPerRadSq,
Real  nominalAngleInDeg 
) [inherited]

Define a harmonic bond bending term applying to a triple of atom classes using the built-in functional form.

You are only allowed to specify a built-in term once per class triple, regardless of the order they are listed. It is permitted to redefine the same term as long as the parameters are the same both times, otherwise this method will throw an exception if (class1,class2, class3) or (class3,class2,class1) has already been assigned. Stiffness (energy per angle^2) must be given in (kJ/mol)/rad^2, and nominal bond angle is in degrees(!). Note that the nominal angle is in degrees while the stiffness is in radians. Odd, I know, but that seems to be how it's done! CAUTION: Energy is ka^2 using this definition, while torque is -2ka; note factor of 2 in torque -- conventions vary.

Parameters:
[in]class1,class2,class3The atom class triple to which this term applies; 1-2-3 or 3-2-1 order doesn't matter.
[in]stiffnessInKJPerRadSqThe bond stiffness in kilojoules per radian^2.
[in]nominalAngleInDegBond angle at which energy and force are zero, in degrees(!).
See also:
defineBondBend_KA() for kcal-angstrom units
defineCustomBondBend() to define your own functional form
void SimTK::DuMMForceFieldSubsystem::defineBondBend_KA ( DuMM::AtomClassIndex  class1,
DuMM::AtomClassIndex  class2,
DuMM::AtomClassIndex  class3,
Real  stiffnessInKcalPerRadSq,
Real  nominalAngleInDeg 
) [inline, inherited]

Same as defineBondBend() except that for convenience this takes stiffness in (kcal/mol)/rad^2 (nominal angle is still in degrees).

Note that the stiffness is immediately converted to our standard MD units internally.

See also:
defineBondBend() for a complete description
Parameters:
[in]class1,class2,class3The atom class triple to which this term applies; 1-2-3 or 3-2-1 order doesn't matter.
[in]stiffnessInKcalPerRadSqThe bond stiffness in kilocalories per radian^2.
[in]nominalAngleInDegBond angle at which energy and force are zero, in degrees(!).
void SimTK::DuMMForceFieldSubsystem::defineBondBend_KA ( int  class1,
int  class2,
int  class3,
Real  stiffnessInKcalPerRadSq,
Real  nominalAngleInDeg 
) [inline, inherited]

Same as defineBondBend_KA() but takes integer class arguments for backwards compatibility.

void SimTK::DuMMForceFieldSubsystem::defineCustomBondBend ( DuMM::AtomClassIndex  class1,
DuMM::AtomClassIndex  class2,
DuMM::AtomClassIndex  class3,
DuMM::CustomBondBend bondBendTerm 
) [inherited]

Define a custom bond bend term to be applied to the indicated triple of atom classes whenever they are found in a 1-2-3 bonded sequence.

Pass in a pointer to a newly- allocated CustomBondBend object; the subsystem takes over ownership of the object so it should NOT be deleted by the caller.

Parameters:
[in]class1,class2,class3The atom class triple to which this term applies; 1-2-3 or 3-2-1 order doesn't matter.
[in]bondBendTermPointer to a heap-allocated object of a concrete type derived from DuMM::CustomBondBend.
See also:
DuMM::CustomBondBend
defineBondBend() to use the built-in harmonic functional form
void SimTK::DuMMForceFieldSubsystem::defineBondTorsion ( DuMM::AtomClassIndex  class1,
DuMM::AtomClassIndex  class2,
DuMM::AtomClassIndex  class3,
DuMM::AtomClassIndex  class4,
int  periodicity,
Real  ampInKJ,
Real  phaseInDegrees 
) [inherited]

Define bond torsion terms applying to a quadruple of atom classes using the built-in functional form, a combination of sinusoids with different periods (see full discussion in the group header for bond torsion terms).

You are only allowed to specify a built-in term of a given period once per class quadruple, whether specified 1-2-3-4 or 4-3-2-1. It is permitted to redefine the same term as long as the parameters are the same both times, otherwise this method will throw an exception if (class1,class2,class3,class4) or (class4,class3,class2, class1) has already been assigned. The amplitude of the energy functions must be given in (kJ/mol), the period is an integer which will multiply the dihedral angle, and the phase (offset) angle is given in degrees (not radians). CAUTION: various conventions exist for specifying the parameters for bond torsion terms. Pay particular attention to the amplitude -- in our convention it is really a half-amplitude since the sinusoids range from -1 to 1 making the full energy and torque "excursion" twice the amplitude.

Parameters:
[in]class1,class2,class3,class4The atom class quadruple to which this term applies; 1-2-3-4 or 4-3-2-1 order doesn't matter.
[in]periodicityThe periodicity (dependence on input angle theta) of the sinusoidal term.
[in]ampInKJThe amplitude (half-range) of the energy function, in kilojoules.
[in]phaseInDegreesAngle at which periodicity*theta yields maximum energy (2*ampInKJ).
See also:
defineBondTorsion_KA() for kcal-angstrom units
defineCustomBondTorsion() to define your own functional form
void SimTK::DuMMForceFieldSubsystem::defineBondTorsion ( DuMM::AtomClassIndex  class1,
DuMM::AtomClassIndex  class2,
DuMM::AtomClassIndex  class3,
DuMM::AtomClassIndex  class4,
int  periodicity1,
Real  amp1InKJ,
Real  phase1InDegrees,
int  periodicity2,
Real  amp2InKJ,
Real  phase2InDegrees 
) [inherited]

Same as defineBondTorsion() but permits two torsion terms (with different periods) to be specified simultaneously.

void SimTK::DuMMForceFieldSubsystem::defineBondTorsion ( DuMM::AtomClassIndex  class1,
DuMM::AtomClassIndex  class2,
DuMM::AtomClassIndex  class3,
DuMM::AtomClassIndex  class4,
int  periodicity1,
Real  amp1InKJ,
Real  phase1InDegrees,
int  periodicity2,
Real  amp2InKJ,
Real  phase2InDegrees,
int  periodicity3,
Real  amp3InKJ,
Real  phase3InDegrees 
) [inherited]

Same as defineBondTorsion() but permits three torsion terms (with different periods) to be specified simultaneously.

void SimTK::DuMMForceFieldSubsystem::defineBondTorsion_KA ( DuMM::AtomClassIndex  class1,
DuMM::AtomClassIndex  class2,
DuMM::AtomClassIndex  class3,
DuMM::AtomClassIndex  class4,
int  periodicity1,
Real  amp1InKcal,
Real  phase1InDegrees 
) [inline, inherited]

Same as defineBondTorsion() but permits takes the amplitude in kcal/mol (but note that this is converted immediately to our MD unit system of kJ/mol).

void SimTK::DuMMForceFieldSubsystem::defineBondTorsion_KA ( DuMM::AtomClassIndex  class1,
DuMM::AtomClassIndex  class2,
DuMM::AtomClassIndex  class3,
DuMM::AtomClassIndex  class4,
int  periodicity1,
Real  amp1InKcal,
Real  phase1InDegrees,
int  periodicity2,
Real  amp2InKcal,
Real  phase2InDegrees 
) [inline, inherited]

Same as defineBondTorsion_KA() but permits two torsion terms (with different periods) to be specified simultaneously.

void SimTK::DuMMForceFieldSubsystem::defineBondTorsion_KA ( DuMM::AtomClassIndex  class1,
DuMM::AtomClassIndex  class2,
DuMM::AtomClassIndex  class3,
DuMM::AtomClassIndex  class4,
int  periodicity1,
Real  amp1InKcal,
Real  phase1InDegrees,
int  periodicity2,
Real  amp2InKcal,
Real  phase2InDegrees,
int  periodicity3,
Real  amp3InKcal,
Real  phase3InDegrees 
) [inline, inherited]

Same as defineBondTorsion_KA() but permits three torsion terms (with different periods) to be specified simultaneously.

void SimTK::DuMMForceFieldSubsystem::defineBondTorsion_KA ( int  class1,
int  class2,
int  class3,
int  class4,
int  periodicity1,
Real  amp1InKcal,
Real  phase1InDegrees 
) [inline, inherited]

Same as defineBondTorsion_KA() but takes integer class arguments for backwards compatibility.

void SimTK::DuMMForceFieldSubsystem::defineBondTorsion_KA ( int  class1,
int  class2,
int  class3,
int  class4,
int  periodicity1,
Real  amp1InKcal,
Real  phase1InDegrees,
int  periodicity2,
Real  amp2InKcal,
Real  phase2InDegrees 
) [inline, inherited]

Same as defineBondTorsion_KA() but takes integer class arguments for backwards compatibility.

void SimTK::DuMMForceFieldSubsystem::defineBondTorsion_KA ( int  class1,
int  class2,
int  class3,
int  class4,
int  periodicity1,
Real  amp1InKcal,
Real  phase1InDegrees,
int  periodicity2,
Real  amp2InKcal,
Real  phase2InDegrees,
int  periodicity3,
Real  amp3InKcal,
Real  phase3InDegrees 
) [inline, inherited]

Same as defineBondTorsion_KA() but takes integer class arguments for backwards compatibility.

void SimTK::DuMMForceFieldSubsystem::defineCustomBondTorsion ( DuMM::AtomClassIndex  class1,
DuMM::AtomClassIndex  class2,
DuMM::AtomClassIndex  class3,
DuMM::AtomClassIndex  class4,
DuMM::CustomBondTorsion bondTorsionTerm 
) [inherited]

Define a custom bond torsion term to be applied to the indicated quadruple of atom classes whenever they are found in a 1-2-3-4 bonded sequence.

Pass in a pointer to a newly-allocated CustomBondTorsion object; the subsystem takes over ownership of the object so it should NOT be deleted by the caller.

Parameters:
[in]class1,class2,class3,class4The atom class quadruple to which this term applies; 1-2-3-4 or 4-3-2-1 order doesn't matter.
[in]bondTorsionTermPointer to a heap-allocated object of a concrete type derived from DuMM::CustomBondTorsion.
See also:
DuMM::CustomBondTorsion
defineBondTorsion() to use the built-in sinusoidal functional form
void SimTK::DuMMForceFieldSubsystem::defineAmberImproperTorsion ( DuMM::AtomClassIndex  class1,
DuMM::AtomClassIndex  class2,
DuMM::AtomClassIndex  class3,
DuMM::AtomClassIndex  class4,
int  periodicity,
Real  ampInKJ,
Real  phaseInDegrees 
) [inherited]

Provide one torsion term in MD units, using kilojoules for amplitude.

void SimTK::DuMMForceFieldSubsystem::defineAmberImproperTorsion ( DuMM::AtomClassIndex  class1,
DuMM::AtomClassIndex  class2,
DuMM::AtomClassIndex  class3,
DuMM::AtomClassIndex  class4,
int  periodicity1,
Real  amp1InKJ,
Real  phase1InDegrees,
int  periodicity2,
Real  amp2InKJ,
Real  phase2InDegrees 
) [inherited]

Provide two torsion terms in MD units, using kilojoules/mole for amplitude.

void SimTK::DuMMForceFieldSubsystem::defineAmberImproperTorsion ( DuMM::AtomClassIndex  class1,
DuMM::AtomClassIndex  class2,
DuMM::AtomClassIndex  class3,
DuMM::AtomClassIndex  class4,
int  periodicity1,
Real  amp1InKJ,
Real  phase1InDegrees,
int  periodicity2,
Real  amp2InKJ,
Real  phase2InDegrees,
int  periodicity3,
Real  amp3InKJ,
Real  phase3InDegrees 
) [inherited]

Provide three torsion terms in MD units, using kilojoules/mole for amplitude.

void SimTK::DuMMForceFieldSubsystem::defineAmberImproperTorsion_KA ( DuMM::AtomClassIndex  class1,
DuMM::AtomClassIndex  class2,
DuMM::AtomClassIndex  class3,
DuMM::AtomClassIndex  class4,
int  periodicity1,
Real  amp1InKcal,
Real  phase1InDegrees 
) [inline, inherited]

Provide one torsion term in KA units, using kilocalories/mole for amplitude.

void SimTK::DuMMForceFieldSubsystem::defineAmberImproperTorsion_KA ( DuMM::AtomClassIndex  class1,
DuMM::AtomClassIndex  class2,
DuMM::AtomClassIndex  class3,
DuMM::AtomClassIndex  class4,
int  periodicity1,
Real  amp1InKcal,
Real  phase1InDegrees,
int  periodicity2,
Real  amp2InKcal,
Real  phase2InDegrees 
) [inline, inherited]

Provide two torsion terms in KA units, using kilocalories/mole for amplitude.

void SimTK::DuMMForceFieldSubsystem::defineAmberImproperTorsion_KA ( DuMM::AtomClassIndex  class1,
DuMM::AtomClassIndex  class2,
DuMM::AtomClassIndex  class3,
DuMM::AtomClassIndex  class4,
int  periodicity1,
Real  amp1InKcal,
Real  phase1InDegrees,
int  periodicity2,
Real  amp2InKcal,
Real  phase2InDegrees,
int  periodicity3,
Real  amp3InKcal,
Real  phase3InDegrees 
) [inline, inherited]

Provide three torsion terms in KA units, using kilocalories/mole for amplitude.

void SimTK::DuMMForceFieldSubsystem::setSolventDielectric ( Real  ) [inherited]
void SimTK::DuMMForceFieldSubsystem::setSoluteDielectric ( Real  ) [inherited]
Real SimTK::DuMMForceFieldSubsystem::getSolventDielectric ( ) const [inherited]
Real SimTK::DuMMForceFieldSubsystem::getSoluteDielectric ( ) const [inherited]
void SimTK::DuMMForceFieldSubsystem::setGbsaIncludeAceApproximation ( bool  ) [inherited]
void SimTK::DuMMForceFieldSubsystem::setGbsaIncludeAceApproximationOn ( ) [inline, inherited]
void SimTK::DuMMForceFieldSubsystem::setGbsaIncludeAceApproximationOff ( ) [inline, inherited]
void SimTK::DuMMForceFieldSubsystem::setVdwGlobalScaleFactor ( Real  ) [inherited]

scale all van der Waals terms

void SimTK::DuMMForceFieldSubsystem::setCoulombGlobalScaleFactor ( Real  ) [inherited]

scale all Coulomb terms

void SimTK::DuMMForceFieldSubsystem::setGbsaGlobalScaleFactor ( Real  ) [inherited]

scale all GBSA terms

void SimTK::DuMMForceFieldSubsystem::setBondStretchGlobalScaleFactor ( Real  ) [inherited]

scale all built-in bond stretch terms

void SimTK::DuMMForceFieldSubsystem::setBondBendGlobalScaleFactor ( Real  ) [inherited]

scale all built-in bond bending terms

void SimTK::DuMMForceFieldSubsystem::setBondTorsionGlobalScaleFactor ( Real  ) [inherited]

scale all built-in bond torsion terms

void SimTK::DuMMForceFieldSubsystem::setAmberImproperTorsionGlobalScaleFactor ( Real  ) [inherited]

scale all improper torsion terms

void SimTK::DuMMForceFieldSubsystem::setCustomBondStretchGlobalScaleFactor ( Real  ) [inherited]

scale all custom bond stretch terms

void SimTK::DuMMForceFieldSubsystem::setCustomBondBendGlobalScaleFactor ( Real  ) [inherited]

scale all custom bond bending terms

void SimTK::DuMMForceFieldSubsystem::setCustomBondTorsionGlobalScaleFactor ( Real  ) [inherited]

scale all custom bond torsion terms

Real SimTK::DuMMForceFieldSubsystem::getVdwGlobalScaleFactor ( ) const [inherited]

get current scale factor for van der Waals terms

Real SimTK::DuMMForceFieldSubsystem::getCoulombGlobalScaleFactor ( ) const [inherited]

get current scale factor for Coulomb terms

Real SimTK::DuMMForceFieldSubsystem::getGbsaGlobalScaleFactor ( ) const [inherited]

get current scale factor for GBSA terms

Real SimTK::DuMMForceFieldSubsystem::getBondStretchGlobalScaleFactor ( ) const [inherited]

get current scale factor for built-in bond stretch terms

Real SimTK::DuMMForceFieldSubsystem::getBondBendGlobalScaleFactor ( ) const [inherited]

get current scale factor for built-in bond bending terms

Real SimTK::DuMMForceFieldSubsystem::getBondTorsionGlobalScaleFactor ( ) const [inherited]

get current scale factor for built-in bond torsion terms

Real SimTK::DuMMForceFieldSubsystem::getAmberImproperTorsionGlobalScaleFactor ( ) const [inherited]

get current scale factor for improper torsion terms

Real SimTK::DuMMForceFieldSubsystem::getCustomBondStretchGlobalScaleFactor ( ) const [inherited]

get current scale factor for custom bond stretch terms

Real SimTK::DuMMForceFieldSubsystem::getCustomBondBendGlobalScaleFactor ( ) const [inherited]

get current scale factor for custom bond bending terms

Real SimTK::DuMMForceFieldSubsystem::getCustomBondTorsionGlobalScaleFactor ( ) const [inherited]

get current scale factor for custom bond torsion terms

void SimTK::DuMMForceFieldSubsystem::setAllGlobalScaleFactors ( Real  s) [inline, inherited]

Set all the global scale factors to the same value.

This is commonly used to turn everything on or off, followed by selectively disabling or enabling individual terms.

void SimTK::DuMMForceFieldSubsystem::setTracing ( bool  ) [inherited]

For debugging, you can ask DuMM to dump some information to std::clog about its attempts to figure out the best way to use the available hardware.

void SimTK::DuMMForceFieldSubsystem::setUseMultithreadedComputation ( bool  ) [inherited]

Enable or disable multithreaded computation (enabled by default).

bool SimTK::DuMMForceFieldSubsystem::getUseMultithreadedComputation ( ) const [inherited]

Is multithreaded computation enabled?

void SimTK::DuMMForceFieldSubsystem::setNumThreadsRequested ( int  ) [inherited]

Request how many threads we want to use if multithreading is enabled; zero (the default) means let DuMM choose, which will likely be one thread per processor.

DuMM will choose to use single threaded code if there is only one processor, but requesting 1 thread here explicitly will cause it to use multithreaded code but with a single thread (useful for debugging sometimes).

int SimTK::DuMMForceFieldSubsystem::getNumThreadsRequested ( ) const [inherited]

What was the last value passed to setNumThreadsRequested()? The default is zero meaning DuMM chooses the number of threads.

This doesn't tell you how many threads are actually in use -- see getNumThreadsInUse() for that.

bool SimTK::DuMMForceFieldSubsystem::isUsingMultithreadedComputation ( ) const [inherited]

Is DuMM using the multithreaded code? This could return true even if there is just one thread, if you forced it with setNumThreadsToUse().

int SimTK::DuMMForceFieldSubsystem::getNumThreadsInUse ( ) const [inherited]

Find out how many threads DuMM is actually using; will be zero until after realizeTopology().

void SimTK::DuMMForceFieldSubsystem::setUseOpenMMAcceleration ( bool  ) [inherited]

This determines whether we use OpenMM GPU acceleration if it is available.

By default, this is set false because OpenMM will compute only to single precision. Note that even if you set this flag, we won't use OpenMM unless (a) it is installed correctly on your machine, and (b) it can run with GPU acceleration. If you want to allow use of the non-accelerated Reference Platform provided by OpenMM, use setAllowOpenMMReference().

bool SimTK::DuMMForceFieldSubsystem::getUseOpenMMAcceleration ( ) const [inherited]

Return the current setting of the flag set by setUseOpenMMAcceleration().

void SimTK::DuMMForceFieldSubsystem::setAllowOpenMMReference ( bool  ) [inherited]

This allows us to use OpenMM even if only the Reference platform is available.

This is for testing/debugging; one should never use the Reference platform in production since it will likely be slower than the CPU implementation.

bool SimTK::DuMMForceFieldSubsystem::getAllowOpenMMReference ( ) const [inherited]

Return the current setting of the flag set by setAllowOpenMMReference().

bool SimTK::DuMMForceFieldSubsystem::isUsingOpenMM ( ) const [inherited]

Return true if DuMM is currently using OpenMM for its computations.

std::string SimTK::DuMMForceFieldSubsystem::getOpenMMPlatformInUse ( ) const [inherited]

Return the OpenMM Platform currently in use, or the empty string if we're not using OpenMM.

long long SimTK::DuMMForceFieldSubsystem::getForceEvaluationCount ( ) const [inherited]

How many times has the forcefield been evaluated?

void SimTK::DuMMForceFieldSubsystem::dump ( ) const [inherited]

Produce an ugly but comprehensive dump of the contents of DuMM's internal data structures, sent to std::cout (stdout).

void SimTK::DuMMForceFieldSubsystem::dumpCForceFieldParameters ( std::ostream &  os,
const String methodName = "loadParameters" 
) const [inherited]

Generate C++ code to reproduce forceField parameters presently in memory.

void SimTK::DuMMForceFieldSubsystem::loadTestMoleculeParameters ( ) [inherited]

Load test parameters.

void SimTK::DuMMForceFieldSubsystem::setTraceOpenMM ( bool  shouldTrace) [inline, inherited]

OBSOLETE NAME -- just use setTracing().

void SimTK::DuMMForceFieldSubsystem::defineIncompleteAtomClass ( DuMM::AtomClassIndex  classIx,
const char *  name,
int  elementNumber,
int  valence 
) [protected, inherited]
void SimTK::DuMMForceFieldSubsystem::defineIncompleteAtomClass_KA ( DuMM::AtomClassIndex  classIx,
const char *  name,
int  elementNumber,
int  valence 
) [inline, protected, inherited]
void SimTK::DuMMForceFieldSubsystem::setAtomClassVdwParameters ( DuMM::AtomClassIndex  atomClassIx,
Real  vdwRadiusInNm,
Real  vdwWellDepthInKJPerMol 
) [protected, inherited]
void SimTK::DuMMForceFieldSubsystem::setAtomClassVdwParameters_KA ( DuMM::AtomClassIndex  atomClassIx,
Real  radiusInAng,
Real  wellDepthInKcal 
) [inline, protected, inherited]
bool SimTK::DuMMForceFieldSubsystem::isValidAtomClass ( DuMM::AtomClassIndex  ) const [protected, inherited]
void SimTK::DuMMForceFieldSubsystem::defineIncompleteChargedAtomType ( DuMM::ChargedAtomTypeIndex  typeIx,
const char *  name,
DuMM::AtomClassIndex  classIx 
) [protected, inherited]
void SimTK::DuMMForceFieldSubsystem::defineIncompleteChargedAtomType_KA ( DuMM::ChargedAtomTypeIndex  typeIx,
const char *  name,
DuMM::AtomClassIndex  classIx 
) [inline, protected, inherited]
void SimTK::DuMMForceFieldSubsystem::setChargedAtomTypeCharge ( DuMM::ChargedAtomTypeIndex  ,
Real  charge 
) [protected, inherited]
void SimTK::DuMMForceFieldSubsystem::setChargedAtomTypeCharge_KA ( DuMM::ChargedAtomTypeIndex  chargedAtomTypeIx,
Real  charge 
) [inline, protected, inherited]

Friends

friend class MolecularMechanicsSystem [friend, inherited]
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines