Open Babel  3.0
mol.h
Go to the documentation of this file.
1 /**********************************************************************
2 mol.h - Handle molecules. Declarations of OBMol, OBAtom, OBBond, OBResidue.
3  (the main header for Open Babel)
4 
5 Copyright (C) 1998-2001 by OpenEye Scientific Software, Inc.
6 Some portions Copyright (C) 2001-2006 by Geoffrey R. Hutchison
7 Some portions Copyright (C) 2003 by Michael Banck
8 
9 This file is part of the Open Babel project.
10 For more information, see <http://openbabel.org/>
11 
12 This program is free software; you can redistribute it and/or modify
13 it under the terms of the GNU General Public License as published by
14 the Free Software Foundation version 2 of the License.
15 
16 This program is distributed in the hope that it will be useful,
17 but WITHOUT ANY WARRANTY; without even the implied warranty of
18 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
19 GNU General Public License for more details.
20 ***********************************************************************/
21 
22 #ifndef OB_MOL_H
23 #define OB_MOL_H
24 
25 #include <openbabel/babelconfig.h>
26 
27 #ifndef EXTERN
28 # define EXTERN extern
29 #endif
30 #ifndef THREAD_LOCAL
31 #ifdef SWIG
32 # define THREAD_LOCAL
33 # elif (__cplusplus >= 201103L)
34 //this is required for correct multi-threading
35 # define THREAD_LOCAL thread_local
36 # else
37 # define THREAD_LOCAL
38 # endif
39 #endif
40 
41 #include <math.h>
42 #include <float.h>
43 
44 #include <vector>
45 #include <string>
46 #include <map>
47 
48 #include <openbabel/base.h>
49 
50 
51 namespace OpenBabel
52 {
53  class OBAtom;
54  class OBBond;
55  class OBResidue;
56  class OBRing;
57  class OBInternalCoord;
58  class OBConversion; //used only as a pointer
59 
60  class vector3;
61  class OBBitVec;
62  class OBMolAtomDFSIter;
63  class OBChainsParser;
64 
65  typedef std::vector<OBAtom*>::iterator OBAtomIterator;
66  typedef std::vector<OBBond*>::iterator OBBondIterator;
67  typedef std::vector<OBResidue*>::iterator OBResidueIterator;
68 
69  // Class OBMol
70  //MOL Property Macros (flags) -- 32+ bits
72 #define OB_SSSR_MOL (1<<1)
73 #define OB_RINGFLAGS_MOL (1<<2)
75 #define OB_AROMATIC_MOL (1<<3)
77 #define OB_ATOMTYPES_MOL (1<<4)
79 #define OB_CHIRALITY_MOL (1<<5)
81 #define OB_PCHARGE_MOL (1<<6)
83 #define OB_HYBRID_MOL (1<<8)
85 #define OB_CLOSURE_MOL (1<<11)
87 #define OB_H_ADDED_MOL (1<<12)
89 #define OB_PH_CORRECTED_MOL (1<<13)
91 #define OB_CHAINS_MOL (1<<15)
93 #define OB_TCHARGE_MOL (1<<16)
95 #define OB_TSPIN_MOL (1<<17)
97 #define OB_RINGTYPES_MOL (1<<18)
99 #define OB_PATTERN_STRUCTURE (1<<19)
101 #define OB_LSSR_MOL (1<<20)
103 #define OB_ATOMSPIN_MOL (1<<21)
105 #define OB_REACTION_MOL (1<<22)
107  // flags 22-32 unspecified
108 
109 #define SET_OR_UNSET_FLAG(X) \
110  if (value) SetFlag(X); \
111  else UnsetFlag(X);
112 
113 #define OB_CURRENT_CONFORMER -1
114 
116 
117  // class introduction in mol.cpp
118  class OBAPI OBMol: public OBBase
119  {
120  protected:
121  int _flags;
124  std::string _title;
125  std::vector<OBAtom*> _vatom;
126  std::vector<OBAtom*> _atomIds;
127  std::vector<OBBond*> _vbond;
128  std::vector<OBBond*> _bondIds;
129  unsigned short int _dimension;
131  unsigned int _totalSpin;
132  double *_c;
133  std::vector<double*> _vconf;
134  double _energy;
135  unsigned int _natoms;
136  unsigned int _nbonds;
137  std::vector<OBResidue*> _residue;
138  std::vector<OBInternalCoord*> _internals;
139  unsigned short int _mod;
140 
141  public:
142 
144 
145  OBMol();
148  OBMol(const OBMol &);
150  virtual ~OBMol();
152  OBMol &operator=(const OBMol &mol);
154  OBMol &operator+=(const OBMol &mol);
155 
158  void ReserveAtoms(int natoms)
159  {
160  if (natoms > 0 && _mod) {
161  _vatom.reserve(natoms);
162  _atomIds.reserve(natoms);
163  }
164  }
165 
168  virtual void DestroyAtom(OBAtom*);
171  virtual void DestroyBond(OBBond*);
174  virtual void DestroyResidue(OBResidue*);
175 
180  bool AddAtom(OBAtom& atom, bool forceNewId = false);
183  bool InsertAtom(OBAtom &);
191  bool AddBond(int beginIdx, int endIdx, int order,
192  int flags=0,int insertpos=-1);
195  bool AddBond(OBBond&);
198  bool AddResidue(OBResidue&);
199 
203  OBAtom *NewAtom();
207  OBAtom *NewAtom(unsigned long id);
211  OBBond *NewBond();
215  OBBond *NewBond(unsigned long id);
217  OBResidue *NewResidue();
222  bool DeleteAtom(OBAtom*, bool destroyAtom = true);
225  bool DeleteBond(OBBond*, bool destroyBond = true);
228  bool DeleteResidue(OBResidue*, bool destroyResidue = true);
230 
232 
233  virtual void BeginModify(void);
240  virtual void EndModify(bool nukePerceivedData=true);
242  int GetMod() { return(_mod); }
245  void IncrementMod() { _mod++; }
248  void DecrementMod() { _mod--; }
250 
252 
253  int GetFlags() { return(_flags); }
257  const char *GetTitle(bool replaceNewlines = true) const;
259  unsigned int NumAtoms() const { return(_natoms); }
261  unsigned int NumBonds() const { return(_nbonds); }
263  unsigned int NumHvyAtoms();
265  unsigned int NumResidues() const { return(static_cast<unsigned int> (_residue.size())); }
267  unsigned int NumRotors(bool sampleRingBonds=false);
268 
271  OBAtom *GetAtom(int idx) const;
273  OBAtom *GetAtomById(unsigned long id) const;
276  OBAtom *GetFirstAtom() const;
279  OBBond *GetBond(int idx) const;
281  OBBond *GetBondById(unsigned long id) const;
284  OBBond *GetBond(int a, int b) const;
285  // The safer version of the above method
287  OBBond *GetBond(OBAtom* bgn, OBAtom* end) const;
290  OBResidue *GetResidue(int idx) const;
291  std::vector<OBInternalCoord*> GetInternalCoord();
296  double GetTorsion(int,int,int,int);
301  double GetTorsion(OBAtom* a,OBAtom* b,OBAtom* c,OBAtom* d);
304  double GetAngle(OBAtom* a, OBAtom* b, OBAtom* c);
307  int AreInSameRing(OBAtom *a, OBAtom *b);
309  std::string GetFormula();
311  std::string GetSpacedFormula(int ones=0, const char* sp=" ", bool implicitH = true);
313  double GetEnergy() const { return _energy; }
315  double GetMolWt(bool implicitH = true);
317  double GetExactMass(bool implicitH = true);
319  int GetTotalCharge();
321  unsigned int GetTotalSpinMultiplicity();
323  unsigned short int GetDimension() const { return _dimension; }
325  double *GetCoordinates() { return(_c); }
327  std::vector<OBRing*> &GetSSSR();
329  std::vector<OBRing*> &GetLSSR();
331  bool AutomaticFormalCharge() { return(_autoFormalCharge); }
333  bool AutomaticPartialCharge() { return(_autoPartialCharge); }
335 
336 
338 
339  void SetTitle(const char *title);
342  void SetTitle(std::string &title);
344  void SetFormula(std::string molFormula);
346  void SetEnergy(double energy) { _energy = energy; }
348  void SetDimension(unsigned short int d) { _dimension = d; }
350  void SetTotalCharge(int charge);
353  void SetTotalSpinMultiplicity(unsigned int spinMultiplicity);
359  void SetInternalCoord(std::vector<OBInternalCoord*> int_coord);
362  { _autoFormalCharge=val; }
365  { _autoPartialCharge=val; }
366 
370  void SetSSSRPerceived(bool value = true) { SET_OR_UNSET_FLAG(OB_SSSR_MOL); }
372  void SetLSSRPerceived(bool value = true) { SET_OR_UNSET_FLAG(OB_LSSR_MOL); }
380  void SetChainsPerceived(bool value = true) { SET_OR_UNSET_FLAG(OB_CHAINS_MOL); }
390  void SetHydrogensAdded(bool value = true) { SET_OR_UNSET_FLAG(OB_H_ADDED_MOL); }
395  void SetIsReaction(bool value = true) { SET_OR_UNSET_FLAG(OB_REACTION_MOL) };
396  bool HasFlag(int flag) { return (_flags & flag) ? true : false; }
397  void SetFlag(int flag) { _flags |= flag; }
398  void UnsetFlag(int flag) { _flags &= (~(flag)); }
399  void SetFlags(int flags) { _flags = flags; }
401 
403 
404  // Description in transform.cpp (command-line transformations to this molecule)
405  virtual OBBase* DoTransformations(const std::map<std::string,std::string>* pOptions,OBConversion* pConv);
406  // Ditto (documentation on transformation options)
407  static const char* ClassDescription();
409  bool Clear();
411  void RenumberAtoms(std::vector<OBAtom*>&);
413  void RenumberAtoms(std::vector<int>);
416  void SetCoordinates(double *c);
418  void ToInertialFrame(int conf, double *rmat);
420  void ToInertialFrame();
422  void Translate(const vector3 &v);
424  void Translate(const vector3 &v, int conf);
426  void Rotate(const double u[3][3]);
428  void Rotate(const double m[9]);
430  void Rotate(const double m[9],int nconf);
432  void Center();
435  bool DeleteHydrogens();
438  bool DeleteHydrogens(OBAtom*);
442  bool DeletePolarHydrogens();
445  bool DeleteNonPolarHydrogens();
448  bool DeleteHydrogen(OBAtom*);
455  bool AddHydrogens(bool polaronly=false,bool correctForPH=false, double pH=7.4);
457  bool AddHydrogens(OBAtom*);
459  bool AddPolarHydrogens();
462  bool AddNonPolarHydrogens();
465  bool AddNewHydrogens(HydrogenType whichHydrogen, bool correctForPH=false, double pH=7.4);
466 
470  bool StripSalts(unsigned int threshold=0);
472  std::vector<OBMol> Separate(int StartIndex=1);
474  bool GetNextFragment( OpenBabel::OBMolAtomDFSIter& iter, OBMol& newMol );
475  // docs in mol.cpp
476  bool CopySubstructure(OBMol& newmol, OBBitVec *includeatoms, OBBitVec *excludebonds = (OBBitVec*)0,
477  unsigned int correctvalence=1,
478  std::vector<unsigned int> *atomorder=(std::vector<unsigned int>*)0,
479  std::vector<unsigned int> *bondorder=(std::vector<unsigned int>*)0);
481  bool ConvertDativeBonds();
485  bool MakeDativeBonds();
491  bool ConvertZeroBonds();
492 
494  bool CorrectForPH(double pH=7.4);
495  // docs in mol.cpp
496  bool AssignSpinMultiplicity(bool NoImplicitH=false);
497 
501  bool AssignTotalChargeToAtoms(int charge);
502 
505  vector3 Center(int nconf);
511  void SetTorsion(OBAtom*,OBAtom*,OBAtom*,OBAtom*,double ang);
513 
515 
516  void FindSSSR();
519  void FindLSSR();
521  void FindRingAtomsAndBonds();
522  // documented in mol.cpp -- locates all atom indexes which can reach 'end'
523  void FindChildren(std::vector<int> & children,int bgnIdx,int endIdx);
524  // documented in mol.cpp -- locates all atoms which can reach 'end'
525  void FindChildren(std::vector<OBAtom*>& children,OBAtom* bgn,OBAtom* end);
530  void FindLargestFragment(OBBitVec &frag);
533  void ContigFragList(std::vector<std::vector<int> >&);
535  void Align(OBAtom*,OBAtom*,vector3&,vector3&);
537  void ConnectTheDots();
539  void PerceiveBondOrders();
541  void FindAngles();
543  void FindTorsions();
544  // documented in mol.cpp: graph-theoretical distance for each atom
545  bool GetGTDVector(std::vector<int> &);
546  // documented in mol.cpp: graph-invariant index for each atom
547  void GetGIVector(std::vector<unsigned int> &);
548  // documented in mol.cpp: calculate symmetry-unique identifiers
549  void GetGIDVector(std::vector<unsigned int> &);
551 
553 
554  bool Has2D(bool Not3D=false);
557  bool Has3D();
559  bool HasNonZeroCoords();
561  bool HasAromaticPerceived() { return(HasFlag(OB_AROMATIC_MOL)); }
563  bool HasSSSRPerceived() { return(HasFlag(OB_SSSR_MOL)); }
565  bool HasLSSRPerceived() { return(HasFlag(OB_LSSR_MOL)); }
569  bool HasAtomTypesPerceived() { return(HasFlag(OB_ATOMTYPES_MOL));}
571  bool HasRingTypesPerceived() { return(HasFlag(OB_RINGTYPES_MOL));}
573  bool HasChiralityPerceived() { return(HasFlag(OB_CHIRALITY_MOL));}
575  bool HasPartialChargesPerceived() { return(HasFlag(OB_PCHARGE_MOL));}
577  bool HasHybridizationPerceived() { return(HasFlag(OB_HYBRID_MOL)); }
579  bool HasClosureBondsPerceived() { return(HasFlag(OB_CLOSURE_MOL)); }
581  bool HasChainsPerceived() { return(HasFlag(OB_CHAINS_MOL)); }
583  bool HasHydrogensAdded() { return(HasFlag(OB_H_ADDED_MOL)); }
585  bool IsCorrectedForPH() { return(HasFlag(OB_PH_CORRECTED_MOL)); }
587  bool HasSpinMultiplicityAssigned() { return(HasFlag(OB_ATOMSPIN_MOL)); }
589  bool IsReaction() { return HasFlag(OB_REACTION_MOL); }
591  bool Empty() { return(_natoms == 0); }
593 
595 
596  int NumConformers() { return((_vconf.empty())?0:static_cast<int> (_vconf.size())); }
599  void SetConformers(std::vector<double*> &v);
601  void AddConformer(double *f) { _vconf.push_back(f); }
604  void SetConformer(unsigned int i);
607  void CopyConformer(double* c,int nconf);
609  void DeleteConformer(int nconf);
611  double *GetConformer(int i) { return(_vconf[i]); }
613  void SetEnergies(std::vector<double> &energies);
615  std::vector<double> GetEnergies();
618  double GetEnergy(int ci);
621  double *BeginConformer(std::vector<double*>::iterator&i)
622  { i = _vconf.begin();
623  return((i == _vconf.end()) ? NULL:*i); }
626  double *NextConformer(std::vector<double*>::iterator&i)
627  { ++i;
628  return((i == _vconf.end()) ? NULL:*i); }
630  std::vector<double*> &GetConformers() { return(_vconf); }
632 
634 
635  OBAtomIterator BeginAtoms() { return _vatom.begin(); }
638  OBAtomIterator EndAtoms() { return _vatom.begin() + NumAtoms() ; }
640  OBBondIterator BeginBonds() { return _vbond.begin(); }
642  OBBondIterator EndBonds() { return _vbond.begin() + NumBonds() ; }
644  OBResidueIterator BeginResidues() { return _residue.begin(); }
646  OBResidueIterator EndResidues() { return _residue.end(); }
647 
650  OBAtom *BeginAtom(OBAtomIterator &i);
653  OBAtom *NextAtom(OBAtomIterator &i);
656  OBBond *BeginBond(OBBondIterator &i);
659  OBBond *NextBond(OBBondIterator &i);
662  OBResidue *BeginResidue(OBResidueIterator &i)
663  {
664  i = _residue.begin();
665  return((i == _residue.end()) ? NULL:*i);
666  }
669  OBResidue *NextResidue(OBResidueIterator &i)
670  {
671  ++i;
672  return((i == _residue.end()) ? NULL:*i);
673  }
677  OBInternalCoord *BeginInternalCoord(std::vector<OBInternalCoord*>::iterator &i)
678  {
679  i = _internals.begin();
680  return((i == _internals.end()) ? NULL:*i);
681  }
685  OBInternalCoord *NextInternalCoord(std::vector<OBInternalCoord*>::iterator &i)
686  {
687  ++i;
688  return((i == _internals.end()) ? NULL:*i);
689  }
691 
692  };
693 
694  // Utility function prototypes
695  //tokenize and Trim declarations moved to base.h
696  // Deprecated -- use OBMessageHandler class instead (docs in obutil.cpp)
697  OBAPI void ThrowError(char *str);
698  // Deprecated -- use OBMessageHandler class instead (docs in obutil.cpp)
699  OBAPI void ThrowError(std::string &str);
701  OBAPI void CartesianToInternal(std::vector<OBInternalCoord*>&,OBMol&);
703  OBAPI void InternalToCartesian(std::vector<OBInternalCoord*>&,OBMol&);
704  // Replace the last extension in str with a new one (docs in obutil.cpp)
705  OBAPI std::string NewExtension(std::string&,char*);
706 
708  namespace detail {
711  template<typename T, int size = sizeof(T)>
712  struct max_value
713  {
714  static const T result = (static_cast<T>(0xFF) << (size-1)*8) + max_value<T, size-1>::result;
715  };
716 
718  template<typename T>
719  struct max_value<T, 0>
720  {
721  static const T result = 0;
722  };
723  }
724 
725  // No unique id
726  static const unsigned long NoId = detail::max_value<unsigned long>::result;
727 
728  //Utility Macros
729 
730 #ifndef BUFF_SIZE
731 #define BUFF_SIZE 32768
732 #endif
733 
734 #ifndef EQ
735 #define EQ(a,b) (!strcmp((a), (b)))
736 #endif
737 
738 #ifndef EQn
739 #define EQn(a,b,n) (!strncmp((a), (b), (n)))
740 #endif
741 
742 #ifndef SQUARE
743 #define SQUARE(x) ((x)*(x))
744 #endif
745 
746 #ifndef IsUnsatType
747 #define IsUnsatType(x) (EQ(x,"Car") || EQ(x,"C2") || EQ(x,"Sox") || EQ(x,"Sac") || EQ(x,"Pac") || EQ(x,"So2"))
748 #endif
749 
750 #ifndef __KCC
751  extern "C"
752  {
753  OBAPI void get_rmat(double*,double*,double*,int);
754  OBAPI void ob_make_rmat(double mat[3][3],double rmat[9]);
755  OBAPI void qtrfit (double *r,double *f,int size,double u[3][3]);
756  OBAPI double superimpose(double*,double*,int);
757  }
758 #else
759  OBAPI void get_rmat(double*,double*,double*,int);
760  OBAPI void ob_make_rmat(double mat[3][3],double rmat[9]);
761  OBAPI void qtrfit (double *r,double *f,int size,double u[3][3]);
762  OBAPI double superimpose(double*,double*,int);
763 #endif // __KCC
764 
765 // extern OBMol* (*CreateMolecule) (void);
766 
767 } // end namespace OpenBabel
768 
769 #endif // OB_MOL_H
770 
Definition: residue.h:336
double * BeginConformer(std::vector< double *>::iterator &i)
Definition: mol.h:621
bool HasClosureBondsPerceived()
Have ring "closure" bonds been assigned? (e.g., OBBond::IsClosure())
Definition: mol.h:579
OBInternalCoord * NextInternalCoord(std::vector< OBInternalCoord *>::iterator &i)
Definition: mol.h:685
void get_rmat(double *, double *, double *, int)
Definition: obutil.cpp:1197
bool HasChiralityPerceived()
Has atom chirality been assigned?
Definition: mol.h:573
#define OB_CHAINS_MOL
Biomolecular chains and residues have been set. See OBChainsParser.
Definition: mol.h:92
#define OB_LSSR_MOL
Largest Set of Smallest Rings (LSSR) done. See OBRing and OBMol::FindLSSR.
Definition: mol.h:102
OBAtomIterator EndAtoms()
Definition: mol.h:638
double superimpose(double *, double *, int)
Definition: obutil.cpp:1103
unsigned int NumAtoms() const
Definition: mol.h:259
void SetPartialChargesPerceived(bool value=true)
Mark that partial charges have been assigned.
Definition: mol.h:384
void SetSpinMultiplicityAssigned(bool value=true)
Definition: mol.h:392
Base classes to build a graph.
#define OB_HYBRID_MOL
Atom hybridizations have been set. See OBAtomTyper.
Definition: mol.h:84
void IncrementMod()
Definition: mol.h:245
Class to convert from one format to another.
Definition: obconversion.h:59
std::vector< OBBond * > _vbond
vector of bonds
Definition: mol.h:127
int _flags
bitfield of flags
Definition: mol.h:121
#define OB_H_ADDED_MOL
Hyrdogen atoms have been added where needed. See OBMol::AddHydrogens.
Definition: mol.h:88
bool HasChainsPerceived()
Have biomolecule chains and residues been assigned by OBChainsParser?
Definition: mol.h:581
unsigned int _totalSpin
Total spin on the molecule (if not specified, assumes lowest possible spin)
Definition: mol.h:131
OBInternalCoord * BeginInternalCoord(std::vector< OBInternalCoord *>::iterator &i)
Definition: mol.h:677
#define OB_PH_CORRECTED_MOL
pH correction for hydrogen addition has been performed.
Definition: mol.h:90
void SetIsReaction(bool value=true)
Definition: mol.h:395
std::string _title
Molecule title.
Definition: mol.h:124
#define OB_RINGTYPES_MOL
Ring typing has been performed. See OBRingTyper.
Definition: mol.h:98
Bond class.
Definition: bond.h:58
void SetRingAtomsAndBondsPerceived(bool value=true)
Mark that rings have been perceived (see OBRing class for details)
Definition: mol.h:374
std::vector< OBResidue * > _residue
Residue information (if applicable)
Definition: mol.h:137
#define OB_REACTION_MOL
Treat as reaction.
Definition: mol.h:106
bool HasAtomTypesPerceived()
Have atom types been assigned by OBAtomTyper?
Definition: mol.h:569
double _energy
heat of formation
Definition: mol.h:134
Molecule Class.
Definition: mol.h:118
double * GetConformer(int i)
Definition: mol.h:611
bool Empty()
Are there any atoms in this molecule?
Definition: mol.h:591
#define SET_OR_UNSET_FLAG(X)
Definition: mol.h:109
void SetLSSRPerceived(bool value=true)
Mark that Largest Set of Smallest Rings has been run (see OBRing class)
Definition: mol.h:372
#define OB_CHIRALITY_MOL
Chirality detection has been performed.
Definition: mol.h:80
bool HasAromaticPerceived()
Has aromatic perception been performed?
Definition: mol.h:561
void SetChainsPerceived(bool value=true)
Mark that chains and residues have been perceived (see OBChainsParser)
Definition: mol.h:380
void SetDimension(unsigned short int d)
Set the dimension of this molecule (i.e., 0, 1 , 2, 3)
Definition: mol.h:348
std::vector< OBInternalCoord * > _internals
Internal Coordinates (if applicable)
Definition: mol.h:138
static const unsigned long NoId
Definition: mol.h:726
void ReserveAtoms(int natoms)
Definition: mol.h:158
OBBondIterator BeginBonds()
Definition: mol.h:640
void qtrfit(double *r, double *f, int size, double u[3][3])
Definition: obutil.cpp:686
void ThrowError(char *str)
Definition: obutil.cpp:55
A speed-optimized vector of bits.
Definition: bitvec.h:57
#define OB_PCHARGE_MOL
Partial charges have been set or percieved.
Definition: mol.h:82
void SetEnergy(double energy)
Set the heat of formation for this molecule (in kcal/mol)
Definition: mol.h:346
bool IsReaction()
Does this OBMol represent a reaction?
Definition: mol.h:589
void SetCorrectedForPH(bool value=true)
Definition: mol.h:391
#define OB_RINGFLAGS_MOL
Ring flags have been set: See OBRing::FindRingAtomsAndBonds.
Definition: mol.h:74
std::vector< double * > & GetConformers()
Definition: mol.h:630
std::vector< OBAtom * >::iterator OBAtomIterator
A standard iterator over a vector of atoms.
Definition: atom.h:48
std::vector< OBBond * > _bondIds
vector of bonds
Definition: mol.h:128
void SetAtomTypesPerceived(bool value=true)
Mark that atom types have been perceived (see OBAtomTyper for details)
Definition: mol.h:376
double GetEnergy() const
Definition: mol.h:313
void SetClosureBondsPerceived(bool value=true)
Mark that ring closure bonds have been assigned by graph traversal.
Definition: mol.h:388
std::vector< OBBond * >::iterator OBBondIterator
A standard iterator over a vector of bonds.
Definition: atom.h:46
std::vector< OBAtom * > _vatom
vector of atoms
Definition: mol.h:125
double * NextConformer(std::vector< double *>::iterator &i)
Definition: mol.h:626
void SetRingTypesPerceived(bool value=true)
Mark that ring types have been perceived (see OBRingTyper for details)
Definition: mol.h:378
void CartesianToInternal(std::vector< OBInternalCoord *> &, OBMol &)
Convert Cartesian XYZ to a set of OBInternalCoord coordinates.
Definition: obutil.cpp:547
void InternalToCartesian(std::vector< OBInternalCoord *> &, OBMol &)
Convert set of OBInternalCoord coordinates into Cartesian XYZ.
Definition: obutil.cpp:444
double * _c
coordinate array
Definition: mol.h:132
Used to transform from z-matrix to cartesian coordinates.
Definition: internalcoord.h:61
OBResidueIterator EndResidues()
Definition: mol.h:646
OBResidue * BeginResidue(OBResidueIterator &i)
Definition: mol.h:662
Represents a vector in 3-dimensional real space.
Definition: vector3.h:44
bool HasSSSRPerceived()
Has the smallest set of smallest rings (FindSSSR) been performed?
Definition: mol.h:563
a C++ template to return the maximum value of a type (e.g., int)
Definition: mol.h:712
bool AutomaticPartialCharge()
Get the current flag for whether partial charges are auto-determined.
Definition: mol.h:333
HydrogenType
Definition: mol.h:115
double GetExactMass(unsigned int atomic_number, unsigned int isotope=0)
Definition: elements.cpp:831
void SetHydrogensAdded(bool value=true)
Mark that explicit hydrogen atoms have been added.
Definition: mol.h:390
unsigned short int GetDimension() const
Definition: mol.h:323
void SetFlags(int flags)
Definition: mol.h:399
bool HasLSSRPerceived()
Has the largest set of smallest rings (FindLSSR) been performed?
Definition: mol.h:565
bool HasRingTypesPerceived()
Have ring types been assigned by OBRingTyper?
Definition: mol.h:571
unsigned short int _dimension
Dimensionality of coordinates.
Definition: mol.h:129
void SetIsPatternStructure(bool value=true)
The OBMol is a pattern, not a complete molecule. Left unchanged by Clear().
Definition: mol.h:394
bool HasPartialChargesPerceived()
Have atomic Gasteiger partial charges been assigned by OBGastChrg?
Definition: mol.h:575
#define OB_AROMATIC_MOL
Aromatic flags have been set for atoms and bonds.
Definition: mol.h:76
void AddConformer(double *f)
Add a new set of coordinates f as a new conformer.
Definition: mol.h:601
#define OB_SSSR_MOL
Smallest Set of Smallest Rings (SSSR) done. See OBRing and OBMol::FindSSSR.
Definition: mol.h:72
void SetHybridizationPerceived(bool value=true)
Mark that hybridization of all atoms has been assigned.
Definition: mol.h:386
bool _autoFormalCharge
Assign formal charges automatically.
Definition: mol.h:123
void SetFlag(int flag)
Definition: mol.h:397
void SetChiralityPerceived(bool value=true)
Mark that chirality has been perceived.
Definition: mol.h:382
unsigned int NumBonds() const
Definition: mol.h:261
std::string NewExtension(std::string &, char *)
Utility function: replace the last extension in string &src with new extension char *ext...
Definition: obutil.cpp:109
void SetAromaticPerceived(bool value=true)
Mark that aromaticity has been perceived for this molecule (see OBAromaticTyper)
Definition: mol.h:368
void SetAutomaticPartialCharge(bool val)
Set the flag for determining partial charges automatically (default=true)
Definition: mol.h:364
Definition: mol.h:115
bool HasHybridizationPerceived()
Has atomic hybridization been assigned by OBAtomTyper?
Definition: mol.h:577
std::vector< double * > _vconf
vector of conformers
Definition: mol.h:133
std::vector< OBAtom * > _atomIds
vector of atoms indexed by id
Definition: mol.h:126
void DecrementMod()
Definition: mol.h:248
unsigned short int _mod
Number of nested calls to BeginModify()
Definition: mol.h:139
unsigned int _natoms
Number of atoms.
Definition: mol.h:135
OBBondIterator EndBonds()
Definition: mol.h:642
unsigned int NumResidues() const
Definition: mol.h:265
Definition: mol.h:115
#define OB_ATOMSPIN_MOL
SpinMultiplicities on atoms have been set in OBMol::AssignSpinMultiplicity()
Definition: mol.h:104
bool IsCorrectedForPH()
Has the molecule been corrected for pH by CorrectForPH?
Definition: mol.h:585
double * GetCoordinates()
Definition: mol.h:325
bool HasRingAtomsAndBondsPerceived()
Have ring atoms and bonds been assigned?
Definition: mol.h:567
#define OB_PATTERN_STRUCTURE
A pattern, not a complete molecule.
Definition: mol.h:100
void SetSSSRPerceived(bool value=true)
Mark that Smallest Set of Smallest Rings has been run (see OBRing class)
Definition: mol.h:370
unsigned int _nbonds
Number of bonds.
Definition: mol.h:136
#define OB_CLOSURE_MOL
Ring "closure" bonds have been set. See OBBond::IsClosure.
Definition: mol.h:86
#define OB_ATOMTYPES_MOL
Atom typing has been performed. See OBAtomTyper.
Definition: mol.h:78
bool HasHydrogensAdded()
Have hydrogens been added to the molecule?
Definition: mol.h:583
OBResidueIterator BeginResidues()
Definition: mol.h:644
void SetAutomaticFormalCharge(bool val)
Set the flag for determining automatic formal charges with pH (default=true)
Definition: mol.h:361
int GetMod()
Definition: mol.h:242
void UnsetFlag(int flag)
Definition: mol.h:398
int _totalCharge
Total charge on the molecule.
Definition: mol.h:130
std::vector< OBResidue * >::iterator OBResidueIterator
Definition: mol.h:67
Base Class.
Definition: base.h:239
bool HasSpinMultiplicityAssigned()
Has total spin multiplicity been assigned?
Definition: mol.h:587
Residue information.
Definition: residue.h:50
bool AutomaticFormalCharge()
Get the current flag for whether formal charges are set with pH correction.
Definition: mol.h:331
Definition: mol.h:115
bool _autoPartialCharge
Assign partial charges automatically.
Definition: mol.h:122
OBResidue * NextResidue(OBResidueIterator &i)
Definition: mol.h:669
bool HasFlag(int flag)
Definition: mol.h:396
void ob_make_rmat(double mat[3][3], double rmat[9])
Definition: obutil.cpp:896
Iterate over all atoms in an OBMol in a depth-first search (DFS)
Definition: obiter.h:67
Global namespace for all Open Babel code.
Definition: alias.h:22
Atom class.
Definition: atom.h:71