OBMol Class Reference

Molecule Class. More...

#include <mol.h>

Inheritance diagram for OBMol:

Inheritance graph
[legend]
List of all members.

Initialization and data (re)size methods

 OBMol ()
 Constructor.
 OBMol (const OBMol &)
virtual ~OBMol ()
 Destructor.
OBMoloperator= (const OBMol &mol)
 Assignment.
OBMoloperator+= (const OBMol &mol)
void ReserveAtoms (int natoms)
virtual OBAtomCreateAtom (void)
virtual OBBondCreateBond (void)
virtual void DestroyAtom (OBNodeBase *)
virtual void DestroyBond (OBEdgeBase *)
bool AddAtom (OBAtom &)
 Add an atom to a molecule.
bool AddBond (int, int, int, int flags=0, int insertpos=-1)
bool AddBond (OBBond &)
bool AddResidue (OBResidue &)
bool InsertAtom (OBAtom &)
bool DeleteAtom (OBAtom *)
bool DeleteBond (OBBond *)
bool DeleteResidue (OBResidue *)
OBAtomNewAtom ()
 Get a new atom to add to a molecule.
OBResidueNewResidue ()

Molecule modification methods

virtual OBBaseDoTransformations (const std::map< std::string, std::string > *pOptions)
bool Clear ()
 Clear all information from a molecule.
void RenumberAtoms (std::vector< OBNodeBase * > &)
 Renumber the atoms of this molecule according to the order in the supplied vector.
void ToInertialFrame (int conf, double *rmat)
 Translate one conformer and rotate by a rotation matrix (which is returned) to the inertial frame-of-reference.
void ToInertialFrame ()
 Translate all conformers to the inertial frame-of-reference.
void Translate (const vector3 &v)
 Translates all conformers in the molecule by the supplied vector.
void Translate (const vector3 &v, int conf)
 Translates one conformer in the molecule by the supplied vector.
void Rotate (const double u[3][3])
void Rotate (const double m[9])
void Rotate (const double m[9], int nconf)
void Center ()
 Translate to the center of all coordinates (for this conformer).
bool Kekulize ()
 Transform to standard Kekule bond structure (presumably from an aromatic form).
bool PerceiveKekuleBonds ()
void NewPerceiveKekuleBonds ()
 Kekulize aromatic rings without using implicit valence.
bool DeleteHydrogen (OBAtom *)
bool DeleteHydrogens ()
bool DeleteHydrogens (OBAtom *)
bool DeleteNonPolarHydrogens ()
bool AddHydrogens (bool polaronly=false, bool correctForPH=true)
bool AddHydrogens (OBAtom *)
bool AddPolarHydrogens ()
bool StripSalts ()
 Deletes all atoms except for the largest contiguous fragment.
bool ConvertDativeBonds ()
 Converts for instance [N+]([O-])=O to N(=O)=O.
bool CorrectForPH ()
bool AssignSpinMultiplicity ()
 set spin multiplicity for H-deficient atoms
vector3 Center (int nconf)
void SetTorsion (OBAtom *, OBAtom *, OBAtom *, OBAtom *, double)
 Set the torsion defined by these atoms, rotating bonded neighbors.
static const char * ClassDescription ()

Public Member Functions

unsigned int NumNodes ()
unsigned int NumEdges ()
void ResetVisitFlags ()
bool SetVisitLock (bool)
bool GetVisitLock ()
OBNodeBaseBegin (std::vector< OBNodeBase * >::iterator &)
OBEdgeBaseBegin (std::vector< OBEdgeBase * >::iterator &)
OBNodeBaseNext (std::vector< OBNodeBase * >::iterator &)
OBEdgeBaseNext (std::vector< OBEdgeBase * >::iterator &)
Molecule modification methods
virtual void BeginModify (void)
 Call when making many modifications -- clears conformer/rotomer data.
virtual void EndModify (bool nukePerceivedData=true)
 Call when done with modificaions -- re-perceive data as needed.
int GetMod ()
void IncrementMod ()
void DecrementMod ()
Generic data handling methods (via OBGenericData)
bool HasData (std::string &)
 
Returns:
whether the generic attribute/value pair exists

bool HasData (const char *)
 
Returns:
whether the generic attribute/value pair exists

bool HasData (unsigned int type)
 
Returns:
whether the generic attribute/value pair exists

void DeleteData (unsigned int type)
void DeleteData (OBGenericData *)
void DeleteData (std::vector< OBGenericData * > &)
void SetData (OBGenericData *d)
unsigned int DataSize ()
 
Returns:
the number of OBGenericData items attached to this molecule.

OBGenericDataGetData (unsigned int type)
OBGenericDataGetData (std::string &)
 Returns the value given an attribute name.
OBGenericDataGetData (const char *)
 Returns the value given an attribute name.
