mol.h

Go to the documentation of this file.
00001 /**********************************************************************
00002 mol.h - Handle molecules. Declarations of OBMol, OBAtom, OBBond, OBResidue.
00003         (the main header for Open Babel)
00004  
00005 Copyright (C) 1998-2001 by OpenEye Scientific Software, Inc.
00006 Some portions Copyright (C) 2001-2005 by Geoffrey R. Hutchison
00007 Some portions Copyright (C) 2003 by Michael Banck
00008  
00009 This file is part of the Open Babel project.
00010 For more information, see <http://openbabel.sourceforge.net/>
00011  
00012 This program is free software; you can redistribute it and/or modify
00013 it under the terms of the GNU General Public License as published by
00014 the Free Software Foundation version 2 of the License.
00015  
00016 This program is distributed in the hope that it will be useful,
00017 but WITHOUT ANY WARRANTY; without even the implied warranty of
00018 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00019 GNU General Public License for more details.
00020 ***********************************************************************/
00021 
00022 #ifndef OB_MOL_H
00023 #define OB_MOL_H
00024 
00025 #include "babelconfig.h"
00026 
00027 #ifndef EXTERN
00028 #  define EXTERN extern
00029 #endif
00030 
00031 #include <math.h>
00032 
00033 #include <algorithm>
00034 #include <vector>
00035 #include <string>
00036 #include <map>
00037 
00038 #if HAVE_IOSTREAM
00039 #include <iostream>
00040 #elif HAVE_IOSTREAM_H
00041 #include <iostream.h>
00042 #endif
00043 
00044 #if HAVE_FSTREAM
00045 #include <fstream>
00046 #elif HAVE_FSTREAM_H
00047 #include <fstream.h>
00048 #endif
00049 
00050 #include "base.h"
00051 #include "data.h"
00052 #include "chains.h"
00053 #include "math/vector3.h"
00054 #include "bitvec.h"
00055 #include "ring.h"
00056 #include "generic.h"
00057 #include "typer.h"
00058 #include "oberror.h"
00059 #include "obiter.h"
00060 #include "reaction.h" //so it gets notices in DLL builds
00061 
00062 namespace OpenBabel
00063 {
00064 
00065 class OBAtom;
00066 class OBBond;
00067 class OBMol;
00068 class OBInternalCoord;
00069 
00070 // Class OBResidue
00071 // class introduction in residue.cpp
00072 class OBAPI OBResidue
00073 {
00074 public:
00075 
00077     OBResidue(void);
00082     OBResidue(const OBResidue &);
00084     virtual ~OBResidue(void);
00085 
00086     OBResidue &operator=(const OBResidue &);
00087 
00088     void    AddAtom(OBAtom *atom);
00089     void    InsertAtom(OBAtom *atom);
00090     void    RemoveAtom(OBAtom *atom);
00091     void    Clear(void);
00092 
00093     void    SetName(const std::string &resname);
00094     void    SetNum(unsigned int resnum);
00095     void    SetChain(char chain);
00096     void    SetChainNum(unsigned int chainnum);
00097     void    SetIdx(unsigned int idx);
00098 
00099     void    SetAtomID(OBAtom *atom, const std::string &id);
00100     void    SetHetAtom(OBAtom *atom, bool hetatm);
00102     void    SetSerialNum(OBAtom *atom, unsigned int sernum);
00103 
00104     std::string    GetName(void)                const;
00105     unsigned int   GetNum(void)                 const;
00106     unsigned int   GetNumAtoms()                const;
00107     char           GetChain(void)               const;
00108     unsigned int   GetChainNum(void)            const;
00109     unsigned int   GetIdx(void)                 const;
00110     unsigned int   GetResKey(void)              const;
00111 
00112     std::vector<OBAtom*> GetAtoms(void)         const;
00113     std::vector<OBBond*> GetBonds(bool = true)  const;
00114 
00115     std::string    GetAtomID(OBAtom *atom)      const;
00117     unsigned       GetSerialNum(OBAtom *atom)   const;
00118 
00119     bool           GetAminoAcidProperty(int)      const;
00120     bool           GetAtomProperty(OBAtom *, int) const;
00121     bool           GetResidueProperty(int)        const;
00122 
00123     bool           IsHetAtom(OBAtom *atom)      const;
00124     bool           IsResidueType(int)           const;
00125 
00127     OBAtom *BeginAtom(std::vector<OBAtom*>::iterator &i);
00129     OBAtom *NextAtom(std::vector<OBAtom*>::iterator &i);
00130 
00132 
00133     bool                              HasData(std::string &);
00134     bool                              HasData(const char *);
00135     bool                              HasData(unsigned int type);
00136     void                              DeleteData(unsigned int type);
00137     void                              DeleteData(OBGenericData*);
00138     void                              DeleteData(std::vector<OBGenericData*>&);
00139     void                              SetData(OBGenericData *d)
00140       { _vdata.push_back(d); }
00142     unsigned int                      DataSize()
00143       { return(_vdata.size()); }
00144     OBGenericData                    *GetData(unsigned int type);
00145     OBGenericData                    *GetData(std::string&);
00146     OBGenericData                    *GetData(const char *);
00147     std::vector<OBGenericData*>      &GetData()
00148       { return(_vdata); }
00149     std::vector<OBGenericData*>::iterator  BeginData()
00150       { return(_vdata.begin()); }
00151     std::vector<OBGenericData*>::iterator  EndData()
00152       { return(_vdata.end()); }
00154 
00155 protected: // members
00156 
00157     unsigned int                _idx;   
00158     char                        _chain; 
00159     unsigned int                _aakey; 
00160     unsigned int                _reskey;
00161     unsigned int                _resnum;
00162     std::string                 _resname;
00163 
00164     std::vector<bool>           _hetatm;
00165     std::vector<std::string>    _atomid;
00166     std::vector<OBAtom*>        _atoms; 
00167     std::vector<unsigned int>   _sernum;
00168     std::vector<OBGenericData*> _vdata; 
00169 }; // OBResidue
00170 
00171 
00172 //ATOM Property Macros (flags)
00174 #define OB_4RING_ATOM     (1<<1)
00176 #define OB_3RING_ATOM     (1<<2)
00178 #define OB_AROMATIC_ATOM  (1<<3)
00180 #define OB_RING_ATOM      (1<<4)
00182 #define OB_CSTEREO_ATOM   (1<<5)
00184 #define OB_ACSTEREO_ATOM  (1<<6)
00186 #define OB_DONOR_ATOM     (1<<7)
00188 #define OB_ACCEPTOR_ATOM  (1<<8)
00190 #define OB_CHIRAL_ATOM    (1<<9)
00192 #define OB_POS_CHIRAL_ATOM (1<<10)
00194 #define OB_NEG_CHIRAL_ATOM (1<<11)
00195 // 12-16 currently unused
00196 
00197 // Class OBAtom
00198 // class introduction in atom.cpp
00199 class OBAPI OBAtom : public OBNodeBase
00200 {
00201 protected:
00202     char                          _ele;         
00203     char                          _impval;      
00204     char                          _type[6];     
00205     short                         _fcharge;     
00206     unsigned short                _isotope;     
00207     short                           _spinmultiplicity;
00208 
00209     //unsigned short int          _idx;         //!< index in parent (inherited)
00210     unsigned short            _cidx;            
00211     unsigned short                _hyb;         
00212     unsigned short                _flags;       
00213     double                         _pcharge;    
00214     double                       **_c;          
00215     vector3                       _v;           
00216     OBResidue                    *_residue;     
00217     //OBMol                      *_parent;      //!< parent molecule (inherited)
00218     //vector<OBBond*>             _bond;        //!< connections (inherited)
00219     std::vector<OBGenericData*>   _vdata;       
00220 
00221     int  GetFlag() const    {  return(_flags);  }
00222     void SetFlag(int flag)  { _flags |= flag;   }
00223     bool HasFlag(int flag)  {  return((_flags & flag) ? true : false); }
00224 
00225 public:
00226 
00228     OBAtom();
00230     virtual ~OBAtom();
00232     OBAtom &operator = (OBAtom &);
00234     void Clear();
00235 
00237 
00238 
00239     void SetIdx(int idx)    { _idx = idx; _cidx = (idx-1)*3; }
00241     void SetHyb(int hyb)    { _hyb = hyb; }
00243     void SetAtomicNum(int atomicnum)    { _ele = (char)atomicnum; }
00245     void SetIsotope(unsigned int iso);
00246     void SetImplicitValence(int val)    { _impval = (char)val; }
00247     void IncrementImplicitValence()     { _impval++; }
00248     void DecrementImplicitValence()     { _impval--; }
00249     void SetFormalCharge(int fcharge)   { _fcharge = fcharge; }
00250     void SetSpinMultiplicity(short spin){ _spinmultiplicity = spin; }
00251     void SetType(char *type);
00252     void SetType(std::string &type);
00253     void SetPartialCharge(double pcharge){ _pcharge = pcharge; }
00254     void SetVector(vector3 &v);
00255     void SetVector(const double x,const double y,const double z);
00257     void SetCoordPtr(double **c)        { _c = c; _cidx = (GetIdx()-1)*3; }
00259     void SetVector();
00260     void SetResidue(OBResidue *res)     { _residue=res; }
00261     //  void SetParent(OBMol *ptr)      { _parent=ptr; } // inherited
00262     void SetAromatic()                  { SetFlag(OB_AROMATIC_ATOM); }
00263     void UnsetAromatic()                { _flags &= (~(OB_AROMATIC_ATOM)); }
00265     void SetClockwiseStereo()           { SetFlag(OB_CSTEREO_ATOM|OB_CHIRAL_ATOM); }
00267     void SetAntiClockwiseStereo()       { SetFlag(OB_ACSTEREO_ATOM|OB_CHIRAL_ATOM); }
00269     void SetPositiveStereo() { SetFlag(OB_POS_CHIRAL_ATOM|OB_CHIRAL_ATOM); }
00271     void SetNegativeStereo() { SetFlag(OB_NEG_CHIRAL_ATOM|OB_CHIRAL_ATOM); }
00273     void UnsetStereo()
00274     {
00275         _flags &= ~(OB_ACSTEREO_ATOM);
00276         _flags &= ~(OB_CSTEREO_ATOM);
00277         _flags &= ~(OB_POS_CHIRAL_ATOM);
00278         _flags &= ~(OB_NEG_CHIRAL_ATOM);
00279         _flags &= ~(OB_CHIRAL_ATOM);
00280     }
00282     void SetInRing()         { SetFlag(OB_RING_ATOM); }
00284     void SetChiral()         { SetFlag(OB_CHIRAL_ATOM); }
00286     void ClearCoordPtr()     { _c = NULL; _cidx=0; }
00288 
00290 
00291     //int        GetStereo()        const { return((int)_stereo);}
00292     int          GetFormalCharge()  const { return(_fcharge);    }
00293     unsigned int GetAtomicNum()     const { return((unsigned int)_ele); }
00294     unsigned short int GetIsotope() const { return(_isotope);    }
00295     int          GetSpinMultiplicity() const { return(_spinmultiplicity); }
00297     double       GetAtomicMass()    const;
00299     double       GetExactMass()     const;
00300     unsigned int GetIdx()           const { return((int)_idx);  }
00301     unsigned int GetCoordinateIdx() const { return((int)_cidx); }
00303     unsigned int GetCIdx()          const { return((int)_cidx); }
00305     unsigned int GetValence()       const
00306     {
00307         return((_vbond.empty()) ? 0 : _vbond.size());
00308     }
00310     unsigned int GetHyb()             const;
00312     unsigned int GetImplicitValence() const;
00314     unsigned int GetHvyValence()      const;
00316     unsigned int GetHeteroValence()   const;
00317     char        *GetType();
00318 
00320     double      GetX()    {        return(x());    }
00322     double      GetY()    {        return(y());    }
00324     double      GetZ()    {        return(z());    }
00325     double      x()
00326     {
00327         if (_c)
00328             return((*_c)[_cidx]);
00329         else
00330             return _v.x();
00331     }
00332     double      y()
00333     {
00334         if (_c)
00335             return((*_c)[_cidx+1]);
00336         else
00337             return _v.y();
00338     }
00339     double      z()
00340     {
00341         if (_c)
00342             return((*_c)[_cidx+2]);
00343         else
00344             return _v.z();
00345     }
00347     double     *GetCoordinate()
00348     {
00349       if (_c)
00350         return(&(*_c)[_cidx]);
00351       else
00352         return NULL;
00353     }
00355     vector3   &GetVector();
00357     double     GetPartialCharge();
00358     OBResidue *GetResidue();
00359     //OBMol   *GetParent()        {return((OBMol*)_parent);}
00361     bool       GetNewBondVector(vector3 &v,double length);
00362     OBBond    *GetBond(OBAtom *);
00363     OBAtom    *GetNextAtom();
00365 
00367 
00368 
00369     std::vector<OBEdgeBase*>::iterator BeginBonds()
00370       { return(_vbond.begin()); }
00372     std::vector<OBEdgeBase*>::iterator EndBonds()
00373       { return(_vbond.end());   }
00375     OBBond *BeginBond(std::vector<OBEdgeBase*>::iterator &i);
00377     OBBond *NextBond(std::vector<OBEdgeBase*>::iterator &i);
00379     OBAtom *BeginNbrAtom(std::vector<OBEdgeBase*>::iterator &);
00381     OBAtom *NextNbrAtom(std::vector<OBEdgeBase*>::iterator &);
00383 
00385     double GetDistance(int index);
00387     double GetDistance(OBAtom*);
00389     double GetAngle(int b, int c);
00391     double GetAngle(OBAtom *b, OBAtom *c);
00392 
00394 
00395     void NewResidue()
00396     {
00397         if (!_residue)
00398             _residue = new OBResidue;
00399     }
00400     void DeleteResidue()
00401     {
00402         if (_residue)
00403             delete _residue;
00404     }
00405     void AddBond(OBBond *bond)
00406     {
00407         _vbond.push_back((OBEdgeBase*)bond);
00408     }
00409     void InsertBond(std::vector<OBEdgeBase*>::iterator &i, OBBond *bond)
00410     {
00411         _vbond.insert(i, (OBEdgeBase*)bond);
00412     }
00413     bool DeleteBond(OBBond*);
00414         void ClearBond() {_vbond.clear();}
00416 
00418 
00419 
00420     unsigned int  CountFreeOxygens()      const;
00422     unsigned int  ImplicitHydrogenCount() const;
00424     unsigned int  ExplicitHydrogenCount() const;
00426     unsigned int  MemberOfRingCount()     const;
00428     unsigned int  MemberOfRingSize()      const;
00430     double        SmallestBondAngle();
00432     double        AverageBondAngle();
00434     unsigned int  BOSum()                 const;
00436     unsigned int  KBOSum()                const;
00438 
00440 
00441 
00442     bool HtoMethyl();
00444     bool SetHybAndGeom(int);
00446 
00448 
00449 
00450     bool HasResidue()    { return(_residue != NULL);    }
00451     bool IsHydrogen()    { return(GetAtomicNum() == 1); }
00452     bool IsCarbon()      { return(GetAtomicNum() == 6); }
00453     bool IsNitrogen()    { return(GetAtomicNum() == 7); }
00454     bool IsOxygen()      { return(GetAtomicNum() == 8); }
00455     bool IsSulfur()      { return(GetAtomicNum() == 16);}
00456     bool IsPhosphorus()  { return(GetAtomicNum() == 15);}
00457     bool IsAromatic()      const;
00458     bool IsInRing()        const;
00459     bool IsInRingSize(int) const;
00461     bool IsHeteroatom();
00463     bool IsNotCorH();
00465     bool IsConnected(OBAtom*);
00467     bool IsOneThree(OBAtom*);
00469     bool IsOneFour(OBAtom*);
00471     bool IsCarboxylOxygen();
00473     bool IsPhosphateOxygen();
00475     bool IsSulfateOxygen();
00477     bool IsNitroOxygen();
00478     bool IsAmideNitrogen();
00479     bool IsPolarHydrogen();
00480     bool IsNonPolarHydrogen();
00481     bool IsAromaticNOxide();
00483     bool IsChiral();
00484     bool IsAxial();
00486     bool IsClockwise()         { return(HasFlag(OB_CSTEREO_ATOM));  }
00488     bool IsAntiClockwise()     { return(HasFlag(OB_ACSTEREO_ATOM)); }
00490     bool IsPositiveStereo() { return(HasFlag(OB_POS_CHIRAL_ATOM)); }
00492     bool IsNegativeStereo() { return(HasFlag(OB_NEG_CHIRAL_ATOM)); }
00494     bool HasChiralitySpecified()
00495       { return(HasFlag(OB_CSTEREO_ATOM|OB_ACSTEREO_ATOM)); }
00497     bool HasChiralVolume()
00498       { return(HasFlag(OB_POS_CHIRAL_ATOM|OB_NEG_CHIRAL_ATOM)); }
00500     bool IsHbondAcceptor();
00502     bool IsHbondDonor();
00504     bool IsHbondDonorH();
00505     bool HasAlphaBetaUnsat(bool includePandS=true);
00506     bool HasBondOfOrder(unsigned int);
00507     int  CountBondsOfOrder(unsigned int);
00508     bool HasNonSingleBond();
00509     bool HasSingleBond()    {        return(HasBondOfOrder(1));    }
00510     bool HasDoubleBond()    {        return(HasBondOfOrder(2));    }
00511     bool HasAromaticBond()  {        return(HasBondOfOrder(5));    }
00513     bool MatchesSMARTS(const char *);
00515 
00517 
00518     bool                              HasData(std::string &);
00519     bool                              HasData(const char *);
00520     bool                              HasData(unsigned int type);
00521     void                              DeleteData(unsigned int type);
00522     void                              DeleteData(OBGenericData*);
00523     void                              DeleteData(std::vector<OBGenericData*>&);
00524     void                              SetData(OBGenericData *d)
00525     {        _vdata.push_back(d);    }
00527     unsigned int                      DataSize()
00528     {        return(_vdata.size());    }
00529     OBGenericData                    *GetData(unsigned int type);
00530     OBGenericData                    *GetData(std::string&);
00531     OBGenericData                    *GetData(const char *);
00532     std::vector<OBGenericData*>      &GetData() { return(_vdata); }
00533     std::vector<OBGenericData*>::iterator  BeginData()
00534     {        return(_vdata.begin());    }
00535     std::vector<OBGenericData*>::iterator  EndData()
00536     {        return(_vdata.end());      }
00538 }; // class OBAtom
00539 
00540 
00541 // Class OBBond
00542 
00543 //BOND Property Macros (flags)
00545 #define OB_AROMATIC_BOND  (1<<1)
00547 #define OB_WEDGE_BOND     (1<<2)
00549 #define OB_HASH_BOND      (1<<3)
00551 #define OB_RING_BOND      (1<<4)
00553 #define OB_TORUP_BOND     (1<<5)
00555 #define OB_TORDOWN_BOND   (1<<6)
00557 #define OB_KSINGLE_BOND   (1<<7)
00559 #define OB_KDOUBLE_BOND   (1<<8)
00561 #define OB_KTRIPLE_BOND   (1<<9)
00562 #define OB_CLOSURE_BOND   (1<<10)
00563 // 11-16 currently unused
00564 
00565 // class introduction in bond.cpp
00566 class OBAPI OBBond : public OBEdgeBase
00567 {
00568 protected:
00569     char                          _order; 
00570     unsigned short int            _flags; 
00571     //OBAtom                     *_bgn;   //!< Not needed, inherited from OBEdgeBase
00572     //OBAtom                     *_end;   //!< Not needed, inherited from OBEdgeBase
00573     //OBMol                      *_parent;//!< Not needed, inherited from OBEdgeBase
00574     //unsigned short int          _idx;   //!< Not needed, inherited from OBEdgeBase
00575     std::vector<OBGenericData*>   _vdata; 
00576 
00577     bool HasFlag(int flag)    { return((_flags & flag) != 0); }
00578     void SetFlag(int flag)    { _flags |= flag;               }
00579 
00580 public:
00582     OBBond();
00584     virtual ~OBBond();
00585 
00587 
00588     void SetIdx(int idx)
00589     {
00590         _idx = idx;
00591     }
00592     void SetBO(int order);
00593     void SetBegin(OBAtom *begin)
00594     {
00595         _bgn = begin;
00596     }
00597     void SetEnd(OBAtom *end)
00598     {
00599         _end = end;
00600     }
00601     // void SetParent(OBMol *ptr)               {_parent=ptr;} // (inherited)
00602     void SetLength(OBAtom*,double);
00603     void Set(int,OBAtom*,OBAtom*,int,int);
00604     void SetKSingle();
00605     void SetKDouble();
00606     void SetKTriple();
00607     void SetAromatic()    { SetFlag(OB_AROMATIC_BOND); }
00608     void SetHash()        { SetFlag(OB_HASH_BOND);     }
00609     void SetWedge()       { SetFlag(OB_WEDGE_BOND);    }
00610     void SetUp()          { SetFlag(OB_TORUP_BOND);    }
00611     void SetDown()        { SetFlag(OB_TORDOWN_BOND);  }
00612     void SetInRing()      { SetFlag(OB_RING_BOND);     }
00613     void SetClosure()     { SetFlag(OB_CLOSURE_BOND);  }
00614 
00615     void UnsetAromatic()  { _flags &= (~(OB_AROMATIC_BOND)); }
00616     void UnsetKekule()
00617     {
00618         _flags &= (~(OB_KSINGLE_BOND|OB_KDOUBLE_BOND|OB_KTRIPLE_BOND));
00619     }
00621 
00623 
00624     unsigned int     GetBO()            const { return((int)_order); }
00625     unsigned int     GetBondOrder()     const { return((int)_order); }
00626     unsigned int     GetFlags()         const { return(_flags);      }
00627     unsigned int     GetBeginAtomIdx()  const { return(_bgn->GetIdx()); }
00628     unsigned int     GetEndAtomIdx()    const { return(_end->GetIdx()); }
00629     OBAtom *GetBeginAtom()    { return((OBAtom*)_bgn);    }
00630     OBAtom *GetEndAtom()      { return((OBAtom*)_end);    }
00631     OBAtom *GetNbrAtom(OBAtom *ptr)
00632     {
00633         return((ptr != _bgn)? (OBAtom*)_bgn : (OBAtom*)_end);
00634     }
00635     // OBMol  *GetParent()                 {return(_parent);}  // (inherited)
00636     double   GetEquibLength();
00637     double   GetLength();
00638     int     GetNbrAtomIdx(OBAtom *ptr)
00639     {
00640         return((ptr!=_bgn)?_bgn->GetIdx():_end->GetIdx());
00641     }
00643 
00645 
00646     bool IsAromatic() const;
00647     bool IsInRing() const;
00652     bool IsRotor();
00653     bool IsAmide();
00654     bool IsPrimaryAmide();
00655     bool IsSecondaryAmide();
00656     bool IsEster();
00657     bool IsCarbonyl();
00658     bool IsSingle();
00659     bool IsDouble();
00660     bool IsTriple();
00661     bool IsKSingle();
00662     bool IsKDouble();
00663     bool IsKTriple();
00664     bool IsClosure();
00667     bool IsUp()    {    return(HasFlag(OB_TORUP_BOND));    }
00670     bool IsDown()  {    return(HasFlag(OB_TORDOWN_BOND));  }
00671     bool IsWedge() {    return(HasFlag(OB_WEDGE_BOND));    }
00672     bool IsHash()  {    return(HasFlag(OB_HASH_BOND));     }
00674     bool IsDoubleBondGeometry();
00676 
00678 
00679     bool                              HasData(std::string &);
00680     bool                              HasData(const char *);
00681     bool                              HasData(unsigned int type);
00682     void                              DeleteData(unsigned int type);
00683     void                              DeleteData(OBGenericData*);
00684     void                              DeleteData(std::vector<OBGenericData*>&);
00685     void                              SetData(OBGenericData *d)
00686     {
00687         _vdata.push_back(d);
00688     }
00690     unsigned int                      DataSize()
00691     {
00692         return(_vdata.size());
00693     }
00694     OBGenericData                    *GetData(unsigned int type);
00695     OBGenericData                    *GetData(std::string&);
00696     OBGenericData                    *GetData(const char *);
00697     std::vector<OBGenericData*>           &GetData()
00698     {
00699         return(_vdata);
00700     }
00701     std::vector<OBGenericData*>::iterator  BeginData()
00702     {
00703         return(_vdata.begin());
00704     }
00705     std::vector<OBGenericData*>::iterator  EndData()
00706     {
00707         return(_vdata.end());
00708     }
00710 }
00711 ; // class OBBond
00712 
00713 
00714 // Class OBMol
00715 
00716 //MOL Property Macros (flags) -- 32+ bits
00717 #define OB_SSSR_MOL              (1<<1)
00718 #define OB_RINGFLAGS_MOL         (1<<2)
00719 #define OB_AROMATIC_MOL          (1<<3)
00720 #define OB_ATOMTYPES_MOL         (1<<4)
00721 #define OB_CHIRALITY_MOL         (1<<5)
00722 #define OB_PCHARGE_MOL           (1<<6)
00723 #define OB_HYBRID_MOL            (1<<8)
00724 #define OB_IMPVAL_MOL            (1<<9)
00725 #define OB_KEKULE_MOL            (1<<10)
00726 #define OB_CLOSURE_MOL           (1<<11)
00727 #define OB_H_ADDED_MOL           (1<<12)
00728 #define OB_PH_CORRECTED_MOL      (1<<13)
00729 #define OB_AROM_CORRECTED_MOL    (1<<14)
00730 #define OB_CHAINS_MOL            (1<<15)
00731 #define OB_TCHARGE_MOL           (1<<16)
00732 #define OB_TSPIN_MOL             (1<<17)
00733 // flags 18-32 unspecified
00734 #define OB_CURRENT_CONFORMER     -1
00735 
00736 // class introduction in mol.cpp
00737 class OBAPI OBMol : public OBGraphBase
00738 {
00739 protected:
00740     int                           _flags;       
00741     bool                          _autoPartialCharge; 
00742     bool                          _autoFormalCharge; 
00743     std::string                   _title;       
00744     //vector<OBAtom*>             _atom;        //!< not needed (inherited)
00745     //vector<OBBond*>             _bond;        //!< not needed (inherited)
00746     unsigned short int            _dimension;   
00747     double                        _energy;      
00748     int                           _totalCharge; 
00749     unsigned int                  _totalSpin;   
00750     double                       *_c;           
00751     std::vector<double*>          _vconf;       
00752     unsigned short int            _natoms;      
00753     unsigned short int            _nbonds;      
00754     std::vector<OBResidue*>       _residue;     
00755     std::vector<OBInternalCoord*> _internals;   
00756     std::vector<OBGenericData*>   _vdata;       
00757     unsigned short int            _mod;         
00758 
00759     bool  HasFlag(int flag)    { return((_flags & flag) ? true : false); }
00760     void  SetFlag(int flag)    { _flags |= flag; }
00761 
00763 
00764     void start_kekulize(std::vector <OBAtom*> &cycle, std::vector<int> &electron);
00765     int expand_kekulize(OBAtom *atom1, OBAtom *atom2, std::vector<int> &currentState, std::vector<int> &initState, std::vector<int> &bcurrentState, std::vector<int> &binitState, std::vector<bool> &mark);
00766     int getorden(OBAtom *atom);
00767     void expandcycle(OBAtom *atom, OBBitVec &avisit);
00769 
00770 public:
00771 
00773 
00774 
00775     OBMol();
00780     OBMol(const OBMol &);
00782     virtual ~OBMol();
00784     OBMol &operator=(const OBMol &mol);
00785     OBMol &operator+=(const OBMol &mol);
00786     void ReserveAtoms(int natoms)
00787     {
00788         if (natoms && _mod)
00789             _vatom.reserve(natoms);
00790     }
00791     virtual OBAtom *CreateAtom(void);
00792     virtual OBBond *CreateBond(void);
00793     virtual void DestroyAtom(OBNodeBase*);
00794     virtual void DestroyBond(OBEdgeBase*);
00795     bool AddAtom(OBAtom&);
00796     bool AddBond(int,int,int,int flags=0,int insertpos=-1);
00797     bool AddBond(OBBond&);
00798     bool AddResidue(OBResidue&);
00799     bool InsertAtom(OBAtom &);
00800     bool DeleteAtom(OBAtom*);
00801     bool DeleteBond(OBBond*);
00802     bool DeleteResidue(OBResidue*);
00803     OBAtom    *NewAtom();
00804     OBResidue *NewResidue();
00806 
00808 
00809 
00810     virtual void BeginModify(void);
00812     virtual void EndModify(bool nukePerceivedData=true);
00813     int GetMod()
00814     {
00815         return(_mod);
00816     }
00817     void IncrementMod()
00818     {
00819         _mod++;
00820     }
00821     void DecrementMod()
00822     {
00823         _mod--;
00824     }
00826 
00828 
00829 
00830     bool                              HasData(std::string &);
00832     bool                              HasData(const char *);
00834     bool                              HasData(unsigned int type);
00835     void                              DeleteData(unsigned int type);
00836     void                              DeleteData(OBGenericData*);
00837     void                              DeleteData(std::vector<OBGenericData*>&);
00838     void                              SetData(OBGenericData *d)
00839     {
00840         _vdata.push_back(d);
00841     }
00843     unsigned int                      DataSize(){ return(_vdata.size()); }
00844     OBGenericData                    *GetData(unsigned int type);
00845     OBGenericData                    *GetData(std::string&);
00846     OBGenericData                    *GetData(const char *);
00847     std::vector<OBGenericData*>      &GetData() { return(_vdata); }
00848     std::vector<OBGenericData*>::iterator  BeginData()
00849     {
00850         return(_vdata.begin());
00851     }
00852     std::vector<OBGenericData*>::iterator  EndData()
00853     {
00854         return(_vdata.end());
00855     }
00857 
00859 
00860     int          GetFlags()               { return(_flags); }
00862     const char  *GetTitle() const         { return(_title.c_str()); }
00864     unsigned int NumAtoms() const         {  return(_natoms); }
00866     unsigned int NumBonds() const         {  return(_nbonds); }
00868     unsigned int NumHvyAtoms();
00870     unsigned int NumResidues() const      { return(_residue.size()); }
00872     unsigned int NumRotors();
00873     
00874     OBAtom      *GetAtom(int);
00875     OBAtom      *GetFirstAtom();
00876     OBBond      *GetBond(int);
00877     OBBond      *GetBond(int, int);
00878     OBBond      *GetBond(OBAtom*,OBAtom*);
00879     OBResidue   *GetResidue(int);
00880     std::vector<OBInternalCoord*> GetInternalCoord();
00882     double       GetTorsion(int,int,int,int);
00884     double       GetTorsion(OBAtom*,OBAtom*,OBAtom*,OBAtom*);
00886     std::string  GetFormula();
00888     double       GetEnergy() const { return(_energy); }
00890     double       GetMolWt();
00892     double       GetExactMass();
00894     int          GetTotalCharge();
00896     unsigned int GetTotalSpinMultiplicity();
00898     unsigned short int GetDimension() const { return _dimension; }
00899     double      *GetCoordinates() { return(_c); }
00901     std::vector<OBRing*> &GetSSSR();
00903     bool AutomaticFormalCharge()   { return(_autoFormalCharge);  }
00905     bool AutomaticPartialCharge()  { return(_autoPartialCharge); }
00907 
00908 
00910 
00911     void   SetTitle(const char *title);
00912     void   SetTitle(std::string &title);
00914     void   SetFormula(std::string molFormula);
00916     void   SetEnergy(double energy) { _energy = energy; }
00918     void   SetDimension(unsigned short int d) { _dimension = d; }
00919     void   SetTotalCharge(int charge);
00920     void   SetTotalSpinMultiplicity(unsigned int spin);
00921     void   SetInternalCoord(std::vector<OBInternalCoord*> int_coord)
00922       { _internals = int_coord; }
00924     void SetAutomaticFormalCharge(bool val)
00925     { _autoFormalCharge=val;  }
00927     void SetAutomaticPartialCharge(bool val)
00928     { _autoPartialCharge=val; }
00929 
00931     void   SetAromaticPerceived()    { SetFlag(OB_AROMATIC_MOL);    }
00933     void   SetSSSRPerceived()        { SetFlag(OB_SSSR_MOL);        }
00935     void   SetRingAtomsAndBondsPerceived(){SetFlag(OB_RINGFLAGS_MOL);}
00937     void   SetAtomTypesPerceived()   { SetFlag(OB_ATOMTYPES_MOL);   }
00939     void   SetChainsPerceived()      { SetFlag(OB_CHAINS_MOL);      }
00941     void   SetChiralityPerceived()   { SetFlag(OB_CHIRALITY_MOL);   }
00943     void   SetPartialChargesPerceived(){ SetFlag(OB_PCHARGE_MOL);   }
00944     void   SetHybridizationPerceived() { SetFlag(OB_HYBRID_MOL);    }
00945     void   SetImplicitValencePerceived(){ SetFlag(OB_IMPVAL_MOL);   }
00946     void   SetKekulePerceived()      { SetFlag(OB_KEKULE_MOL);      }
00947     void   SetClosureBondsPerceived(){ SetFlag(OB_CLOSURE_MOL);     }
00948     void   SetHydrogensAdded()       { SetFlag(OB_H_ADDED_MOL);     }
00949     void   SetCorrectedForPH()       { SetFlag(OB_PH_CORRECTED_MOL);}
00950     void   SetAromaticCorrected()    { SetFlag(OB_AROM_CORRECTED_MOL);}
00951     void   SetSpinMultiplicityAssigned(){ SetFlag(OB_TSPIN_MOL);    }
00952     void   SetFlags(int flags)       { _flags = flags;              }
00953 
00954     void   UnsetAromaticPerceived()  { _flags &= (~(OB_AROMATIC_MOL));   }
00955     void   UnsetPartialChargesPerceived(){ _flags &= (~(OB_PCHARGE_MOL));}
00956     void   UnsetImplicitValencePerceived(){_flags &= (~(OB_IMPVAL_MOL)); }
00957     void   UnsetFlag(int flag)       { _flags &= (~(flag));              }
00958 
00960 
00961     // Description in transform.cpp
00962     virtual OBBase*    DoTransformations(const std::map<std::string,std::string>* pOptions);
00963     static const char* ClassDescription();
00965     bool Clear();
00967     void RenumberAtoms(std::vector<OBNodeBase*>&);
00969     void ToInertialFrame(int conf, double *rmat);
00971     void ToInertialFrame();
00973     void Translate(const vector3 &v);
00975     void Translate(const vector3 &v, int conf);
00976     void Rotate(const double u[3][3]);
00977     void Rotate(const double m[9]);
00978     void Rotate(const double m[9],int nconf);
00980     void Center();
00982     bool Kekulize();
00983     bool PerceiveKekuleBonds();
00984 
00985     void NewPerceiveKekuleBonds();
00986 
00987     bool DeleteHydrogen(OBAtom*);
00988     bool DeleteHydrogens();
00989     bool DeleteHydrogens(OBAtom*);
00990     bool DeleteNonPolarHydrogens();
00991     bool AddHydrogens(bool polaronly=false,bool correctForPH=true);
00992     bool AddHydrogens(OBAtom*);
00993     bool AddPolarHydrogens();
00994 
00996     bool StripSalts();
00998     bool ConvertDativeBonds();
00999 
01000     bool CorrectForPH();
01001     bool AssignSpinMultiplicity();
01002     vector3 Center(int nconf);
01004     void SetTorsion(OBAtom*,OBAtom*,OBAtom*,OBAtom*,double);
01006 
01008 
01009 
01010     void FindSSSR();
01011     void FindRingAtomsAndBonds();
01012     void FindChiralCenters();
01013     void FindChildren(std::vector<int> &,int,int);
01014     void FindChildren(std::vector<OBAtom*>&,OBAtom*,OBAtom*);
01015     void FindLargestFragment(OBBitVec &);
01018     void ContigFragList(std::vector<std::vector<int> >&);
01020     void Align(OBAtom*,OBAtom*,vector3&,vector3&);
01022     void ConnectTheDots();
01024     void PerceiveBondOrders();
01025     void FindTorsions();
01026     // documented in mol.cpp: graph-theoretical distance for each atom
01027     bool         GetGTDVector(std::vector<int> &);
01028     // documented in mol.cpp: graph-invariant index for each atom
01029     void         GetGIVector(std::vector<unsigned int> &);
01030     // documented in mol.cpp: calculate symmetry-unique identifiers
01031     void         GetGIDVector(std::vector<unsigned int> &);
01033 
01035 
01036 
01037     bool Has2D();
01039     bool Has3D();
01041     bool HasNonZeroCoords();
01042     bool HasAromaticPerceived()     { return(HasFlag(OB_AROMATIC_MOL)); }
01043     bool HasSSSRPerceived()         { return(HasFlag(OB_SSSR_MOL));     }
01044     bool HasRingAtomsAndBondsPerceived(){return(HasFlag(OB_RINGFLAGS_MOL));}
01045     bool HasAtomTypesPerceived()    { return(HasFlag(OB_ATOMTYPES_MOL));}
01046     bool HasChiralityPerceived()    { return(HasFlag(OB_CHIRALITY_MOL));}
01047     bool HasPartialChargesPerceived() { return(HasFlag(OB_PCHARGE_MOL));}
01048     bool HasHybridizationPerceived() { return(HasFlag(OB_HYBRID_MOL));  }
01049     bool HasImplicitValencePerceived() { return(HasFlag(OB_IMPVAL_MOL));}
01050     bool HasKekulePerceived() { return(HasFlag(OB_KEKULE_MOL));         }
01051     bool HasClosureBondsPerceived() { return(HasFlag(OB_CLOSURE_MOL));  }
01052     bool HasChainsPerceived() { return(HasFlag(OB_CHAINS_MOL));         }
01053     bool HasHydrogensAdded() { return(HasFlag(OB_H_ADDED_MOL));         }
01054     bool HasAromaticCorrected() { return(HasFlag(OB_AROM_CORRECTED_MOL));}
01055     bool IsCorrectedForPH() { return(HasFlag(OB_PH_CORRECTED_MOL));     }
01056     bool HasSpinMultiplicityAssigned() { return(HasFlag(OB_TSPIN_MOL)); }
01058     bool IsChiral();
01060     bool Empty()                       { return(_natoms == 0);          }
01062 
01064 
01065     int     NumConformers()    { return((_vconf.empty())?0:_vconf.size()); }
01066     void    SetConformers(std::vector<double*> &v);
01067     void    AddConformer(double *f)    {  _vconf.push_back(f);    }
01068     void    SetConformer(int i)        {  _c = _vconf[i];         }
01069     void    CopyConformer(double*,int);
01070     void    DeleteConformer(int);
01071     double  *GetConformer(int i)       {  return(_vconf[i]);      }
01072     double  *BeginConformer(std::vector<double*>::iterator&i)
01073       { i = _vconf.begin();
01074         return((i == _vconf.end()) ? NULL:*i); }
01075     double  *NextConformer(std::vector<double*>::iterator&i)
01076       { i++;
01077         return((i == _vconf.end()) ? NULL:*i); }
01078     std::vector<double*> &GetConformers() {   return(_vconf);     }
01080 
01082 
01083 
01084     OBAtom *BeginAtom(std::vector<OBNodeBase*>::iterator &i);
01086     OBAtom *NextAtom(std::vector<OBNodeBase*>::iterator &i);
01088     OBBond *BeginBond(std::vector<OBEdgeBase*>::iterator &i);
01090     OBBond *NextBond(std::vector<OBEdgeBase*>::iterator &i);
01092     OBResidue *BeginResidue(std::vector<OBResidue*>::iterator &i)
01093     {
01094         i = _residue.begin();
01095         return((i == _residue.end()) ? NULL:*i);
01096     }
01098     OBResidue *NextResidue(std::vector<OBResidue*>::iterator &i)
01099     {
01100         i++;
01101         return((i == _residue.end()) ? NULL:*i);
01102     }
01103     OBInternalCoord *BeginInternalCoord(std::vector<OBInternalCoord*>::iterator &i)
01104     {
01105         i = _internals.begin();
01106         return((i == _internals.end()) ? NULL:*i);
01107     }
01108     OBInternalCoord *NextInternalCoord(std::vector<OBInternalCoord*>::iterator &i)
01109     {
01110         i++;
01111         return((i == _internals.end()) ? NULL:*i);
01112     }
01114 
01115     //  Removed with OBConversion framework -- see OBConversion class instead
01117 
01118     // friend std::ostream&       operator<< ( std::ostream&, OBMol& ) ;
01119     // friend std::istream&       operator>> ( std::istream&, OBMol& ) ;
01121 };
01122 
01124 class OBAPI OBInternalCoord
01125 {
01126 public:
01127     //class members
01128     OBAtom *_a,*_b,*_c;
01129     double   _dst,_ang,_tor;
01131     OBInternalCoord(OBAtom *a=(OBAtom*)NULL,
01132                     OBAtom *b=(OBAtom*)NULL,
01133                     OBAtom *c=(OBAtom*)NULL)
01134     {
01135         _a = a;
01136         _b = b;
01137         _c = c;
01138         _dst = _ang = _tor = 0.0;
01139     }
01140 };
01141 
01142 //function prototypes
01143 
01144 OBAPI bool tokenize(std::vector<std::string>&, const char *buf, const char *delimstr=" \t\n");
01145 OBAPI bool tokenize(std::vector<std::string>&, std::string&, const char *delimstr=" \t\n", int limit=-1);
01147 OBAPI void Trim(std::string& txt);
01149 OBAPI void ThrowError(char *str);
01151 OBAPI void ThrowError(std::string &str);
01152 OBAPI void CartesianToInternal(std::vector<OBInternalCoord*>&,OBMol&);
01153 OBAPI void InternalToCartesian(std::vector<OBInternalCoord*>&,OBMol&);
01154 OBAPI std::string NewExtension(std::string&,char*);
01155 // Now handled by OBConversion class
01156 // OBAPI bool SetInputType(OBMol&,std::string&);
01157 // OBAPI bool SetOutputType(OBMol&,std::string&);
01158 
01159 //global definitions
01161 EXTERN  OBElementTable   etab;
01164 EXTERN  OBTypeTable      ttab;
01166 EXTERN  OBIsotopeTable   isotab;
01168 EXTERN  OBAromaticTyper  aromtyper;
01171 EXTERN  OBAtomTyper      atomtyper;
01173 EXTERN  OBChainsParser   chainsparser;
01175 EXTERN  OBMessageHandler obErrorLog;
01177 EXTERN  OBResidueData    resdat;
01178 
01179 //Utility Macros
01180 
01181 #ifndef BUFF_SIZE
01182 #define BUFF_SIZE 32768
01183 #endif
01184 
01185 #ifndef EQ
01186 #define EQ(a,b) (!strcmp((a), (b)))
01187 #endif
01188 
01189 #ifndef EQn
01190 #define EQn(a,b,n) (!strncmp((a), (b), (n)))
01191 #endif
01192 
01193 #ifndef SQUARE
01194 #define SQUARE(x) ((x)*(x))
01195 #endif
01196 
01197 #ifndef IsUnsatType
01198 #define IsUnsatType(x)  (EQ(x,"Car") || EQ(x,"C2") || EQ(x,"Sox") || EQ(x,"Sac") || EQ(x,"Pac") || EQ(x,"So2"))
01199 #endif
01200 
01201 #ifndef __KCC
01202 extern "C"
01203 {
01204     OBAPI void  get_rmat(double*,double*,double*,int);
01205     OBAPI void  ob_make_rmat(double mat[3][3],double rmat[9]);
01206     OBAPI void  qtrfit (double *r,double *f,int size,double u[3][3]);
01207     OBAPI double superimpose(double*,double*,int);
01208 }
01209 #else
01210 OBAPI void get_rmat(double*,double*,double*,int);
01211 OBAPI void ob_make_rmat(double mat[3][3],double rmat[9]);
01212 OBAPI void qtrfit (double *r,double *f,int size,double u[3][3]);
01213 OBAPI double superimpose(double*,double*,int);
01214 #endif // __KCC
01215 
01216 } // end namespace OpenBabel
01217 
01218 #endif // OB_MOL_H
01219