Compound Class Reference

The base class for atoms, molecules, and chemical groups. More...

#include <Compound.h>

Inheritance diagram for Compound:
PIMPLHandle< Compound, CompoundRep > AlcoholOHGroup AromaticSixMemberedCHGroup BiopolymerResidue CarboxylateGroup SingleAtom LigandDroplet MethyleneGroup MethylGroup Molecule P12 PrimaryAmineGroup PurineBaseCore PyrimidineBaseCore TwoNMethylGuanidineGroup Water WaterDroplet

List of all members.

Classes

class  SingleAtom
 Base class for single-atom Compound building blocks. More...

Public Types

enum  MatchStratagem { Match_Exact, Match_Idealized, Match_TopologyOnly }

Public Member Functions

CompoundsetTopLevelTransform (const Transform &transform)
 Set the orientation and location of this Compound.
const TransformgetTopLevelTransform () const
int getNAtoms () const
size_t getNBondCenters () const
size_t getNBondCenters (Compound::AtomIndex atomIndex) const
size_t getNBonds () const
Molecular topology methods



 Compound ()
 default Compound constructor.
 Compound (const Name &name)
 Construct a Compound with a type name.
int getNumAtoms () const
size_t getNumBondCenters () const
size_t getNumBondCenters (Compound::AtomIndex atomIndex) const
size_t getNumBonds () const
bool hasAtom (const AtomPathName &name) const
bool hasBondCenter (const BondCenterPathName &bondCenter) const
Compound::AtomIndex getBondAtomIndex (Compound::BondIndex bondIndex, int which) const
CompoundsetBaseAtom (const Compound::AtomName &name, const Element &element, const Transform &location=Vec3(0))
 Add the first atom unconnected to anything else (yet).
CompoundsetBaseAtom (const Compound::AtomName &name, const Biotype &biotype, const Transform &location=Vec3(0))
 Add the first atom unconnected to anything else (yet).
CompoundsetBaseAtom (const Compound::SingleAtom &c, const Transform &=Transform())
 Add a first subcompound containing exactly one atom, so the Compound::AtomName can be reused for the Compound::Name.
CompoundsetBaseCompound (const Compound::Name &n, const Compound &c, const Transform &location=Transform())
 Add a first subcompound without attaching it to anything.
CompoundbondAtom (const Compound::SingleAtom &c, const BondCenterPathName &parentBondName, mdunits::Length distance, Angle dihedral=0, BondMobility::Mobility=BondMobility::Default)
 Add a subcompound containing exactly one atom, so the Compound::AtomName can be reused for the Compound::Name.
CompoundbondAtom (const Compound::SingleAtom &c, const BondCenterPathName &parentBondName)
 Bond atom using default bond length and dihedral angle.
CompoundbondCompound (const Compound::Name &n, const Compound &c, const BondCenterPathName &parentBondName, mdunits::Length distance, Angle dihedral=180 *Deg2Rad, BondMobility::Mobility mobility=BondMobility::Default)
 Add a subcompound attached by its inboard bond to an existing bond center.
CompoundbondCompound (const Compound::Name &n, const Compound &c, const BondCenterPathName &parentBondName)
 Add a subcompound attached by its inboard bond to an existing bond center.
CompoundbondCompound (const Compound::Name &n, const Compound &c, const BondCenterPathName &parentBondName, BondMobility::Mobility mobility)
 Add a subcompound attached by its inboard bond to an existing bond center.
CompoundsetInboardBondCenter (const Compound::BondCenterName &centerName)
 setInboardBondCenter assigns special status to a bond center.
CompoundconvertInboardBondCenterToOutboard ()
 Make so that this compound can no longer be a child to the geometry of another compound Raises an error if the inboard bond center is already bonded.
CompoundaddFirstBondCenter (const Compound::BondCenterName &centerName, const Compound::AtomPathName &atomName)
 Assign the very first bond center on a particular atom.
CompoundaddSecondBondCenter (const Compound::BondCenterName &centerName, const Compound::AtomName &atomName, Angle bondAngle1)
 Place a second bond center on an atom, placed a particular angle away from the first bond center.
CompoundaddPlanarBondCenter (const Compound::BondCenterName &centerName, const Compound::AtomName &atomName, Angle bondAngle1, Angle bondAngle2)
 Places a third or later bond center on an atom, in the same plane as the first two bond centers.
CompoundaddRightHandedBondCenter (const Compound::BondCenterName &centerName, const Compound::AtomName &atomName, Angle bondAngle1, Angle bondAngle2)
 Place a third or later bond center on an atom, defined by two angular displacements from the first and second bond centers, respectively.
CompoundaddLeftHandedBondCenter (const Compound::BondCenterName &centerName, const Compound::AtomName &atomName, Angle bondAngle1, Angle bondAngle2)
 Place a third or later bond center on an atom, defined by two angular displacements from the first and second bond centers, respectively.
CompoundaddRingClosingBond (const Compound::BondCenterPathName &centerName1, const Compound::BondCenterPathName &centerName2, mdunits::Length bondLength, Angle dihedral=180 *Deg2Rad, BondMobility::Mobility mobility=BondMobility::Default)
 Adds a covalent bond that is not part of the main tree-topology of bonds in a Compound.
CompoundaddRingClosingBond (const Compound::BondCenterPathName &centerName1, const Compound::BondCenterPathName &centerName2)
 Adds a covalent bond that is not part of the main tree-topology of bonds in a molecule.
const CompoundgetSubcompound (const Compound::Name &subcompoundName) const
CompoundupdSubcompound (const Compound::Name &subcompoundName)
Compound nomenclature methods



CompoundsetCompoundName (const Name &)
 Name a type of Compound.
const NamegetCompoundName () const
CompoundaddCompoundSynonym (const Name &)
 Add an additional name for this class of compound.
const AtomPathName getAtomName (Compound::AtomIndex) const
 Returns the most recently assigned name, if any, given to an atom in this compound.
const ElementgetAtomElement (Compound::AtomIndex) const
const ElementgetAtomElement (const Compound::AtomName &) const
 Name a type of Compound.
CompoundnameAtom (const Compound::AtomName &newName, const AtomPathName &oldName)
CompoundnameAtom (const Compound::AtomName &newName, const AtomPathName &oldName, BiotypeIndex biotype)
CompounddefineDihedralAngle (const Compound::DihedralName &angleName, const Compound::AtomPathName &atom1, const Compound::AtomPathName &atom2, const Compound::AtomPathName &atom3, const Compound::AtomPathName &atom4, Angle offset=0 *Deg2Rad)
 Define a named dihedral angle using four atoms.
CompounddefineDihedralAngle (const Compound::DihedralName &angleName, const Compound::BondCenterPathName &bond1, const Compound::BondCenterPathName &bond2, Angle offset=0 *Deg2Rad)
 Define a named dihedral in terms of two bond centers.
CompoundsetPdbResidueNumber (int resNum)
 Set value to populate "residue number" field for PDB file output.
int getPdbResidueNumber () const
CompoundsetPdbResidueName (const String &resName)
 Set string to populate "residue type" field for PDB file output.
const StringgetPdbResidueName () const
CompoundsetPdbChainId (char chainId)
 Set character to populate "chain id" field for PDB file output.
char getPdbChainId () const
CompoundnameBondCenter (const Compound::BondCenterName &newName, const BondCenterPathName &oldName)
 Define a local name for a particular BondCenter in this Compound.
CompoundinheritAtomNames (const Compound::Name &)
 Convenience method to locally import all of the atom names that a particular subcompound uses.
CompoundinheritCompoundSynonyms (const Compound &otherCompound)
 Name a type of Compound.