std::vector< OBGenericData * > & GetData ()
std::vector< OBGenericData
* >::iterator 
BeginData ()
std::vector< OBGenericData
* >::iterator 
EndData ()
Data retrieval methods
int GetFlags ()
const char * GetTitle () const
 
Returns:
the title of this molecule (often the filename)

unsigned int NumAtoms () const
 
Returns:
the number of atoms (i.e. OBAtom children)

unsigned int NumBonds () const
 
Returns:
the number of bonds (i.e. OBBond children)

unsigned int NumHvyAtoms ()
 
Returns:
the number of non-hydrogen atoms

unsigned int NumResidues () const
 
Returns:
the number of residues (i.e. OBResidue substituents)

unsigned int NumRotors ()
 
Returns:
the number of rotatble bonds. See OBBond::IsRotor() for details

OBAtomGetAtom (int)
OBAtomGetFirstAtom ()
OBBondGetBond (int)
OBBondGetBond (int, int)
OBBondGetBond (OBAtom *, OBAtom *)
OBResidueGetResidue (int)
std::vector< OBInternalCoord * > GetInternalCoord ()
double GetTorsion (int, int, int, int)
 
Returns:
the dihedral angle between the four atoms supplied a1-a2-a3-a4)

double GetTorsion (OBAtom *, OBAtom *, OBAtom *, OBAtom *)
 
Returns:
the dihedral angle between the four atoms supplied a1-a2-a3-a4)

std::string GetFormula ()
 
Returns:
the stochoimetric formula (e.g., C4H6O)

double GetEnergy () const
 
Returns:
the heat of formation for this molecule (in kcal/mol)

double GetMolWt ()
 
Returns:
the standard molar mass given by IUPAC atomic masses (amu)

double GetExactMass ()
 
Returns:
the mass given by isotopes (or most abundant isotope, if not specified)

int GetTotalCharge ()
 
Returns:
the total charge on this molecule (i.e., 0 = neutral, +1, -1...)

unsigned int GetTotalSpinMultiplicity ()
 
Returns:
the total spin on this molecule (i.e., 1 = singlet, 2 = doublet...)

unsigned short int GetDimension () const
 
Returns:
the dimensionality of coordinates (i.e., 0 = unknown or no coord, 2=2D, 3=3D)

double * GetCoordinates ()
std::vector< OBRing * > & GetSSSR ()
 Implements blue-obelisk:findSmallestSetOfSmallestRings.
bool AutomaticFormalCharge ()
 Get the current flag for whether formal charges are set with pH correction.
bool AutomaticPartialCharge ()
 Get the current flag for whether partial charges are auto-determined.
Data modification methods
void SetTitle (const char *title)
void SetTitle (std::string &title)
void SetFormula (std::string molFormula)
 Set the stochiometric formula for this molecule.
void SetEnergy (double energy)
 Set the heat of formation for this molecule (in kcal/mol).
void SetDimension (unsigned short int d)
 Set the dimension of this molecule (i.e., 0, 1 , 2, 3).
void SetTotalCharge (int charge)
void SetTotalSpinMultiplicity (unsigned int spin)
void SetInternalCoord (std::vector< OBInternalCoord * > int_coord)
void SetAutomaticFormalCharge (bool val)
 Set the flag for determining automatic formal charges with pH (default=true).
void SetAutomaticPartialCharge (bool val)
 Set the flag for determining partial charges automatically (default=true).
void SetAromaticPerceived ()
 Mark that aromaticity has been perceived for this molecule (see OBAromaticTyper).
void SetSSSRPerceived ()
 Mark that Smallest Set of Smallest Rings has been run (see OBRing class).
void SetRingAtomsAndBondsPerceived ()
 Mark that rings have been perceived (see OBRing class for details).
void SetAtomTypesPerceived ()
 Mark that atom types have been perceived (see OBAtomTyper for details).
void SetChainsPerceived ()
 Mark that chains and residues have been perceived (see OBChainsParser).
void SetChiralityPerceived ()
 Mark that chirality has been perceived.
void SetPartialChargesPerceived ()
 Mark that partial charges have been assigned.
void SetHybridizationPerceived ()
void SetImplicitValencePerceived ()
void SetKekulePerceived ()
void SetClosureBondsPerceived ()
void SetHydrogensAdded ()
void SetCorrectedForPH ()
void SetAromaticCorrected ()
void SetSpinMultiplicityAssigned ()
void SetFlags (int flags)
void UnsetAromaticPerceived ()
void UnsetPartialChargesPerceived ()
void UnsetImplicitValencePerceived ()
void UnsetFlag (int flag)
Molecule utilities and perception methods
void FindSSSR ()
 Find Smallest Set of Smallest Rings (see OBRing class for more details).
void FindRingAtomsAndBonds ()
void FindChiralCenters ()
 Sets atom->IsChiral() to true for chiral atoms.
