Open Babel  3.0
Public Member Functions | Protected Attributes | List of all members
OBMol Class Reference

#include <openbabel/mol.h>

Inheritance diagram for OBMol:
OBBase

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 void DestroyAtom (OBAtom *)
 
virtual void DestroyBond (OBBond *)
 
virtual void DestroyResidue (OBResidue *)
 
bool AddAtom (OBAtom &atom, bool forceNewId=false)
 
bool InsertAtom (OBAtom &)
 
bool AddBond (int beginIdx, int endIdx, int order, int flags=0, int insertpos=-1)
 
bool AddBond (OBBond &)
 
bool AddResidue (OBResidue &)
 
OBAtomNewAtom ()
 
OBAtomNewAtom (unsigned long id)
 
OBBondNewBond ()
 
OBBondNewBond (unsigned long id)
 
OBResidueNewResidue ()
 
bool DeleteAtom (OBAtom *, bool destroyAtom=true)
 
bool DeleteBond (OBBond *, bool destroyBond=true)
 
bool DeleteResidue (OBResidue *, bool destroyResidue=true)
 
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 (bool sampleRingBonds=false)
 
OBAtomGetAtom (int idx) const
 
OBAtomGetAtomById (unsigned long id) const
 
OBAtomGetFirstAtom () const
 
OBBondGetBond (int idx) const
 
OBBondGetBondById (unsigned long id) 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)
 
int AreInSameRing (OBAtom *a, OBAtom *b)
 
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 ()
 
std::vector< OBRing * > & GetLSSR ()
 
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 (bool value=true)
 
void SetSSSRPerceived (bool value=true)
 
void SetLSSRPerceived (bool value=true)
 
void SetRingAtomsAndBondsPerceived (bool value=true)
 
void SetAtomTypesPerceived (bool value=true)
 
void SetRingTypesPerceived (bool value=true)
 
void SetChainsPerceived (bool value=true)
 
void SetChiralityPerceived (bool value=true)
 
void SetPartialChargesPerceived (bool value=true)
 
void SetHybridizationPerceived (bool value=true)
 
void SetClosureBondsPerceived (bool value=true)
 
void SetHydrogensAdded (bool value=true)
 
void SetCorrectedForPH (bool value=true)
 
void SetSpinMultiplicityAssigned (bool value=true)
 
void SetIsPatternStructure (bool value=true)
 
void SetIsReaction (bool value=true)
 
bool HasFlag (int flag)
 
void SetFlag (int flag)
 
void UnsetFlag (int flag)
 
void SetFlags (int flags)
 
Molecule utilities and perception methods
void FindSSSR ()
 
void FindLSSR ()
 
void FindRingAtomsAndBonds ()
 
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 Not3D=false)
 
bool Has3D ()
 
bool HasNonZeroCoords ()
 
bool HasAromaticPerceived ()
 
bool HasSSSRPerceived ()
 
bool HasLSSRPerceived ()
 
bool HasRingAtomsAndBondsPerceived ()
 
bool HasAtomTypesPerceived ()
 
bool HasRingTypesPerceived ()
 
bool HasChiralityPerceived ()
 
bool HasPartialChargesPerceived ()
 
bool HasHybridizationPerceived ()
 
bool HasClosureBondsPerceived ()
 
bool HasChainsPerceived ()
 
bool HasHydrogensAdded ()
 
bool IsCorrectedForPH ()
 
bool HasSpinMultiplicityAssigned ()
 
bool IsReaction ()
 
bool Empty ()
 
Multiple conformer member functions
int NumConformers ()
 
void SetConformers (std::vector< double *> &v)
 
void AddConformer (double *f)
 
void SetConformer (unsigned 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)
 
size_t 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 Attributes

int _flags
 
bool _autoPartialCharge
 
bool _autoFormalCharge
 
std::string _title
 
std::vector< OBAtom * > _vatom
 
std::vector< OBAtom * > _atomIds
 
std::vector< OBBond * > _vbond
 
std::vector< OBBond * > _bondIds
 
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
 

Molecule modification methods

virtual void BeginModify (void)
 
virtual void EndModify (bool nukePerceivedData=true)
 
int GetMod ()
 
void IncrementMod ()
 
void DecrementMod ()
 
virtual OBBaseDoTransformations (const std::map< std::string, std::string > *pOptions, OBConversion *pConv)
 
bool Clear ()
 
void RenumberAtoms (std::vector< OBAtom *> &)
 
void RenumberAtoms (std::vector< int >)
 
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 DeleteHydrogens ()
 
bool DeleteHydrogens (OBAtom *)
 
bool DeletePolarHydrogens ()
 
bool DeleteNonPolarHydrogens ()
 
bool DeleteHydrogen (OBAtom *)
 
bool AddHydrogens (bool polaronly=false, bool correctForPH=false, double pH=7.4)
 
bool AddHydrogens (OBAtom *)
 
bool AddPolarHydrogens ()
 
bool AddNonPolarHydrogens ()
 
bool AddNewHydrogens (HydrogenType whichHydrogen, bool correctForPH=false, double pH=7.4)
 
bool StripSalts (unsigned int threshold=0)
 
std::vector< OBMolSeparate (int StartIndex=1)
 
bool GetNextFragment (OpenBabel::OBMolAtomDFSIter &iter, OBMol &newMol)
 
bool CopySubstructure (OBMol &newmol, OBBitVec *includeatoms, OBBitVec *excludebonds=(OBBitVec *) 0, unsigned int correctvalence=1, std::vector< unsigned int > *atomorder=(std::vector< unsigned int > *) 0, std::vector< unsigned int > *bondorder=(std::vector< unsigned int > *) 0)
 
bool ConvertDativeBonds ()
 
bool MakeDativeBonds ()
 
bool ConvertZeroBonds ()
 
bool CorrectForPH (double pH=7.4)
 
bool AssignSpinMultiplicity (bool NoImplicitH=false)
 
bool AssignTotalChargeToAtoms (int charge)
 
vector3 Center (int nconf)
 
