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 SetCoordinates (double *c)
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=false, double pH=7.4)
bool AddHydrogens (OBAtom *)
bool AddPolarHydrogens ()
bool StripSalts (int threshold)
std::vector< OBMolSeparate (int StartIndex=1)
bool GetNextFragment (OpenBabel::OBMolAtomDFSIter &iter, OBMol &newMol)
bool ConvertDativeBonds ()
bool CorrectForPH (double pH=7.4)
bool AssignSpinMultiplicity (bool NoImplicitH=false)
vector3 Center (int nconf)
void SetTorsion (OBAtom *, OBAtom *, OBAtom *, OBAtom *, double ang)
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 &)
OBAtomNewAtom ()
OBBondNewBond ()
OBResidueNewResidue ()
bool DeleteAtom (OBAtom *, bool destroyAtom=true)
bool DeleteBond (OBBond *, bool destroyBond=true)
bool DeleteResidue (OBResidue *, bool destroyResidue=true)
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 (bool replaceNewlines=true) const
unsigned int NumAtoms () const
unsigned int NumBonds () const
unsigned int NumHvyAtoms ()
unsigned int NumResidues () const
unsigned int NumRotors ()
OBAtomGetAtom (int idx) const
OBAtomGetFirstAtom () const
OBBondGetBond (int idx) const
OBBondGetBond (int a, int b) const
OBBondGetBond (OBAtom *bgn, OBAtom *end) const
OBResidueGetResidue (int idx) const
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=" ", bool implicitH=true)
double GetEnergy () const
double GetMolWt (bool implicitH=true)
double GetExactMass (bool implicitH=true)
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 spinMultiplicity)
void SetInternalCoord (std::vector< OBInternalCoord * > int_coord)
void SetAutomaticFormalCharge (bool val)
void SetAutomaticPartialCharge (bool val)
void SetAromaticPerceived ()
void SetSSSRPerceived ()
void SetRingAtomsAndBondsPerceived ()
void SetAtomTypesPerceived ()
void SetRingTypesPerceived ()
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 UnsetSSSRPerceived ()
void UnsetRingTypesPerceived ()
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 HasRingTypesPerceived ()
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)
void SetEnergies (std::vector< double > &energies)
std::vector< double > GetEnergies ()
double GetEnergy (int ci)
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 * > &)
bool DeleteData (const std::string &s)
void SetData (OBGenericData *d)
void CloneData (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)
std::vector< OBGenericData * > GetAllData (const unsigned int type)
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)
bool expandcycle (OBAtom *atom, OBBitVec &avisit, OBAtom *first=NULL, int depth=0)

Protected Attributes

int _flags
bool _autoPartialCharge
bool _autoFormalCharge
std::string _title
std::vector< OBAtom * > _vatom
std::vector< OBBond * > _vbond
unsigned short int _dimension
int _totalCharge
unsigned int _totalSpin
double * _c
std::vector< double * > _vconf
double _energy
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 mol  ) 

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.

Referenced by OBMol::NewPerceiveKekuleBonds().

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

Referenced by OBMol::start_kekulize().

int getorden ( OBAtom atom  )  [protected]

Give the priority to give two electrons instead of 1.

Referenced by OBMol::NewPerceiveKekuleBonds().

bool expandcycle ( OBAtom atom,
OBBitVec avisit,
OBAtom first = NULL,
int  depth = 0 
) [protected]

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

Referenced by OBMol::NewPerceiveKekuleBonds().

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

Referenced by OBMol::AddAtom(), and OBMol::NewAtom().

OBBond * CreateBond ( void   )  [virtual]

Create a new OBBond pointer. Does no bookkeeping

Deprecated:
Use NewBond instead, which ensures internal connections

Referenced by OBMol::AddBond(), and OBMol::NewBond().

OBResidue * CreateResidue ( void   )  [virtual]

Create a new OBResidue pointer. Does no bookkeeping

Deprecated:
Use NewResidue instead, which ensures internal connections

Referenced by OBMol::AddResidue(), and OBMol::NewResidue().

void DestroyAtom ( OBAtom atom  )  [virtual]

Free an OBAtom pointer if defined. Does no bookkeeping

See also:
DeleteAtom which ensures internal connections

Referenced by OBMol::Clear(), OBMol::DeleteAtom(), OBMol::DeleteHydrogen(), and OBMol::~OBMol().

void DestroyBond ( OBBond bond  )  [virtual]

Free an OBBond pointer if defined. Does no bookkeeping

See also:
DeleteBond which ensures internal connections

Referenced by OBMol::Clear(), OBMol::DeleteBond(), and OBMol::~OBMol().

void DestroyResidue ( OBResidue residue  )  [virtual]

Free an OBResidue pointer if defined. Does no bookkeeping

See also:
DeleteResidue which ensures internal connections