void FindChildren (std::vector< int > &, int, int)
void FindChildren (std::vector< OBAtom * > &, OBAtom *, OBAtom *)
void FindLargestFragment (OBBitVec &)
void ContigFragList (std::vector< std::vector< int > > &)
void Align (OBAtom *, OBAtom *, vector3 &, vector3 &)
 Aligns atom a on p1 and atom b along p1->p2 vector.
void ConnectTheDots ()
 Adds single bonds based on atom proximity.
void PerceiveBondOrders ()
 Attempts to perceive multiple bonds based on geometries.
void FindTorsions ()
 Fills the OBGeneric OBTorsionData with torsions from the mol.
bool GetGTDVector (std::vector< int > &)
 Calculates the graph theoretical distance of each atom. Vector is indexed from zero.
void GetGIVector (std::vector< unsigned int > &)
 Calculates a set of graph invariant indexes using the graph theoretical distance, number of connected heavy atoms, aromatic boolean, ring boolean, atomic number, and summation of bond orders connected to the atom. Vector is indexed from zero.
void GetGIDVector (std::vector< unsigned int > &)
 Calculates a set of symmetry identifiers for a molecule. Atoms with the same symmetry ID are symmetrically equivalent. Vector is indexed from zero.
Methods to check for existence of properties
bool Has2D ()
 Are there non-zero coordinates in two dimensions (i.e. X and Y)?
bool Has3D ()
 Are there non-zero coordinates in all three dimensions (i.e. X, Y, Z)?
bool HasNonZeroCoords ()
 Are there any non-zero coordinates?
bool HasAromaticPerceived ()
bool HasSSSRPerceived ()
bool HasRingAtomsAndBondsPerceived ()
bool HasAtomTypesPerceived ()
bool HasChiralityPerceived ()
bool HasPartialChargesPerceived ()
bool HasHybridizationPerceived ()
bool HasImplicitValencePerceived ()
bool HasKekulePerceived ()
bool HasClosureBondsPerceived ()
bool HasChainsPerceived ()
bool HasHydrogensAdded ()
bool HasAromaticCorrected ()
bool IsCorrectedForPH ()
bool HasSpinMultiplicityAssigned ()
bool IsChiral ()
 Is this molecule chiral?
bool Empty ()
 Are there any atoms in this molecule?
Multiple conformer member functions
int NumConformers ()
void SetConformers (std::vector< double * > &v)
void AddConformer (double *f)
void SetConformer (int i)
void CopyConformer (double *, int)
void DeleteConformer (int)
double * GetConformer (int i)
double * BeginConformer (std::vector< double * >::iterator &i)
double * NextConformer (std::vector< double * >::iterator &i)
std::vector< double * > & GetConformers ()
Iterator methods
OBAtomBeginAtom (std::vector< OBNodeBase * >::iterator &i)
OBAtomNextAtom (std::vector< OBNodeBase * >::iterator &i)
OBBondBeginBond (std::vector< OBEdgeBase * >::iterator &i)
OBBondNextBond (std::vector< OBEdgeBase * >::iterator &i)
OBResidueBeginResidue (std::vector< OBResidue * >::iterator &i)
OBResidueNextResidue (std::vector< OBResidue * >::iterator &i)
OBInternalCoordBeginInternalCoord (std::vector< OBInternalCoord * >::iterator &i)
OBInternalCoordNextInternalCoord (std::vector< OBInternalCoord * >::iterator &i)

Protected Member Functions

bool HasFlag (int flag)
void SetFlag (int flag)
Internal Kekulization routines -- see kekulize.cpp and NewPerceiveKekuleBonds()
void start_kekulize (std::vector< OBAtom * > &cycle, std::vector< int > &electron)
 Start kekulizing one or a fused set of aromatic ring(s).
int expand_kekulize (OBAtom *atom1, OBAtom *atom2, std::vector< int > &currentState, std::vector< int > &initState, std::vector< int > &bcurrentState, std::vector< int > &binitState, std::vector< bool > &mark)
 recursively assign single and double bonds according to the electronical state of the atoms of the current bond
int getorden (OBAtom *atom)
 Give the priority to give two electrons instead of 1.
void expandcycle (OBAtom *atom, OBBitVec &avisit)
 Recursively find the aromatic atoms with an aromatic bond to the current atom.

Protected Attributes

int _flags
 bitfield of flags
bool _autoPartialCharge
 Assign partial charges automatically.
bool _autoFormalCharge
 Assign formal charges automatically.
std::string _title
 Molecule title.
unsigned short int _dimension
 Dimensionality of coordinates.
double _energy
 Molecular heat of formation (if applicable).
int _totalCharge
 Total charge on the molecule.
unsigned int _totalSpin
 Total spin on the molecule (if not specified, assumes lowest possible spin).
double * _c
 coordinate array
std::vector< double * > _vconf
 vector of conformers
unsigned short int _natoms
 Number of atoms.
unsigned short int _nbonds
 Number of bonds.
std::vector< OBResidue * > _residue
 Residue information (if applicable).
