OBMol Class Reference

Molecule Class. More...

#include <openbabel/mol.h>

Inheritance diagram for OBMol:

Inheritance graph
[legend]
List of all members.

Molecule modification methods

virtual OBBaseDoTransformations (const std::map< std::string, std::string > *pOptions)
bool Clear ()
void RenumberAtoms (std::vector< OBAtom * > &)
void ToInertialFrame (int conf, double *rmat)
void ToInertialFrame ()
void Translate (const vector3 &v)
void Translate (const vector3 &v, int conf)
void Rotate (const double u[3][3])
void Rotate (const double m[9])
void Rotate (const double m[9], int nconf)
void Center ()
bool Kekulize ()
bool PerceiveKekuleBonds ()
void NewPerceiveKekuleBonds ()
bool DeleteHydrogens ()
bool DeleteHydrogens (OBAtom *)
bool DeleteNonPolarHydrogens ()
bool DeleteHydrogen (OBAtom *)
bool AddHydrogens (bool polaronly=false, bool correctForPH=true)
bool AddHydrogens (OBAtom *)
bool AddPolarHydrogens ()
bool StripSalts ()
std::vector< OBMolSeparate (int StartIndex=1)
bool ConvertDativeBonds ()
bool CorrectForPH ()
bool AssignSpinMultiplicity ()
vector3 Center (int nconf)
void SetTorsion (OBAtom *, OBAtom *, OBAtom *, OBAtom *, double)
static const char * ClassDescription ()

Public Member Functions

template<class T>
T * CastAndClear (bool clear=true)
Initialization and data (re)size methods
 OBMol ()
 OBMol (const OBMol &)