CompoundinheritBondCenterNames (const Compound::Name &)
 Convenience method to locally import all of the BondCenter names that a particular subcompound uses.
Compound::AtomIndex getAtomIndex (const Compound::AtomPathName &) const
Compound simulation methods



void setMultibodySystem (MultibodySystem &system)
 Add this Compound, including all of its subcompounds, to a CompoundSystem, for simulation etc.
CompoundsetBondMobility (BondMobility::Mobility mobility, const AtomPathName &atom1, const AtomPathName &atom2)
 Override the default rotatability of a bond.
CompoundsetBondMobility (BondMobility::Mobility mobility, Compound::BondIndex bondIndex)
 Override the default rotatability of a bond.
CompoundsetCompoundBondMobility (BondMobility::Mobility mobility)
 Set BondMobility for every bond in the Compound.
MobilizedBodyIndex getAtomMobilizedBodyIndex (Compound::AtomIndex) const
 get the simbody MobilizedBody to which a particular atom is attached
Vec3 getAtomLocationInMobilizedBodyFrame (Compound::AtomIndex) const
 get the location of an Atom in the frame of the simbody MobilizedBody to which the Atom is attached
CompoundsetBiotypeIndex (const Compound::AtomPathName &atomName, BiotypeIndex biotype)
 define the Biotype for an Atom in this Compound
CompoundsetAtomBiotype (const Compound::AtomPathName &atomName, const String &biotypeResidueName, const String &biotypeAtomName, SimTK::Ordinality::Residue ordinality=SimTK::Ordinality::Any)
 Add this Compound, including all of its subcompounds, to a CompoundSystem, for simulation etc.
BiotypeIndex getAtomBiotypeIndex (Compound::AtomIndex) const
DuMM::AtomIndex getDuMMAtomIndex (Compound::AtomIndex) const
 get DuMMForceFieldSubsystem atom index corresponding to an atom in this Compound

Protected Member Functions

void setDuMMAtomIndex (Compound::AtomIndex, DuMM::AtomIndex)
 Stores relationship between a Compound Atom and an Atom defined in a DuMMForcefieldSubsystem.
 Compound (CompoundRep *ip)

Friends

class CompoundSystem

Molecular geometry methods



enum  PlanarBondMatchingPolicy { KeepPlanarBonds, DistortPlanarBonds }
 

Compute atom location in local Compound frame.

More...
Vec3 calcDefaultAtomLocationInGroundFrame (const AtomPathName &atomName) const
 Compute atom location in local Compound frame.
Vec3 calcDefaultAtomLocationInCompoundFrame (const AtomPathName &atomName) const
 Compute atom location in local Compound frame.
CompoundsetDefaultInboardBondLength (mdunits::Length)
 Sets default (initial) bond length of current or future Bond using this Compound's inboard BondCenter.
CompoundsetDefaultInboardDihedralAngle (Angle)
 Stores a default (initial) dihedral angle in the inboard BondCenter of this Compound.
CompoundsetDefaultBondAngle (Angle angle, const AtomPathName &atom1, const AtomPathName &atom2, const AtomPathName &atom3)
 Sets a default(initial) bond angle defined by three atoms.
CompoundsetDefaultBondLength (mdunits::Length length, const AtomPathName &atom1, const AtomPathName &atom2)
 Sets a default (inital) bond length defined by two atoms.
CompoundsetDefaultDihedralAngle (const DihedralName &dihedralName, Angle angleInRadians)
 Sets a default (initial) dihedral angle of a previously named dihedral.
CompoundsetDefaultDihedralAngle (Angle angle, Compound::AtomIndex atom1, Compound::AtomIndex atom2, Compound::AtomIndex atom3, Compound::AtomIndex atom4)
 Compute atom location in local Compound frame.
CompoundsetDefaultDihedralAngle (Angle angle, const Compound::AtomName &atom1, const Compound::AtomName &atom2, const Compound::AtomName &atom3, const Compound::AtomName &atom4)
 Compute atom location in local Compound frame.
Angle calcDefaultDihedralAngle (const DihedralName &dihedralName) const
 Computes default (initial) dihedral angle of a previously named dihedral.
CompoundsetDihedralAngle (State &state, const DihedralName &dihedralName, Angle)
 Sets dynamic dihedral angle of a previously named dihedral.
Angle calcDihedralAngle (const State &state, const DihedralName &dihedralName) const
 Computes dynamic dihedral angle of a previously named dihedral.
Transform calcDefaultAtomFrameInCompoundFrame (Compound::AtomIndex) const
Vec3 calcAtomLocationInGroundFrame (const State &state, Compound::AtomIndex atomId) const
Vec3 calcAtomLocationInCompoundFrame (const State &state, Compound::AtomIndex atomId) const
 Compute atom location in local Compound frame.
Vec3 calcAtomVelocityInGroundFrame (const State &state, Compound::AtomIndex atomId) const
Vec3 calcAtomAccelerationInGroundFrame (const State &state, Compound::AtomIndex atomId) const
Transform getSubcompoundFrameInParentFrame (const Compound::Name &subcompoundName) const
virtual AtomTargetLocations createAtomTargets (const class PdbStructure &targetStructure) const
 Create a mapping between this Compound and atom locations in a PdbStructure.
virtual AtomTargetLocations createAtomTargets (const class PdbChain &targetChain) const
 Compute atom location in local Compound frame.
virtual AtomTargetLocations createAtomTargets (const class PdbResidue &targetResidue) const
 Compute atom location in local Compound frame.
CompoundmatchDefaultAtomChirality (const AtomTargetLocations &atomTargets, Angle planarityTolerance=90.0 *Deg2Rad)
 Adjust stereochemistry about chiral atoms to match that seen in a set of atomic locations.
CompoundmatchDefaultBondLengths (const AtomTargetLocations &atomTargets)
 Compute atom location in local Compound frame.
CompoundmatchDefaultBondAngles (const AtomTargetLocations &atomTargets)
 Compute atom location in local Compound frame.
CompoundmatchDefaultDihedralAngles (const AtomTargetLocations &atomTargets, PlanarBondMatchingPolicy policy=KeepPlanarBonds)
 Compute atom location in local Compound frame.
CompoundmatchDefaultTopLevelTransform (const AtomTargetLocations &atomTargets)
 Compute atom location in local Compound frame.
TransformAndResidual getTransformAndResidual (const Compound::AtomTargetLocations &atomTargets) const
 Compute atom location in local Compound frame.
CompoundmatchDefaultConfiguration (const AtomTargetLocations &atomTargets, MatchStratagem matchStratagem=Match_Exact)
 Compute atom location in local Compound frame.
CompoundfitDefaultConfiguration (const AtomTargetLocations &atomTargets, SimTK::Real targetRms)
 Optimize adjustable degrees of freedom to best match atom targets.
const CompoundpopulateDefaultPdbChain (class PdbChain &, int &defaultNextResidueNumber, const Transform &transform=Transform()) const
 Write current default(initial) Compound configuration into a PdbChain object.
const CompoundpopulatePdbChain (const State &state, class PdbChain &, int &defaultNextResidueNumber, const Transform &transform=Transform()) const
 Write dynamic Compound configuration into a PdbChain object.
std::ostream & writeDefaultPdb (std::ostream &os, const Transform &transform=Transform()) const
 Write the default (initial) configuration in Protein Data Bank (PDB) format.
void writeDefaultPdb (const char *outFileName, const Transform &transform) const
 /brief This polymorphism takes a char* file name rather than ostream, to save the user a couple of lines of code.
std::ostream & writeDefaultPdb (std::ostream &os, int &nextAtomSerialNumber, const Transform &transform=Transform()) const
 Write the default (initial) configuration in Protein Data Bank (PDB) format.
