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 OB_EXTERN
28 #  define OB_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
71   //! Smallest Set of Smallest Rings (SSSR) done. See OBRing and OBMol::FindSSSR
72 #define OB_SSSR_MOL              (1<<1)
73   //! Ring flags have been set: See OBRing::FindRingAtomsAndBonds
74 #define OB_RINGFLAGS_MOL         (1<<2)
75   //! Aromatic flags have been set for atoms and bonds
76 #define OB_AROMATIC_MOL          (1<<3)
77   //! Atom typing has been performed. See OBAtomTyper
78 #define OB_ATOMTYPES_MOL         (1<<4)
79   //! Chirality detection has been performed.
80 #define OB_CHIRALITY_MOL         (1<<5)
81   //! Partial charges have been set or percieved
82 #define OB_PCHARGE_MOL           (1<<6)
83   //! Atom hybridizations have been set. See OBAtomTyper
84 #define OB_HYBRID_MOL            (1<<8)
85   //! Ring "closure" bonds have been set. See OBBond::IsClosure
86 #define OB_CLOSURE_MOL           (1<<11)
87   //! Hyrdogen atoms have been added where needed. See OBMol::AddHydrogens
88 #define OB_H_ADDED_MOL           (1<<12)
89   //! pH correction for hydrogen addition has been performed.
90 #define OB_PH_CORRECTED_MOL      (1<<13)
91   //! Biomolecular chains and residues have been set. See OBChainsParser
92 #define OB_CHAINS_MOL            (1<<15)
93   //! Total charge on this molecule has been set. See OBMol::SetTotalCharge
94 #define OB_TCHARGE_MOL		       (1<<16)
95   //! Total spin on this molecule has been set. See OBMol::SetTotalSpinMultiplicity
96 #define OB_TSPIN_MOL             (1<<17)
97   //! Ring typing has been performed. See OBRingTyper
98 #define OB_RINGTYPES_MOL         (1<<18)
99   //! A pattern, not a complete molecule.
100 #define OB_PATTERN_STRUCTURE     (1<<19)
101   //! Largest Set of Smallest Rings (LSSR) done. See OBRing and OBMol::FindLSSR
102 #define OB_LSSR_MOL              (1<<20)
103   //! SpinMultiplicities on atoms have been set in OBMol::AssignSpinMultiplicity()
104 #define OB_ATOMSPIN_MOL          (1<<21)
105   //! Treat as reaction
106 #define OB_REACTION_MOL          (1<<22)
107   //! Molecule is repeating in a periodic unit cell
108 #define OB_PERIODIC_MOL          (1<<23)
109   // flags 24-32 unspecified
110 
111 #define SET_OR_UNSET_FLAG(X) \
112   if (value) SetFlag(X); \
113   else     UnsetFlag(X);
114 
115 #define OB_CURRENT_CONFORMER	 -1
116 
117 enum HydrogenType { AllHydrogen, PolarHydrogen, NonPolarHydrogen };
118 
119   // class introduction in mol.cpp
120  class OBAPI OBMol: public OBBase
121   {
122   protected:
123     int                           _flags;	//!< bitfield of flags
124     bool                          _autoPartialCharge;//!< Assign partial charges automatically
125     bool                          _autoFormalCharge;//!< Assign formal charges automatically
126     std::string                   _title;     	//!< Molecule title
127     std::vector<OBAtom*>          _vatom;      	//!< vector of atoms
128     std::vector<OBAtom*>          _atomIds;    	//!< vector of atoms indexed by id
129     std::vector<OBBond*>          _vbond;      	//!< vector of bonds
130     std::vector<OBBond*>          _bondIds;     //!< vector of bonds
131     unsigned short int            _dimension;   //!< Dimensionality of coordinates
132     int				  _totalCharge; //!< Total charge on the molecule
133     unsigned int                  _totalSpin;   //!< Total spin on the molecule (if not specified, assumes lowest possible spin)
134     double                        *_c;	        //!< coordinate array
135     std::vector<double*>          _vconf;       //!< vector of conformers
136     double                        _energy;      //!< heat of formation
137     unsigned int                  _natoms;      //!< Number of atoms
138     unsigned int                  _nbonds;      //!< Number of bonds
139     std::vector<OBResidue*>       _residue;     //!< Residue information (if applicable)
140     std::vector<OBInternalCoord*> _internals;   //!< Internal Coordinates (if applicable)
141     unsigned short int            _mod;	        //!< Number of nested calls to BeginModify()
142 
143   public:
144 
145     //! \name Initialization and data (re)size methods
146     //@{
147     //! Constructor
148     OBMol();
149     //! Copy constructor, copies atoms,bonds and OBGenericData
150     OBMol(const OBMol &);
151     //! Destructor
152     virtual ~OBMol();
153     //! Assignment, copies atoms,bonds and OBGenericData
154     OBMol &operator=(const OBMol &mol);
155     //! Copies atoms and bonds but not OBGenericData
156     OBMol &operator+=(const OBMol &mol);
157 
158     //! Reserve a minimum number of atoms for internal storage
159     //! This improves performance since the internal atom vector does not grow.
ReserveAtoms(int natoms)160     void ReserveAtoms(int natoms)
161     {
162       if (natoms > 0 && _mod) {
163         _vatom.reserve(natoms);
164         _atomIds.reserve(natoms);
165       }
166     }
167 
168     //! Free an OBAtom pointer if defined. Does no bookkeeping
169     //! \see DeleteAtom which ensures internal connections
170     virtual void DestroyAtom(OBAtom*);
171     //! Free an OBBond pointer if defined. Does no bookkeeping
172     //! \see DeleteBond which ensures internal connections
173     virtual void DestroyBond(OBBond*);
174     //! Free an OBResidue pointer if defined. Does no bookkeeping
175     //! \see DeleteResidue which ensures internal connections
176     virtual void DestroyResidue(OBResidue*);
177 
178     //! Add the specified atom to this molecule
179     //! \param atom        the atom to add
180     //! \param forceNewId  whether to make a new atom Id even if the atom already has one (default is false)
181     //! \return Whether the method was successful
182     bool AddAtom(OBAtom& atom, bool forceNewId = false);
183     //! Add a new atom to this molecule (like AddAtom)
184     //! Calls BeginModify() before insertion and EndModify() after insertion
185     bool InsertAtom(OBAtom &);
186     //! Add a new bond to the molecule with the specified parameters
187     //! \param beginIdx  the atom index of the "start" atom
188     //! \param endIdx    the atom index of the "end" atom
189     //! \param order     the bond order (see OBBond::GetBondOrder())
190     //! \param flags     any bond flags such as stereochemistry (default = none)
191     //! \param insertpos the position index to insert the bond (default = none)
192     //! \return Whether the new bond creation was successful
193     bool AddBond(int beginIdx, int endIdx, int order,
194                  int flags=0,int insertpos=-1);
195     //! Add the specified residue to this molecule and update connections
196     //! \return Whether the method was successful
197     bool AddBond(OBBond&);
198     //! Add the specified residue to this molecule and update connections
199     //! \return Whether the method was successful
200     bool AddResidue(OBResidue&);
201 
202     //! Create a new OBAtom in this molecule and ensure connections
203     //! (e.g. OBAtom::GetParent(). A new unique id will be assigned
204     //! to this atom.
205     OBAtom    *NewAtom();
206     //! Create a new OBAtom in this molecule and ensure connections.
207     //! (e.g. OBAtom::GetParent(). The @p id will be assigned to this
208     //! atom.
209     OBAtom    *NewAtom(unsigned long id);
210     //! Create a new OBBond in this molecule and ensure connections
211     //! (e.g. OBBond::GetParent(). A new unique id will be assigned
212     //! to this bond.
213     OBBond    *NewBond();
214     //! Create a new OBBond in this molecule and ensure connections
215     //! (e.g. OBBond::GetParent(). The @p id will be assigned to this
216     //! bond.
217     OBBond    *NewBond(unsigned long id);
218     //! Create a new OBResidue in this molecule and ensure connections.
219     OBResidue *NewResidue();
220     //! Deletes an atom from this molecule and all appropriate bonds.
221     //! Updates the molecule and atom and bond indexes accordingly.
222     //! \warning Does not update any residues which may contain this atom
223     //! \return Whether deletion was successful
224     bool DeleteAtom(OBAtom*, bool destroyAtom = true);
225     //! Deletes an bond from this molecule and updates accordingly
226     //! \return Whether deletion was successful
227     bool DeleteBond(OBBond*, bool destroyBond = true);
228     //! Deletes a residue from this molecule and updates accordingly.
229     //! \return Whether deletion was successful
230     bool DeleteResidue(OBResidue*, bool destroyResidue = true);
231     //@}
232 
233     //! \name Molecule modification methods
234     //@{
235     //! Call when making many modifications -- clears conformer/rotomer data.
236     //! The method "turns off" perception routines, improving performance.
237     //! Changes in molecular structure will be re-considered after modifications.
238     virtual void BeginModify(void);
239     //! Call when done with modificaions -- re-perceive data as needed.
240     //! This method "turns on" perception routines and re-evaluates molecular
241     //! structure.
242     virtual void EndModify(bool nukePerceivedData=true);
243     //! \return The number of nested BeginModify() calls. Used internally.
GetMod()244     int GetMod()           {      return(_mod);    }
245     //! Increase the number of nested BeginModify calls. Dangerous!
246     //! Instead, properly use BeginModify as needed.
IncrementMod()247     void IncrementMod()    {      _mod++;          }
248     //! Decrease the number of nested BeginModify calls. Dangerous!
249     //! Instead, properly use EndModify as needed.
DecrementMod()250     void DecrementMod()    {      _mod--;          }
251     //@}
252 
253     //! \name Data retrieval methods
254     //@{
255     //! \return the entire set of flags. (Internal use, mainly.)
GetFlags()256     int          GetFlags()               { return(_flags); }
257     //! \return the title of this molecule (often the filename)
258     //! \param replaceNewlines whether to replace any newline characters with spaces
259     const char  *GetTitle(bool replaceNewlines = true) const;
260     //! \return the number of atoms (i.e. OBAtom children)
NumAtoms()261     unsigned int NumAtoms() const         {  return(_natoms); }
262     //! \return the number of bonds (i.e. OBBond children)
NumBonds()263     unsigned int NumBonds() const         {  return(_nbonds); }
264     //! \return the number of non-hydrogen atoms
265     unsigned int NumHvyAtoms();
266     //! \return the number of residues (i.e. OBResidue substituents)
NumResidues()267     unsigned int NumResidues() const      { return(static_cast<unsigned int> (_residue.size())); }
268     //! \return the number of rotatable bonds. If sampleRingBonds is true, will include rotors within rings (see OBBond::IsRotor() for details)
269     unsigned int NumRotors(bool sampleRingBonds=false);
270 
271     //! \return the atom at index @p idx or NULL if it does not exist.
272     //! \warning Atom indexing will change. Use iterator methods instead.
273     OBAtom      *GetAtom(int idx) const;
274     //! \return the atom with @p id or NULL if it does not exist.
275     OBAtom      *GetAtomById(unsigned long id) const;
276     //! \return the first atom in this molecule, or NULL if none exist.
277     //! \deprecated Will be removed in favor of more standard iterator methods
278     OBAtom      *GetFirstAtom() const;
279     //! \return the bond at index @p idx or NULL if it does not exist.
280     //! \warning Bond indexing may change. Use iterator methods instead.
281     OBBond      *GetBond(int idx) const;
282     //! \return the bond with @p id or NULL if it does not exist.
283     OBBond      *GetBondById(unsigned long id) const;
284     //! \return the bond connecting the atom indexed by @p a and @p b or NULL if none exists.
285     //! \warning Atom indexing will change. Use atom objects and iterators instead.
286     OBBond      *GetBond(int a, int b) const;
287     // The safer version of the above method
288     //! \return the bond between the atoms @p bgn and @p end or NULL if none exists
289     OBBond      *GetBond(OBAtom* bgn, OBAtom* end) const;
290     //! \return the residue indexed by @p idx, or NULL if none exists
291     //! \warning Residue indexing may change. Use iterator methods instead.
292     OBResidue   *GetResidue(int idx) const;
293     std::vector<OBInternalCoord*> GetInternalCoord();
294     /*! \return the dihedral angle (in degrees) between the four atoms supplied a1-a2-a3-a4)
295      *  WARNING: SetTorsion takes an angle in radians while GetTorsion returns it
296      *  in degrees
297      */
298     double       GetTorsion(int,int,int,int);
299     /*! \return the dihedral angle (in degrees) between the four atoms @p a, @p b, @p c, and @p d)
300      *  WARNING: SetTorsion takes an angle in radians while GetTorsion returns it
301      *  in degrees
302      */
303     double       GetTorsion(OBAtom* a,OBAtom* b,OBAtom* c,OBAtom* d);
304     //! \return the angle (in degrees) between the three atoms @p a, @p b and @p c
305     //!  (where  a-> b (vertex) -> c )
306     double GetAngle(OBAtom* a, OBAtom* b, OBAtom* c);
307     //! \return the size of the smallest ring if a and b are in the same ring, 0 otherwise
308     //! \since version 2.4
309     int AreInSameRing(OBAtom *a, OBAtom *b);
310     //! \return the stochoimetric formula (e.g., C4H6O)
311     std::string  GetFormula();
312     //! \return the stochoimetric formula in spaced format e.g. C 4 H 6 O 1
313     std::string  GetSpacedFormula(int ones=0, const char* sp=" ", bool implicitH = true);
314     //! \return the heat of formation for this molecule (in kcal/mol)
GetEnergy()315     double       GetEnergy() const { return _energy; }
316     //! \return the standard molar mass given by IUPAC atomic masses (amu)
317     double       GetMolWt(bool implicitH = true);
318     //! \return the mass given by isotopes (or most abundant isotope, if not specified)
319     double	 GetExactMass(bool implicitH = true);
320     //! \return the total charge on this molecule (i.e., 0 = neutral, +1, -1...)
321     int		 GetTotalCharge();
322     //! \return the total spin on this molecule (i.e., 1 = singlet, 2 = doublet...)
323     unsigned int GetTotalSpinMultiplicity();
324     //! \return the dimensionality of coordinates (i.e., 0 = unknown or no coord, 2=2D, 3=3D)
GetDimension()325     unsigned short int GetDimension() const { return _dimension; }
326     //! \return the set of all atomic coordinates. See OBAtom::GetCoordPtr for more
GetCoordinates()327     double      *GetCoordinates() { return(_c); }
328     //! \return the Smallest Set of Smallest Rings has been run (see OBRing class)
329     std::vector<OBRing*> &GetSSSR();
330     //! \return the Largest Set of Smallest Rings has been run (see OBRing class)
331     std::vector<OBRing*> &GetLSSR();
332     //! Get the current flag for whether formal charges are set with pH correction
AutomaticFormalCharge()333     bool AutomaticFormalCharge()   { return(_autoFormalCharge);  }
334     //! Get the current flag for whether partial charges are auto-determined
AutomaticPartialCharge()335     bool AutomaticPartialCharge()  { return(_autoPartialCharge); }
336     //@}
337 
338 
339     //! \name Data modification methods
340     //@{
341     //! Set the title of this molecule to @p title
342     void   SetTitle(const char *title);
343     //! Set the title of this molecule to @p title
344     void   SetTitle(std::string &title);
345     //! Set the stochiometric formula for this molecule
346     void   SetFormula(std::string molFormula);
347     //! Set the heat of formation for this molecule (in kcal/mol)
SetEnergy(double energy)348     void   SetEnergy(double energy) { _energy = energy; }
349     //! Set the dimension of this molecule (i.e., 0, 1 , 2, 3)
SetDimension(unsigned short int d)350     void   SetDimension(unsigned short int d) { _dimension = d; }
351     //! Set the total charge of this molecule to @p charge
352     void   SetTotalCharge(int charge);
353     //! Set the total spin multiplicity of this molecule to @p spinMultiplicity
354     //! Overrides the calculation from spin multiplicity of OBAtoms
355     void   SetTotalSpinMultiplicity(unsigned int spinMultiplicity);
356     //! Set the internal coordinates to @p int_coord
357     //! (Does not call InternalToCartesian to update the 3D cartesian
358     //! coordinates).
359     //! The size of the @p int_coord has to be the same as the number of atoms
360     //! in molecule (+ NULL at the beginning).
361     void SetInternalCoord(std::vector<OBInternalCoord*> int_coord);
362     //! Set the flag for determining automatic formal charges with pH (default=true)
SetAutomaticFormalCharge(bool val)363     void SetAutomaticFormalCharge(bool val)
364     { _autoFormalCharge=val;  }
365     //! Set the flag for determining partial charges automatically (default=true)
SetAutomaticPartialCharge(bool val)366     void SetAutomaticPartialCharge(bool val)
367     { _autoPartialCharge=val; }
368 
369     //! Mark that aromaticity has been perceived for this molecule (see OBAromaticTyper)
370     void   SetAromaticPerceived(bool value = true)    { SET_OR_UNSET_FLAG(OB_AROMATIC_MOL);    }
371     //! Mark that Smallest Set of Smallest Rings has been run (see OBRing class)
372     void   SetSSSRPerceived(bool value = true)        { SET_OR_UNSET_FLAG(OB_SSSR_MOL);        }
373     //! Mark that Largest Set of Smallest Rings has been run (see OBRing class)
374     void   SetLSSRPerceived(bool value = true)        { SET_OR_UNSET_FLAG(OB_LSSR_MOL);        }
375     //! Mark that rings have been perceived (see OBRing class for details)
376     void   SetRingAtomsAndBondsPerceived(bool value = true) { SET_OR_UNSET_FLAG(OB_RINGFLAGS_MOL); }
377     //! Mark that atom types have been perceived (see OBAtomTyper for details)
378     void   SetAtomTypesPerceived(bool value = true)   { SET_OR_UNSET_FLAG(OB_ATOMTYPES_MOL);   }
379     //! Mark that ring types have been perceived (see OBRingTyper for details)
380     void   SetRingTypesPerceived(bool value = true)   { SET_OR_UNSET_FLAG(OB_RINGTYPES_MOL);   }
381     //! Mark that chains and residues have been perceived (see OBChainsParser)
382     void   SetChainsPerceived(bool value = true)      { SET_OR_UNSET_FLAG(OB_CHAINS_MOL);      }
383     //! Mark that chirality has been perceived
384     void   SetChiralityPerceived(bool value = true)   { SET_OR_UNSET_FLAG(OB_CHIRALITY_MOL);   }
385     //! Mark that partial charges have been assigned
386     void   SetPartialChargesPerceived(bool value = true) { SET_OR_UNSET_FLAG(OB_PCHARGE_MOL);  }
387     //! Mark that hybridization of all atoms has been assigned
388     void   SetHybridizationPerceived(bool value = true)  { SET_OR_UNSET_FLAG(OB_HYBRID_MOL);   }
389     //! Mark that ring closure bonds have been assigned by graph traversal
390     void   SetClosureBondsPerceived(bool value = true)   { SET_OR_UNSET_FLAG(OB_CLOSURE_MOL);  }
391     //! Mark that explicit hydrogen atoms have been added
392 
393     void   SetHydrogensAdded(bool value = true) { SET_OR_UNSET_FLAG(OB_H_ADDED_MOL); }
394     void   SetCorrectedForPH(bool value = true) { SET_OR_UNSET_FLAG(OB_PH_CORRECTED_MOL); }
395     void   SetSpinMultiplicityAssigned(bool value = true) { SET_OR_UNSET_FLAG(OB_ATOMSPIN_MOL); }
396     //! The OBMol is a pattern, not a complete molecule. Left unchanged by Clear().
397     void   SetIsPatternStructure(bool value = true) { SET_OR_UNSET_FLAG(OB_PATTERN_STRUCTURE); }
398     void   SetIsReaction(bool value = true)               { SET_OR_UNSET_FLAG(OB_REACTION_MOL); }
399     //! Mark that distance calculations, etc., should apply periodic boundary conditions through the minimimum image convention.
400     //! Does not automatically recalculate bonding.
401     void   SetPeriodicMol(bool value = true){ SET_OR_UNSET_FLAG(OB_PERIODIC_MOL); }
HasFlag(int flag)402     bool   HasFlag(int flag)   { return (_flags & flag) ? true : false; }
SetFlag(int flag)403     void   SetFlag(int flag)   { _flags |= flag; }
UnsetFlag(int flag)404     void   UnsetFlag(int flag) { _flags &= (~(flag)); }
SetFlags(int flags)405     void   SetFlags(int flags) { _flags = flags; }
406 
407     //@}
408 
409     //! \name Molecule modification methods
410     //@{
411     // Description in transform.cpp (command-line transformations to this molecule)
412     virtual OBBase*    DoTransformations(const std::map<std::string,std::string>* pOptions,OBConversion* pConv);
413     // Ditto (documentation on transformation options)
414     static const char* ClassDescription();
415     //! Clear all information from a molecule except OB_PATTERN_STRUCTURE left unchanged
416     bool Clear();
417     //! Renumber the atoms of this molecule according to the order in the supplied vector
418     void RenumberAtoms(std::vector<OBAtom*>&);
419     //! Renumber the atoms of this molecule using the initial indexes in the supplied vector
420     void RenumberAtoms(std::vector<int>);
421     //! Set the coordinates for all atoms in this conformer.
422     //! \sa OBMol::GetCoordinates()
423     void SetCoordinates(double *c);
424     //! Translate one conformer and rotate by a rotation matrix (which is returned) to the inertial frame-of-reference
425     void ToInertialFrame(int conf, double *rmat);
426     //! Translate all conformers to the inertial frame-of-reference
427     void ToInertialFrame();
428     //! Translates all conformers in the molecule by the supplied vector
429     void Translate(const vector3 &v);
430     //! Translates one conformer in the molecule by the supplied vector
431     void Translate(const vector3 &v, int conf);
432     //! Rotate all conformers using the supplied matrix @p u (a 3x3 array of double)
433     void Rotate(const double u[3][3]);
434     //! Rotate all conformers using the supplied matrix @p m (a linear 3x3 row-major array of double)
435     void Rotate(const double m[9]);
436     //! Rotate a specific conformer @p nconf using the supplied rotation matrix @p m
437     void Rotate(const double m[9],int nconf);
438     //! Translate to the center of all coordinates (for this conformer)
439     void Center();
440     //! Suppress hydrogens by converting explicit hydrogen atoms to implicit
441     //! \return Success
442     bool DeleteHydrogens();
443     //! Suppress explicit hydrogen atoms on the supplied atom
444     //! \return Success
445     bool DeleteHydrogens(OBAtom*);
446     //! Suppress explicit hydrogen atoms connected to a polar atom
447     //! \see OBAtom::IsPolarHydrogen
448     //! \since version 2.4
449     bool DeletePolarHydrogens();
450     //! Suppress explicit hydrogen atoms connected to a non-polar atom
451     //! \see OBAtom::IsNonPolarHydrogen
452     bool DeleteNonPolarHydrogens();
453     //! Suppress the supplied atom if it is a hydrogen
454     //! (Helper function for DeleteHydrogens)
455     bool DeleteHydrogen(OBAtom*);
456     //! Convert implicit hydrogens to explicit atoms in the molecular graph
457     //! \param polaronly    Whether to add hydrogens only to polar atoms
458     //! (i.e., not to C atoms)
459     //! \param correctForPH Whether to call CorrectForPH() first
460     //! \param pH The pH to use for CorrectForPH() modification
461     //! \return Whether any hydrogens were added
462     bool AddHydrogens(bool polaronly=false,bool correctForPH=false, double pH=7.4);
463     //! For a particular atom, convert implicit hydrogens to explicit atoms in the molecular graph
464     bool AddHydrogens(OBAtom*);
465     //! For polar atoms only, convert implicit hydrogens to explicit atoms in the molecular graph
466     bool AddPolarHydrogens();
467     //! For non-polar atoms only, convert implicit hydrogens to explicit atoms in the molecular graph
468     //! \since version 2.4
469     bool AddNonPolarHydrogens();
470     //! For polar and/or non-polar atoms, convert implicit hydrogens to explicit atoms in the molecular graph
471     //! \since verison 2.4
472     bool AddNewHydrogens(HydrogenType whichHydrogen, bool correctForPH=false, double pH=7.4);
473 
474     //! If @p threshold is not specified or is zero, remove all but the largest
475     //! contiguous fragment. If @p threshold is non-zero, remove any fragments with fewer
476     //! than @p threshold atoms.
477     bool StripSalts(unsigned int threshold=0);
478     //! Copies each disconnected fragment as a separate OBMol
479     std::vector<OBMol> Separate(int StartIndex=1);
480     //! Iterative component of Separate to copy one fragment at a time
481     bool GetNextFragment( OpenBabel::OBMolAtomDFSIter& iter, OBMol& newMol );
482     // docs in mol.cpp
483     bool CopySubstructure(OBMol& newmol, OBBitVec *includeatoms, OBBitVec *excludebonds = (OBBitVec*)nullptr,
484       unsigned int correctvalence=1,
485       std::vector<unsigned int> *atomorder=(std::vector<unsigned int>*)nullptr,
486       std::vector<unsigned int> *bondorder=(std::vector<unsigned int>*)nullptr);
487     //! Converts the charged form of coordinate bonds, e.g.[N+]([O-])=O to N(=O)=O
488     bool ConvertDativeBonds();
489     //! Converts 5-valent N and P only. Return true if conversion occurred.
490     //! \return has charged form of dative bonds(e.g.[N+]([O-])=O from N(=O)=O).
491     //! \since version 2.4
492     bool MakeDativeBonds();
493     /** Convert zero-order bonds to single or double bonds and adjust adjacent atom
494      *  charges in an attempt to achieve the correct valence state.
495      *  @return Whether any modifications were made
496      *  @since version 2.4
497      */
498     bool ConvertZeroBonds();
499 
500     //! Correct for pH by applying the OBPhModel transformations
501     bool CorrectForPH(double pH=7.4);
502     // docs in mol.cpp
503     bool AssignSpinMultiplicity(bool NoImplicitH=false);
504 
505     //! Put the specified molecular charge on appropriate atoms.
506     //! Assumes all the hydrogen is explicitly included in the molecule.
507     //! \since version 2.4
508     bool AssignTotalChargeToAtoms(int charge);
509 
510     //! \return the center of the supplied conformer @p nconf
511     //! \see Center() to actually center all conformers at the origin
512     vector3 Center(int nconf);
513     /*! Set the torsion defined by these atoms, rotating bonded neighbors
514      *  \par ang The torsion angle in radians
515      *  WARNING: SetTorsion takes an angle in radians while GetTorsion returns it
516      *  in degrees
517      */
518     void SetTorsion(OBAtom*,OBAtom*,OBAtom*,OBAtom*,double ang);
519     //@}
520 
521     //! \name Molecule utilities and perception methods
522     //@{
523     //! Find Smallest Set of Smallest Rings (see OBRing class for more details)
524     void FindSSSR();
525     //! Find Largest Set of Smallest Rings
526     void FindLSSR();
527     //! Find all ring atoms and bonds. Does not need to call FindSSSR().
528     void FindRingAtomsAndBonds();
529     // documented in mol.cpp -- locates all atom indexes which can reach 'end'
530     void FindChildren(std::vector<int> & children,int bgnIdx,int endIdx);
531     // documented in mol.cpp -- locates all atoms which can reach 'end'
532     void FindChildren(std::vector<OBAtom*>& children,OBAtom* bgn,OBAtom* end);
533     //! Find the largest fragment in OBMol
534     //! (which may include multiple non-connected fragments)
535     //! \param frag   Return (by reference) a bit vector indicating the atoms
536     //! in the largest fragment
537     void FindLargestFragment(OBBitVec &frag);
538     //! Sort a list of contig fragments by size from largest to smallest
539     //! Each vector<int> contains the atom numbers of a contig fragment
540     void ContigFragList(std::vector<std::vector<int> >&);
541     //! Aligns atom a on p1 and atom b along p1->p2 vector
542     void Align(OBAtom*,OBAtom*,vector3&,vector3&);
543     //! Adds single bonds based on atom proximity
544     void ConnectTheDots();
545     //! Attempts to perceive multiple bonds based on geometries
546     void PerceiveBondOrders();
547     //! Fills out an OBAngleData with angles from the molecule
548     void FindAngles();
549     //! Fills out an OBTorsionData with angles from the molecule
550     void FindTorsions();
551     // documented in mol.cpp: graph-theoretical distance for each atom
552     bool         GetGTDVector(std::vector<int> &);
553     // documented in mol.cpp: graph-invariant index for each atom
554     void         GetGIVector(std::vector<unsigned int> &);
555     // documented in mol.cpp: calculate symmetry-unique identifiers
556     void         GetGIDVector(std::vector<unsigned int> &);
557     //@}
558 
559     //! \name Methods to check for existence of properties
560     //@{
561     //! Are there non-zero coordinates in two dimensions (i.e. X and Y)- and, if Not3D is true, no Z coordinates?
562     bool Has2D(bool Not3D=false);
563     //! Are there non-zero coordinates in all three dimensions (i.e. X, Y, Z)?
564     bool Has3D();
565     //! Are there any non-zero coordinates?
566     bool HasNonZeroCoords();
567     //! Has aromatic perception been performed?
HasAromaticPerceived()568     bool HasAromaticPerceived()     { return(HasFlag(OB_AROMATIC_MOL)); }
569     //! Has the smallest set of smallest rings (FindSSSR) been performed?
HasSSSRPerceived()570     bool HasSSSRPerceived()         { return(HasFlag(OB_SSSR_MOL));     }
571     //! Has the largest set of smallest rings (FindLSSR) been performed?
HasLSSRPerceived()572     bool HasLSSRPerceived()         { return(HasFlag(OB_LSSR_MOL));     }
573     //! Have ring atoms and bonds been assigned?
HasRingAtomsAndBondsPerceived()574     bool HasRingAtomsAndBondsPerceived(){return(HasFlag(OB_RINGFLAGS_MOL));}
575     //! Have atom types been assigned by OBAtomTyper?
HasAtomTypesPerceived()576     bool HasAtomTypesPerceived()    { return(HasFlag(OB_ATOMTYPES_MOL));}
577     //! Have ring types been assigned by OBRingTyper?
HasRingTypesPerceived()578     bool HasRingTypesPerceived()    { return(HasFlag(OB_RINGTYPES_MOL));}
579     //! Has atom chirality been assigned?
HasChiralityPerceived()580     bool HasChiralityPerceived()    { return(HasFlag(OB_CHIRALITY_MOL));}
581     //! Have atomic Gasteiger partial charges been assigned by OBGastChrg?
HasPartialChargesPerceived()582     bool HasPartialChargesPerceived() { return(HasFlag(OB_PCHARGE_MOL));}
583     //! Has atomic hybridization been assigned by OBAtomTyper?
HasHybridizationPerceived()584     bool HasHybridizationPerceived() { return(HasFlag(OB_HYBRID_MOL));  }
585     //! Have ring "closure" bonds been assigned? (e.g., OBBond::IsClosure())
HasClosureBondsPerceived()586     bool HasClosureBondsPerceived() { return(HasFlag(OB_CLOSURE_MOL));  }
587     //! Have biomolecule chains and residues been assigned by OBChainsParser?
HasChainsPerceived()588     bool HasChainsPerceived() { return(HasFlag(OB_CHAINS_MOL));         }
589     //! Have hydrogens been added to the molecule?
HasHydrogensAdded()590     bool HasHydrogensAdded() { return(HasFlag(OB_H_ADDED_MOL));         }
591     //! Has the molecule been corrected for pH by CorrectForPH?
IsCorrectedForPH()592     bool IsCorrectedForPH() { return(HasFlag(OB_PH_CORRECTED_MOL));     }
593     //! Has total spin multiplicity been assigned?
HasSpinMultiplicityAssigned()594     bool HasSpinMultiplicityAssigned() { return(HasFlag(OB_ATOMSPIN_MOL)); }
595     //! Does this OBMol represent a reaction?
IsReaction()596     bool IsReaction()                  { return HasFlag(OB_REACTION_MOL); }
597     //! Is this molecule periodic? Should periodic boundary conditions be applied?
IsPeriodic()598     bool IsPeriodic() { return(HasFlag(OB_PERIODIC_MOL)); }
599 
600     //! Are there any atoms in this molecule?
Empty()601     bool Empty()                       { return(_natoms == 0);          }
602     //@}
603 
604     //! \name Multiple conformer member functions
605     //@{
606     //! \return the number of conformers in this molecule
NumConformers()607     int     NumConformers()    { return((_vconf.empty())?0:static_cast<int> (_vconf.size())); }
608     //! Set the entire set of conformers for this molecule to @p v
609     void    SetConformers(std::vector<double*> &v);
610     //! Add a new set of coordinates @p f as a new conformer
AddConformer(double * f)611     void    AddConformer(double *f)    {  _vconf.push_back(f);    }
612     //! Set the molecule's current conformer to @p i
613     //! Does nothing if @p i is larger than NumConformers()
614     void    SetConformer(unsigned int i);
615     //! Copy the conformer @p nconf into the array @p c
616     //! \warning Does no checking to see if @p c is large enough
617     void    CopyConformer(double* c,int nconf);
618     //! Delete the conformer @p nconf
619     void    DeleteConformer(int nconf);
620     //! \return the coordinates to conformer @p i
GetConformer(int i)621     double  *GetConformer(int i)       {  return(_vconf[i]);      }
622     //! Set the entire set of conformer energies
623     void    SetEnergies(std::vector<double> &energies);
624     //! Set the entire set of conformer energies
625     std::vector<double> GetEnergies();
626     //! Get the energy for conformer ci
627     //! \par ci conformer index
628     double  GetEnergy(int ci);
629     //! Set the iterator to the beginning of the conformer list
630     //! \return the array of coordinates for the first conformer
BeginConformer(std::vector<double * >::iterator & i)631     double  *BeginConformer(std::vector<double*>::iterator&i)
632     { i = _vconf.begin();
633       return((i == _vconf.end()) ? nullptr:*i); }
634     //! Advance the iterator to the next confomer, if possible
635     //! \return The array of coordinates for the next conformer, or NULL if none exist
NextConformer(std::vector<double * >::iterator & i)636     double  *NextConformer(std::vector<double*>::iterator&i)
637     { ++i;
638       return((i == _vconf.end()) ? nullptr:*i); }
639     //! \return the entire set of conformers for this molecule as a vector of floating point arrays
GetConformers()640     std::vector<double*> &GetConformers() {   return(_vconf);     }
641     //@}
642 
643     //! \name Iterator methods
644     //@{
645     //! \return An atom iterator pointing to the beginning of the atom list
BeginAtoms()646     OBAtomIterator BeginAtoms()   { return _vatom.begin(); }
647     //! \return An atom iterator pointing to the end of the atom list
EndAtoms()648     OBAtomIterator EndAtoms()    { return _vatom.begin() + NumAtoms() ; }
649     //! \return A bond iterator pointing to the beginning of the bond list
BeginBonds()650     OBBondIterator BeginBonds()   { return _vbond.begin(); }
651     //! \return A bond iterator pointing to the end of the bond list
EndBonds()652     OBBondIterator EndBonds()     { return _vbond.begin() + NumBonds() ; }
653     //! \return A residue iterator pointing to the beginning of the residue list
BeginResidues()654     OBResidueIterator BeginResidues() { return _residue.begin(); }
655     //! \return A residue iterator pointing to the end of the residue list
EndResidues()656     OBResidueIterator EndResidues()   { return _residue.end();   }
657 
658     //! Set the iterator @p i to the beginning of the atom list
659     //! \return the first atom (or NULL if none exist)
660     OBAtom *BeginAtom(OBAtomIterator &i);
661     //! Advance the iterator @p i to the next atom in the molecule
662     //! \return the next atom (if any, or NULL if none exist)
663     OBAtom *NextAtom(OBAtomIterator &i);
664     //! Set the iterator @p i to the beginning of the bond list
665     //! \return the first bond (or NULL if none exist)
666     OBBond *BeginBond(OBBondIterator &i);
667     //! Advance the iterator @p i to the next bond in the molecule
668     //! \return the next bond (if any, or NULL if none exist)
669     OBBond *NextBond(OBBondIterator &i);
670     //! Set the iterator @p i to the beginning of the resdiue list
671     //! \return the first residue (or NULL if none exist)
BeginResidue(OBResidueIterator & i)672     OBResidue *BeginResidue(OBResidueIterator &i)
673     {
674       i = _residue.begin();
675       return((i == _residue.end()) ? nullptr:*i);
676     }
677     //! Advance the iterator @p i to the next residue in the molecule
678     //! \return the next residue (if any, or NULL if not possible)
NextResidue(OBResidueIterator & i)679     OBResidue *NextResidue(OBResidueIterator &i)
680     {
681       ++i;
682       return((i == _residue.end()) ? nullptr:*i);
683     }
684     //! Set the iterator to the beginning of the internal coordinate list
685     //! \return the first internal coordinate record, or NULL if none exist
686     //! \see SetInternalCoord
BeginInternalCoord(std::vector<OBInternalCoord * >::iterator & i)687     OBInternalCoord *BeginInternalCoord(std::vector<OBInternalCoord*>::iterator &i)
688     {
689       i = _internals.begin();
690       return((i == _internals.end()) ? nullptr:*i);
691     }
692     //! Advance the iterator to the next internal coordinate record
693     //! \return the next first internal coordinate record, or NULL if none exist
694     //! \see SetInternalCoord
NextInternalCoord(std::vector<OBInternalCoord * >::iterator & i)695     OBInternalCoord *NextInternalCoord(std::vector<OBInternalCoord*>::iterator &i)
696     {
697       ++i;
698       return((i == _internals.end()) ? nullptr:*i);
699     }
700     //@}
701 
702   };
703 
704   // Utility function prototypes
705   //tokenize and Trim declarations moved to base.h
706   // Deprecated -- use OBMessageHandler class instead (docs in obutil.cpp)
707   OBAPI void ThrowError(char *str);
708   // Deprecated -- use OBMessageHandler class instead (docs in obutil.cpp)
709   OBAPI void ThrowError(std::string &str);
710   //! Convert Cartesian XYZ to a set of OBInternalCoord coordinates
711   OBAPI void CartesianToInternal(std::vector<OBInternalCoord*>&,OBMol&);
712   //! Convert set of OBInternalCoord coordinates into Cartesian XYZ
713   OBAPI void InternalToCartesian(std::vector<OBInternalCoord*>&,OBMol&);
714   // Replace the last extension in str with a new one (docs in obutil.cpp)
715   OBAPI std::string NewExtension(std::string&,char*);
716 
717   //! \brief Nested namespace for max_value templates
718   namespace detail {
719     //! \struct max_value mol.h <openbabel/mol.h>
720     //! \brief a C++ template to return the maximum value of a type (e.g., int)
721     template<typename T, int size = sizeof(T)>
722     struct max_value
723     {
724       static const T result = (static_cast<T>(0xFF) << (size-1)*8) + max_value<T, size-1>::result;
725     };
726 
727     //! \brief a C++ template to return the maximum value of a type (e.g., int)
728     template<typename T>
729     struct max_value<T, 0>
730     {
731       static const T result = 0;
732     };
733   }
734 
735   // No unique id
736   static const unsigned long NoId = detail::max_value<unsigned long>::result;
737 
738   //Utility Macros
739 
740 #ifndef BUFF_SIZE
741 #define BUFF_SIZE 32768
742 #endif
743 
744 #ifndef EQ
745 #define EQ(a,b) (!strcmp((a), (b)))
746 #endif
747 
748 #ifndef EQn
749 #define EQn(a,b,n) (!strncmp((a), (b), (n)))
750 #endif
751 
752 #ifndef SQUARE
753 #define SQUARE(x) ((x)*(x))
754 #endif
755 
756 #ifndef IsUnsatType
757 #define IsUnsatType(x)  (EQ(x,"Car") || EQ(x,"C2") || EQ(x,"Sox") || EQ(x,"Sac") || EQ(x,"Pac") || EQ(x,"So2"))
758 #endif
759 
760 #ifndef __KCC
761   extern "C"
762   {
763     OBAPI void  get_rmat(double*,double*,double*,int);
764     OBAPI void  ob_make_rmat(double mat[3][3],double rmat[9]);
765     OBAPI void  qtrfit (double *r,double *f,int size,double u[3][3]);
766     OBAPI double superimpose(double*,double*,int);
767   }
768 #else
769   OBAPI void get_rmat(double*,double*,double*,int);
770   OBAPI void ob_make_rmat(double mat[3][3],double rmat[9]);
771   OBAPI void qtrfit (double *r,double *f,int size,double u[3][3]);
772   OBAPI double superimpose(double*,double*,int);
773 #endif // __KCC
774 
775 //  extern OBMol* (*CreateMolecule) (void);
776 
777 } // end namespace OpenBabel
778 
779 #endif // OB_MOL_H
780 
781 //! \file mol.h
782 //! \brief Handle molecules. Declarations of OBMol, OBAtom, OBBond, OBResidue.
783 //!        (the main header for Open Babel)
784