virtual ~OBMol ()
OBMoloperator= (const OBMol &mol)
OBMoloperator+= (const OBMol &mol)
void ReserveAtoms (int natoms)
virtual OBAtomCreateAtom (void)
virtual OBBondCreateBond (void)
virtual OBResidueCreateResidue (void)
virtual void DestroyAtom (OBAtom *)
virtual void DestroyBond (OBBond *)
virtual void DestroyResidue (OBResidue *)
bool AddAtom (OBAtom &)
bool InsertAtom (OBAtom &)
bool AddBond (int beginIdx, int endIdx, int order, int flags=0, int insertpos=-1)
bool AddBond (OBBond &)
bool AddResidue (OBResidue &)
virtual OBAtomNewAtom ()
virtual OBBondNewBond ()
virtual OBResidueNewResidue ()
bool DeleteAtom (OBAtom *)
bool DeleteBond (OBBond *)
bool DeleteResidue (OBResidue *)
Molecule modification methods
virtual void BeginModify (void)
virtual void EndModify (bool nukePerceivedData=true)
int GetMod ()
void IncrementMod ()
void DecrementMod ()
Data retrieval methods
int GetFlags ()
const char * GetTitle () const
unsigned int NumAtoms () const
unsigned int NumBonds () const
unsigned int NumHvyAtoms ()
unsigned int NumResidues () const
unsigned int NumRotors ()
OBAtomGetAtom (int idx)
OBAtomGetFirstAtom ()
OBBondGetBond (int idx)
OBBondGetBond (int a, int b)
OBBondGetBond (OBAtom *bgn, OBAtom *end)
OBResidueGetResidue (int)
std::vector< OBInternalCoord * > GetInternalCoord ()
double GetTorsion (int, int, int, int)
double GetTorsion (OBAtom *a, OBAtom *b, OBAtom *c, OBAtom *d)
double GetAngle (OBAtom *a, OBAtom *b, OBAtom *c)
std::string GetFormula ()
std::string GetSpacedFormula (int ones=0, const char *sp=" ")
double GetEnergy () const
double GetMolWt ()
double GetExactMass ()
int GetTotalCharge ()
unsigned int GetTotalSpinMultiplicity ()
unsigned short int GetDimension () const
double * GetCoordinates ()
std::vector< OBRing * > & GetSSSR ()
bool AutomaticFormalCharge ()
bool AutomaticPartialCharge ()
Data modification methods
void SetTitle (const char *title)
void SetTitle (std::string &title)
void SetFormula (std::string molFormula)
void SetEnergy (double energy)
void SetDimension (unsigned short int d)
void SetTotalCharge (int charge)
void SetTotalSpinMultiplicity (unsigned int spin)
void SetInternalCoord (std::vector< OBInternalCoord * > int_coord)
void SetAutomaticFormalCharge (bool val)
void SetAutomaticPartialCharge (bool val)
void SetAromaticPerceived ()
void SetSSSRPerceived ()
void SetRingAtomsAndBondsPerceived ()
void SetAtomTypesPerceived ()
void SetChainsPerceived ()
void SetChiralityPerceived ()
void SetPartialChargesPerceived ()
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 UnsetHydrogensAdded ()
void UnsetFlag (int flag)
Molecule utilities and perception methods
void FindSSSR ()
void FindRingAtomsAndBonds ()
void FindChiralCenters ()
void FindChildren (std::vector< int > &children, int bgnIdx, int endIdx)
void FindChildren (std::vector< OBAtom * > &children, OBAtom *bgn, OBAtom *end)
void FindLargestFragment (OBBitVec &frag)
void ContigFragList (std::vector< std::vector< int > > &)
void Align (OBAtom *, OBAtom *, vector3 &, vector3 &)
void ConnectTheDots ()
void PerceiveBondOrders ()
void FindAngles ()
void FindTorsions ()
bool GetGTDVector (std::vector< int > &)
void GetGIVector (std::vector< unsigned int > &)
void GetGIDVector (std::vector< unsigned int > &)
Methods to check for existence of properties
bool Has2D ()
bool Has3D ()
bool HasNonZeroCoords ()
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 ()
bool Empty ()
Multiple conformer member functions
int NumConformers ()
void SetConformers (std::vector< double * > &v)
void AddConformer (double *f)
void SetConformer (int i)
void CopyConformer (double *c, int nconf)
void DeleteConformer (int nconf)
double * GetConformer (int i)
double * BeginConformer (std::vector< double * >::iterator &i)
double * NextConformer (std::vector< double * >::iterator &i)
std::vector< double * > & GetConformers ()
Iterator methods
OBAtomIterator BeginAtoms ()
OBAtomIterator EndAtoms ()
OBBondIterator BeginBonds ()
OBBondIterator EndBonds ()
OBResidueIterator BeginResidues ()
OBResidueIterator EndResidues ()
OBAtomBeginAtom (OBAtomIterator &i)
OBAtomNextAtom (OBAtomIterator &i)
OBBondBeginBond (OBBondIterator &i)
OBBondNextBond (OBBondIterator &i)
OBResidueBeginResidue (OBResidueIterator &i)
OBResidueNextResidue (OBResidueIterator &i)
OBInternalCoordBeginInternalCoord (std::vector< OBInternalCoord * >::iterator &i)
OBInternalCoordNextInternalCoord (std::vector< OBInternalCoord * >::iterator &i)
Generic data handling methods (via OBGenericData)
bool HasData (const std::string &)
bool HasData (const char *)
bool HasData (const unsigned int type)
void DeleteData (unsigned int type)
void DeleteData (OBGenericData *)
void DeleteData (std::vector< OBGenericData * > &)
void SetData (OBGenericData *d)
unsigned int DataSize () const
OBGenericDataGetData (const unsigned int type)
OBGenericDataGetData (const std::string &)
OBGenericDataGetData (const char *)
std::vector< OBGenericData * > & GetData ()
std::vector< OBGenericData * > GetData (DataOrigin source)
OBDataIterator BeginData ()
OBDataIterator EndData ()

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)
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)
int getorden (OBAtom *atom)
void expandcycle (OBAtom *atom, OBBitVec &avisit)

Protected Attributes

int _flags
bool _autoPartialCharge
bool _autoFormalCharge
std::string _title
std::vector< OBAtom * > _vatom
std::vector< OBBond * > _vbond
unsigned short int _dimension
double _energy
int _totalCharge
unsigned int _totalSpin
double * _c
std::vector< double * > _vconf
unsigned int _natoms
unsigned int _nbonds
std::vector< OBResidue * > _residue
std::vector< OBInternalCoord * > _internals
unsigned short int _mod
std::vector< OBGenericData * > _vdata

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 <openbabel/mol.h>
      #include <openbabel/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_OF_MOL(atom, mol) // iterator access (see OBMolAtomIter)