std::ostream & writePdb (const State &state, std::ostream &os, const Transform &transform=Transform()) const
 Write the dynamic Compound configuration in Protein Data Bank (PDB) format.
std::ostream & writePdb (const State &state, std::ostream &os, int &nextAtomSerialNumber, const Transform &transform=Transform()) const
 Write the dynamic Compound configuration in Protein Data Bank (PDB) format.

Data types that identify subcomponents of molecular compounds



typedef String Name
 Type for name of a particular subcompound within a Compound.
typedef String AtomName
 Type for name of a particular atom within a Compound.
typedef String BondCenterName
 Type for name of a particular bond center within a Compound or atom.
typedef String DihedralName
 Type for name of a particular named dihedral angle within a Compound.
typedef String BondCenterPathName
 Type for name of a particular bond center within a Compound, possibly including subcompound indirection.
typedef String AtomPathName
 Type for name of a particular atom within a Compound, possibly including subcompound path.
typedef std::map< AtomIndex, Vec3AtomTargetLocations
 Type for set of target atom locations to be used for structure matching.
 SimTK_DEFINE_UNIQUE_LOCAL_INDEX_TYPE (Compound, AtomIndex)
 Compound::Index type is an integer index into subcompounds of a Compound.
 SimTK_DEFINE_UNIQUE_LOCAL_INDEX_TYPE (Compound, LocalAtomIndex)
 Compound::LocalAtomIndex type is an integer index into atoms directly attached to a Compound, that is to say that the atom does not belong to any subcompounds of the Compound.
 SimTK_DEFINE_UNIQUE_LOCAL_INDEX_TYPE (Compound, BondCenterIndex)
 Compound::BondCenterIndex type is an integer index into BondCenters of a Compound.
 SimTK_DEFINE_UNIQUE_LOCAL_INDEX_TYPE (Compound, BondIndex)
 Compound::BondIndex type is an integer index into Bonds of a Compound.

Detailed Description

The base class for atoms, molecules, and chemical groups.

The Compound class is the base for all molecular entities in the SimTK Molmodel API.


Member Typedef Documentation

typedef String AtomName

Type for name of a particular atom within a Compound.

Compound::AtomName is not intrinsic to the atom itself, but rather the relationship between an atom and a particular Compound. A particular atom may have different names in a compound and its subcompounds. For example, a particular atom might have the name "methyl H1" in a Compound, and have the name "H1" within a subcompound of that Compound. Unlike AtomPathNames, AtomNames must not be qualified by subcompound indirection. e.g. "H1" is a valid AtomName and a valid AtomPathName, while "methyl/H1" (because of the "/" character) is not a valid AtomName, but is a valid AtomPathName.

Type for name of a particular atom within a Compound, possibly including subcompound path.

Compound::AtomPathName is not intrinsic to the atom itself, but rather the relationship between an atom and a particular Compound. A particular atom may have different names in a compound and its subcompounds. For example, a particular atom might have the name "methyl H1" in a Compound, and have the name "H1" within a subcompound of that Compound. AtomPathNames may be qualified by subcompound indirection, in contrast to AtomNames. e.g. "H1" is a valid AtomName and a valid AtomPathName, while "methyl/H1" (because of the "/" character) is not a valid AtomName, but is a valid AtomPathName.

typedef std::map<AtomIndex, Vec3> AtomTargetLocations

Type for set of target atom locations to be used for structure matching.

Type for name of a particular bond center within a Compound or atom.

Compound::BondCenterName is not intrinsic to the bond center itself, but rather the relationship between a bond center and a particular atom or Compound. BondCenterNames are local to a particular compound. In contrast, BondCenterPathNames may be qualified by subcompound indirection. e.g. "bond1" is a valid BondCenterName and a valid BondCenterPathName, while "methyl gamma/bond1" is not a valid BondCenterName, but is a valid BondCenterPathName.

Type for name of a particular bond center within a Compound, possibly including subcompound indirection.

Compound::BondCenterPathName is not intrinsic to the bond center itself, but rather the relationship between a bond center and a particular Compound. BondCenterPathNames may be qualified by subcompound indirection, in contrast to BondCenterNames. e.g. "bond1" is a valid BondCenterName and a valid BondCenterPathName, while "methyl gamma/bond1" is not a valid BondCenterName, but is a valid BondCenterPathName.

Type for name of a particular named dihedral angle within a Compound.

Compound::DihedralName is not intrinsic to the dihedral angle itself, but rather the relationship between a dihedral angle and a particular Compound.

typedef String Name

Type for name of a particular subcompound within a Compound.

Compound::Name is used as an index to subcompounds of a parent compound, e.g. "methyl group gamma"


Member Enumeration Documentation

Enumerator:
Match_Exact 
Match_Idealized 
Match_TopologyOnly 

Compute atom location in local Compound frame.

Returns:
default (initial) location of atom in orthogonal nanometers with respect to Ground frame
Enumerator:
KeepPlanarBonds 
DistortPlanarBonds 

Constructor & Destructor Documentation

Compound (  ) 

default Compound constructor.

Create an empty compound object representing a simulatable molecular structure.

Compound ( const Name name  )  [explicit]

Construct a Compound with a type name.

Create an empty compound object representing a simulatable molecular structure. The name is intended to represent the class of molecule being represented. e.g. "ethane", not "ethane number 3"

Parameters:
name name of the compound type, not the compound instance
Compound ( CompoundRep *  ip  )  [explicit, protected]

Member Function Documentation

Compound& addCompoundSynonym ( const Name  ) 

Add an additional name for this class of compound.

Adding synonyms can be useful for automatically resolving biotypes, if the compound name in the force field parameter file differs from that used to define the compound. See also setCompoundName()

Returns:
a reference to this compound object

Referenced by Adenylate::Adenylate(), Aspartate::Aspartate(), Cysteine::Cysteine(), Cytidylate::Cytidylate(), Glutamate::Glutamate(), Guanylate::Guanylate(), Histidine::Histidine(), NMethylAmideResidue::NMethylAmideResidue(), TwoNMethylGuanineBaseGroup::TwoNMethylGuanineBaseGroup(), TwoNMethylGuanylate::TwoNMethylGuanylate(), and Uridylate::Uridylate().

Compound& addFirstBondCenter ( const Compound::BondCenterName centerName,
const Compound::AtomPathName atomName 
)

Assign the very first bond center on a particular atom.

Returns:
a reference to this Compound object
a reference to this compound object
Parameters:
centerName name for the new BondCenter, must be unique within the Compound
atomName name of the Atom to attach the bond center to

Referenced by BivalentAtom::BivalentAtom(), QuadrivalentAtom::QuadrivalentAtom(), TrivalentAtom::TrivalentAtom(), and UnivalentAtom::UnivalentAtom().

Compound& addLeftHandedBondCenter ( const Compound::BondCenterName centerName,
const Compound::AtomName atomName,
Angle  bondAngle1,
Angle  bondAngle2 
)

Place a third or later bond center on an atom, defined by two angular displacements from the first and second bond centers, respectively.

Placed the new bond center such that the directions of the first, second, and new bond centers form a left-handed geometry.

Returns:
a reference to this compound object
Parameters:
centerName name for the new BondCenter, must be unique within the Compound
atomName name of the Atom to attach the bond center to
bondAngle1 default bond bend angle relating the first bond center to the new bond center
bondAngle2 default bond bend angle relating the second bond center to the new bond center

Referenced by QuadrivalentAtom::QuadrivalentAtom().

Compound& addPlanarBondCenter ( const Compound::BondCenterName centerName,
const Compound::AtomName atomName,
Angle  bondAngle1,
Angle  bondAngle2 
)