Referenced by OBMol::Clear(), OBMol::DeleteResidue(), and OBMol::~OBMol().

bool AddAtom ( OBAtom atom  ) 

Add an atom to a molecule.

Add the specified atom to this molecule

Returns:
Whether the method was successful
Also checks bond_queue for any bonds that should be made to the new atom

Referenced by OBMol::GetNextFragment(), OBMol::InsertAtom(), OBMol::operator+=(), and OBMol::operator=().

bool InsertAtom ( OBAtom atom  ) 

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

Referenced by OBMol::AddAtom(), OBMol::AddBond(), OBMol::AddHydrogens(), OBResidueData::AssignBonds(), OBBuilder::Build(), OBMol::ConnectTheDots(), AliasData::Expand(), OBMol::GetNextFragment(), OBAtom::HtoMethyl(), OBMol::NewAtom(), OBMol::operator+=(), OBMol::operator=(), and OBAtom::SetHybAndGeom().

bool AddBond ( OBBond bond  ) 

Add the specified residue to this molecule and update connections

Returns:
Whether the method was successful

bool AddResidue ( OBResidue residue  ) 

Add the specified residue to this molecule and update connections

Returns:
Whether the method was successful

Referenced by OBMol::operator+=().

OBAtom * NewAtom (  ) 

Instantiate a New Atom and add it to the molecule.

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

Checks bond_queue for any bonds that should be made to the new atom and updates atom indexes.

Referenced by OBMol::AddHydrogens(), AliasData::Expand(), OBUnitCell::FillUnitCell(), OBForceField::GetGrid(), OBAtom::HtoMethyl(), and OBAtom::SetHybAndGeom().

OBBond * NewBond (  ) 

Instantiate a New Bond and add it to the molecule.

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

Since:
version 2.1 Sets the proper Bond index and insures this molecule is set as the parent.

Referenced by OBBuilder::Connect().

OBResidue * NewResidue (  ) 

Create a new OBResidue in this molecule and ensure connections.

Referenced by OBMol::operator=().

bool DeleteAtom ( OBAtom atom,
bool  destroyAtom = true 
)

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

Referenced by OBChemTsfm::Apply(), OBMol::DeleteHydrogens(), OBMol::DeleteNonPolarHydrogens(), OBForceField::GetGrid(), OpenBabel::InternalToCartesian(), OBAtom::SetHybAndGeom(), and OBMol::StripSalts().

bool DeleteBond ( OBBond bond,
bool  destroyBond = true 
)

Deletes an bond from this molecule and updates accordingly

Returns:
Whether deletion was successful

Referenced by OBBuilder::Build(), OBMol::ConnectTheDots(), OBMol::DeleteAtom(), OBMol::DeleteHydrogen(), and OBBuilder::Swap().

bool DeleteResidue ( OBResidue residue,
bool  destroyResidue = true 
)

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.

Referenced by OBMol::AddHydrogens(), OBMol::AddResidue(), OBMol::DeleteAtom(), OBMol::DeleteBond(), OBForceField::GetGrid(), OBAtom::HtoMethyl(), OBMol::InsertAtom(), OBMol::operator+=(), OBMol::operator=(), OBAtom::SetHybAndGeom(), and OBMol::StripSalts().

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.

Referenced by OBMol::AddHydrogens(), OBMol::AddResidue(), OBMol::DeleteAtom(), OBMol::DeleteBond(), OBForceField::GetGrid(), OBAtom::HtoMethyl(), OBMol::InsertAtom(), OBMol::operator+=(), OBMol::operator=(), OBAtom::SetHybAndGeom(), and OBMol::StripSalts().

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.

Referenced by OBMol::AddHydrogens(), OBMol::DeleteHydrogen(), OBMol::DeleteHydrogens(), OBMol::DeleteNonPolarHydrogens(), and OBAtom::SetHybAndGeom().

void DecrementMod (  )  [inline]

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

Referenced by OBMol::AddHydrogens(), OBMol::DeleteHydrogen(), OBMol::DeleteHydrogens(), OBMol::DeleteNonPolarHydrogens(), and OBAtom::SetHybAndGeom().

int GetFlags (  )  [inline]

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

Referenced by OBAtomTyper::AssignImplicitValence().

const char * GetTitle ( bool  replaceNewlines = true  )  const

unsigned int NumAtoms (  )  const [inline]

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