or
      FOR_BONDS_OF_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<OBAtom*> and vector<OBBond*> 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 <openbabel/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 <openbabel/obiter.h>
      #include <openbabel/mol.h>

      OBMol mol;
      double exactMass = 0.0;
      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, copies atoms,bonds and OBGenericData.

~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, copies atoms,bonds and OBGenericData.

OBMol & operator+= ( const OBMol mol  ) 

Copies atoms and bonds but not OBGenericData.

void ReserveAtoms ( int  natoms  )  [inline]

Reserve a minimum number of atoms for internal storage This improves performance since the internal atom vector does not grow.

OBAtom * CreateAtom ( void   )  [virtual]

Create a new OBAtom pointer. Does no bookkeeping

Deprecated:
Use NewAtom instead, which ensures internal connections

OBBond * CreateBond ( void   )  [virtual]

Create a new OBBond pointer. Does no bookkeeping

Deprecated:
Use NewBond instead, which ensures internal connections

OBResidue * CreateResidue ( void   )  [virtual]

Create a new OBResidue pointer. Does no bookkeeping

Deprecated:
Use NewResidue instead, which ensures internal connections

void DestroyAtom ( OBAtom  )  [virtual]

Free an OBAtom pointer if defined. Does no bookkeeping

See also:
DeleteAtom which ensures internal connections

void DestroyBond ( OBBond  )  [virtual]

Free an OBBond pointer if defined. Does no bookkeeping

See also:
DeleteBond which ensures internal connections

void DestroyResidue ( OBResidue  )  [virtual]

Free an OBResidue pointer if defined. Does no bookkeeping

See also:
DeleteResidue which ensures internal connections

bool AddAtom ( OBAtom  ) 

Add an atom to a molecule.

Add the specified atom to this molecule

Returns:
Whether the method was successful

bool InsertAtom ( OBAtom  ) 

Add a new atom to this molecule (like AddAtom) Calls BeginModify() before insertion and EndModify() after insertion

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

Add a new bond to the molecule with the specified parameters

Parameters:
beginIdx the atom index of the "start" atom
endIdx the atom index of the "end" atom
order the bond order (see OBBond::GetBO())
flags any bond flags such as stereochemistry (default = none)
insertpos the position index to insert the bond (default = none)
Returns:
Whether the new bond creation was successful

bool AddBond ( OBBond  ) 

Add the specified residue to this molecule and update connections

Returns:
Whether the method was successful

bool AddResidue ( OBResidue  ) 

Add the specified residue to this molecule and update connections

Returns:
Whether the method was successful

OBAtom * NewAtom (  )  [virtual]

Instantiate a New Atom and add it to the molecule.

Create a new OBAtom in this molecule and ensure connections. (e.g. OBAtom::GetParent()

OBBond * NewBond (  )  [virtual]

Instantiate a New Bond and add it to the molecule.

Create a new OBBond in this molecule and ensure connections. (e.g. OBBond::GetParent()

OBResidue * NewResidue (  )  [virtual]

Create a new OBResidue in this molecule and ensure connections.

bool DeleteAtom ( OBAtom  ) 

Deletes an atom from this molecule and all appropriate bonds. Updates the molecule and atom and bond indexes accordingly.

Warning:
Does not update any residues which may contain this atom
Returns:
Whether deletion was successful

bool DeleteBond ( OBBond  ) 

Deletes an bond from this molecule and updates accordingly

Returns:
Whether deletion was successful

bool DeleteResidue ( OBResidue  ) 

Deletes a residue from this molecule and updates accordingly.

Returns:
Whether deletion was successful

void BeginModify ( void   )  [virtual]

Call when making many modifications -- clears conformer/rotomer data. The method "turns off" perception routines, improving performance. Changes in molecular structure will be re-considered after modifications.

void EndModify ( bool  nukePerceivedData = true  )  [virtual]

Call when done with modificaions -- re-perceive data as needed. This method "turns on" perception routines and re-evaluates molecular structure.

int GetMod (  )  [inline]

Returns:
The number of nested BeginModify() calls. Used internally.

void IncrementMod (  )  [inline]

Increase the number of nested BeginModify calls. Dangerous! Instead, properly use BeginModify as needed.

void DecrementMod (  )  [inline]

Decrease the number of nested BeginModify calls. Dangerous! Instead, properly use EndModify as needed.

int GetFlags (  )  [inline]

Returns:
the entire set of flags. (Internal use, mainly.)

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:
the atom at index idx or NULL if it does not exist.
Warning:
Atom indexing will change. Use iterator methods instead.

OBAtom * GetFirstAtom (  ) 

Returns:
the first atom in this molecule, or NULL if none exist.
Deprecated:
Will be removed in favor of more standard iterator methods

OBBond * GetBond ( int  idx  ) 

Returns:
the bond at index idx or NULL if it does not exist.
Warning:
Bond indexing may change. Use iterator methods instead.

OBBond * GetBond ( int  a,
int  b 
)

Returns:
the bond connecting the atom indexed by a and b or NULL if none exists.
Warning:
Atom indexing will change. Use atom objects and iterators instead.

OBBond * GetBond ( OBAtom bgn,
OBAtom end 
)

Returns:
the bond between the atoms bgn and end or NULL if none exists

OBResidue * GetResidue ( int   ) 

Returns:
the residue indexed by idx, or NULL if none exists
Warning:
Residue indexing may change. Use iterator methods instead.

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 a,
OBAtom b,
OBAtom c,
OBAtom d 
)

Returns:
the dihedral angle between the four atoms a, b, c, and d)