Places a third or later bond center on an atom, in the same plane as the first two bond centers.

The angle represents the angular distance from the first bond center on the atom. The second angle is used to determine whether the new bond center is on the same side of the plane as the second bond center.

Returns:
a reference to this compound object
Parameters:
centerName name for the new BondCenter, must be unique within the Compound
atomName name of the Atom to attach the bond center to
bondAngle1 default bond bend angle relating the first bond center to the new bond center
bondAngle2 APPROXIMATE default bond bend angle relating the second bond center to the new bond center

Referenced by TrivalentAtom::TrivalentAtom().

Compound& addRightHandedBondCenter ( const Compound::BondCenterName centerName,
const Compound::AtomName atomName,
Angle  bondAngle1,
Angle  bondAngle2 
)

Place a third or later bond center on an atom, defined by two angular displacements from the first and second bond centers, respectively.

Placed the new bond center such that the directions of the first, second, and new bond centers form a right-handed geometry.

Returns:
a reference to this compound object
Parameters:
centerName name for the new BondCenter, must be unique within the Compound
atomName name of the Atom to attach the bond center to
bondAngle1 default bond bend angle relating the first bond center to the new bond center
bondAngle2 default bond bend angle relating the second bond center to the new bond center

Referenced by QuadrivalentAtom::QuadrivalentAtom().

Compound& addRingClosingBond ( const Compound::BondCenterPathName centerName1,
const Compound::BondCenterPathName centerName2 
)

Adds a covalent bond that is not part of the main tree-topology of bonds in a molecule.

Such bonds are necessary when there are ring structures. These bonds will not correspond to MobilizedBodies when the molecule is realized as a simbody multibody model, and thus the bond lengths, bond angles, and dihedral angles for a ring closing bond might not be enforced in the simulation.

Returns:
a reference to this compound object
Parameters:
centerName1 name of the first existing bond center in the new bond
centerName2 name of the other existing bond center in the new bond
Compound& addRingClosingBond ( const Compound::BondCenterPathName centerName1,
const Compound::BondCenterPathName centerName2,
mdunits::Length  bondLength,
Angle  dihedral = 180 *Deg2Rad,
BondMobility::Mobility  mobility = BondMobility::Default 
)

Adds a covalent bond that is not part of the main tree-topology of bonds in a Compound.

Such bonds are necessary when there are ring structures. These bonds will not correspond to MobilizedBodies when the molecule is realized as a simbody multibody model, and thus the bond lengths, bond angles, dihedral angles, and mobility for a ring closing bond might not be enforced in the simulation.

Returns:
a reference to this compound object
Parameters:
centerName1 name of the first existing bond center in the new bond
centerName2 name of the other existing bond center in the new bond
bondLength default bond length in nanometers of the new bond (might have no effect)
dihedral default dihedral angle about new bond, with respect to the first other bond center on each atom
mobility type of motion permitted in the bond connecting parent to new subcompound

Referenced by AdenineBase::AdenineBase(), CytosineBase::CytosineBase(), GuanineBase::GuanineBase(), Histidine::Histidine(), P12::P12(), Phenylalanine::Phenylalanine(), Proline::Proline(), PurineBaseCore::PurineBaseCore(), Tryptophan::Tryptophan(), TwoNMethylGuanidineGroup::TwoNMethylGuanidineGroup(), Tyrosine::Tyrosine(), and UracilBase::UracilBase().

Compound& addSecondBondCenter ( const Compound::BondCenterName centerName,
const Compound::AtomName atomName,
Angle  bondAngle1 
)

Place a second bond center on an atom, placed a particular angle away from the first bond center.

Returns:
a reference to this compound object
Parameters:
centerName name for the new BondCenter, must be unique within the Compound
atomName name of the Atom to attach the bond center to
bondAngle1 default bond bend angle relating the first two bond centers of the atom

Referenced by BivalentAtom::BivalentAtom(), QuadrivalentAtom::QuadrivalentAtom(), and TrivalentAtom::TrivalentAtom().

Compound& bondAtom ( const Compound::SingleAtom c,
const BondCenterPathName parentBondName 
)

Bond atom using default bond length and dihedral angle.

Bond length and dihedral angle must have already been predefined in the parent *or* the child BondCenter, but not both

Returns:
a reference to this compound object
Parameters:
c the new subcompound to attach to this compound (a copy of the subcompound will be attached)
Compound& bondAtom ( const Compound::SingleAtom c,
const BondCenterPathName parentBondName,
mdunits::Length  distance,
Angle  dihedral = 0,
BondMobility::Mobility  = BondMobility::Default 
)

Add a subcompound containing exactly one atom, so the Compound::AtomName can be reused for the Compound::Name.

This atom is connected to existing material.

Returns:
a reference to this compound object
Parameters:
c the new subcompound to attach to this compound (a copy of the subcompound will be attached)
parentBondName name of the bond center on the parent Compound to which the new subcompound will be attached
distance default bond length in nanometers of new bond connecting new subcompound to parent Compound
dihedral default dihedral angle about new bond, with respect to the first bond center on each atom

Referenced by AcetylResidue::AcetylResidue(), AdenineBase::AdenineBase(), Alanine::Alanine(), AlcoholOHGroup::AlcoholOHGroup(), Arginine::Arginine(), AromaticSixMemberedCHGroup::AromaticSixMemberedCHGroup(), Asparagine::Asparagine(), BetaBranchedAminoAcidResidue::BetaBranchedAminoAcidResidue(), CarboxylateGroup::CarboxylateGroup(), Cysteine::Cysteine(), CytosineBase::CytosineBase(), Glutamate::Glutamate(), Glutamine::Glutamine(), Glycine::Glycine(), GuanineBase::GuanineBase(), Histidine::Histidine(), HNAminoAcidResidue::HNAminoAcidResidue(), Leucine::Leucine(), Lysine::Lysine(), Methane::Methane(), Methionine::Methionine(), MethyleneGroup::MethyleneGroup(), MethylGroup::MethylGroup(), NMethylAmideResidue::NMethylAmideResidue(), P12::P12(), Phenylalanine::Phenylalanine(), PrimaryAmineGroup::PrimaryAmineGroup(), Proline::Proline(), PurineBaseCore::PurineBaseCore(), Tryptophan::Tryptophan(), TwoNMethylGuanidineGroup::TwoNMethylGuanidineGroup(), Tyrosine::Tyrosine(), UracilBase::UracilBase(), Valine::Valine(), and Water::Water().

Compound& bondCompound ( const Compound::Name n,
const Compound c,
const BondCenterPathName parentBondName,
BondMobility::Mobility  mobility 
)

Add a subcompound attached by its inboard bond to an existing bond center.

Shorter version uses default bond length and dihedral angle, which must have been predefined in the parent *or* the child BondCenter, but not both.

Returns:
a reference to this compound object
Parameters:
n name for the new subcompound from the viewpoint of the parent Compound
c the new subcompound to attach to this compound (a copy of the subcompound will be attached)
parentBondName name of the bond center on the parent Compound to which the new subcompound will be attached
mobility type of motion permitted in the bond connecting parent to new subcompound
Compound& bondCompound ( const Compound::Name n,
const Compound c,
const BondCenterPathName parentBondName 
)

Add a subcompound attached by its inboard bond to an existing bond center.

Shorter version uses default bond length and dihedral angle, which must have been predefined in the parent *or* the child BondCenter, but not both.

Returns:
a reference to this compound object
Parameters:
n name for the new subcompound from the viewpoint of the parent Compound
c the new subcompound to attach to this compound (a copy of the subcompound will be attached)
parentBondName name of the bond center on the parent Compound to which the new subcompound will be attached
Compound& bondCompound ( const Compound::Name n,
const Compound c,
const BondCenterPathName parentBondName,
mdunits::Length  distance,
Angle  dihedral = 180 *Deg2Rad,
BondMobility::Mobility  mobility = BondMobility::Default 
)