Referenced by OBMol::AddBond(), OBMol::AddHydrogens(), OBRingSearch::AddRingFromClosure(), OBChemTsfm::Apply(), patty::assign_types(), OBAromaticTyper::AssignAromaticFlags(), OBGastChrg::AssignPartialCharges(), OBMol::AssignSpinMultiplicity(), OBRotorList::AssignTorVals(), OpenBabel::BreakChiralTies(), OpenBabel::CalcSignedVolume(), OpenBabel::CalculateSymmetry(), OpenBabel::CanonicalLabels(), OpenBabel::CartesianToInternal(), OBMol::Center(), OBForceField::ConjugateGradientsInitialize(), OBForceField::ConjugateGradientsTakeNSteps(), OBBuilder::Connect(), OBMol::ConnectTheDots(), OpenBabel::construct_c_matrix(), OpenBabel::construct_g_matrix(), OBMol::ContigFragList(), OBMol::CopyConformer(), OBForceField::CorrectVelocities(), OBRotamerList::CreateConformerList(), OBMol::DeleteHydrogen(), OpenBabel::DetermineFRJ(), OBForceField::DistanceGeometry(), OBMol::EndModify(), AliasData::Expand(), OpenBabel::FastSingleMatch(), OBMol::FindLargestFragment(), OpenBabel::FindRingAtoms(), OBMol::FindRingAtomsAndBonds(), OBRotorList::FindRotors(), OBMol::FindSSSR(), OpenBabel::FixCisTransBonds(), OBForceField::GenerateVelocities(), OBMol::GetAtom(), OBForceField::GetAtomTypes(), OpenBabel::GetChirality(), OBForceField::GetConformers(), OBForceField::GetCoordinates(), OpenBabel::GetDFFVector(), OBMol::GetExactMass(), OBMol::GetGIDVector(), OBMol::GetGIVector(), OpenBabel::GetGIVector(), OBForceField::GetGrid(), OBMol::GetGTDVector(), OpenBabel::GetGTDVector(), OBMol::GetInternalCoord(), OBMol::GetMolWt(), OBAtom::GetNextAtom(), OBMol::GetNextFragment(), OBForceField::GetPartialCharges(), OBMol::GetSpacedFormula(), OpenBabel::GraphPotentials(), OBAtom::HtoMethyl(), OBBond::IsClosure(), OBForceField::IsSetupNeeded(), OBForceField::LineSearch(), OBMoleculeFormat::MakeCombinedMolecule(), OBSSMatch::Match(), OpenBabel::MinimumPairRMS(), OBMol::NewPerceiveKekuleBonds(), OBMolAtomBFSIter::OBMolAtomBFSIter(), OBMolAtomDFSIter::OBMolAtomDFSIter(), OBSSMatch::OBSSMatch(), OBMol::operator+=(), OBMol::operator=(), OBMol::PerceiveKekuleBonds(), OBMoleculeFormat::ReadChemObjectImpl(), OBRotorList::RemoveSymVals(), OBMol::RenumberAtoms(), OBMol::Rotate(), OBMol::Separate(), OBRotamerList::SetBaseCoordinateSets(), OBForceField::SetConformers(), OBMol::SetCoordinates(), OBForceField::SetCoordinates(), OBRotorList::SetRotAtoms(), OBProxGrid::Setup(), OBForceField::Setup(), OpenBabel::SetupAtomMatchTable(), OBMol::start_kekulize(), OBForceField::SteepestDescentTakeNSteps(), OBMol::ToInertialFrame(), OBMol::Translate(), OBForceField::WeightedRotorSearch(), and OBMoleculeFormat::WriteChemObjectImpl().

unsigned int NumBonds (  )  const [inline]

unsigned int NumHvyAtoms (  ) 

Returns:
the number of non-hydrogen atoms

Referenced by OBMol::DoTransformations(), OBMol::GetExactMass(), OBMol::GetMolWt(), and OBMol::GetSpacedFormula().

unsigned int NumResidues (  )  const [inline]

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

Referenced by OBMol::GetResidue(), and OBMol::operator=().

unsigned int NumRotors (  ) 

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

OBAtom * GetAtom ( int  idx  )  const

Returns:
the atom at index idx or NULL if it does not exist.
Warning:
Atom indexing will change. Use iterator methods instead.
Returns a pointer to the atom after a safety check 0 < idx <= NumAtoms