double GetAngle ( OBAtom a,
OBAtom b,
OBAtom c 
)

Returns:
the angle between the three atoms a, b and c (where a-> b (vertex) -> c )

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.

string GetSpacedFormula ( int  ones = 0,
const char *  sp = " " 
)

Returns:
the stochoimetric formula in spaced format e.g. C 4 H 6 O 1

Stochoimetric formula in spaced format e.g. C 4 H 6 O 1 No pair data is stored. Normally use without parameters: GetSpacedFormula()

Since:
version 2.1

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]

Returns:
the set of all atomic coordinates. See OBAtom::GetCoordPtr for more

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  ) 

Set the title of this molecule to title.

void SetTitle ( std::string &  title  ) 

Set the title of this molecule to 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  ) 

Set the total charge of this molecule to charge.

void SetTotalSpinMultiplicity ( unsigned int  spin  ) 

Set the total spin multiplicity of this molecule to spin (i.e., 0 = singlet (default), 1 = doublet, 2 = triplet, etc.)

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

Set the internal coordinates to int_coord (Does not call InternalToCartesian to update the 3D cartesian coordinates)

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]

Mark that hybridization of all atoms has been assigned.

void SetImplicitValencePerceived (  )  [inline]

Mark that the implicit hydrogen valence of all atoms has been assigned.

void SetKekulePerceived (  )  [inline]

Mark that Kekule forms have been assigned by Kekulize().

void SetClosureBondsPerceived (  )  [inline]

Mark that ring closure bonds have been assigned by graph traversal.

void SetHydrogensAdded (  )  [inline]

Mark that explicit hydrogen atoms have been added.

void SetCorrectedForPH (  )  [inline]

void SetAromaticCorrected (  )  [inline]

void SetSpinMultiplicityAssigned (  )  [inline]

void SetFlags ( int  flags  )  [inline]

void UnsetAromaticPerceived (  )  [inline]

void UnsetPartialChargesPerceived (  )  [inline]

void UnsetImplicitValencePerceived (  )  [inline]

void UnsetHydrogensAdded (  )  [inline]

void UnsetFlag ( int  flag  )  [inline]

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

Perform a set of transformations specified by the user

Typically these are program options to filter or modify an object For example, see OBMol::DoTransformations() and OBMol::ClassDescription()

Reimplemented from OBBase.

const char * ClassDescription (  )  [static]

Returns:
A list of descriptions of command-line options for DoTransformations()

Reimplemented from OBBase.

bool Clear (  )  [virtual]

Clear all information from a molecule.

Reimplemented from OBBase.

void RenumberAtoms ( std::vector< OBAtom * > &  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]  ) 

Rotate all conformers using the supplied matrix u (a 3x3 array of double).

void Rotate ( const double  m[9]  ) 

Rotate all conformers using the supplied matrix m (a linear 3x3 row-major array of double).

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

Rotate a specific conformer nconf using the supplied rotation matrix m.

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 DeleteHydrogens (  ) 

Delete all hydrogens from the molecule

Returns:
Success

bool DeleteHydrogens ( OBAtom  ) 

Delete all hydrogens from the supplied atom

Returns:
Success

bool DeleteNonPolarHydrogens (  ) 