std::vector< OBInternalCoord * > _internals
 Internal Coordinates (if applicable).
std::vector< OBGenericData * > _vdata
 Custom data -- see OBGenericData class for more.
unsigned short int _mod
 Number of nested calls to BeginModify().
bool _vlock
std::vector< OBNodeBase * > _vatom
std::vector< OBEdgeBase * > _vbond

Detailed Description

Molecule Class.

The most important class in Open Babel is OBMol, or the molecule class. The OBMol class is designed to store all the basic information associated with a molecule, to make manipulations on the connection table of a molecule facile, and to provide member functions which automatically perceive information about a molecule. A guided tour of the OBMol class is a good place to start.

An OBMol class can be declared:

      OBMol mol;

For example:

      #include <iostream.h>
      #include "mol.h"
      #include "obconversion.h"
      int main(int argc,char **argv)
      {
      OBConversion conv(&cin,&cout);
      if(conv.SetInAndOutFormats("SDF","MOL2"))
      { 
      OBMol mol;
      if(conv.Read(&mol))
      ...manipulate molecule 
                
      conv->Write(&mol);
      }
      return(1);
      }

will read in a molecule in SD file format from stdin (or the C++ equivalent cin) and write a MOL2 format file out to standard out. Additionally, The input and output formats can be altered using the OBConversion class

Once a molecule has been read into an OBMol (or created via other methods) the atoms and bonds can be accessed by the following methods:

      OBAtom *atom;
      atom = mol.GetAtom(5); //random access of an atom
or
      OBBond *bond;
      bond = mol.GetBond(14); //random access of a bond
or
      FOR_ATOMS_IN_MOL(atom, mol) // iterator access (see OBMolAtomIter)
or
      FOR_BONDS_IN_MOL(bond, mol) // iterator access (see OBMolBondIter)