Referenced by OBMol::AddBond(), OBMol::Align(), OBChemTsfm::Apply(), OpenBabel::ApplyRotMatToBond(), OBBondTyper::AssignFunctionalGroupBonds(), OBAtomTyper::AssignHyb(), OBAtomTyper::AssignImplicitValence(), OBPhModel::AssignSeedPartialCharge(), OBRingTyper::AssignTypes(), OBAtomTyper::AssignTypes(), OBBuilder::Build(), OpenBabel::BuildOBRTreeVector(), OpenBabel::CalcSignedVolume(), OpenBabel::CartesianToInternal(), OBBuilder::Connect(), OBMol::ConnectTheDots(), OBMol::ContigFragList(), OpenBabel::CorrectBadResonanceForm(), OBBuilder::CorrectStereoAtoms(), OBBuilder::CorrectStereoBonds(), OBForceField::DistanceGeometry(), AliasData::Expand(), OpenBabel::FastSingleMatch(), OBRing::findCenterAndNormal(), OBMol::FindChildren(), OBMol::FindLargestFragment(), OpenBabel::FindRings(), OBAtom::GetAngle(), OBMol::GetBond(), OBForceField::GetCoordinates(), OpenBabel::GetDFFVector(), OBAtom::GetDistance(), OBForceField::GetGrid(), OBMol::GetGTDVector(), OpenBabel::GetGTDVector(), OBAtom::GetNextAtom(), OBMol::GetNextFragment(), OBRing::GetRootAtom(), OBRotorRules::GetRotorIncrements(), OBRing::IsAromatic(), OBForceField::IsSetupNeeded(), OBSSMatch::Match(), OpenBabel::match(), OBAtom::MatchesSMARTS(), OpenBabel::MinimumPairRMS(), OBMol::NewPerceiveKekuleBonds(), OBMolAtomBFSIter::operator++(), OBMolAtomDFSIter::operator++(), OBMol::operator+=(), OBMol::operator=(), OBMol::PerceiveBondOrders(), OBAromaticTyper::SelectRootAtoms(), OBForceField::SetConformers(), OBForceField::SetCoordinates(), OBRotorList::SetEvalAtoms(), OBBond::SetLength(), OBRotamerList::Setup(), OBFFConstraints::Setup(), OBMol::StripSalts(), OBBuilder::Swap(), and OBForceField::UpdatePairsSimple().

OBAtom * GetFirstAtom (  )  const

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

Returns:
the bond at index idx or NULL if it does not exist.
Warning:
Bond indexing may change. Use iterator methods instead.
Returns a pointer to the bond after a safety check 0 <= idx < NumBonds

Referenced by OBMol::AddBond(), OBChemTsfm::Apply(), OBBuilder::Build(), OpenBabel::FastSingleMatch(), OpenBabel::FindRings(), OBMol::GetBond(), OBSSMatch::Match(), OBMol::PerceiveBondOrders(), OBMol::start_kekulize(), and OBBuilder::Swap().

OBBond * GetBond ( int  a,
int  b 
) const

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

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

OBResidue * GetResidue ( int  idx  )  const

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

Referenced by OBMol::operator=().

std::vector< OBInternalCoord * > GetInternalCoord (  ) 

double GetTorsion ( int  a,
int  b,
int  c,
int  d 
)

Returns:
the dihedral angle (in degrees) between the four atoms supplied a1-a2-a3-a4) WARNING: SetTorsion takes an angle in radians while GetTorsion returns it in degrees

Referenced by OBBuilder::CorrectStereoAtoms(), OBRotorRules::GetRotorIncrements(), and OBMol::PerceiveBondOrders().

double GetTorsion ( OBAtom a,
OBAtom b,
OBAtom c,
OBAtom d 
)

Returns:
the dihedral angle (in degrees) between the four atoms a, b, c, and d) WARNING: SetTorsion takes an angle in radians while GetTorsion returns it in degrees

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

Returns:
the angle (in degrees) 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 = " ",
bool  implicitH = true 
)

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

Referenced by OBMol::DoTransformations(), OBMol::GetFormula(), and OBMoleculeFormat::MakeCombinedMolecule().

double GetEnergy (  )  const [inline]

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

Referenced by OBMol::operator=().

double GetMolWt ( bool  implicitH = true  ) 

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

double GetExactMass ( bool  implicitH = true  ) 

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 number of unpaired electrons assuming no further pairing of spins.

unsigned short int GetDimension (  )  const [inline]

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

Referenced by OBPhModel::CorrectForPH(), OBMol::GetNextFragment(), OBMoleculeFormat::MakeCombinedMolecule(), OBMol::operator+=(), and OBMol::operator=().

double* GetCoordinates (  )  [inline]

vector< OBRing * > & GetSSSR (  ) 

bool AutomaticFormalCharge (  )  [inline]

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

Referenced by OBPhModel::CorrectForPH().

bool AutomaticPartialCharge (  )  [inline]

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

Referenced by OBPhModel::AssignSeedPartialCharge().

void SetTitle ( const char *  title  ) 

Set the title of this molecule to title.

Referenced by OBMol::DoTransformations(), and OBMoleculeFormat::MakeCombinedMolecule().

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

Referenced by OBBuilder::Build(), and OBMol::GetNextFragment().

void SetTotalCharge ( int  charge  ) 

Set the total charge of this molecule to charge.

void SetTotalSpinMultiplicity ( unsigned int  spinMultiplicity  ) 

Set the total spin multiplicity of this molecule to spinMultiplicity Overrides the calculation from spin multiplicity of OBAtoms

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

Referenced by OBAromaticTyper::AssignAromaticFlags(), OBAtomTyper::AssignImplicitValence(), and OBMol::EndModify().