Delete all hydrogen atoms connected to a non-polar atom

See also:
OBAtom::IsNonPolarHydrogen

bool DeleteHydrogen ( OBAtom  ) 

Delete the supplied atom if it is a hydrogen (Helper function for DeleteHydrogens)

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

Add hydrogens to the entire molecule to fill out implicit valence spots

Parameters:
polaronly Whether to add hydrogens only to polar atoms (i.e., not to C atoms)
correctForPH Whether to call CorrectForPH() first
Returns:
Whether any hydrogens were added

bool AddHydrogens ( OBAtom  ) 

Add hydrogens only to the supplied atom to fill out implicit valence.

bool AddPolarHydrogens (  ) 

Add only polar hydrogens (i.e., attached to polar atoms, not C).

bool StripSalts (  ) 

Deletes all atoms except for the largest contiguous fragment.

vector< OBMol > Separate ( int  StartIndex = 1  ) 

Copies each disconnected fragment as a separate OBMol.

bool ConvertDativeBonds (  ) 

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

bool CorrectForPH (  ) 

Correct for pH by applying the OBPhModel transformations.

bool AssignSpinMultiplicity (  ) 

set spin multiplicity for H-deficient atoms

vector3 Center ( int  nconf  ) 

Returns:
the center of the supplied conformer nconf
See also:
Center() to actually center all conformers at the origin

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 (  ) 

Find all ring atoms and bonds. Does not need to call FindSSSR().

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 frag  ) 

Find the largest fragment in OBMol (which may include multiple non-connected fragments)

Parameters:
frag Return (by reference) a bit vector indicating the atoms in the largest fragment

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 FindAngles (  ) 

Fills out an OBAngleData with angles from the molecule.

void FindTorsions (  ) 

Fills out an OBTorsionData with angles from the molecule.

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]

Has aromatic perception been performed?

bool HasSSSRPerceived (  )  [inline]

Has the smallest set of smallest rings (FindSSSR) been performed?

bool HasRingAtomsAndBondsPerceived (  )  [inline]

Have ring atoms and bonds been assigned?

bool HasAtomTypesPerceived (  )  [inline]

Have atom types been assigned by OBAtomTyper?

bool HasChiralityPerceived (  )  [inline]

Has atom chirality been assigned?

bool HasPartialChargesPerceived (  )  [inline]

Have atomic Gasteiger partial charges been assigned by OBGastChrg?

bool HasHybridizationPerceived (  )  [inline]

Has atomic hybridization been assigned by OBAtomTyper?

bool HasImplicitValencePerceived (  )  [inline]

Has implicit hydrogen valence been assigned by OBAtomTyper?

bool HasKekulePerceived (  )  [inline]

Has aromaticity and Kekule forms been assigned by Kekulize?

bool HasClosureBondsPerceived (  )  [inline]

Have ring "closure" bonds been assigned? (e.g., OBBond::IsClosure()).

bool HasChainsPerceived (  )  [inline]

Have biomolecule chains and residues been assigned by OBChainsParser?

bool HasHydrogensAdded (  )  [inline]

Have hydrogens been added to the molecule?

bool HasAromaticCorrected (  )  [inline]

Have aromatic nitrogens been "corrected?" (deprecated).

bool IsCorrectedForPH (  )  [inline]

Has the molecule been corrected for pH by CorrectForPH?

bool HasSpinMultiplicityAssigned (  )  [inline]

Has total spin multiplicity been assigned?

bool IsChiral (  ) 

Is this molecule chiral?

bool Empty (  )  [inline]

Are there any atoms in this molecule?

int NumConformers (  )  [inline]

Returns:
the number of conformers in this molecule

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

Set the entire set of conformers for this molecule to v.

void AddConformer ( double *  f  )  [inline]

Add a new set of coordinates f as a new conformer.

void SetConformer ( int  i  ) 

Set the molecule's current conformer to i Does nothing if i is less than 0 or i is larger than NumConformers()

void CopyConformer ( double *  c,
int  nconf 
)

Copy the conformer nconf into the array c

Warning:
Does no checking to see if c is large enough

void DeleteConformer ( int  nconf  ) 

Delete the conformer nconf.

double* GetConformer ( int  i  )  [inline]

Returns:
the coordinates to conformer i

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

Set the iterator to the beginning of the conformer list

Returns:
the array of coordinates for the first conformer

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