Add a subcompound attached by its inboard bond to an existing bond center.

Returns:
a reference to this compound object
Parameters:
n name for the new subcompound from the viewpoint of the parent Compound
c the new subcompound to attach to this compound (a copy of the subcompound will be attached)
parentBondName name of the bond center on the parent Compound to which the new subcompound will be attached
distance default bond length in nanometers of new bond connecting new subcompound to parent Compound
dihedral default dihedral angle about new bond, with respect to the first other bond center on each atom
mobility allowed motion in new bond

Referenced by Adenylate::Adenylate(), Aspartate::Aspartate(), BetaUnbranchedAminoAcidResidue::BetaUnbranchedAminoAcidResidue(), Cytidylate::Cytidylate(), Ethane::Ethane(), Glutamate::Glutamate(), Guanylate::Guanylate(), Isoleucine::Isoleucine(), Lysine::Lysine(), P12::P12(), RNA::RNA(), Serine::Serine(), Threonine::Threonine(), TwoNMethylGuanylate::TwoNMethylGuanylate(), Tyrosine::Tyrosine(), Uridylate::Uridylate(), and Valine::Valine().

Vec3 calcAtomAccelerationInGroundFrame ( const State state,
Compound::AtomIndex  atomId 
) const
Returns:
vector acceleration in nanometers per picosecond squared
Parameters:
state simbody State representing current configuration
atomId integer index of Atom with respect to this Compound
Vec3 calcAtomLocationInCompoundFrame ( const State state,
Compound::AtomIndex  atomId 
) const

Compute atom location in local Compound frame.

Returns:
default (initial) location of atom in orthogonal nanometers with respect to Ground frame
Parameters:
state simbody State representing current configuration
atomId integer index of Atom with respect to this Compound
Vec3 calcAtomLocationInGroundFrame ( const State state,
Compound::AtomIndex  atomId 
) const
Returns:
location in orthogonal nanometers
Parameters:
state simbody State representing current configuration
atomId integer index of Atom with respect to this Compound
Vec3 calcAtomVelocityInGroundFrame ( const State state,
Compound::AtomIndex  atomId 
) const
Returns:
vector velocity in nanometers per picosecond
Parameters:
state simbody State representing current configuration
atomId integer index of Atom with respect to this Compound
Transform calcDefaultAtomFrameInCompoundFrame ( Compound::AtomIndex   )  const
Returns:
location and orientation of Atom with respect to this Compound
Vec3 calcDefaultAtomLocationInCompoundFrame ( const AtomPathName atomName  )  const

Compute atom location in local Compound frame.

Returns:
default (initial) location of atom in orthogonal nanometers with respect to Ground frame
Parameters:
atomName name of the atom from viewpoint of Compound
Vec3 calcDefaultAtomLocationInGroundFrame ( const AtomPathName atomName  )  const

Compute atom location in local Compound frame.

Returns:
default (initial) location of atom in orthogonal nanometers with respect to Ground frame
Parameters:
atomName name of the atom from viewpoint of Compound

Referenced by LigandDroplet::LigandDroplet(), and WaterDroplet::WaterDroplet().

Angle calcDefaultDihedralAngle ( const DihedralName dihedralName  )  const

Computes default (initial) dihedral angle of a previously named dihedral.

Returns:
dihedral angle in radians with respect to first other bond centers of bonded atoms
Parameters:
dihedralName name of predefined dihedral angle with respect to this Compound

Referenced by Ethane::calcDefaultTorsionAngle().

Angle calcDihedralAngle ( const State state,
const DihedralName dihedralName 
) const

Computes dynamic dihedral angle of a previously named dihedral.

Returns:
dihedral angle in radians with respect to first other bond centers of bonded atoms
Parameters:
state simbody State object representing the current configuration
dihedralName name of predefined dihedral angle with respect to this Compound
Compound& convertInboardBondCenterToOutboard (  ) 

Make so that this compound can no longer be a child to the geometry of another compound Raises an error if the inboard bond center is already bonded.

Returns:
a reference to this compound object

Referenced by AcetylResidue::AcetylResidue(), Ethane::Ethane(), Methane::Methane(), and RNA::RNA().

virtual AtomTargetLocations createAtomTargets ( const class PdbResidue targetResidue  )  const [virtual]

Compute atom location in local Compound frame.

Returns:
default (initial) location of atom in orthogonal nanometers with respect to Ground frame
virtual AtomTargetLocations createAtomTargets ( const class PdbChain targetChain  )  const [virtual]

Compute atom location in local Compound frame.

Returns:
default (initial) location of atom in orthogonal nanometers with respect to Ground frame

Reimplemented in Biopolymer.

virtual AtomTargetLocations createAtomTargets ( const class PdbStructure targetStructure  )  const [virtual]

Create a mapping between this Compound and atom locations in a PdbStructure.

For use in the process of importing a configuration from a PDB file

Reimplemented in Biopolymer.

Referenced by P12::P12().

Compound& defineDihedralAngle ( const Compound::DihedralName angleName,
const Compound::BondCenterPathName bond1,
const Compound::BondCenterPathName bond2,
Angle  offset = 0 *Deg2Rad 
)

Define a named dihedral in terms of two bond centers.

Useful in cases where not all four atoms have been placed yet. The bond centers must be from two atoms that are bonded (atoms 2 and 3 of 4) But are NOT the bond centers that connect the two middle atoms.

The offset parameter is used in cases where the dihedral angle definition differs from the IUPAC definition of dihedral angles using four atoms. For example the definition of protein phi and psi angles is 180 degrees offset from the standard four-atom dihedral angle definition.

Returns:
a reference to this compound object
Parameters:
angleName unique name for new dihedral angle definition
bond1 first bond center, connecting atom 2 to atom 1
bond2 second bond center, connecting atom 3 to atom 4
offset nomenclature offset
Compound& defineDihedralAngle ( const Compound::DihedralName angleName,
const Compound::AtomPathName atom1,
const Compound::AtomPathName atom2,
const Compound::AtomPathName atom3,
const Compound::AtomPathName atom4,
Angle  offset = 0 *Deg2Rad 
)

Define a named dihedral angle using four atoms.

The offset parameter is used in cases where the dihedral angle definition differs from the IUPAC definition of dihedral angles using four atoms. For example the definition of protein phi and psi angles is 180 degrees offset from the standard four-atom dihedral angle definition.

Returns:
a reference to this compound object
Parameters:
angleName unique name for new dihedral angle definition
atom1 first atom name
atom2 second atom name
atom3 third atom name
atom4 fourth atom name
offset nomenclature offset

Referenced by Adenylate::Adenylate(), Arginine::Arginine(), Asparagine::Asparagine(), Aspartate::Aspartate(), BetaBranchedAminoAcidResidue::BetaBranchedAminoAcidResidue(), BetaUnbranchedAminoAcidResidue::BetaUnbranchedAminoAcidResidue(), Cysteine::Cysteine(), Cytidylate::Cytidylate(), Ethane::Ethane(), Glutamate::Glutamate(), Glutamine::Glutamine(), Guanylate::Guanylate(), Histidine::Histidine(), Isoleucine::Isoleucine(), Leucine::Leucine(), Lysine::Lysine(), Methionine::Methionine(), NMethylAmideResidue::NMethylAmideResidue(), P12::P12(), Phenylalanine::Phenylalanine(), Proline::Proline(), RNA::RNA(), Serine::Serine(), Tryptophan::Tryptophan(), TwoNMethylGuanylate::TwoNMethylGuanylate(), Tyrosine::Tyrosine(), Uridylate::Uridylate(), and Valine::Valine().