void SetSSSRPerceived (  )  [inline]

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

Referenced by OBMol::FindSSSR().

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

Referenced by OBAtomTyper::AssignTypes().

void SetRingTypesPerceived (  )  [inline]

Mark that ring types have been perceived (see OBRingTyper for details).

Referenced by OBRingTyper::AssignTypes().

void SetChainsPerceived (  )  [inline]

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

void SetChiralityPerceived (  )  [inline]

Mark that chirality has been perceived.

Referenced by OBMol::FindChiralCenters().

void SetPartialChargesPerceived (  )  [inline]

Mark that partial charges have been assigned.

Referenced by OBPhModel::AssignSeedPartialCharge().

void SetHybridizationPerceived (  )  [inline]

Mark that hybridization of all atoms has been assigned.

Referenced by OBAtomTyper::AssignHyb(), OBBuilder::Build(), and OBMol::PerceiveBondOrders().

void SetImplicitValencePerceived (  )  [inline]

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

Referenced by OBAtomTyper::AssignImplicitValence().

void SetKekulePerceived (  )  [inline]

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

Referenced by OBMol::NewPerceiveKekuleBonds(), and OBMol::PerceiveKekuleBonds().

void SetClosureBondsPerceived (  )  [inline]

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

Referenced by OBBond::IsClosure().

void SetHydrogensAdded (  )  [inline]

Mark that explicit hydrogen atoms have been added.

Referenced by OBMol::AddHydrogens().

void SetCorrectedForPH (  )  [inline]

Referenced by OBPhModel::CorrectForPH().

void SetAromaticCorrected (  )  [inline]

void SetSpinMultiplicityAssigned (  )  [inline]

void SetFlags ( int  flags  )  [inline]

void UnsetAromaticPerceived (  )  [inline]

Referenced by OBMol::EndModify().

void UnsetSSSRPerceived (  )  [inline]

void UnsetRingTypesPerceived (  )  [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.

Referenced by OBMoleculeFormat::ReadChemObjectImpl().

const char * ClassDescription (  )  [static]

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

Reimplemented from OBBase.

bool Clear ( void   )  [virtual]

Clear all information from a molecule.

Reimplemented from OBBase.

Referenced by OBMol::operator=(), OBMoleculeFormat::ReadNameIndex(), and OBMol::Separate().

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 SetCoordinates ( double *  c  ) 

Set the coordinates for all atoms in this conformer.

See also:
OBMol::GetCoordinates()

Referenced by OBForceField::WeightedRotorSearch().

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

Referenced by OBMol::Rotate().

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

Referenced by OBMol::DoTransformations().

bool Kekulize (  ) 

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

Referenced by OBMol::EndModify(), and OBMol::PerceiveBondOrders().

bool PerceiveKekuleBonds (  ) 

void NewPerceiveKekuleBonds (  ) 

Kekulize aromatic rings without using implicit valence.

This new perceive kekule bonds function has been designed to handle molecule files without explicit hydrogens such as pdb or xyz. (It can, of course, easily handle explicit hydrogens too.) 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

Referenced by OBPhModel::CorrectForPH(), and OBMol::DoTransformations().

bool DeleteHydrogens ( OBAtom atom  ) 

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

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

Referenced by OBMol::DeleteAtom(), and OBMol::DeleteHydrogens().

bool AddHydrogens ( bool  polaronly = false,
bool  correctForPH = false,
double  pH = 7.4 
)

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
pH The pH to use for CorrectForPH() modification
Returns:
Whether any hydrogens were added

Referenced by OBMol::AddPolarHydrogens(), OBMol::DoTransformations(), and OBSmartsPattern::Match().

bool AddHydrogens ( OBAtom atom  ) 

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

Deletes all atoms except for the largest contiguous fragment.

vector< OBMol > Separate ( int  StartIndex = 1  ) 

Copies each disconnected fragment as a separate OBMol.

Referenced by OBMoleculeFormat::ReadChemObjectImpl().

bool GetNextFragment ( OpenBabel::OBMolAtomDFSIter iter,
OBMol newMol 
)

Iterative component of Separate to copy one fragment at a time.

Referenced by OBMol::Separate().

bool ConvertDativeBonds (  ) 

Converts the charged form of coordinate bonds, e.g.[N+]([O-])=O to N(=O)=O.

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

Referenced by OBMol::DoTransformations().

bool CorrectForPH ( double  pH = 7.4  ) 

Correct for pH by applying the OBPhModel transformations.

Referenced by OBMol::AddHydrogens().

bool AssignSpinMultiplicity ( bool  NoImplicitH = false  ) 

set spin multiplicity for H-deficient atoms

If NoImplicitH is true then the molecule has no implicit hydrogens. Individual atoms on which ForceNoH() has been called also have no implicit hydrogens. If NoImplicitH is false (the default), then if there are any explicit hydrogens on an atom then they constitute all the hydrogen on that atom. However, a hydrogen atom with its _isotope!=0 is not considered explicit hydrogen for this purpose. In addition, an atom which has had ForceImplH()called for it is never considered hydrogen deficient, e.g. unbracketed atoms in SMILES. Any discrepancy with the expected atom valency is interpreted as the atom being a radical of some sort and iits _spinMultiplicity is set to 2 when it is one hydrogen short and 3 when it is two hydrogens short and similarly for greater hydrogen deficiency.

So SMILES C[CH] is interpreted as methyl carbene, CC[H][H] as ethane, and CC[2H] as CH3CH2D.

Referenced by OBMol::PerceiveBondOrders().

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

Set the torsion defined by these atoms, rotating bonded neighbors

ang The torsion angle in radians
WARNING: SetTorsion takes an angle in radians while GetTorsion returns it in degrees

Referenced by OBBuilder::CorrectStereoBonds().

void FindSSSR (  ) 

void FindRingAtomsAndBonds (  ) 

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

Referenced by OBRotorList::FindRotors(), OBMol::FindSSSR(), OBBond::IsInRing(), and OBAtom::IsInRing().

void FindChiralCenters (  ) 

Find all chiral atom centers. See OBAtom::IsChiral() for more details.

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'

Referenced by OBMol::Align(), OpenBabel::ApplyRotMatToBond(), OBRotorList::AssignTorVals(), OBBond::SetLength(), OBRotorList::SetRotAtoms(), OBRotorList::SetRotAtomsByFix(), OBMol::SetTorsion(), and OBRotamerList::Setup().

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 > > &  cfl  ) 

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