void SetTorsion (OBAtom *, OBAtom *, OBAtom *, OBAtom *, double ang)
 
static const char * ClassDescription ()
 

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

...
#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/mol.h>
OBMol mol;
double exactMass = 0.0;
{
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.

Examples:
obconformersearch_default.cpp, obconversion_readstring.cpp, and obforcefield_energy.cpp.

Constructor & Destructor Documentation

◆ OBMol() [1/2]

OBMol ( )

Constructor.

◆ OBMol() [2/2]

OBMol ( const OBMol mol)

Copy constructor, copies atoms,bonds and OBGenericData.

◆ ~OBMol()

~OBMol ( )
virtual

Destructor.

Member Function Documentation

◆ operator=()

OBMol & operator= ( const OBMol mol)

Assignment, copies atoms,bonds and OBGenericData.

◆ operator+=()

OBMol & operator+= ( const OBMol mol)

Copies atoms and bonds but not OBGenericData.

◆ ReserveAtoms()

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.

◆ DestroyAtom()

void DestroyAtom ( OBAtom atom)
virtual

Free an OBAtom pointer if defined. Does no bookkeeping

See also
DeleteAtom which ensures internal connections

Referenced by OBMol::~OBMol().

◆ DestroyBond()

void DestroyBond ( OBBond bond)
virtual

Free an OBBond pointer if defined. Does no bookkeeping

See also
DeleteBond which ensures internal connections

Referenced by OBMol::~OBMol().

◆ DestroyResidue()

void DestroyResidue ( OBResidue residue)
virtual

Free an OBResidue pointer if defined. Does no bookkeeping

See also
DeleteResidue which ensures internal connections

Referenced by OBMol::~OBMol().

◆ AddAtom()

bool AddAtom ( OBAtom atom,
bool  forceNewId = false 
)

Add an atom to a molecule.

Add the specified atom to this molecule

Parameters
atomthe atom to add
forceNewIdwhether to make a new atom Id even if the atom already has one (default is false)
Returns
Whether the method was successful

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

Referenced by OpenBabel::addFragment(), and OBMol::CopySubstructure().

◆ InsertAtom()

bool InsertAtom ( OBAtom atom)

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

◆ AddBond() [1/2]

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
beginIdxthe atom index of the "start" atom
endIdxthe atom index of the "end" atom
orderthe bond order (see OBBond::GetBondOrder())
flagsany bond flags such as stereochemistry (default = none)
insertposthe position index to insert the bond (default = none)
Returns
Whether the new bond creation was successful

Referenced by OpenBabel::addFragment(), OBResidueData::AssignBonds(), OBBuilder::Build(), OBMol::ConnectTheDots(), OBMol::CopySubstructure(), AliasData::Expand(), OBAtom::HtoMethyl(), and OBAtom::SetHybAndGeom().

◆ AddBond() [2/2]

bool AddBond ( OBBond bond)

Add the specified residue to this molecule and update connections

Returns
Whether the method was successful

◆ AddResidue()

bool AddResidue ( OBResidue residue)

Add the specified residue to this molecule and update connections

Returns
Whether the method was successful

◆ NewAtom() [1/2]

OBAtom * NewAtom ( )

Create a new OBAtom in this molecule and ensure connections (e.g. OBAtom::GetParent(). A new unique id will be assigned to this atom.

Referenced by OBMol::CopySubstructure(), OBUnitCell::FillUnitCell(), OBAtom::HtoMethyl(), and OBAtom::SetHybAndGeom().

◆ NewAtom() [2/2]

OBAtom * NewAtom ( unsigned long  id)

Instantiate a New Atom and add it to the molecule.

Create a new OBAtom in this molecule and ensure connections. (e.g. OBAtom::GetParent(). The id will be assigned to this atom.

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

◆ NewBond() [1/2]

OBBond * NewBond ( )

Create a new OBBond in this molecule and ensure connections (e.g. OBBond::GetParent(). A new unique id will be assigned to this bond.

Referenced by OBBuilder::Connect().

◆ NewBond() [2/2]

OBBond * NewBond ( unsigned long  id)

Instantiate a New Bond and add it to the molecule.

Create a new OBBond in this molecule and ensure connections (e.g. OBBond::GetParent(). The id will be assigned to this bond.

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

◆ NewResidue()

OBResidue * NewResidue ( )

Create a new OBResidue in this molecule and ensure connections.

Referenced by OBMol::CopySubstructure(), and OBChainsParser::~OBChainsParser().

◆ DeleteAtom()

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(), AliasData::Expand(), OBUnitCell::FillUnitCell(), OpenBabel::InternalToCartesian(), and OBAtom::SetHybAndGeom().

◆ DeleteBond()

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(), OBBuilder::Connect(), OBMol::ConnectTheDots(), OBBuilder::CorrectStereoAtoms(), OBBuilder::IsSpiroAtom(), and OBBuilder::Swap().

◆ DeleteResidue()

bool DeleteResidue ( OBResidue residue,
bool  destroyResidue = true 
)

Deletes a residue from this molecule and updates accordingly.

Returns
Whether deletion was successful

Referenced by OBChainsParser::~OBChainsParser().

◆ BeginModify()

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 OBChemTsfm::Apply(), OBMol::ConnectTheDots(), OBAtom::HtoMethyl(), OBMol::MakeDativeBonds(), and OBAtom::SetHybAndGeom().

◆ EndModify()

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 OpenBabel::addFragment(), OBChemTsfm::Apply(), OBMol::ConnectTheDots(), OBAtom::HtoMethyl(), OBMol::MakeDativeBonds(), and OBAtom::SetHybAndGeom().

◆ GetMod()

int GetMod ( )
inline
Returns
The number of nested BeginModify() calls. Used internally.

Referenced by OpenBabel::intToStr().

◆ IncrementMod()

void IncrementMod ( )
inline

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

Referenced by OBAtom::SetHybAndGeom().

◆ DecrementMod()

void DecrementMod ( )
inline

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

Referenced by OBAtom::SetHybAndGeom().

◆ GetFlags()

int GetFlags ( )
inline
Returns
the entire set of flags. (Internal use, mainly.)

◆ GetTitle()

const char * GetTitle ( bool  replaceNewlines = true) const
virtual

◆ NumAtoms()

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

Referenced by OBRingSearch::AddRingFromClosure(), OBAlign::Align(), OpenBabel::alternate(), OBChemTsfm::Apply(), patty::assign_types(), OpenBabel::AssignOBAromaticityModel(), OBGastChrg::AssignPartialCharges(), OBRotorList::AssignTorVals(), OpenBabel::CartesianToInternal(), OBMol::Center(), OpenBabel::ComparePairSecond(), OpenBabel::CompileMoleculeQuery(), OBBuilder::Connect(), OBMol::ConnectTheDots(), OBMol::CopyConformer(), OBMol::CopySubstructure(), OBRotamerList::CreateConformerList(), OBDepict::DrawMolecule(), AliasData::Expand(), OBSmartsMatcher::FastSingleMatch(), OBChargeModel::FillChargeVectors(), OpenBabel::FindAutomorphisms(), OpenBabel::findMetalloceneBonds(), OpenBabel::FindRingAtomsAndBonds2(), OBRotorList::FindRotors(), OpenBabel::generateDiagram(), OBAlign::GetAlignment(), OBForceField::GetAtomTypes(), OBForceField::GetConformers(), OBForceField::GetCoordinates(), OpenBabel::GetDFFVector(), OpenBabel::GetHeavyAtomCoords(), OpenBabel::GetMaxAtomIdx(), OBMol::GetNextFragment(), OBForceField::GetPartialCharges(), OBSpectrophore::GetSpectrophore(), OBAtom::HtoMethyl(), OpenBabel::InternalToCartesian(), OpenBabel::intToStr(), OBStericConformerFilter::IsGood(), OBForceField::IsSetupNeeded(), OBBuilder::IsSpiroAtom(), OBMoleculeFormat::MakeCombinedMolecule(), OpenBabel::OBBondGetSmallestRingSize(), OBMolAtomBFSIter::OBMolAtomBFSIter(), OBMolAtomDFSIter::OBMolAtomDFSIter(), OBSSMatch::OBSSMatch(), OBMol::operator=(), OBChainsParser::PerceiveChains(), OBMoleculeFormat::ReadChemObjectImpl(), OBMol::RenumberAtoms(), OBMol::Rotate(), OBRMSDConformerScore::Score(), OBEnergyConformerScore::Score(), OBMinimizingEnergyConformerScore::Score(), OBMinimizingRMSDConformerScore::Score(), OBMol::Separate(), OBRotamerList::SetBaseCoordinateSets(), OBForceField::SetConformers(), OBMol::SetCoordinates(), OBForceField::SetCoordinates(), OBAlign::SetRefMol(), OBRotorList::SetRotAtoms(), OBAlign::SetTargetMol(), OBProxGrid::Setup(), OBForceField::Setup(), OBSmartsMatcher::SetupAtomMatchTable(), OBMol::Translate(), OpenBabel::UpdateConformersFromTree(), OBAlign::UpdateCoords(), OpenBabel::visitRing(), OBMoleculeFormat::WriteChemObjectImpl(), and OBChainsParser::~OBChainsParser().

◆ NumBonds()

unsigned int NumBonds ( ) const
inline

◆ NumHvyAtoms()

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

◆ NumResidues()

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

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

◆ NumRotors()

unsigned int NumRotors ( bool  sampleRingBonds = false)
Returns
the number of rotatable bonds. If sampleRingBonds is true, will include rotors within rings (see OBBond::IsRotor() for details)

◆ GetAtom()

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 OpenBabel::alternate(), OBChemTsfm::Apply(), OpenBabel::ApplyRotMatToBond(), OBBondTyper::AssignFunctionalGroupBonds(), OBAtomTyper::AssignHyb(), OBPhModel::AssignSeedPartialCharge(), OBAtomTyper::AssignTypes(), OBRingTyper::AssignTypes(), OBBuilder::Build(), OpenBabel::BuildOBRTreeVector(), OpenBabel::CartesianToInternal(), OBBuilder::Connect(), OBMol::ConnectTheDots(), OBMol::ConvertZeroBonds(), OBMol::CopySubstructure(), OBBuilder::CorrectStereoAtoms(), OBDepict::DrawMolecule(), AliasData::Expand(), OBSmartsMatcher::FastSingleMatch(), OpenBabel::FindAutomorphisms(), OBRing::findCenterAndNormal(), OpenBabel::findMetalloceneBonds(), OpenBabel::FindRingAtomsAndBonds2(), OpenBabel::generateDiagram(), OBAlign::GetAlignment(), OBAtom::GetAngle(), OBForceField::GetAtomTypes(), OBForceField::GetCoordinates(), OpenBabel::GetDFFVector(), OBAtom::GetDistance(), OpenBabel::GetHeavyAtomCoords(), OBForceField::GetPartialCharges(), OBRing::GetRootAtom(), OBRotorRules::GetRotorIncrements(), OpenBabel::groupRedraw(), OpenBabel::intToStr(), OBRing::IsAromatic(), OBStericConformerFilter::IsGood(), OBForceField::IsSetupNeeded(), OBBuilder::IsSpiroAtom(), OBSmartsMatcher::match(), OBAtom::MatchesSMARTS(), OBMolAtomDFSIter::operator++(), OBMolAtomBFSIter::operator++(), OBMol::PerceiveBondOrders(), OBChainsParser::PerceiveChains(), OBMol::RenumberAtoms(), AliasData::RevertToAliasForm(), OBRotorList::SetEvalAtoms(), OBBond::SetLength(), OBAlign::SetRefMol(), OBAlign::SetTargetMol(), OBRotamerList::Setup(), OBFFConstraints::Setup(), OBBuilder::Swap(), OBPointGroup::Symmetrize(), and OBChainsParser::~OBChainsParser().

◆ GetAtomById()

OBAtom * GetAtomById ( unsigned long  id) const

◆ GetFirstAtom()

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

Referenced by OBBuilder::Build().

◆ GetBond() [1/3]

OBBond * GetBond ( int  idx) const

◆ GetBondById()

OBBond * GetBondById ( unsigned long  id) const
Returns
the bond with id or NULL if it does not exist.

Referenced by OpenBabel::findMetalloceneBonds().

◆ GetBond() [2/3]

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.

◆ GetBond() [3/3]

OBBond * GetBond ( OBAtom bgn,
OBAtom end 
) const
Returns
the bond between the atoms bgn and end or NULL if none exists

◆ GetResidue()

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=(), and OBChainsParser::~OBChainsParser().

◆ GetInternalCoord()

std::vector< OBInternalCoord * > GetInternalCoord ( )

◆ GetTorsion() [1/2]

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::CorrectStereoBonds(), OBRotorRules::GetRotorIncrements(), OBMol::PerceiveBondOrders(), and OBRotamerList::Setup().

◆ GetTorsion() [2/2]

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

◆ GetAngle()

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 )

◆ AreInSameRing()

int AreInSameRing ( OBAtom a,
OBAtom b 
)
Returns
the size of the smallest ring if a and b are in the same ring, 0 otherwise
Since
version 2.4
version 2.4

Referenced by OBBuilder::IsSpiroAtom().

◆ GetFormula()

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.

◆ GetSpacedFormula()

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 OBMoleculeFormat::MakeCombinedMolecule().

◆ GetEnergy() [1/2]

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

Referenced by OBMol::operator=().

◆ GetMolWt()

double GetMolWt ( bool  implicitH = true)
Returns
the standard molar mass given by IUPAC atomic masses (amu)

◆ GetExactMass()

double GetExactMass ( bool  implicitH = true)
Returns
the mass given by isotopes (or most abundant isotope, if not specified)

◆ GetTotalCharge()

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

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

◆ GetTotalSpinMultiplicity()

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. if it fails (gives singlet for odd number of electronic systems), then assign wrt parity of the total electrons.

Referenced by OBMol::operator=().

◆ GetDimension()

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

Referenced by OBBuilder::Build(), OpenBabel::CanonicalLabels(), OBMol::CopySubstructure(), OBPhModel::CorrectForPH(), OBMoleculeFormat::MakeCombinedMolecule(), OBMol::operator+=(), and OBMol::operator=().

◆ GetCoordinates()

double* GetCoordinates ( )
inline

◆ GetSSSR()

vector< OBRing * > & GetSSSR ( )

◆ GetLSSR()

vector< OBRing * > & GetLSSR ( )
Returns
the Largest Set of Smallest Rings has been run (see OBRing class)

Referenced by OBMol::AreInSameRing().

◆ AutomaticFormalCharge()

bool AutomaticFormalCharge ( )
inline

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

Referenced by OBPhModel::CorrectForPH().

◆ AutomaticPartialCharge()

bool AutomaticPartialCharge ( )
inline

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

Referenced by OBPhModel::AssignSeedPartialCharge().

◆ SetTitle() [1/2]

void SetTitle ( const char *  title)
virtual

Set the title of this molecule to title.

Reimplemented from OBBase.

Referenced by OBMoleculeFormat::MakeCombinedMolecule().

◆ SetTitle() [2/2]

void SetTitle ( std::string &  title)

Set the title of this molecule to title.

◆ SetFormula()

void SetFormula ( std::string  molFormula)

Set the stochiometric formula for this molecule.

◆ SetEnergy()

void SetEnergy ( double  energy)
inline

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

◆ SetDimension()

void SetDimension ( unsigned short int  d)
inline

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

Referenced by OBBuilder::Build(), OBMol::CopySubstructure(), and AliasData::Expand().

◆ SetTotalCharge()

void SetTotalCharge ( int  charge)

Set the total charge of this molecule to charge.

◆ SetTotalSpinMultiplicity()

void SetTotalSpinMultiplicity ( unsigned int  spinMultiplicity)

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

◆ SetInternalCoord()

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

Set the internal coordinates to int_coord (Does not call InternalToCartesian to update the 3D cartesian coordinates). The size of the int_coord has to be the same as the number of atoms in molecule (+ NULL at the beginning).

◆ SetAutomaticFormalCharge()

void SetAutomaticFormalCharge ( bool  val)
inline

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

◆ SetAutomaticPartialCharge()

void SetAutomaticPartialCharge ( bool  val)
inline

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

◆ SetAromaticPerceived()

void SetAromaticPerceived ( bool  value = true)
inline

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

Referenced by OBAromaticTyper::AssignAromaticFlags(), and OBMol::PerceiveBondOrders().

◆ SetSSSRPerceived()

void SetSSSRPerceived ( bool  value = true)
inline

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

◆ SetLSSRPerceived()

void SetLSSRPerceived ( bool  value = true)
inline

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

◆ SetRingAtomsAndBondsPerceived()

void SetRingAtomsAndBondsPerceived ( bool  value = true)
inline

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

Referenced by OpenBabel::FindRingAtomsAndBonds2().

◆ SetAtomTypesPerceived()

void SetAtomTypesPerceived ( bool  value = true)
inline

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

Referenced by OBAtomTyper::AssignTypes().

◆ SetRingTypesPerceived()

void SetRingTypesPerceived ( bool  value = true)
inline

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

Referenced by OBRingTyper::AssignTypes().

◆ SetChainsPerceived()

void SetChainsPerceived ( bool  value = true)
inline

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

Referenced by OBChainsParser::PerceiveChains().

◆ SetChiralityPerceived()

void SetChiralityPerceived ( bool  value = true)
inline

Mark that chirality has been perceived.

Referenced by OBBuilder::Build().

◆ SetPartialChargesPerceived()

void SetPartialChargesPerceived ( bool  value = true)
inline

Mark that partial charges have been assigned.

Referenced by OBPhModel::AssignSeedPartialCharge().

◆ SetHybridizationPerceived()

void SetHybridizationPerceived ( bool  value = true)
inline

Mark that hybridization of all atoms has been assigned.

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

◆ SetClosureBondsPerceived()

void SetClosureBondsPerceived ( bool  value = true)
inline

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

Referenced by OpenBabel::FindRingAtomsAndBonds2().

◆ SetHydrogensAdded()

void SetHydrogensAdded ( bool  value = true)
inline

Mark that explicit hydrogen atoms have been added.

◆ SetCorrectedForPH()

void SetCorrectedForPH ( bool  value = true)
inline

Referenced by OBPhModel::CorrectForPH().

◆ SetSpinMultiplicityAssigned()

void SetSpinMultiplicityAssigned ( bool  value = true)
inline

◆ SetIsPatternStructure()

void SetIsPatternStructure ( bool  value = true)
inline

The OBMol is a pattern, not a complete molecule. Left unchanged by Clear().

Referenced by OpenBabel::alternate(), and AliasData::Expand().

◆ SetIsReaction()

void SetIsReaction ( bool  value = true)
inline

◆ HasFlag()

bool HasFlag ( int  flag)
inline

Referenced by OBMol::operator=().

◆ SetFlag()

void SetFlag ( int  flag)
inline

Referenced by OBMol::CopySubstructure().

◆ UnsetFlag()

void UnsetFlag ( int  flag)
inline

Referenced by OBMol::RenumberAtoms().

◆ SetFlags()

void SetFlags ( int  flags)
inline

◆ DoTransformations()

OBBase * DoTransformations ( const std::map< std::string, std::string > *  ,
OBConversion  
)
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() Base type does nothing

Reimplemented from OBBase.

Referenced by OBMoleculeFormat::ReadChemObjectImpl().

◆ ClassDescription()

const char * ClassDescription ( )
static

◆ Clear()

bool Clear ( void  )
virtual

Clear all information from a molecule except OB_PATTERN_STRUCTURE left unchanged.

Reimplemented from OBBase.

Referenced by OpenBabel::addFragment(), OpenBabel::alternate(), OBMoleculeFormat::ReadNameIndex(), and OBMol::Separate().

◆ RenumberAtoms() [1/2]

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.

Referenced by OBMol::RenumberAtoms().

◆ RenumberAtoms() [2/2]

void RenumberAtoms ( std::vector< int >  v)

Renumber the atoms of this molecule using the initial indexes in the supplied vector.

Renumber the atoms according to the order of indexes in the supplied vector This with assemble an atom vector and call RenumberAtoms(vector<OBAtom*>) It will return without action if the supplied vector is empty or does not have the same number of atoms as the molecule.

Since
version 2.3

◆ SetCoordinates()

void SetCoordinates ( double *  c)

Set the coordinates for all atoms in this conformer.

See also
OBMol::GetCoordinates()

◆ ToInertialFrame() [1/2]

void ToInertialFrame ( int  conf,
double *  rmat 
)

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

◆ ToInertialFrame() [2/2]

void ToInertialFrame ( )

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

◆ Translate() [1/2]

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

◆ Translate() [2/2]

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.

◆ Rotate() [1/3]

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

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

Referenced by OBMol::Rotate().

◆ Rotate() [2/3]

void Rotate ( const double  m[9])

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

◆ Rotate() [3/3]

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

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

◆ Center() [1/2]

void Center ( )

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

Referenced by OBPointGroup::Setup().

◆ DeleteHydrogens() [1/2]

bool DeleteHydrogens ( )

Delete all hydrogens from the molecule

Returns
Success

Referenced by OBPhModel::CorrectForPH(), and AliasData::Expand().

◆ DeleteHydrogens() [2/2]

bool DeleteHydrogens ( OBAtom atom)

Delete all hydrogens from the supplied atom

Returns
Success

◆ DeletePolarHydrogens()

bool DeletePolarHydrogens ( )

Delete all hydrogen atoms connected to a polar atom

See also
OBAtom::IsPolarHydrogen
Since
version 2.4

◆ DeleteNonPolarHydrogens()

bool DeleteNonPolarHydrogens ( )

Delete all hydrogen atoms connected to a non-polar atom

See also
OBAtom::IsNonPolarHydrogen

◆ DeleteHydrogen()

bool DeleteHydrogen ( OBAtom atom)

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

◆ AddHydrogens() [1/2]

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
polaronlyWhether to add hydrogens only to polar atoms (i.e., not to C atoms)
correctForPHWhether to call CorrectForPH() first
pHThe pH to use for CorrectForPH() modification
Returns
Whether any hydrogens were added

Referenced by OBSmartsPattern::Match().

◆ AddHydrogens() [2/2]

bool AddHydrogens ( OBAtom atom)

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

◆ AddPolarHydrogens()

bool AddPolarHydrogens ( )

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

◆ AddNonPolarHydrogens()

bool AddNonPolarHydrogens ( )

Add only nonpolar hydrogens (i.e., attached to C)

Since
version 2.4

◆ AddNewHydrogens()

bool AddNewHydrogens ( HydrogenType  whichHydrogen,
bool  correctForPH = false,
double  pH = 7.4 
)

Add polar and/or nonpolar hydrogens

Since
verison 2.4

◆ StripSalts()

bool StripSalts ( unsigned int  threshold = 0)

If threshold is not specified or is zero, remove all but the largest contiguous fragment. If threshold is non-zero, remove any fragments with fewer than threshold atoms.

◆ Separate()

vector< OBMol > Separate ( int  StartIndex = 1)

Copies each disconnected fragment as a separate OBMol.

Referenced by OBBuilder::Build(), and OBMoleculeFormat::ReadChemObjectImpl().

◆ GetNextFragment()

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

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

Referenced by OBMol::Separate().

◆ CopySubstructure()

bool CopySubstructure ( OBMol newmol,
OBBitVec atoms,
OBBitVec excludebonds = (OBBitVec*)0,
unsigned int  correctvalence = 1,
std::vector< unsigned int > *  atomorder = (std::vector<unsigned int>*)0,
std::vector< unsigned int > *  bondorder = (std::vector<unsigned int>*)0 
)

Copy part of a molecule to another molecule.

This function copies a substructure of a molecule to another molecule. The key information needed is an OBBitVec indicating which atoms to include and (optionally) an OBBitVec indicating which bonds to exclude. By default, only bonds joining included atoms are copied.

When an atom is copied, but not all of its bonds are, by default hydrogen counts are adjusted to account for the missing bonds. That is, given the SMILES "CF", if we copy the two atoms but exclude the bond, we will end up with "C.F". This behavior can be changed by specifiying a value other than 1 for the correctvalence parameter. A value of 0 will yield "[C].[F]" while 2 will yield "C*.F*" (see correctvalence below for more information).

Aromaticity is preserved as present in the original OBMol. If this is not desired, the user should call OBMol::SetAromaticPerceived(false) on the new OBMol.

Stereochemistry is only preserved if the corresponding elements are wholly present in the substructure. For example, all four atoms and bonds of a tetrahedral stereocenter must be copied.

Residue information is preserved if the original OBMol is marked as having its residues perceived. If this is not desired, either call OBMol::SetChainsPerceived(false) in advance on the original OBMol to avoid copying the residues (and then reset it afterwards), or else call it on the new OBMol so that residue information will be reperceived (when requested).

Here is an example of using this method to copy ring systems to a new molecule. Given the molecule represented by the SMILES string, "FC1CC1c2ccccc2I", we will end up with a new molecule represented by the SMILES string, "C1CC1.c2ccccc2".

OBBitVec atoms(mol.NumAtoms() + 1); // the maximum size needed
FOR_ATOMS_OF_MOL(atom, mol) {
if(atom->IsInRing())
atoms.SetBitOn(atom->Idx());
}
OBBitVec excludebonds(mol.NumBonds()); // the maximum size needed
FOR_BONDS_OF_MOL(bond, mol) {
if(!bond->IsInRing())
excludebonds.SetBitOn(bond->Idx());
}
OBMol newmol;
mol.CopySubstructure(&newmol, &atoms, &excludebonds);

When used from Python, note that "None" may be used to specify an empty value for the excludebonds parameter.

Remarks
Some alternatives to using this function, which may be preferred in some instances due to efficiency or convenience are:
  1. Copying the entire OBMol, and then deleting the unwanted parts
  2. Modifiying the original OBMol, and then restoring it
  3. Using the SMILES writer option -xf to specify fragment atom idxs
Returns
A boolean indicating success or failure. Currently failure is only reported if one of the specified atoms is not present, or atoms is a NULL pointer.
Parameters
newmolThe molecule to which to add the substructure. Note that atoms are appended to this molecule.
atomsAn OBBitVec, indexed by atom Idx, specifying which atoms to copy
excludebondsAn OBBitVec, indexed by bond Idx, specifying a list of bonds to exclude. By default, all bonds between the specified atoms are included - this parameter overrides that.
correctvalenceA value of 0, 1 (default) or 2 that indicates how atoms with missing bonds are handled: 0 - Leave the implicit hydrogen count unchanged; 1 - Adjust the implicit hydrogen count to correct for the missing bonds; 2 - Replace the missing bonds with bonds to dummy atoms
atomorderRecord the Idxs of the original atoms. That is, the first element in this vector will be the Idx of the atom in the original OBMol that corresponds to the first atom in the new OBMol. Note that the information is appended to this vector.
bondorderRecord the Idxs of the original bonds. See atomorder above.

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

◆ ConvertDativeBonds()

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.

◆ MakeDativeBonds()

bool MakeDativeBonds ( )

Converts 5-valent N and P only. Return true if conversion occurred.

Returns
has charged form of dative bonds(e.g.[N+]([O-])=O from N(=O)=O).
Since
version 2.4

Converts 5-valent N to charged form of dative bonds, e.g. -N(=O)=O converted to -[N+]([O-])=O. Returns true if conversion occurs.

◆ ConvertZeroBonds()

bool ConvertZeroBonds ( )

Convert zero-order bonds to single or double bonds and adjust adjacent atom charges in an attempt to achieve the correct valence state.

Returns
Whether any modifications were made
Since
version 2.4

This function is useful when writing to legacy formats (such as MDL MOL) that do not support zero-order bonds. It is worth noting that some compounds cannot be well represented using just single, double and triple bonds, even with adjustments to adjacent charges. In these cases, simply converting zero-order bonds to single bonds is all that can be done.

Algorithm from:
Clark, A. M. Accurate Specification of Molecular Structures: The Case for
Zero-Order Bonds and Explicit Hydrogen Counting. Journal of Chemical Information
and Modeling, 51, 3149-3157 (2011). http://pubs.acs.org/doi/abs/10.1021/ci200488k

◆ CorrectForPH()

bool CorrectForPH ( double  pH = 7.4)

Correct for pH by applying the OBPhModel transformations.

◆ AssignSpinMultiplicity()

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 OpenBabel::alternate().

◆ AssignTotalChargeToAtoms()

bool AssignTotalChargeToAtoms ( int  charge)

Put the specified molecular charge on appropriate atoms. Assumes all the hydrogen is explicitly included in the molecule.

Since
version 2.4

◆ Center() [2/2]

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

◆ SetTorsion()

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::Connect(), and OBBuilder::CorrectStereoBonds().

◆ FindSSSR()

void FindSSSR ( )

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

Referenced by OBAtom::IsInRingSize(), OBAtom::MemberOfRingCount(), OBAtom::MemberOfRingSize(), OBMolRingIter::OBMolRingIter(), and OBMol::PerceiveBondOrders().

◆ FindLSSR()

void FindLSSR ( )

Find Largest Set of Smallest Rings.

◆ FindRingAtomsAndBonds()

void FindRingAtomsAndBonds ( )

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

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

◆ FindChildren() [1/2]

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 OpenBabel::ApplyRotMatToBond(), OBRotorList::AssignTorVals(), OBBuilder::IsSpiroAtom(), OBBond::SetLength(), OBRotorList::SetRotAtoms(), OBRotorList::SetRotAtomsByFix(), and OBRotamerList::Setup().

◆ FindChildren() [2/2]

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'

◆ FindLargestFragment()

void FindLargestFragment ( OBBitVec frag)

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

Parameters
fragReturn (by reference) a bit vector indicating the atoms in the largest fragment

◆ ContigFragList()

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 OBMol::ConvertZeroBonds().

◆ Align()

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

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

◆ ConnectTheDots()

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.

◆ PerceiveBondOrders()

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.

◆ FindAngles()

void FindAngles ( )

Fills out an OBAngleData with angles from the molecule.

Referenced by OBMolAngleIter::OBMolAngleIter().

◆ FindTorsions()

void FindTorsions ( )

Fills out an OBTorsionData with angles from the molecule.

Referenced by OBMolTorsionIter::OBMolTorsionIter().

◆ GetGTDVector()

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

Calculates the graph theoretical distance (GTD) of each atom.

Creates a vector (indexed from zero) containing, for each atom in the molecule, the number of bonds plus one to the most distant non-H atom.

For example, for the molecule H3CC(=O)Cl the GTD value for C1 would be 3, as the most distant non-H atom (either Cl or O) is 2 bonds away.

Since the GTD measures the distance to non-H atoms, the GTD values for terminal H atoms tend to be larger than for non-H terminal atoms. In the example above, the GTD values for the H atoms are all 4.

Referenced by OBRotorList::FindRotors(), and OpenBabel::intToStr().

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

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

◆ Has2D()

bool Has2D ( bool  Not3D = false)

Are there non-zero coordinates in two dimensions (i.e. X and Y)- and, if Not3D is true, no Z coordinates?

Referenced by AliasData::Expand().

◆ Has3D()

bool Has3D ( )

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

Referenced by AliasData::Expand().

◆ HasNonZeroCoords()

bool HasNonZeroCoords ( )

Are there any non-zero coordinates?

◆ HasAromaticPerceived()

bool HasAromaticPerceived ( )
inline

Has aromatic perception been performed?

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

◆ HasSSSRPerceived()

bool HasSSSRPerceived ( )
inline

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

Referenced by OBAtom::IsInRingSize(), OBAtom::MemberOfRingCount(), OBAtom::MemberOfRingSize(), OBMolRingIter::OBMolRingIter(), and OBMol::PerceiveBondOrders().

◆ HasLSSRPerceived()

bool HasLSSRPerceived ( )
inline

Has the largest set of smallest rings (FindLSSR) been performed?

◆ HasRingAtomsAndBondsPerceived()

bool HasRingAtomsAndBondsPerceived ( )
inline

Have ring atoms and bonds been assigned?

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

◆ HasAtomTypesPerceived()

bool HasAtomTypesPerceived ( )
inline

Have atom types been assigned by OBAtomTyper?

Referenced by OBAtom::GetType().

◆ HasRingTypesPerceived()

bool HasRingTypesPerceived ( )
inline

Have ring types been assigned by OBRingTyper?

Referenced by OBRing::GetType().

◆ HasChiralityPerceived()

bool HasChiralityPerceived ( )
inline

Has atom chirality been assigned?

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

◆ HasPartialChargesPerceived()

bool HasPartialChargesPerceived ( )
inline

Have atomic Gasteiger partial charges been assigned by OBGastChrg?

◆ HasHybridizationPerceived()

bool HasHybridizationPerceived ( )
inline

Has atomic hybridization been assigned by OBAtomTyper?

Referenced by OBAtom::GetHyb().

◆ HasClosureBondsPerceived()

bool HasClosureBondsPerceived ( )
inline

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

Referenced by OpenBabel::DetermineFRJ(), and OBBond::IsClosure().

◆ HasChainsPerceived()

bool HasChainsPerceived ( )
inline

Have biomolecule chains and residues been assigned by OBChainsParser?

Referenced by OBMol::CopySubstructure(), and OBAtom::GetResidue().

◆ HasHydrogensAdded()

bool HasHydrogensAdded ( )
inline

Have hydrogens been added to the molecule?

◆ IsCorrectedForPH()

bool IsCorrectedForPH ( )
inline

Has the molecule been corrected for pH by CorrectForPH?

Referenced by OBPhModel::CorrectForPH().

◆ HasSpinMultiplicityAssigned()

bool HasSpinMultiplicityAssigned ( )
inline

Has total spin multiplicity been assigned?

◆ IsReaction()

bool IsReaction ( )
inline

Does this OBMol represent a reaction?

Referenced by OBMoleculeFormat::ReadChemObjectImpl().

◆ Empty()

bool Empty ( )
inline

Are there any atoms in this molecule?

Referenced by OBMol::ConnectTheDots(), OBSSMatch::OBSSMatch(), OBMol::PerceiveBondOrders(), and OBMol::RenumberAtoms().

◆ NumConformers()

int NumConformers ( )
inline

◆ SetConformers()

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

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

Referenced by OBConformerSearch::GetConformers(), and OBForceField::GetConformers().

◆ AddConformer()

void AddConformer ( double *  f)
inline

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

Referenced by OpenBabel::UpdateConformersFromTree().

◆ SetConformer()

void SetConformer ( unsigned int  i)

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

Referenced by OBMol::Center(), OBMoleculeFormat::DoOutputOptions(), and OBForceField::GetConformers().

◆ CopyConformer()

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

◆ DeleteConformer()

void DeleteConformer ( int  nconf)

Delete the conformer nconf.

◆ GetConformer()

double* GetConformer ( int  i)
inline

◆ SetEnergies()

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

Set the entire set of conformer energies.

◆ GetEnergies()

vector< double > GetEnergies ( )

Set the entire set of conformer energies.

◆ GetEnergy() [2/2]

double GetEnergy ( int  ci)

Get the energy for conformer ci

ci conformer index

◆ BeginConformer()

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

◆ NextConformer()

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

◆ GetConformers()

std::vector<double*>& GetConformers ( )
inline
Returns
the entire set of conformers for this molecule as a vector of floating point arrays

Referenced by OBRotamerList::SetBaseCoordinateSets().

◆ BeginAtoms()

OBAtomIterator BeginAtoms ( )
inline
Returns
An atom iterator pointing to the beginning of the atom list

◆ EndAtoms()

OBAtomIterator EndAtoms ( )
inline
Returns
An atom iterator pointing to the end of the atom list

◆ BeginBonds()

OBBondIterator BeginBonds ( )
inline
Returns
A bond iterator pointing to the beginning of the bond list

◆ EndBonds()

OBBondIterator EndBonds ( )
inline
Returns
A bond iterator pointing to the end of the bond list

◆ BeginResidues()

OBResidueIterator BeginResidues ( )
inline
Returns
A residue iterator pointing to the beginning of the residue list

◆ EndResidues()

OBResidueIterator EndResidues ( )
inline
Returns
A residue iterator pointing to the end of the residue list

◆ BeginAtom()

OBAtom * BeginAtom ( OBAtomIterator i)

◆ NextAtom()

OBAtom * NextAtom ( OBAtomIterator i)

◆ BeginBond()

OBBond * BeginBond ( OBBondIterator i)

◆ NextBond()

OBBond * NextBond ( OBBondIterator i)

◆ BeginResidue()

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+=(), OBChainsParser::~OBChainsParser(), and OBMol::~OBMol().

◆ NextResidue()

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+=(), OBChainsParser::~OBChainsParser(), and OBMol::~OBMol().

◆ BeginInternalCoord()

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

◆ NextInternalCoord()

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

◆ CastAndClear()

T* CastAndClear ( bool  clear = true)
inlineinherited

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

◆ HasData() [1/3]

bool HasData ( const std::string &  s)
inherited

◆ HasData() [2/3]

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

◆ HasData() [3/3]

bool HasData ( const unsigned int  type)
inherited
Returns
whether the generic attribute/value pair exists, for a given OBGenericDataType

◆ DeleteData() [1/4]

void DeleteData ( unsigned int  type)
inherited

◆ DeleteData() [2/4]

void DeleteData ( OBGenericData gd)
inherited

Delete the given generic data from this object.

◆ DeleteData() [3/4]

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

Delete all of the given generic data from this object.

◆ DeleteData() [4/4]

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

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

◆ SetData()

void SetData ( OBGenericData d)
inlineinherited

◆ CloneData()

void CloneData ( OBGenericData d)
inherited

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

Since
version 2.2

Referenced by AliasData::Expand().

◆ DataSize()

size_t DataSize ( ) const
inlineinherited
Returns
the number of OBGenericData items attached to this molecule.

◆ GetData() [1/5]

OBGenericData * GetData ( const unsigned int  type)
inherited

◆ GetData() [2/5]

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

◆ GetData() [3/5]

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

◆ GetData() [4/5]

std::vector<OBGenericData*>& GetData ( )
inlineinherited
Returns
all data, suitable for iterating

Referenced by OBMol::GetEnergies(), OBMol::GetEnergy(), and OBMol::SetEnergies().

◆ GetData() [5/5]

std::vector< OBGenericData * > GetData ( DataOrigin  source)
inherited
Returns
all data with a specific origin, suitable for iterating

◆ GetAllData()

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

Referenced by OpenBabel::CanonicalLabels(), OBMol::CopySubstructure(), OBBuilder::CorrectStereoAtoms(), OBBuilder::CorrectStereoBonds(), OpenBabel::DeleteStereoOnAtom(), and OBMol::operator+=().

◆ BeginData()

OBDataIterator BeginData ( )
inlineinherited
Returns
An iterator pointing to the beginning of the data

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

◆ EndData()

OBDataIterator EndData ( )
inlineinherited
Returns
An iterator pointing to the end of the data

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

Member Data Documentation

◆ _flags

int _flags
protected

◆ _autoPartialCharge

bool _autoPartialCharge
protected

Assign partial charges automatically.

Referenced by OBMol::OBMol().

◆ _autoFormalCharge

bool _autoFormalCharge
protected

Assign formal charges automatically.

Referenced by OBMol::OBMol().

◆ _title

std::string _title
protected

Molecule title.

Referenced by OBMol::OBMol().

◆ _vatom

std::vector<OBAtom*> _vatom
protected

◆ _atomIds

std::vector<OBAtom*> _atomIds
protected

vector of atoms indexed by id

Referenced by OBMol::OBMol().

◆ _vbond

std::vector<OBBond*> _vbond
protected

vector of bonds

Referenced by OBMol::BeginBond(), OBMol::NextBond(), and OBMol::OBMol().

◆ _bondIds

std::vector<OBBond*> _bondIds
protected

vector of bonds

Referenced by OBMol::OBMol().

◆ _dimension

unsigned short int _dimension
protected

Dimensionality of coordinates.

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

◆ _totalCharge

int _totalCharge
protected

Total charge on the molecule.

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

◆ _totalSpin

unsigned int _totalSpin
protected

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

◆ _c

double* _c
protected

◆ _vconf

std::vector<double*> _vconf
protected

◆ _energy

double _energy
protected

heat of formation

Referenced by OBMol::OBMol().

◆ _natoms

unsigned int _natoms
protected

Number of atoms.

Referenced by OBMol::ConnectTheDots(), and OBMol::OBMol().

◆ _nbonds

unsigned int _nbonds
protected

Number of bonds.

Referenced by OBMol::OBMol().

◆ _residue

std::vector<OBResidue*> _residue
protected

Residue information (if applicable)

◆ _internals

std::vector<OBInternalCoord*> _internals
protected

Internal Coordinates (if applicable)

◆ _mod

unsigned short int _mod
protected

Number of nested calls to BeginModify()

Referenced by OBMol::OBMol().

◆ _vdata

std::vector<OBGenericData*> _vdata
protectedinherited

Custom data.

Referenced by OBMol::OBMol().


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