It is important to note that atom arrays currently begin at 1 and bond arrays begin at 0. Requesting atom 0 (
      OBAtom *atom = mol.GetAtom(0); 
will result in an error, but
      OBBond *bond = mol.GetBond(0);
is perfectly valid. Note that this is expected to change in the near future to simplify coding and improve efficiency.

The ambiguity of numbering issues and off-by-one errors led to the use of iterators in Open Babel. An iterator is essentially just a pointer, but when used in conjunction with Standard Template Library (STL) vectors it provides an unambiguous way to loop over arrays. OBMols store their atom and bond information in STL vectors. Since vectors are template based, a vector of any user defined type can be declared. OBMols declare vector<OBNodeBase*> and vector<OBEdgeBase*> to store atom and bond information. Iterators are then a natural way to loop over the vectors of atoms and bonds.

A variety of predefined iterators have been created to simplify common looping requests (e.g., looping over all atoms in a molecule, bonds to a given atom, etc.)

      #include "obiter.h"
      ...
      #define FOR_ATOMS_OF_MOL(a,m)     for( OBMolAtomIter     a(m); a; a++ )
      #define FOR_BONDS_OF_MOL(b,m)     for( OBMolBondIter     b(m); b; b++ )
      #define FOR_NBORS_OF_ATOM(a,p)    for( OBAtomAtomIter    a(p); a; a++ )
      #define FOR_BONDS_OF_ATOM(b,p)    for( OBAtomBondIter    b(p); b; b++ )
      #define FOR_RESIDUES_OF_MOL(r,m)  for( OBResidueIter     r(m); r; r++ )
      #define FOR_ATOMS_OF_RESIDUE(a,r) for( OBResidueAtomIter a(r); a; a++ )
      ...

These convenience functions can be used like so:

      #include "obiter.h"
      #include "mol.h"

      OBMol mol;
      double exactMass = 0.0f;
      FOR_ATOMS_OF_MOL(a, mol)
      {
      exactMass +=  a->GetExactMass();
      }

Note that with these convenience macros, the iterator "a" (or whichever name you pick) is declared for you -- you do not need to do it beforehand.


Constructor & Destructor Documentation

OBMol (  ) 

Constructor.

OBMol ( const OBMol  ) 

Copy constructor

Warning:
Currently does not copy all associated OBGenericData This requires a (minor) API change, and will thus only be fixed in 2.1 or later releases.

~OBMol (  )  [virtual]

Destructor.


Member Function Documentation

bool HasFlag ( int  flag  )  [inline, protected]

void SetFlag ( int  flag  )  [inline, protected]

void start_kekulize ( std::vector< OBAtom * > &  cycle,
std::vector< int > &  electron 
) [protected]

Start kekulizing one or a fused set of aromatic ring(s).

The initial electronic state indicates if an atoms must make a double bond or not Kekulizing is attempted recursively for all the atoms bonded to the first atom of the cycle.

int expand_kekulize ( OBAtom atom1,
OBAtom atom2,
std::vector< int > &  currentState,
std::vector< int > &  initState,
std::vector< int > &  bcurrentState,
std::vector< int > &  binitState,
std::vector< bool > &  mark 
) [protected]

recursively assign single and double bonds according to the electronical state of the atoms of the current bond

int getorden ( OBAtom atom  )  [protected]

Give the priority to give two electrons instead of 1.

void expandcycle ( OBAtom atom,
OBBitVec avisit 
) [protected]

Recursively find the aromatic atoms with an aromatic bond to the current atom.

OBMol & operator= ( const OBMol mol  ) 

Assignment.

OBMol & operator+= ( const OBMol mol  ) 

void ReserveAtoms ( int  natoms  )  [inline]

OBAtom * CreateAtom ( void   )  [virtual]

OBBond * CreateBond ( void   )  [virtual]

void DestroyAtom ( OBNodeBase  )  [virtual]

void DestroyBond ( OBEdgeBase  )  [virtual]

bool AddAtom ( OBAtom atom  ) 

Add an atom to a molecule.

Also checks bond_queue for any bonds that should be made to the new atom

bool AddBond ( int  ,
int  ,
int  ,
int  flags = 0,
int  insertpos = -1 
)

bool AddBond ( OBBond  ) 

bool AddResidue ( OBResidue  ) 

bool InsertAtom ( OBAtom  ) 

bool DeleteAtom ( OBAtom  ) 

bool DeleteBond ( OBBond  ) 

bool DeleteResidue ( OBResidue  ) 

OBAtom * NewAtom (  ) 

Get a new atom to add to a molecule.

Also checks bond_queue for any bonds that should be made to the new atom

OBResidue * NewResidue (  ) 

void BeginModify ( void   )  [virtual]

Call when making many modifications -- clears conformer/rotomer data.

void EndModify ( bool  nukePerceivedData = true  )  [virtual]

Call when done with modificaions -- re-perceive data as needed.

int GetMod (  )  [inline]

void IncrementMod (  )  [inline]

void DecrementMod (  )  [inline]

bool HasData ( std::string &   ) 

Returns:
whether the generic attribute/value pair exists

bool HasData ( const char *   ) 

Returns:
whether the generic attribute/value pair exists

bool HasData ( unsigned int  type  ) 

Returns:
whether the generic attribute/value pair exists

void DeleteData ( unsigned int  type  ) 

void DeleteData ( OBGenericData  ) 

void DeleteData ( std::vector< OBGenericData * > &   ) 

void SetData ( OBGenericData d  )  [inline]

unsigned int DataSize (  )  [inline]

Returns:
the number of OBGenericData items attached to this molecule.

OBGenericData * GetData ( unsigned int  type  ) 

OBGenericData * GetData ( std::string &   ) 

Returns the value given an attribute name.

OBGenericData * GetData ( const char *   ) 

Returns the value given an attribute name.

std::vector<OBGenericData*>& GetData (  )  [inline]

std::vector<OBGenericData*>::iterator BeginData (  )  [inline]

std::vector<OBGenericData*>::iterator EndData (  )  [inline]

int GetFlags (  )  [inline]

const char* GetTitle (  )  const [inline]

Returns:
the title of this molecule (often the filename)

unsigned int NumAtoms (  )  const [inline]

Returns:
the number of atoms (i.e. OBAtom children)

unsigned int NumBonds (  )  const [inline]

Returns:
the number of bonds (i.e. OBBond children)

unsigned int NumHvyAtoms (  ) 

Returns:
the number of non-hydrogen atoms

unsigned int NumResidues (  )  const [inline]

Returns:
the number of residues (i.e. OBResidue substituents)

unsigned int NumRotors (  ) 

Returns:
the number of rotatble bonds. See OBBond::IsRotor() for details

OBAtom * GetAtom ( int  idx  ) 

Returns a pointer to the atom after a safety check 0 < idx <= NumAtoms

OBAtom * GetFirstAtom (  ) 

OBBond * GetBond ( int  idx  ) 

Returns a pointer to the bond after a safety check 0 <= idx < NumBonds

OBBond * GetBond ( int  ,
int   
)

OBBond * GetBond ( OBAtom ,
OBAtom  
)

OBResidue * GetResidue ( int   ) 

std::vector< OBInternalCoord * > GetInternalCoord (  ) 

double GetTorsion ( int  ,
int  ,
int  ,
int   
)

Returns:
the dihedral angle between the four atoms supplied a1-a2-a3-a4)

double GetTorsion ( OBAtom ,
OBAtom ,
OBAtom ,
OBAtom  
)

Returns:
the dihedral angle between the four atoms supplied a1-a2-a3-a4)

string GetFormula (  ) 

Returns:
the stochoimetric formula (e.g., C4H6O)

Stochoimetric formula (e.g., C4H6O). This is either set by OBMol::SetFormula() or generated on-the-fly using the "Hill order" -- i.e., C first if present, then H if present all other elements in alphabetical order.

double GetEnergy (  )  const [inline]

Returns:
the heat of formation for this molecule (in kcal/mol)