Referenced by OpenBabel::DetermineFRJ(), and OBMol::StripSalts().

void Align ( OBAtom a1,
OBAtom a2,
vector3 p1,
vector3 p2 
)

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." AssignSpinMultiplicity(true) is called at the end of the function. The true states that there are no implict hydrogens in the molecule.

void FindAngles (  ) 

Fills out an OBAngleData with angles from the molecule.

Referenced by OBMolAngleIter::OBMolAngleIter().

void FindTorsions (  ) 

Fills out an OBTorsionData with angles from the molecule.

Referenced by OBMolTorsionIter::OBMolTorsionIter().

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

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

Referenced by OBRotorList::FindRotors(), and OBMol::GetGIVector().

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

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.

Referenced by OBMol::GetGIDVector().

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

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

Referenced by OBMol::FindChiralCenters().

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

Referenced by OpenBabel::CalcSignedVolume(), and OpenBabel::match().

bool HasNonZeroCoords (  ) 

Are there any non-zero coordinates?

Referenced by OBMol::AddHydrogens().

bool HasAromaticPerceived (  )  [inline]

Has aromatic perception been performed?

Referenced by OBAromaticTyper::AssignAromaticFlags(), OBBond::IsAromatic(), and OBAtom::IsAromatic().

bool HasSSSRPerceived (  )  [inline]

bool HasRingAtomsAndBondsPerceived (  )  [inline]

Have ring atoms and bonds been assigned?

Referenced by OBBond::IsInRing(), and OBAtom::IsInRing().

bool HasAtomTypesPerceived (  )  [inline]

Have atom types been assigned by OBAtomTyper?

Referenced by OBAtom::GetType().

bool HasRingTypesPerceived (  )  [inline]

Have ring types been assigned by OBRingTyper?

Referenced by OBRing::GetType().

bool HasChiralityPerceived (  )  [inline]

Has atom chirality been assigned?

Referenced by OBMol::FindChiralCenters().

bool HasPartialChargesPerceived (  )  [inline]

Have atomic Gasteiger partial charges been assigned by OBGastChrg?

bool HasHybridizationPerceived (  )  [inline]

Has atomic hybridization been assigned by OBAtomTyper?

Referenced by OBAtom::GetHyb().

bool HasImplicitValencePerceived (  )  [inline]

Has implicit hydrogen valence been assigned by OBAtomTyper?

Referenced by OBAtomTyper::AssignImplicitValence(), OBAtom::GetImplicitValence(), and OBAtom::ImplicitHydrogenCount().

bool HasKekulePerceived (  )  [inline]

Has aromaticity and Kekule forms been assigned by Kekulize?

Referenced by OBMol::NewPerceiveKekuleBonds(), and OBMol::PerceiveKekuleBonds().

bool HasClosureBondsPerceived (  )  [inline]

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

Referenced by OBBond::IsClosure().

bool HasChainsPerceived (  )  [inline]

Have biomolecule chains and residues been assigned by OBChainsParser?

bool HasHydrogensAdded (  )  [inline]

Have hydrogens been added to the molecule?

Referenced by OBMol::AddHydrogens().