Compound& fitDefaultConfiguration ( const AtomTargetLocations atomTargets,
SimTK::Real  targetRms 
)

Optimize adjustable degrees of freedom to best match atom targets.

BiotypeIndex getAtomBiotypeIndex ( Compound::AtomIndex   )  const
Returns:
the biotype assigned to an Atom in this Compound
const Element& getAtomElement ( const Compound::AtomName  )  const

Name a type of Compound.

Sets the name for this class of compound, e.g. "ethane", not "ethane number 3". The compound name can be important for automatic resolution of atom types for during force field resolution.

Returns:
a reference to this compound object
const Element& getAtomElement ( Compound::AtomIndex   )  const
Returns:
chemical element of atom
Compound::AtomIndex getAtomIndex ( const Compound::AtomPathName  )  const
Returns:
integer index of named atom, with repect to this Compound
Vec3 getAtomLocationInMobilizedBodyFrame ( Compound::AtomIndex   )  const

get the location of an Atom in the frame of the simbody MobilizedBody to which the Atom is attached

Requires that this Compound has already been modeled in a CompoundSystem

Returns:
location in orthogonal nanometers with respect to the rigid body to which the atom is attached

Referenced by RNA::calcRNABaseNormal().

MobilizedBodyIndex getAtomMobilizedBodyIndex ( Compound::AtomIndex   )  const

get the simbody MobilizedBody to which a particular atom is attached

Requires that this Compound has already been modeled in a CompoundSystem

Returns:
the integer index of the MobilizedBody in the CompoundSystem

Referenced by RNA::calcRNABaseNormal().

const AtomPathName getAtomName ( Compound::AtomIndex   )  const

Returns the most recently assigned name, if any, given to an atom in this compound.

Referenced by LigandDroplet::LigandDroplet(), and WaterDroplet::WaterDroplet().

Compound::AtomIndex getBondAtomIndex ( Compound::BondIndex  bondIndex,
int  which 
) const
Returns:
true if Compound contains a subcompound by that name (relative to this Compound)
integer index of one of the two atoms involved in a particular covalent bond
Parameters:
bondIndex integer index of a bond in this Compound
which zero(0) for parent-side (rootward) atom of bond, one(1) for child-side (leafward) atom
const Name& getCompoundName (  )  const
Returns:
name of Compound type

Referenced by RibonucleosideResidue::withPhosphodiester().

DuMM::AtomIndex getDuMMAtomIndex ( Compound::AtomIndex   )  const

get DuMMForceFieldSubsystem atom index corresponding to an atom in this Compound

Requires that this Compound has already been modeled in a CompoundSystem

Returns:
an AtomIndex in the DuMMForceFieldSubystem
int getNAtoms (  )  const [inline]
size_t getNBondCenters ( Compound::AtomIndex  atomIndex  )  const [inline]
size_t getNBondCenters (  )  const [inline]
size_t getNBonds (  )  const [inline]
int getNumAtoms (  )  const
Returns:
the number of atoms in the Compound, including those within subcompounds.

Referenced by LigandDroplet::LigandDroplet(), and WaterDroplet::WaterDroplet().

size_t getNumBondCenters ( Compound::AtomIndex  atomIndex  )  const
Returns:
total number of BondCenters in a particular atom of a Compound
size_t getNumBondCenters (  )  const
Returns:
total number of BondCenters in a Compound, including its subcompounds
size_t getNumBonds (  )  const
Returns:
total number of Bonds in a Compound, including its subcompounds

Referenced by P12::P12().

char getPdbChainId (  )  const
Returns:
character that would populate "chain id" field for PDB file output for this Compound

Referenced by RibonucleosideResidue::withPhosphodiester().

const String& getPdbResidueName (  )  const
Returns:
string that would be used to populate "residue type" field for this Compound in a PDB file

Referenced by RNA::RNA(), and RibonucleosideResidue::withPhosphodiester().

int getPdbResidueNumber (  )  const
Returns:
integer that would appear in "residue number" field for this Compound in a PDB file

Referenced by RibonucleosideResidue::withPhosphodiester().

const Compound& getSubcompound ( const Compound::Name subcompoundName  )  const
Returns:
read-only reference to a subcompound of this Compound
mutable reference to a subcompound of this Compound
read-only reference to a subcompound of this Compound
Parameters:
subcompoundName name of subcompound
Transform getSubcompoundFrameInParentFrame ( const Compound::Name subcompoundName  )  const
Returns:
default (initial) location and orientation of a subcompound with respect to this Compound
Parameters:
subcompoundName name of subcompound with respect to this Compound
const Transform& getTopLevelTransform (  )  const
Returns:
the orientation and location of this Compound

Referenced by CompoundSystem::adoptCompound().

TransformAndResidual getTransformAndResidual ( const Compound::AtomTargetLocations atomTargets  )  const

Compute atom location in local Compound frame.

Returns:
default (initial) location of atom in orthogonal nanometers with respect to Ground frame
bool hasAtom ( const AtomPathName name  )  const
Returns:
total number of subcompounds in a compoud, including sub-subcompounds etc.
true if Compound contains an atom by that name (relative to this Compound)
Parameters:
name name of the atom relative to this compound.

Referenced by RNA::RNA().

bool hasBondCenter ( const BondCenterPathName bondCenter  )  const
Returns:
true if Compound contains a BondCenter (half-bond) by that name (relative to this Compound)
Parameters:
bondCenter name of the BondCenter relative to this compound.
Compound& inheritAtomNames ( const Compound::Name  ) 

Convenience method to locally import all of the atom names that a particular subcompound uses.

Returns:
a reference to this compound object

Referenced by Adenylate::Adenylate(), Cytidylate::Cytidylate(), Guanylate::Guanylate(), Methane::Methane(), RNA::RNA(), TwoNMethylGuanylate::TwoNMethylGuanylate(), and Uridylate::Uridylate().

Compound& inheritBondCenterNames ( const Compound::Name  ) 

Convenience method to locally import all of the BondCenter names that a particular subcompound uses.

Returns:
a reference to this compound object
Compound& inheritCompoundSynonyms ( const Compound otherCompound  ) 

Name a type of Compound.

Sets the name for this class of compound, e.g. "ethane", not "ethane number 3". The compound name can be important for automatic resolution of atom types for during force field resolution.

Returns:
a reference to this compound object
Compound& matchDefaultAtomChirality ( const AtomTargetLocations atomTargets,
Angle  planarityTolerance = 90.0 *Deg2Rad 
)

Adjust stereochemistry about chiral atoms to match that seen in a set of atomic locations.

Choose a small value of planarityTolerance parameter to break the planarity of out-of-plane atoms from the target set

Returns:
a reference to this compound
Parameters:
atomTargets another set of atom locations to match chirality of
planarityTolerance break planarity of BondCenters more than this angle out of plane (radians)

Referenced by P12::P12().

Compound& matchDefaultBondAngles ( const AtomTargetLocations atomTargets  ) 

Compute atom location in local Compound frame.

Returns:
default (initial) location of atom in orthogonal nanometers with respect to Ground frame

Referenced by P12::P12().

Compound& matchDefaultBondLengths ( const AtomTargetLocations atomTargets  ) 

Compute atom location in local Compound frame.

Returns:
default (initial) location of atom in orthogonal nanometers with respect to Ground frame

Referenced by P12::P12().

Compound& matchDefaultConfiguration ( const AtomTargetLocations atomTargets,
MatchStratagem  matchStratagem = Match_Exact 
) [inline]