double GetMolWt (  ) 

Returns:
the standard molar mass given by IUPAC atomic masses (amu)

double GetExactMass (  ) 

Returns:
the mass given by isotopes (or most abundant isotope, if not specified)

int GetTotalCharge (  ) 

Returns:
the total charge on this molecule (i.e., 0 = neutral, +1, -1...)

Returns the total molecular charge -- if it has not previously been set it is calculated from the atomic formal charge information. (This may or may not be correct!) If you set atomic charges with OBAtom::SetFormalCharge() you really should set the molecular charge with OBMol::SetTotalCharge()

unsigned int GetTotalSpinMultiplicity (  ) 

Returns:
the total spin on this molecule (i.e., 1 = singlet, 2 = doublet...)

Returns the total spin multiplicity -- if it has not previously been set it is calculated from the atomic spin multiplicity information assuming the high-spin case (i.e. it simply sums the atomic spins, making no attempt to pair spins). However, if you set atomic spins with OBAtom::SetSpinMultiplicity() you really should set the molecular spin with OBMol::SetTotalSpinMultiplicity()

unsigned short int GetDimension (  )  const [inline]

Returns:
the dimensionality of coordinates (i.e., 0 = unknown or no coord, 2=2D, 3=3D)

double* GetCoordinates (  )  [inline]

vector< OBRing * > & GetSSSR (  ) 

Implements blue-obelisk:findSmallestSetOfSmallestRings.

bool AutomaticFormalCharge (  )  [inline]

Get the current flag for whether formal charges are set with pH correction.

bool AutomaticPartialCharge (  )  [inline]

Get the current flag for whether partial charges are auto-determined.

void SetTitle ( const char *  title  ) 

void SetTitle ( std::string &  title  ) 

void SetFormula ( std::string  molFormula  ) 

Set the stochiometric formula for this molecule.

void SetEnergy ( double  energy  )  [inline]

Set the heat of formation for this molecule (in kcal/mol).

void SetDimension ( unsigned short int  d  )  [inline]

Set the dimension of this molecule (i.e., 0, 1 , 2, 3).

void SetTotalCharge ( int  charge  ) 

void SetTotalSpinMultiplicity ( unsigned int  spin  ) 

void SetInternalCoord ( std::vector< OBInternalCoord * >  int_coord  )  [inline]

void SetAutomaticFormalCharge ( bool  val  )  [inline]

Set the flag for determining automatic formal charges with pH (default=true).

void SetAutomaticPartialCharge ( bool  val  )  [inline]

Set the flag for determining partial charges automatically (default=true).

void SetAromaticPerceived (  )  [inline]

Mark that aromaticity has been perceived for this molecule (see OBAromaticTyper).

void SetSSSRPerceived (  )  [inline]

Mark that Smallest Set of Smallest Rings has been run (see OBRing class).

void SetRingAtomsAndBondsPerceived (  )  [inline]

Mark that rings have been perceived (see OBRing class for details).

void SetAtomTypesPerceived (  )  [inline]

Mark that atom types have been perceived (see OBAtomTyper for details).

void SetChainsPerceived (  )  [inline]

Mark that chains and residues have been perceived (see OBChainsParser).

void SetChiralityPerceived (  )  [inline]

Mark that chirality has been perceived.

void SetPartialChargesPerceived (  )  [inline]

Mark that partial charges have been assigned.

void SetHybridizationPerceived (  )  [inline]

void SetImplicitValencePerceived (  )  [inline]

void SetKekulePerceived (  )  [inline]

void SetClosureBondsPerceived (  )  [inline]

void SetHydrogensAdded (  )  [inline]

void SetCorrectedForPH (  )  [inline]

void SetAromaticCorrected (  )  [inline]

void SetSpinMultiplicityAssigned (  )  [inline]

void SetFlags ( int  flags  )  [inline]

void UnsetAromaticPerceived (  )  [inline]

void UnsetPartialChargesPerceived (  )  [inline]

void UnsetImplicitValencePerceived (  )  [inline]

void UnsetFlag ( int  flag  )  [inline]

OBBase * DoTransformations ( const std::map< std::string, std::string > *  pOptions  )  [virtual]

Reimplemented from OBBase.

const char * ClassDescription (  )  [static]

Reimplemented from OBBase.

bool Clear (  ) 

Clear all information from a molecule.

void RenumberAtoms ( std::vector< OBNodeBase * > &  v  ) 

Renumber the atoms of this molecule according to the order in the supplied vector.

Renumber the atoms in this molecule according to the order in the supplied vector. This will return without action if the supplied vector is empty or does not have the same number of atoms as the molecule.

void ToInertialFrame ( int  conf,
double *  rmat 
)

Translate one conformer and rotate by a rotation matrix (which is returned) to the inertial frame-of-reference.

void ToInertialFrame (  ) 

Translate all conformers to the inertial frame-of-reference.

void Translate ( const vector3 v  ) 

