Open Babel
3.0
|
#include <openbabel/mol.h>
Public Member Functions | |
template<class T > | |
T * | CastAndClear (bool clear=true) |
Initialization and data (re)size methods | |
OBMol () | |
OBMol (const OBMol &) | |
virtual | ~OBMol () |
OBMol & | operator= (const OBMol &mol) |
OBMol & | operator+= (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 &) |
OBAtom * | NewAtom () |
OBAtom * | NewAtom (unsigned long id) |
OBBond * | NewBond () |
OBBond * | NewBond (unsigned long id) |
OBResidue * | NewResidue () |
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) |
OBAtom * | GetAtom (int idx) const |
OBAtom * | GetAtomById (unsigned long id) const |
OBAtom * | GetFirstAtom () const |
OBBond * | GetBond (int idx) const |
OBBond * | GetBondById (unsigned long id) const |
OBBond * | GetBond (int a, int b) const |
OBBond * | GetBond (OBAtom *bgn, OBAtom *end) const |
OBResidue * | GetResidue (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 () |
OBAtom * | BeginAtom (OBAtomIterator &i) |
OBAtom * | NextAtom (OBAtomIterator &i) |
OBBond * | BeginBond (OBBondIterator &i) |
OBBond * | NextBond (OBBondIterator &i) |
OBResidue * | BeginResidue (OBResidueIterator &i) |
OBResidue * | NextResidue (OBResidueIterator &i) |
OBInternalCoord * | BeginInternalCoord (std::vector< OBInternalCoord *>::iterator &i) |
OBInternalCoord * | NextInternalCoord (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 |
OBGenericData * | GetData (const unsigned int type) |
OBGenericData * | GetData (const std::string &) |
OBGenericData * | GetData (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 OBBase * | DoTransformations (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< OBMol > | Separate (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 () |
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:
For example:
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:
or
or
or
It is important to note that atom arrays currently begin at 1 and bond arrays begin at 0. Requesting atom 0 (
will result in an error, but
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.)
These convenience functions can be used like so:
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.
OBMol | ( | ) |
Constructor.
Copy constructor, copies atoms,bonds and OBGenericData.
|
virtual |
Destructor.
Assignment, copies atoms,bonds and OBGenericData.
Copies atoms and bonds but not OBGenericData.
|
inline |
Reserve a minimum number of atoms for internal storage This improves performance since the internal atom vector does not grow.
|
virtual |
Free an OBAtom pointer if defined. Does no bookkeeping
Referenced by OBMol::~OBMol().
|
virtual |
Free an OBBond pointer if defined. Does no bookkeeping
Referenced by OBMol::~OBMol().
|
virtual |
Free an OBResidue pointer if defined. Does no bookkeeping
Referenced by OBMol::~OBMol().
bool AddAtom | ( | OBAtom & | atom, |
bool | forceNewId = false |
||
) |
Add an atom to a molecule.
Add the specified atom to this molecule
atom | the atom to add |
forceNewId | whether to make a new atom Id even if the atom already has one (default is false) |
Also checks bond_queue for any bonds that should be made to the new atom
Referenced by OpenBabel::addFragment(), and OBMol::CopySubstructure().
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
beginIdx | the atom index of the "start" atom |
endIdx | the atom index of the "end" atom |
order | the bond order (see OBBond::GetBondOrder()) |
flags | any bond flags such as stereochemistry (default = none) |
insertpos | the position index to insert the bond (default = none) |
Referenced by OpenBabel::addFragment(), OBResidueData::AssignBonds(), OBBuilder::Build(), OBMol::ConnectTheDots(), OBMol::CopySubstructure(), AliasData::Expand(), OBAtom::HtoMethyl(), and OBAtom::SetHybAndGeom().
bool AddBond | ( | OBBond & | bond | ) |
Add the specified residue to this molecule and update connections
bool AddResidue | ( | OBResidue & | residue | ) |
Add the specified residue to this molecule and update connections
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().
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.
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().
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.
OBResidue * NewResidue | ( | ) |
Create a new OBResidue in this molecule and ensure connections.
Referenced by OBMol::CopySubstructure(), and OBChainsParser::~OBChainsParser().
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.
Referenced by OBChemTsfm::Apply(), AliasData::Expand(), OBUnitCell::FillUnitCell(), OpenBabel::InternalToCartesian(), and OBAtom::SetHybAndGeom().
bool DeleteBond | ( | OBBond * | bond, |
bool | destroyBond = true |
||
) |
Deletes an bond from this molecule and updates accordingly
Referenced by OBBuilder::Build(), OBBuilder::Connect(), OBMol::ConnectTheDots(), OBBuilder::CorrectStereoAtoms(), OBBuilder::IsSpiroAtom(), and OBBuilder::Swap().
bool DeleteResidue | ( | OBResidue * | residue, |
bool | destroyResidue = true |
||
) |
Deletes a residue from this molecule and updates accordingly.
Referenced by OBChainsParser::~OBChainsParser().
|
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().
|
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().
|
inline |
Referenced by OpenBabel::intToStr().
|
inline |
Increase the number of nested BeginModify calls. Dangerous! Instead, properly use BeginModify as needed.
Referenced by OBAtom::SetHybAndGeom().
|
inline |
Decrease the number of nested BeginModify calls. Dangerous! Instead, properly use EndModify as needed.
Referenced by OBAtom::SetHybAndGeom().
|
inline |
|
virtual |
replaceNewlines | whether to replace any newline characters with spaces |
Reimplemented from OBBase.
Referenced by OBBuilder::CorrectStereoAtoms(), OBMoleculeFormat::DeferMolOutput(), OBMoleculeFormat::MakeCombinedMolecule(), OBMol::operator+=(), OBMol::operator=(), OBMol::PerceiveBondOrders(), OBMoleculeFormat::ReadChemObjectImpl(), OBMoleculeFormat::ReadNameIndex(), OBMoleculeFormat::WriteChemObjectImpl(), and OpenBabel::WriteTitles().
|
inline |
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().
|
inline |
Referenced by OpenBabel::alternate(), OBBuilder::Build(), OBDepict::DrawMolecule(), AliasData::Expand(), OpenBabel::FindRingAtomsAndBonds2(), OpenBabel::GetMaxBondIdx(), OpenBabel::intToStr(), OBForceField::IsSetupNeeded(), OBMoleculeFormat::MakeCombinedMolecule(), OBMolBondBFSIter::OBMolBondBFSIter(), OBMol::operator=(), OpenBabel::visitRing(), and OBChainsParser::~OBChainsParser().
unsigned int NumHvyAtoms | ( | ) |
|
inline |
Referenced by OBMol::operator=(), and OBChainsParser::~OBChainsParser().
unsigned int NumRotors | ( | bool | sampleRingBonds = false | ) |
OBAtom * GetAtom | ( | int | idx | ) | const |
idx
or NULL if it does not exist. 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().
OBAtom * GetAtomById | ( | unsigned long | id | ) | const |
id
or NULL if it does not exist. Referenced by OpenBabel::CanonicalLabels(), OBMol::CopySubstructure(), OBBuilder::CorrectStereoAtoms(), OBBuilder::CorrectStereoBonds(), AliasData::Expand(), OpenBabel::findMetalloceneBonds(), and OBBuilder::IsSpiroAtom().
OBAtom * GetFirstAtom | ( | ) | const |
Referenced by OBBuilder::Build().
OBBond * GetBond | ( | int | idx | ) | const |
idx
or NULL if it does not exist. Returns a pointer to the bond after a safety check 0 <= idx < NumBonds
Referenced by OpenBabel::addNbrs(), OpenBabel::alternate(), OBChemTsfm::Apply(), OpenBabel::atomRingToBondRing(), OBBuilder::Build(), OpenBabel::CanonicalLabels(), OBBuilder::Connect(), OBBuilder::CorrectStereoAtoms(), OBBuilder::CorrectStereoBonds(), OBDepict::DrawMolecule(), AliasData::Expand(), OBSmartsMatcher::FastSingleMatch(), OpenBabel::findMetalloceneBonds(), OpenBabel::intToStr(), OBForceField::IsSetupNeeded(), OBBuilder::IsSpiroAtom(), OBMolBondBFSIter::OBMolBondBFSIter(), OBMolBondBFSIter::operator++(), OBMol::PerceiveBondOrders(), OBRotamerList::Setup(), and OBBuilder::Swap().
OBBond * GetBondById | ( | unsigned long | id | ) | const |
id
or NULL if it does not exist. Referenced by OpenBabel::findMetalloceneBonds().
OBBond * GetBond | ( | int | a, |
int | b | ||
) | const |
a
and b
or NULL if none exists. bgn
and end
or NULL if none exists OBResidue * GetResidue | ( | int | idx | ) | const |
idx
, or NULL if none exists Referenced by OBMol::operator=(), and OBChainsParser::~OBChainsParser().
std::vector< OBInternalCoord * > GetInternalCoord | ( | ) |
double GetTorsion | ( | int | a, |
int | b, | ||
int | c, | ||
int | d | ||
) |
Referenced by OBBuilder::CorrectStereoBonds(), OBRotorRules::GetRotorIncrements(), OBMol::PerceiveBondOrders(), and OBRotamerList::Setup().
a
, b
, c
, and d
) WARNING: SetTorsion takes an angle in radians while GetTorsion returns it in degrees a
, b
and c
(where a-> b (vertex) -> c ) Referenced by OBBuilder::IsSpiroAtom().
string GetFormula | ( | ) |
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 |
||
) |
Stochoimetric formula in spaced format e.g. C 4 H 6 O 1 No pair data is stored. Normally use without parameters: GetSpacedFormula()
Referenced by OBMoleculeFormat::MakeCombinedMolecule().
|
inline |
Referenced by OBMol::operator=().
double GetMolWt | ( | bool | implicitH = true | ) |
double GetExactMass | ( | bool | implicitH = true | ) |
int GetTotalCharge | ( | ) |
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=().
unsigned int GetTotalSpinMultiplicity | ( | ) |
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=().
|
inline |
Referenced by OBBuilder::Build(), OpenBabel::CanonicalLabels(), OBMol::CopySubstructure(), OBPhModel::CorrectForPH(), OBMoleculeFormat::MakeCombinedMolecule(), OBMol::operator+=(), and OBMol::operator=().
|
inline |
Referenced by OBScoreGrid::Eval(), OBEnergyConformerScore::Score(), OBMinimizingEnergyConformerScore::Score(), OBMinimizingRMSDConformerScore::Score(), OBRotamerList::SetCurrentCoordinates(), and OBProxGrid::Setup().
vector< OBRing * > & GetSSSR | ( | ) |
Implements blue-obelisk:findSmallestSetOfSmallestRings.
Referenced by OBRingTyper::AssignTypes(), OBDepict::DrawMolecule(), OBBond::FindSmallestRing(), OBAtom::IsInRingSize(), OBAtom::MemberOfRingCount(), OBAtom::MemberOfRingSize(), OBMol::PerceiveBondOrders(), and OBRotor::SetRings().
vector< OBRing * > & GetLSSR | ( | ) |
Referenced by OBMol::AreInSameRing().
|
inline |
Get the current flag for whether formal charges are set with pH correction.
Referenced by OBPhModel::CorrectForPH().
|
inline |
Get the current flag for whether partial charges are auto-determined.
Referenced by OBPhModel::AssignSeedPartialCharge().
|
virtual |
Set the title of this molecule to title
.
Reimplemented from OBBase.
Referenced by 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.
|
inline |
Set the heat of formation for this molecule (in kcal/mol)
|
inline |
Set the dimension of this molecule (i.e., 0, 1 , 2, 3)
Referenced by OBBuilder::Build(), OBMol::CopySubstructure(), and AliasData::Expand().
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 | ) |
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).
|
inline |
Set the flag for determining automatic formal charges with pH (default=true)
|
inline |
Set the flag for determining partial charges automatically (default=true)
|
inline |
Mark that aromaticity has been perceived for this molecule (see OBAromaticTyper)
Referenced by OBAromaticTyper::AssignAromaticFlags(), and OBMol::PerceiveBondOrders().
|
inline |
Mark that Smallest Set of Smallest Rings has been run (see OBRing class)
|
inline |
Mark that Largest Set of Smallest Rings has been run (see OBRing class)
|
inline |
Mark that rings have been perceived (see OBRing class for details)
Referenced by OpenBabel::FindRingAtomsAndBonds2().
|
inline |
Mark that atom types have been perceived (see OBAtomTyper for details)
Referenced by OBAtomTyper::AssignTypes().
|
inline |
Mark that ring types have been perceived (see OBRingTyper for details)
Referenced by OBRingTyper::AssignTypes().
|
inline |
Mark that chains and residues have been perceived (see OBChainsParser)
Referenced by OBChainsParser::PerceiveChains().
|
inline |
Mark that chirality has been perceived.
Referenced by OBBuilder::Build().
|
inline |
Mark that partial charges have been assigned.
Referenced by OBPhModel::AssignSeedPartialCharge().
|
inline |
Mark that hybridization of all atoms has been assigned.
Referenced by OBAtomTyper::AssignHyb(), OBBuilder::Build(), and OBMol::PerceiveBondOrders().
|
inline |
Mark that ring closure bonds have been assigned by graph traversal.
Referenced by OpenBabel::FindRingAtomsAndBonds2().
|
inline |
Mark that explicit hydrogen atoms have been added.
|
inline |
Referenced by OBPhModel::CorrectForPH().
|
inline |
|
inline |
The OBMol is a pattern, not a complete molecule. Left unchanged by Clear().
Referenced by OpenBabel::alternate(), and AliasData::Expand().
|
inline |
|
inline |
Referenced by OBMol::operator=().
|
inline |
Referenced by OBMol::CopySubstructure().
|
inline |
Referenced by OBMol::RenumberAtoms().
|
inline |
|
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().
|
static |
|
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().
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().
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.
void SetCoordinates | ( | double * | c | ) |
Set the coordinates for all atoms in this conformer.
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 OBPointGroup::Setup().
bool DeleteHydrogens | ( | ) |
Delete all hydrogens from the molecule
Referenced by OBPhModel::CorrectForPH(), and AliasData::Expand().
bool DeleteHydrogens | ( | OBAtom * | atom | ) |
Delete all hydrogens from the supplied atom
bool DeletePolarHydrogens | ( | ) |
Delete all hydrogen atoms connected to a polar atom
bool DeleteNonPolarHydrogens | ( | ) |
Delete all hydrogen atoms connected to a non-polar atom
bool DeleteHydrogen | ( | OBAtom * | atom | ) |
Delete the supplied atom if it is a hydrogen (Helper function for 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
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 |
Referenced by 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 AddNonPolarHydrogens | ( | ) |
Add only nonpolar hydrogens (i.e., attached to C)
bool AddNewHydrogens | ( | HydrogenType | whichHydrogen, |
bool | correctForPH = false , |
||
double | pH = 7.4 |
||
) |
Add polar and/or nonpolar hydrogens
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.
vector< OBMol > Separate | ( | int | StartIndex = 1 | ) |
Copies each disconnected fragment as a separate OBMol.
Referenced by OBBuilder::Build(), and 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 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".
When used from Python, note that "None" may be used to specify an empty value for the excludebonds
parameter.
atoms
is a NULL pointer.newmol | The molecule to which to add the substructure. Note that atoms are appended to this molecule. |
atoms | An OBBitVec, indexed by atom Idx, specifying which atoms to copy |
excludebonds | An 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. |
correctvalence | A 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 |
atomorder | Record 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. |
bondorder | Record the Idxs of the original bonds. See atomorder above. |
Referenced by OBBuilder::Build(), and OBMol::GetNextFragment().
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.
bool MakeDativeBonds | ( | ) |
Converts 5-valent N and P only. Return true if conversion occurred.
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.
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.
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
bool CorrectForPH | ( | double | pH = 7.4 | ) |
Correct for pH by applying the OBPhModel transformations.
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().
bool AssignTotalChargeToAtoms | ( | int | charge | ) |
Put the specified molecular charge on appropriate atoms. Assumes all the hydrogen is explicitly included in the molecule.
vector3 Center | ( | int | nconf | ) |
nconf
Set the torsion defined by these atoms, rotating bonded neighbors
Referenced by OBBuilder::Connect(), and OBBuilder::CorrectStereoBonds().
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().
void FindLSSR | ( | ) |
Find Largest Set of Smallest Rings.
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().
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().
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)
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 OBMol::ConvertZeroBonds().
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 (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().
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.
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.
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().
bool Has3D | ( | ) |
Are there non-zero coordinates in all three dimensions (i.e. X, Y, Z)?
Referenced by AliasData::Expand().
bool HasNonZeroCoords | ( | ) |
Are there any non-zero coordinates?
|
inline |
Has aromatic perception been performed?
Referenced by OBAromaticTyper::AssignAromaticFlags(), OBBond::IsAromatic(), and OBAtom::IsAromatic().
|
inline |
Has the smallest set of smallest rings (FindSSSR) been performed?
Referenced by OBAtom::IsInRingSize(), OBAtom::MemberOfRingCount(), OBAtom::MemberOfRingSize(), OBMolRingIter::OBMolRingIter(), and OBMol::PerceiveBondOrders().
|
inline |
Has the largest set of smallest rings (FindLSSR) been performed?
|
inline |
Have ring atoms and bonds been assigned?
Referenced by OBBond::IsInRing(), and OBAtom::IsInRing().
|
inline |
Have atom types been assigned by OBAtomTyper?
Referenced by OBAtom::GetType().
|
inline |
Have ring types been assigned by OBRingTyper?
Referenced by OBRing::GetType().
|
inline |
Has atom chirality been assigned?
Referenced by OpenBabel::CanonicalLabels(), and OBMol::operator=().
|
inline |
Have atomic Gasteiger partial charges been assigned by OBGastChrg?
|
inline |
Has atomic hybridization been assigned by OBAtomTyper?
Referenced by OBAtom::GetHyb().
|
inline |
Have ring "closure" bonds been assigned? (e.g., OBBond::IsClosure())
Referenced by OpenBabel::DetermineFRJ(), and OBBond::IsClosure().
|
inline |
Have biomolecule chains and residues been assigned by OBChainsParser?
Referenced by OBMol::CopySubstructure(), and OBAtom::GetResidue().
|
inline |
Have hydrogens been added to the molecule?
|
inline |
Has the molecule been corrected for pH by CorrectForPH?
Referenced by OBPhModel::CorrectForPH().
|
inline |
Has total spin multiplicity been assigned?
|
inline |
Does this OBMol represent a reaction?
Referenced by OBMoleculeFormat::ReadChemObjectImpl().
|
inline |
Are there any atoms in this molecule?
Referenced by OBMol::ConnectTheDots(), OBSSMatch::OBSSMatch(), OBMol::PerceiveBondOrders(), and OBMol::RenumberAtoms().
|
inline |
Referenced by OBMol::Center(), OBMoleculeFormat::DoOutputOptions(), OBMol::operator=(), OBMol::RenumberAtoms(), OBMol::Rotate(), OBForceField::SetConformers(), and OBMol::Translate().
void SetConformers | ( | std::vector< double *> & | v | ) |
Set the entire set of conformers for this molecule to v
.
Referenced by OBConformerSearch::GetConformers(), and OBForceField::GetConformers().
|
inline |
Add a new set of coordinates f
as a new conformer.
Referenced by OpenBabel::UpdateConformersFromTree().
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().
void CopyConformer | ( | double * | c, |
int | nconf | ||
) |
Copy the conformer nconf
into the array c
c
is large enough void DeleteConformer | ( | int | nconf | ) |
Delete the conformer nconf
.
|
inline |
i
Referenced by OBMol::operator=(), OBMol::RenumberAtoms(), OBMol::Rotate(), OBForceField::SetConformers(), and OBMol::Translate().
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
|
inline |
Set the iterator to the beginning of the conformer list
|
inline |
Advance the iterator to the next confomer, if possible
|
inline |
Referenced by OBRotamerList::SetBaseCoordinateSets().
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
OBAtom * BeginAtom | ( | OBAtomIterator & | i | ) |
Set the iterator i
to the beginning of the atom list
Referenced by OBResidueData::AssignBonds(), OBAtomTyper::AssignHyb(), OBGastChrg::AssignPartialCharges(), OBAtomTyper::AssignTypes(), OpenBabel::CartesianToInternal(), OBMol::Center(), OBMol::ConnectTheDots(), OBMol::ConvertDativeBonds(), OBSmartsMatcher::FastSingleMatch(), OBChargeModel::FillChargeVectors(), OpenBabel::GetDFFVector(), OBAtom::GetPartialCharge(), OBMol::Has2D(), OBMol::Has3D(), OBMol::HasNonZeroCoords(), OBGrid::Init(), OpenBabel::InternalToCartesian(), OBMol::MakeDativeBonds(), OBMolAtomIter::OBMolAtomIter(), OBMolPairIter::OBMolPairIter(), OBMolPairIter::operator++(), OBMol::operator+=(), OBMol::operator=(), OBMol::PerceiveBondOrders(), OBChainsParser::PerceiveChains(), OBMol::RenumberAtoms(), OBMol::SetCoordinates(), OBProxGrid::Setup(), OBSmartsMatcher::SetupAtomMatchTable(), and OBMol::~OBMol().
OBAtom * NextAtom | ( | OBAtomIterator & | i | ) |
Advance the iterator i
to the next atom in the molecule
Referenced by OBResidueData::AssignBonds(), OBAtomTyper::AssignHyb(), OBGastChrg::AssignPartialCharges(), OBAtomTyper::AssignTypes(), OpenBabel::CartesianToInternal(), OBMol::Center(), OBMol::ConnectTheDots(), OBMol::ConvertDativeBonds(), OBSmartsMatcher::FastSingleMatch(), OBChargeModel::FillChargeVectors(), OpenBabel::GetDFFVector(), OBAtom::GetPartialCharge(), OBMol::Has2D(), OBMol::Has3D(), OBMol::HasNonZeroCoords(), OBGrid::Init(), OpenBabel::InternalToCartesian(), OBMol::MakeDativeBonds(), OBMolPairIter::OBMolPairIter(), OBMolPairIter::operator++(), OBMol::operator+=(), OBMol::operator=(), OBMol::PerceiveBondOrders(), OBChainsParser::PerceiveChains(), OBMol::RenumberAtoms(), OBMol::SetCoordinates(), OBProxGrid::Setup(), OBSmartsMatcher::SetupAtomMatchTable(), and OBMol::~OBMol().
OBBond * BeginBond | ( | OBBondIterator & | i | ) |
Set the iterator i
to the beginning of the bond list
Referenced by OBGastChrg::AssignPartialCharges(), OpenBabel::DetermineFRJ(), OBDepict::DrawMolecule(), OBRotorList::FindRotors(), OBAtomBondIter::OBAtomBondIter(), OBMolBondIter::OBMolBondIter(), OBMol::operator+=(), OBMol::operator=(), OBChainsParser::PerceiveChains(), and OBMol::~OBMol().
OBBond * NextBond | ( | OBBondIterator & | i | ) |
Advance the iterator i
to the next bond in the molecule
Referenced by OBGastChrg::AssignPartialCharges(), OpenBabel::DetermineFRJ(), OBDepict::DrawMolecule(), OBRotorList::FindRotors(), OBMolBondIter::operator++(), OBAtomBondIter::operator++(), OBMol::operator+=(), OBMol::operator=(), OBChainsParser::PerceiveChains(), and OBMol::~OBMol().
|
inline |
Set the iterator i
to the beginning of the resdiue list
Referenced by OBResidueIter::OBResidueIter(), OBMol::operator+=(), OBChainsParser::~OBChainsParser(), and OBMol::~OBMol().
|
inline |
Advance the iterator i
to the next residue in the molecule
Referenced by OBResidueIter::operator++(), OBMol::operator+=(), OBChainsParser::~OBChainsParser(), and OBMol::~OBMol().
|
inline |
Set the iterator to the beginning of the internal coordinate list
|
inline |
Advance the iterator to the next internal coordinate record
|
inlineinherited |
By default clears the object. Called from ReadMolecule of most format classes.
|
inherited |
Referenced by OBDepict::AddAtomLabels(), OBDepict::DrawMolecule(), AliasData::Expand(), OBForceField::GetAtomTypes(), OBForceField::GetConformers(), OBForceField::GetCoordinates(), OBMol::GetEnergies(), OBMol::GetEnergy(), OBForceField::GetPartialCharges(), OBDescriptor::MatchPairData(), OBMoleculeFormat::ReadChemObjectImpl(), and OBMol::SetEnergies().
|
inherited |
|
inherited |
|
inherited |
Delete any data matching the given OBGenericDataType.
Referenced by OpenBabel::CanonicalLabels(), OBDescriptor::DeleteProperties(), OpenBabel::DeleteStereoOnAtom(), and OBMol::RenumberAtoms().
|
inherited |
Delete the given generic data from this object.
|
inherited |
Delete all of the given generic data from this object.
|
inherited |
Deletes the generic data with the specified attribute, returning false if not found.
|
inlineinherited |
Adds a data object; does nothing if d==NULL.
Referenced by OBGastChrg::AssignPartialCharges(), OBMol::CopySubstructure(), AliasData::Expand(), OBForceField::GetAtomTypes(), OBForceField::GetConformers(), OBForceField::GetCoordinates(), OBMol::GetEnergies(), OBMol::GetEnergy(), OBForceField::GetPartialCharges(), OBMoleculeFormat::MakeCombinedMolecule(), OBDescriptor::PredictAndSave(), AliasData::RevertToAliasForm(), and OBMol::SetEnergies().
|
inherited |
Adds a copy of a data object; does nothing if d == NULL
Referenced by AliasData::Expand().
|
inlineinherited |
|
inherited |
Referenced by OBDepict::AddAtomLabels(), OBDepict::DrawMolecule(), OpenBabel::extract_thermochemistry(), OBDescriptor::FilterCompare(), OpenBabel::GetAtomSymClass(), OBForceField::GetAtomTypes(), OBForceField::GetConformers(), OBForceField::GetCoordinates(), OBForceField::GetPartialCharges(), OBDescriptor::GetValues(), OpenBabel::IsSuppressibleHydrogen(), OBMoleculeFormat::MakeCombinedMolecule(), OBMolAngleIter::OBMolAngleIter(), OBMolRingIter::OBMolRingIter(), OBMolTorsionIter::OBMolTorsionIter(), and OBDescriptor::PredictAndSave().
|
inherited |
|
inherited |
|
inlineinherited |
Referenced by OBMol::GetEnergies(), OBMol::GetEnergy(), and OBMol::SetEnergies().
|
inherited |
|
inherited |
Referenced by OpenBabel::CanonicalLabels(), OBMol::CopySubstructure(), OBBuilder::CorrectStereoAtoms(), OBBuilder::CorrectStereoBonds(), OpenBabel::DeleteStereoOnAtom(), and OBMol::operator+=().
|
inlineinherited |
Referenced by OBMol::AddBond(), OBAtom::Duplicate(), OBMoleculeFormat::MakeCombinedMolecule(), and OBMol::operator=().
|
inlineinherited |
Referenced by OBMol::AddBond(), OBAtom::Duplicate(), OBMoleculeFormat::MakeCombinedMolecule(), and OBMol::operator=().
|
protected |
bitfield of flags
Referenced by OBMol::CopySubstructure(), OBMol::OBMol(), and OBMol::PerceiveBondOrders().
|
protected |
Assign partial charges automatically.
Referenced by OBMol::OBMol().
|
protected |
Assign formal charges automatically.
Referenced by OBMol::OBMol().
|
protected |
Molecule title.
Referenced by OBMol::OBMol().
|
protected |
vector of atoms
Referenced by OBMol::BeginAtom(), OBMol::NextAtom(), OBMol::OBMol(), and OBMol::RenumberAtoms().
|
protected |
vector of atoms indexed by id
Referenced by OBMol::OBMol().
|
protected |
vector of bonds
Referenced by OBMol::BeginBond(), OBMol::NextBond(), and OBMol::OBMol().
|
protected |
vector of bonds
Referenced by OBMol::OBMol().
|
protected |
Dimensionality of coordinates.
Referenced by OBMol::ConnectTheDots(), OBMol::OBMol(), and OBMol::PerceiveBondOrders().
|
protected |
Total charge on the molecule.
Referenced by OBMol::OBMol(), and OBMol::PerceiveBondOrders().
|
protected |
Total spin on the molecule (if not specified, assumes lowest possible spin)
|
protected |
coordinate array
Referenced by OBMol::ConnectTheDots(), OBMol::OBMol(), OBMol::operator=(), OBMol::Rotate(), OBMol::SetConformer(), OBMol::SetConformers(), OBMol::SetCoordinates(), and OBMol::Translate().
|
protected |
vector of conformers
Referenced by OBMol::ConnectTheDots(), OBMol::CopyConformer(), OBMol::DeleteConformer(), OBMol::OBMol(), OBMol::SetConformer(), OBMol::SetConformers(), OBMol::SetCoordinates(), and OBMol::~OBMol().
|
protected |
heat of formation
Referenced by OBMol::OBMol().
|
protected |
Number of atoms.
Referenced by OBMol::ConnectTheDots(), and OBMol::OBMol().
|
protected |
Number of bonds.
Referenced by OBMol::OBMol().
|
protected |
Residue information (if applicable)
|
protected |
Internal Coordinates (if applicable)
|
protected |
Number of nested calls to BeginModify()
Referenced by OBMol::OBMol().
|
protectedinherited |
Custom data.
Referenced by OBMol::OBMol().