Compute atom location in local Compound frame.

Returns:
default (initial) location of atom in orthogonal nanometers with respect to Ground frame

References Compound::DistortPlanarBonds, Compound::KeepPlanarBonds, Compound::Match_Exact, Compound::Match_Idealized, and Compound::Match_TopologyOnly.

Compound& matchDefaultDihedralAngles ( const AtomTargetLocations atomTargets,
PlanarBondMatchingPolicy  policy = KeepPlanarBonds 
)

Compute atom location in local Compound frame.

Returns:
default (initial) location of atom in orthogonal nanometers with respect to Ground frame
Parameters:
policy whether to keep ideal torsions on bonds between planar atoms

Referenced by P12::P12().

Compound& matchDefaultTopLevelTransform ( const AtomTargetLocations atomTargets  ) 

Compute atom location in local Compound frame.

Returns:
default (initial) location of atom in orthogonal nanometers with respect to Ground frame

Referenced by P12::P12().

Compound& nameAtom ( const Compound::AtomName newName,
const AtomPathName oldName,
BiotypeIndex  biotype 
)
Returns:
a reference to this compound object
Parameters:
newName new name for this atom; name is local to this Compound
oldName previous name; can be a subcompound-qualified Atom name
biotype new Biotype for the atom, with respect to this Compound
Compound& nameAtom ( const Compound::AtomName newName,
const AtomPathName oldName 
)
Compound& nameBondCenter ( const Compound::BondCenterName newName,
const BondCenterPathName oldName 
)

Define a local name for a particular BondCenter in this Compound.

Returns:
a reference to this compound object
Parameters:
newName new local name for bond center
oldName previous name, can be a subcompound-qualified name

Referenced by AcetylResidue::AcetylResidue(), AlcoholOHGroup::AlcoholOHGroup(), BetaBranchedAminoAcidResidue::BetaBranchedAminoAcidResidue(), BetaUnbranchedAminoAcidResidue::BetaUnbranchedAminoAcidResidue(), CarboxylateGroup::CarboxylateGroup(), MethyleneGroup::MethyleneGroup(), MethylGroup::MethylGroup(), NMethylAmideResidue::NMethylAmideResidue(), and PurineBaseCore::PurineBaseCore().

const Compound& populateDefaultPdbChain ( class PdbChain ,
int &  defaultNextResidueNumber,
const Transform transform = Transform() 
) const

Write current default(initial) Compound configuration into a PdbChain object.

const Compound& populatePdbChain ( const State state,
class PdbChain ,
int &  defaultNextResidueNumber,
const Transform transform = Transform() 
) const

Write dynamic Compound configuration into a PdbChain object.

Compound& setAtomBiotype ( const Compound::AtomPathName atomName,
const String biotypeResidueName,
const String biotypeAtomName,
SimTK::Ordinality::Residue  ordinality = SimTK::Ordinality::Any 
) [inline]

Add this Compound, including all of its subcompounds, to a CompoundSystem, for simulation etc.

This Compound will become a top-level object with the CompoundSystem.

References Biotype::defineBiotype(), Biotype::exists(), Biotype::get(), Biotype::getElement(), Biotype::getIndex(), and Biotype::getValence().

Referenced by AcetylResidue::AcetylResidue(), Alanine::Alanine(), NMethylAmideResidue::NMethylAmideResidue(), and P12::P12().

Compound& setBaseAtom ( const Compound::SingleAtom c,
const Transform = Transform() 
)

Add a first subcompound containing exactly one atom, so the Compound::AtomName can be reused for the Compound::Name.

This atom is not connected to anything else (yet).

Returns:
a reference to this compound object
Parameters:
c single-atom Compound to add as first subcompound of this Compound
Compound& setBaseAtom ( const Compound::AtomName name,
const Biotype biotype,
const Transform location = Vec3(0) 
)

Add the first atom unconnected to anything else (yet).

Returns:
a reference to this compound object
Parameters:
name name of the new atom being created
biotype Biotype of the new atom being created
location default location of the new atom being created in orthogonal nanometers. Defaults to (0,0,0)
Compound& setBaseAtom ( const Compound::AtomName name,
const Element element,
const Transform location = Vec3(0) 
)
Compound& setBaseCompound ( const Compound::Name n,
const Compound c,
const Transform location = Transform() 
)

Add a first subcompound without attaching it to anything.

Returns:
a reference to this compound object
Parameters:
n name for new subcompound, from viewpoint of parent
c new subcompound to be copied and attached to parent Compound

Referenced by AcetylResidue::AcetylResidue(), Ethane::Ethane(), Methane::Methane(), and P12::P12().

Compound& setBiotypeIndex ( const Compound::AtomPathName atomName,
BiotypeIndex  biotype 
)

define the Biotype for an Atom in this Compound

Returns:
a reference to this compound object
Parameters:
atomName name of the atom in the context of this Compound
biotype integer index of an existing Biotype known to the Biotype class

Referenced by Methane::Methane(), RNA::RNA(), Serine::Serine(), TwoNMethylGuanidineGroup::TwoNMethylGuanidineGroup(), TwoNMethylGuanylate::TwoNMethylGuanylate(), and Water::Water().

Compound& setBondMobility ( BondMobility::Mobility  mobility,
Compound::BondIndex  bondIndex 
)

Override the default rotatability of a bond.

Returns:
a reference to this compound object
Parameters:
mobility the new allowed motion of the bond
bondIndex the integer index of the Bond in the context of this Compound
Compound& setBondMobility ( BondMobility::Mobility  mobility,
const AtomPathName atom1,
const AtomPathName atom2 
)

Override the default rotatability of a bond.

Returns:
a reference to this compound object
Parameters:
mobility the new allowed motion of the bond
atom1 the name of the first atom defining the Bond
atom2 the name of the second atom defining the Bond

Referenced by Arginine::Arginine(), Asparagine::Asparagine(), Glutamine::Glutamine(), Histidine::Histidine(), P12::P12(), Phenylalanine::Phenylalanine(), Proline::Proline(), RNA::setRNABondMobility(), Tryptophan::Tryptophan(), and Tyrosine::Tyrosine().

Compound& setCompoundBondMobility ( BondMobility::Mobility  mobility  ) 

Set BondMobility for every bond in the Compound.

Returns:
A reference to this Compound
Compound& setCompoundName ( const Name  ) 
Compound& setDefaultBondAngle ( Angle  angle,
const AtomPathName atom1,
const AtomPathName atom2,
const AtomPathName atom3 
)

Sets a default(initial) bond angle defined by three atoms.

Warning:
setDefaultBondAngle only works if one of the two bond centers is bond center 1 or 2
Todo:
  • remove this restriction
Returns:
a reference to this compound object
Parameters:
angle new default bond angle in radians
atom1 name of first atom defining bond angle
atom2 name of middle atom defining bond angle
atom3 name of third atom defining bond angle

Referenced by Cysteine::Cysteine(), Proline::Proline(), and Water::Water().

Compound& setDefaultBondLength ( mdunits::Length  length,
const AtomPathName atom1,
const AtomPathName atom2 
)

Sets a default (inital) bond length defined by two atoms.

Returns:
a reference to this compound object
Parameters:
length default bond length in nanometers
atom1 name of first atom in bond
atom2 name of second atom in bond
Compound& setDefaultDihedralAngle ( Angle  angle,
const Compound::AtomName atom1,
const Compound::AtomName atom2,
const Compound::AtomName atom3,
const Compound::AtomName atom4 
)

Compute atom location in local Compound frame.