bool HasAromaticCorrected (  )  [inline]

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

Referenced by OBAtomTyper::AssignImplicitValence(), and OBAtomTyper::CorrectAromaticNitrogens().

bool IsCorrectedForPH (  )  [inline]

Has the molecule been corrected for pH by CorrectForPH?

Referenced by OBMol::AddHydrogens(), OBPhModel::CorrectForPH(), and OBMol::CorrectForPH().

bool HasSpinMultiplicityAssigned (  )  [inline]

Has total spin multiplicity been assigned?

Referenced by OBMol::AssignSpinMultiplicity().

bool IsChiral (  ) 

Is this molecule chiral?

bool Empty (  )  [inline]

int NumConformers (  )  [inline]

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

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

Referenced by OBForceField::GetConformers(), OBMol::operator=(), and OBForceField::SetConformers().

void AddConformer ( double *  f  )  [inline]

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

Referenced by OBForceField::WeightedRotorSearch().

void SetConformer ( int  i  ) 

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]

void SetEnergies ( std::vector< double > &  energies  ) 

Set the entire set of conformer energies.

vector< double > GetEnergies (  ) 

Set the entire set of conformer energies.

double GetEnergy ( int  ci  ) 

Get the energy for conformer ci

ci conformer index

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

Referenced by OBForceField::RandomRotorSearchInitialize(), OBRotamerList::SetBaseCoordinateSets(), and OBForceField::SystematicRotorSearchInitialize().

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)

Referenced by OBMol::AddHydrogens(), OBAromaticTyper::AssignAromaticFlags(), OBResidueData::AssignBonds(), OBAtomTyper::AssignHyb(), OBAtomTyper::AssignImplicitValence(), OBGastChrg::AssignPartialCharges(), OBMol::AssignSpinMultiplicity(), OBAtomTyper::AssignTypes(), OBMol::BeginModify(), OpenBabel::BreakChiralTies(), OpenBabel::CalcSignedVolume(), OpenBabel::CalculateSymmetry(), OpenBabel::CanonicalLabels(), OpenBabel::CartesianToInternal(), OBMol::Center(), OBMol::ConnectTheDots(), OpenBabel::construct_c_matrix(), OpenBabel::construct_g_matrix(), OBMol::ContigFragList(), OBMol::ConvertDativeBonds(), OBMol::DeleteAtom(), OBMol::DeleteHydrogen(), OBMol::DeleteHydrogens(), OBMol::DeleteNonPolarHydrogens(), OBMol::EndModify(), OBAromaticTyper::ExcludeSmallRing(), OpenBabel::FastSingleMatch(), OBMol::FindChiralCenters(), OBMol::FindLargestFragment(), OpenBabel::GetChirality(), OpenBabel::GetDFFVector(), OBMol::GetExactMass(), OBMol::GetGIDVector(), OBMol::GetGIVector(), OpenBabel::GetGIVector(), OBMol::GetGTDVector(), OpenBabel::GetGTDVector(), OBMol::GetMolWt(), OBAtom::GetPartialCharge(), OBMol::GetTotalCharge(), OBMol::GetTotalSpinMultiplicity(), OBMol::Has2D(), OBMol::Has3D(), OBMol::HasNonZeroCoords(), OBGrid::Init(), OpenBabel::InternalToCartesian(), OBMol::IsChiral(), OBBond::IsClosure(), OBSSMatch::Match(), OBMol::NumHvyAtoms(), OBMolAtomIter::OBMolAtomIter(), OBMolPairIter::OBMolPairIter(), OBMolPairIter::operator++(), OBMol::operator+=(), OBMol::operator=(), OBMol::PerceiveBondOrders(), OBMol::PerceiveKekuleBonds(), OBMol::RenumberAtoms(), OBMol::SetCoordinates(), OBProxGrid::Setup(), OpenBabel::SetupAtomMatchTable(), OBMol::ToInertialFrame(), and OBMol::~OBMol().

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)