Advance the iterator to the next confomer, if possible

Returns:
The array of coordinates for the next conformer, or NULL if none exist

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

Returns:
the entire set of conformers for this molecule as a vector of floating point arrays

OBAtomIterator BeginAtoms (  )  [inline]

Returns:
An atom iterator pointing to the beginning of the atom list

OBAtomIterator EndAtoms (  )  [inline]

Returns:
An atom iterator pointing to the end of the atom list

OBBondIterator BeginBonds (  )  [inline]

Returns:
A bond iterator pointing to the beginning of the bond list

OBBondIterator EndBonds (  )  [inline]

Returns:
A bond iterator pointing to the end of the bond list

OBResidueIterator BeginResidues (  )  [inline]

Returns:
A residue iterator pointing to the beginning of the residue list

OBResidueIterator EndResidues (  )  [inline]

Returns:
A residue iterator pointing to the end of the residue list

OBAtom * BeginAtom ( OBAtomIterator i  ) 

Set the iterator i to the beginning of the atom list

Returns:
the first atom (or NULL if none exist)

OBAtom * NextAtom ( OBAtomIterator i  ) 

Advance the iterator i to the next atom in the molecule

Returns:
the next atom (if any, or NULL if none exist)

OBBond * BeginBond ( OBBondIterator i  ) 

Set the iterator i to the beginning of the bond list

Returns:
the first bond (or NULL if none exist)

OBBond * NextBond ( OBBondIterator i  ) 

Advance the iterator i to the next bond in the molecule

Returns:
the next bond (if any, or NULL if none exist)

OBResidue* BeginResidue ( OBResidueIterator i  )  [inline]

Set the iterator i to the beginning of the resdiue list

Returns:
the first residue (or NULL if none exist)

OBResidue* NextResidue ( OBResidueIterator i  )  [inline]

Advance the iterator i to the next residue in the molecule

Returns:
the next residue (if any, or NULL if not possible)

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

Set the iterator to the beginning of the internal coordinate list

Returns:
the first internal coordinate record, or NULL if none exist
See also:
SetInternalCoord

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

Advance the iterator to the next internal coordinate record

Returns:
the next first internal coordinate record, or NULL if none exist
See also:
SetInternalCoord

T* CastAndClear ( bool  clear = true  )  [inline, inherited]

By default clears the object. Called from ReadMolecule of most format classes.

bool HasData ( const std::string &   )  [inherited]

Returns:
whether the generic attribute/value pair exists

bool HasData ( const char *   )  [inherited]

Returns:
whether the generic attribute/value pair exists

bool HasData ( const unsigned int  type  )  [inherited]

Returns:
whether the generic attribute/value pair exists, for a given OBGenericDataType

void DeleteData ( unsigned int  type  )  [inherited]

Delete any data matching the given OBGenericDataType.

void DeleteData ( OBGenericData  )  [inherited]

Delete the given generic data from this object.

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

Delete all of the given generic data from this object.

void SetData ( OBGenericData d  )  [inline, inherited]

Adds a data object; does nothing if d==NULL.

unsigned int DataSize (  )  const [inline, inherited]

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

OBGenericData * GetData ( const unsigned int  type  )  [inherited]

Returns:
the first matching data for a given type from OBGenericDataType or NULL if nothing matches

OBGenericData * GetData ( const std::string &   )  [inherited]

Returns:
the value given an attribute name

OBGenericData * GetData ( const char *   )  [inherited]

Returns:
the value given an attribute name

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

Returns:
all data, suitable for iterating

std::vector< OBGenericData * > GetData ( DataOrigin  source  )  [inherited]

Returns:
all data with a specific origin, suitable for iterating

OBDataIterator BeginData (  )  [inline, inherited]

Returns:
An iterator pointing to the beginning of the data

OBDataIterator EndData (  )  [inline, inherited]

Returns:
An iterator pointing to the end of the data


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.

std::vector<OBAtom*> _vatom [protected]

vector of atoms

std::vector<OBBond*> _vbond [protected]

vector of bonds

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 int _natoms [protected]

Number of atoms.

unsigned int _nbonds [protected]

Number of bonds.

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

Residue information (if applicable).

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

Internal Coordinates (if applicable).

unsigned short int _mod [protected]

Number of nested calls to BeginModify().

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

Custom data.


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