Returns:
default (initial) location of atom in orthogonal nanometers with respect to Ground frame
Compound& setDefaultDihedralAngle ( Angle  angle,
Compound::AtomIndex  atom1,
Compound::AtomIndex  atom2,
Compound::AtomIndex  atom3,
Compound::AtomIndex  atom4 
)

Compute atom location in local Compound frame.

Returns:
default (initial) location of atom in orthogonal nanometers with respect to Ground frame
Compound& setDefaultDihedralAngle ( const DihedralName dihedralName,
Angle  angleInRadians 
)
Compound& setDefaultInboardBondLength ( mdunits::Length   ) 

Sets default (initial) bond length of current or future Bond using this Compound's inboard BondCenter.

Default inboard bond length and dihedral should only be set when the future bonding partners of this compound are restricted to a particular atom type

Returns:
a reference to this compound object

Referenced by AlcoholOHGroup::AlcoholOHGroup(), AliphaticCarbon::AliphaticCarbon(), AliphaticHydrogen::AliphaticHydrogen(), AromaticSixMemberedCHGroup::AromaticSixMemberedCHGroup(), CarboxylateGroup::CarboxylateGroup(), MethyleneGroup::MethyleneGroup(), MethylGroup::MethylGroup(), NMethylAmideResidue::NMethylAmideResidue(), and PrimaryAmineGroup::PrimaryAmineGroup().

Compound& setDefaultInboardDihedralAngle ( Angle   ) 

Stores a default (initial) dihedral angle in the inboard BondCenter of this Compound.

If exactly one of the two BondCenters that are joined to form a Bond has default geometry set, then the bondCompound() lacking explicit geometry arguments may be used to create the Bond.

Returns:
a reference to this compound object
Compound& setDihedralAngle ( State state,
const DihedralName dihedralName,
Angle   
)

Sets dynamic dihedral angle of a previously named dihedral.

The dihedral angle must be an adjustable parameter of the dynamic model.

Returns:
a reference to this compound object
Parameters:
state simbody State object representing the current configuration
dihedralName name of predefined dihedral angle with respect to this Compound
void setDuMMAtomIndex ( Compound::AtomIndex  ,
DuMM::AtomIndex   
) [protected]

Stores relationship between a Compound Atom and an Atom defined in a DuMMForcefieldSubsystem.

Compound& setInboardBondCenter ( const Compound::BondCenterName centerName  ) 

setInboardBondCenter assigns special status to a bond center.

There can be at most one inboard bond center in a Compound. Only an inboard bond center can be used to bond to a parent compound

Returns:
a reference to this compound object
Parameters:
centerName name of existing BondCenter to use as inboard BondCenter

Referenced by BivalentAtom::BivalentAtom(), QuadrivalentAtom::QuadrivalentAtom(), TrivalentAtom::TrivalentAtom(), and UnivalentAtom::UnivalentAtom().

void setMultibodySystem ( MultibodySystem system  ) 

Add this Compound, including all of its subcompounds, to a CompoundSystem, for simulation etc.

This Compound will become a top-level object with the CompoundSystem.

Referenced by CompoundSystem::adoptCompound().

Compound& setPdbChainId ( char  chainId  ) 

Set character to populate "chain id" field for PDB file output.

Returns:
a reference to this compound object
Parameters:
chainId Protein Data Bank "chain Id" for this Compound

Referenced by LigandDroplet::LigandDroplet(), P12::P12(), and WaterDroplet::WaterDroplet().

Compound& setPdbResidueName ( const String resName  ) 
Compound& setPdbResidueNumber ( int  resNum  ) 

Set value to populate "residue number" field for PDB file output.

Returns:
a reference to this compound object

Referenced by Protein::initialize(), LigandDroplet::LigandDroplet(), P12::P12(), RNA::RNA(), and WaterDroplet::WaterDroplet().

Compound& setTopLevelTransform ( const Transform transform  ) 

Set the orientation and location of this Compound.

Referenced by CompoundSystem::adoptCompound().

SimTK_DEFINE_UNIQUE_LOCAL_INDEX_TYPE ( Compound  ,
BondIndex   
)

Compound::BondIndex type is an integer index into Bonds of a Compound.

It is NOT instrinsic to the Bond, but represents the relationship between a Bond and precisely one of its parent compounds.

SimTK_DEFINE_UNIQUE_LOCAL_INDEX_TYPE ( Compound  ,
BondCenterIndex   
)

Compound::BondCenterIndex type is an integer index into BondCenters of a Compound.

It is NOT instrinsic to the BondCenter, but represents the relationship between a BondCenter and precisely one of its parent compounds.

SimTK_DEFINE_UNIQUE_LOCAL_INDEX_TYPE ( Compound  ,
LocalAtomIndex   
)

Compound::LocalAtomIndex type is an integer index into atoms directly attached to a Compound, that is to say that the atom does not belong to any subcompounds of the Compound.

This index type should ordinarily not be used by clients of the Molmodel API

SimTK_DEFINE_UNIQUE_LOCAL_INDEX_TYPE ( Compound  ,
AtomIndex   
)

Compound::Index type is an integer index into subcompounds of a Compound.

It is NOT instrinsic to the subcompound, but represents the relationship between a subcompound and precisely one of its parent compounds. Compound::AtomIndex type is an integer index into atoms of a Compound. It is NOT instrinsic to the atom, but represents the relationship between an atom and precisely one of its parent compounds.

Compound& updSubcompound ( const Compound::Name subcompoundName  ) 
Returns:
mutable reference to a subcompound of this Compound
Parameters:
subcompoundName name of subcompound
std::ostream& writeDefaultPdb ( std::ostream &  os,
int &  nextAtomSerialNumber,
const Transform transform = Transform() 
) const

Write the default (initial) configuration in Protein Data Bank (PDB) format.

integer nextAtomSerialNumber reference is incremented within writeDefaultPdb method, so that subsequent calls to this method using the same integer variable will continue numbering atoms where the previous call left off.

Parameters:
os output stream to write PDB coordinates to
nextAtomSerialNumber mutable integer reference containing the next desired atom serial number
transform optional change to location and orientation of molecule
void writeDefaultPdb ( const char *  outFileName,
const Transform transform 
) const

/brief This polymorphism takes a char* file name rather than ostream, to save the user a couple of lines of code.

std::ostream& writeDefaultPdb ( std::ostream &  os,
const Transform transform = Transform() 
) const

Write the default (initial) configuration in Protein Data Bank (PDB) format.

Starting with the first atom's serial number as one(1).

Parameters:
os output stream to write PDB coordinates to
transform optional change to location and orientation of molecule
std::ostream& writePdb ( const State state,
std::ostream &  os,
int &  nextAtomSerialNumber,
const Transform transform = Transform() 
) const

Write the dynamic Compound configuration in Protein Data Bank (PDB) format.

integer nextAtomSerialNumber reference is incremented within writePdb method, so that subsequent calls to this method using the same integer variable will continue numbering atoms where the previous call left off.

Parameters:
state simbody state representing the current configuration of the molecule
os output stream to write PDB coordinates to
nextAtomSerialNumber mutable integer reference containing the next desired atom serial number
transform optional change to location and orientation of molecule
std::ostream& writePdb ( const State state,
std::ostream &  os,
const Transform transform = Transform() 
) const

Write the dynamic Compound configuration in Protein Data Bank (PDB) format.

Starting with the first atom's serial number as one(1).

Parameters:
state simbody state representing the current configuration of the molecule
os output stream to write PDB coordinates to
transform optional change to location and orientation of molecule

Referenced by PeriodicPdbWriter::handleEvent().


Friends And Related Function Documentation

friend class CompoundSystem [friend]

The documentation for this class was generated from the following file:

Generated on Wed Dec 30 11:05:11 2009 for SimTKcore by  doxygen 1.6.1