Translates all conformers in the molecule by the supplied vector.

this method adds the vector v to all atom positions in all conformers

void Translate ( const vector3 v,
int  nconf 
)

Translates one conformer in the molecule by the supplied vector.

this method adds the vector v to all atom positions in the conformer nconf. If nconf == OB_CURRENT_CONFORMER, then the atom positions in the current conformer are translated.

void Rotate ( const double  u[3][3]  ) 

void Rotate ( const double  m[9]  ) 

void Rotate ( const double  m[9],
int  nconf 
)

void Center (  ) 

Translate to the center of all coordinates (for this conformer).

bool Kekulize (  ) 

Transform to standard Kekule bond structure (presumably from an aromatic form).

bool PerceiveKekuleBonds (  ) 

void NewPerceiveKekuleBonds (  ) 

Kekulize aromatic rings without using implicit valence.

This new perceive kekule bonds function has been especifically designed to handle molecule files without explicit hydrogens such as pdb or xyz. The function does not rely on GetImplicitValence function The function looks for groups of aromatic cycle For each group it tries to guess the number of electrons given by each atom in order to satisfy the huckel (4n+2) rule If the huckel rule cannot be satisfied the algorithm try with its best alternative guess Then it recursively walk on the atoms of the cycle and assign single and double bonds

bool DeleteHydrogen ( OBAtom  ) 

bool DeleteHydrogens (  ) 

bool DeleteHydrogens ( OBAtom  ) 

bool DeleteNonPolarHydrogens (  ) 

bool AddHydrogens ( bool  polaronly = false,
bool  correctForPH = true 
)

bool AddHydrogens ( OBAtom  ) 

bool AddPolarHydrogens (  ) 

bool StripSalts (  ) 

Deletes all atoms except for the largest contiguous fragment.

bool ConvertDativeBonds (  ) 

Converts for instance [N+]([O-])=O to N(=O)=O.

bool CorrectForPH (  ) 

bool AssignSpinMultiplicity (  ) 

set spin multiplicity for H-deficient atoms

vector3 Center ( int  nconf  ) 

void SetTorsion ( OBAtom ,
OBAtom ,
OBAtom ,
OBAtom ,
double   
)

Set the torsion defined by these atoms, rotating bonded neighbors.

void FindSSSR (  ) 

Find Smallest Set of Smallest Rings (see OBRing class for more details).

void FindRingAtomsAndBonds (  ) 

void FindChiralCenters (  ) 

Sets atom->IsChiral() to true for chiral atoms.

void FindChildren ( std::vector< int > &  children,
int  first,
int  second 
)

locates all atoms for which there exists a path to 'second' without going through 'first' children must not include 'second'

void FindChildren ( std::vector< OBAtom * > &  children,
OBAtom bgn,
OBAtom end 
)

locates all atoms for which there exists a path to 'end' without going through 'bgn' children must not include 'end'

void FindLargestFragment ( OBBitVec  ) 

void ContigFragList ( std::vector< std::vector< int > > &   ) 

Sort a list of contig fragments by size from largest to smallest Each vector<int> contains the atom numbers of a contig fragment

void Align ( OBAtom ,
OBAtom ,
vector3 ,
vector3  
)

Aligns atom a on p1 and atom b along p1->p2 vector.

void ConnectTheDots ( void   ) 

Adds single bonds based on atom proximity.

This method adds single bonds between all atoms closer than their combined atomic covalent radii, then "cleans up" making sure bonded atoms are not closer than 0.4A and the atom does not exceed its valence. It implements blue-obelisk:rebondFrom3DCoordinates.

void PerceiveBondOrders (  ) 

Attempts to perceive multiple bonds based on geometries.

This method uses bond angles and geometries from current connectivity to guess atom types and then filling empty valences with multiple bonds. It currently has a pass to detect some frequent functional groups. It still needs a pass to detect aromatic rings to "clean up."

void FindTorsions (  ) 

Fills the OBGeneric OBTorsionData with torsions from the mol.

bool GetGTDVector ( std::vector< int > &   ) 

Calculates the graph theoretical distance of each atom. Vector is indexed from zero.

void GetGIVector ( std::vector< unsigned int > &   ) 

Calculates a set of graph invariant indexes using the graph theoretical distance, number of connected heavy atoms, aromatic boolean, ring boolean, atomic number, and summation of bond orders connected to the atom. Vector is indexed from zero.

void GetGIDVector ( std::vector< unsigned int > &   ) 

Calculates a set of symmetry identifiers for a molecule. Atoms with the same symmetry ID are symmetrically equivalent. Vector is indexed from zero.

bool Has2D (  ) 

Are there non-zero coordinates in two dimensions (i.e. X and Y)?

bool Has3D (  ) 

Are there non-zero coordinates in all three dimensions (i.e. X, Y, Z)?

bool HasNonZeroCoords (  ) 

Are there any non-zero coordinates?