Referenced by OBMol::AddHydrogens(), OBAromaticTyper::AssignAromaticFlags(), OBResidueData::AssignBonds(), OBAtomTyper::AssignHyb(), OBAtomTyper::AssignImplicitValence(), OBGastChrg::AssignPartialCharges(), OBMol::AssignSpinMultiplicity(), OBAtomTyper::AssignTypes(), OBMol::BeginModify(), OpenBabel::BreakChiralTies(), OpenBabel::CalcSignedVolume(), OpenBabel::CalculateSymmetry(), OpenBabel::CanonicalLabels(), OpenBabel::CartesianToInternal(), OBMol::Center(), OBMol::ConnectTheDots(), OpenBabel::construct_c_matrix(), OpenBabel::construct_g_matrix(), OBMol::ContigFragList(), OBMol::ConvertDativeBonds(), OBMol::DeleteAtom(), OBMol::DeleteHydrogen(), OBMol::DeleteHydrogens(), OBMol::DeleteNonPolarHydrogens(), OBMol::EndModify(), OBAromaticTyper::ExcludeSmallRing(), OpenBabel::FastSingleMatch(), OBMol::FindChiralCenters(), OBMol::FindLargestFragment(), OpenBabel::GetChirality(), OpenBabel::GetDFFVector(), OBMol::GetExactMass(), OBMol::GetGIDVector(), OBMol::GetGIVector(), OpenBabel::GetGIVector(), OBMol::GetGTDVector(), OpenBabel::GetGTDVector(), OBMol::GetMolWt(), OBAtom::GetPartialCharge(), OBMol::GetTotalCharge(), OBMol::GetTotalSpinMultiplicity(), OBMol::Has2D(), OBMol::Has3D(), OBMol::HasNonZeroCoords(), OBGrid::Init(), OpenBabel::InternalToCartesian(), OBMol::IsChiral(), OBBond::IsClosure(), OBSSMatch::Match(), OBMol::NumHvyAtoms(), OBMolPairIter::OBMolPairIter(), OBMolPairIter::operator++(), OBMolAtomIter::operator++(), OBMol::operator+=(), OBMol::operator=(), OBMol::PerceiveBondOrders(), OBMol::PerceiveKekuleBonds(), OBMol::RenumberAtoms(), OBMol::SetCoordinates(), OBProxGrid::Setup(), OpenBabel::SetupAtomMatchTable(), OBMol::ToInertialFrame(), and OBMol::~OBMol().

OBBond * BeginBond ( OBBondIterator i  ) 

OBBond * NextBond ( OBBondIterator i  ) 

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)

Referenced by OBResidueIter::OBResidueIter(), OBMol::operator+=(), and OBMol::~OBMol().

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)

Referenced by OBResidueIter::operator++(), OBMol::operator+=(), and OBMol::~OBMol().

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 &  s  )  [inherited]

bool HasData ( const char *  s  )  [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]

void DeleteData ( OBGenericData gd  )  [inherited]

Delete the given generic data from this object.

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

Delete all of the given generic data from this object.

bool DeleteData ( const std::string &  s  )  [inherited]

Deletes the generic data with the specified attribute, returning false if not found.

void SetData ( OBGenericData d  )  [inline, inherited]

void CloneData ( OBGenericData d  )  [inherited]

Adds a copy of a data object; does nothing if d == NULL

Since:
version 2.2

unsigned int DataSize (  )  const [inline, inherited]

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

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

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

Returns:
any data matching the given attribute name or NULL if nothing matches

the value given an attribute name

OBGenericData * GetData ( const char *  s  )  [inherited]

Returns:
any data matching the given attribute name or NULL if nothing matches

the value given an attribute name

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

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

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

std::vector< OBGenericData * > GetAllData ( const unsigned int  type  )  [inherited]

Returns:
the all matching data for a given type from OBGenericDataType or an empty vector if nothing matches
Since:
version 2.2

OBDataIterator BeginData (  )  [inline, inherited]

Returns:
An iterator pointing to the beginning of the data

Referenced by OBMol::AddAtom(), OBAtom::Duplicate(), OBMoleculeFormat::MakeCombinedMolecule(), OBMol::NewAtom(), and OBMol::operator=().

OBDataIterator EndData (  )  [inline, inherited]

Returns:
An iterator pointing to the end of the data

Referenced by OBMol::AddAtom(), OBAtom::Duplicate(), OBMoleculeFormat::MakeCombinedMolecule(), OBMol::NewAtom(), and OBMol::operator=().


Member Data Documentation

int _flags [protected]

bool _autoPartialCharge [protected]

Assign partial charges automatically.

Referenced by OBMol::OBMol().

bool _autoFormalCharge [protected]

Assign formal charges automatically.

Referenced by OBMol::OBMol().

std::string _title [protected]

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

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

unsigned short int _dimension [protected]

int _totalCharge [protected]

Total charge on the molecule.

Referenced by OBMol::GetTotalCharge(), OBMol::OBMol(), and OBMol::SetTotalCharge().

unsigned int _totalSpin [protected]

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

Referenced by OBMol::GetTotalSpinMultiplicity(), and OBMol::SetTotalSpinMultiplicity().

double* _c [protected]

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

double _energy [protected]

heat of formation

Referenced by OBMol::OBMol(), and OBMol::operator=().

unsigned int _natoms [protected]

unsigned int _nbonds [protected]

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

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

Internal Coordinates (if applicable).

Referenced by OBMol::GetInternalCoord().

unsigned short int _mod [protected]

Number of nested calls to BeginModify().

Referenced by OBMol::BeginModify(), OBMol::Clear(), OBMol::EndModify(), and OBMol::OBMol().

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


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