bool HasAromaticPerceived (  )  [inline]

bool HasSSSRPerceived (  )  [inline]

bool HasRingAtomsAndBondsPerceived (  )  [inline]

bool HasAtomTypesPerceived (  )  [inline]

bool HasChiralityPerceived (  )  [inline]

bool HasPartialChargesPerceived (  )  [inline]

bool HasHybridizationPerceived (  )  [inline]

bool HasImplicitValencePerceived (  )  [inline]

bool HasKekulePerceived (  )  [inline]

bool HasClosureBondsPerceived (  )  [inline]

bool HasChainsPerceived (  )  [inline]

bool HasHydrogensAdded (  )  [inline]

bool HasAromaticCorrected (  )  [inline]

bool IsCorrectedForPH (  )  [inline]

bool HasSpinMultiplicityAssigned (  )  [inline]

bool IsChiral (  ) 

Is this molecule chiral?

bool Empty (  )  [inline]

Are there any atoms in this molecule?

int NumConformers (  )  [inline]

void SetConformers ( std::vector< double * > &  v  ) 

void AddConformer ( double *  f  )  [inline]

void SetConformer ( int  i  )  [inline]

void CopyConformer ( double *  ,
int   
)

void DeleteConformer ( int   ) 

double* GetConformer ( int  i  )  [inline]

double* BeginConformer ( std::vector< double * >::iterator &  i  )  [inline]

double* NextConformer ( std::vector< double * >::iterator &  i  )  [inline]

std::vector<double*>& GetConformers (  )  [inline]

OBAtom * BeginAtom ( std::vector< OBNodeBase * >::iterator &  i  ) 

Deprecated:
Use FOR_ATOMS_OF_MOL and OBMolAtomIter instead

OBAtom * NextAtom ( std::vector< OBNodeBase * >::iterator &  i  ) 

Deprecated:
Use FOR_ATOMS_OF_MOL and OBMolAtomIter instead

OBBond * BeginBond ( std::vector< OBEdgeBase * >::iterator &  i  ) 

Deprecated:
Use FOR_BONDS_OF_MOL and OBMolBondIter instead

OBBond * NextBond ( std::vector< OBEdgeBase * >::iterator &  i  ) 

Deprecated:
Use FOR_BONDS_OF_MOL and OBMolBondIter instead

OBResidue* BeginResidue ( std::vector< OBResidue * >::iterator &  i  )  [inline]

Deprecated:
Use FOR_RESIDUES_OF_MOL and OBResidueIter instead

OBResidue* NextResidue ( std::vector< OBResidue * >::iterator &  i  )  [inline]

Deprecated:
Use FOR_RESIDUES_OF_MOL and OBResidueIter instead

OBInternalCoord* BeginInternalCoord ( std::vector< OBInternalCoord * >::iterator &  i  )  [inline]

OBInternalCoord* NextInternalCoord ( std::vector< OBInternalCoord * >::iterator &  i  )  [inline]

unsigned int NumNodes (  )  [inline, inherited]

unsigned int NumEdges (  )  [inline, inherited]

void ResetVisitFlags (  )  [inherited]

bool SetVisitLock ( bool   )  [inherited]

bool GetVisitLock (  )  [inline, inherited]

OBNodeBase * Begin ( std::vector< OBNodeBase * >::iterator &   )  [inherited]

OBEdgeBase * Begin ( std::vector< OBEdgeBase * >::iterator &   )  [inherited]

OBNodeBase * Next ( std::vector< OBNodeBase * >::iterator &   )  [inherited]

OBEdgeBase * Next ( std::vector< OBEdgeBase * >::iterator &   )  [inherited]


Member Data Documentation

int _flags [protected]

bitfield of flags

bool _autoPartialCharge [protected]

Assign partial charges automatically.

bool _autoFormalCharge [protected]

Assign formal charges automatically.

std::string _title [protected]

Molecule title.

unsigned short int _dimension [protected]

Dimensionality of coordinates.

double _energy [protected]

Molecular heat of formation (if applicable).

int _totalCharge [protected]

Total charge on the molecule.

unsigned int _totalSpin [protected]

Total spin on the molecule (if not specified, assumes lowest possible spin).

double* _c [protected]

coordinate array

std::vector<double*> _vconf [protected]

vector of conformers

unsigned short int _natoms [protected]

Number of atoms.

unsigned short int _nbonds [protected]

Number of bonds.

std::vector<OBResidue*> _residue [protected]

Residue information (if applicable).

std::vector<OBInternalCoord*> _internals [protected]

Internal Coordinates (if applicable).

std::vector<OBGenericData*> _vdata [protected]

Custom data -- see OBGenericData class for more.

unsigned short int _mod [protected]

Number of nested calls to BeginModify().

bool _vlock [protected, inherited]

std::vector<OBNodeBase*> _vatom [protected, inherited]

std::vector<OBEdgeBase*> _vbond [protected, inherited]


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