ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/openbabel/mol.hpp
Revision: 1081
Committed: Thu Oct 19 20:49:05 2006 UTC (18 years, 6 months ago) by gezelter
File size: 49810 byte(s)
Log Message:
updated OpenBabel to version 2.0.2

File Contents

# Content
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-2005 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.sourceforge.net/>
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 "config.h"
26
27 #ifndef EXTERN
28 # define EXTERN extern
29 #endif
30
31 #include <math.h>
32
33 #include <algorithm>
34 #include <vector>
35 #include <string>
36 #include <map>
37
38 #if HAVE_IOSTREAM
39 #include <iostream>
40 #elif HAVE_IOSTREAM_H
41 #include <iostream.h>
42 #endif
43
44 #if HAVE_FSTREAM
45 #include <fstream>
46 #elif HAVE_FSTREAM_H
47 #include <fstream.h>
48 #endif
49
50 #include "base.hpp"
51 #include "data.hpp"
52 #include "chains.hpp"
53 #include "vector3.hpp"
54 #include "bitvec.hpp"
55 #include "ring.hpp"
56 #include "generic.hpp"
57 #include "typer.hpp"
58 #include "oberror.hpp"
59 #include "obiter.hpp"
60 #include "reaction.hpp" //so it gets notices in DLL builds
61
62 namespace OpenBabel
63 {
64
65 class OBAtom;
66 class OBBond;
67 class OBMol;
68 class OBInternalCoord;
69
70 // Class OBResidue
71 // class introduction in residue.cpp
72 class OBAPI OBResidue
73 {
74 public:
75
76 //! Constructor
77 OBResidue(void);
78 //! Copy constructor
79 //! \warning Currently does not copy all associated OBGenericData
80 //! This requires a (minor) API change, and will thus only be fixed in 2.1
81 //! or later releases.
82 OBResidue(const OBResidue &);
83 //! Destructor
84 virtual ~OBResidue(void);
85
86 OBResidue &operator=(const OBResidue &);
87
88 void AddAtom(OBAtom *atom);
89 void InsertAtom(OBAtom *atom);
90 void RemoveAtom(OBAtom *atom);
91 void Clear(void);
92
93 void SetName(const std::string &resname);
94 void SetNum(unsigned int resnum);
95 void SetChain(char chain);
96 void SetChainNum(unsigned int chainnum);
97 void SetIdx(unsigned int idx);
98
99 void SetAtomID(OBAtom *atom, const std::string &id);
100 void SetHetAtom(OBAtom *atom, bool hetatm);
101 //! Set the atomic serial number for a given atom (see OBSerialNums)
102 void SetSerialNum(OBAtom *atom, unsigned int sernum);
103
104 std::string GetName(void) const;
105 unsigned int GetNum(void) const;
106 unsigned int GetNumAtoms() const;
107 char GetChain(void) const;
108 unsigned int GetChainNum(void) const;
109 unsigned int GetIdx(void) const;
110 unsigned int GetResKey(void) const;
111
112 std::vector<OBAtom*> GetAtoms(void) const;
113 std::vector<OBBond*> GetBonds(bool = true) const;
114
115 std::string GetAtomID(OBAtom *atom) const;
116 //! \return the serial number of the supplied atom (uses OBSerialNums)
117 unsigned GetSerialNum(OBAtom *atom) const;
118
119 bool GetAminoAcidProperty(int) const;
120 bool GetAtomProperty(OBAtom *, int) const;
121 bool GetResidueProperty(int) const;
122
123 bool IsHetAtom(OBAtom *atom) const;
124 bool IsResidueType(int) const;
125
126 //! \deprecated Use FOR_ATOMS_OF_RESIDUE and OBResidueAtomIter instead
127 OBAtom *BeginAtom(std::vector<OBAtom*>::iterator &i);
128 //! \deprecated Use FOR_ATOMS_OF_RESIDUE and OBResidueAtomIter instead
129 OBAtom *NextAtom(std::vector<OBAtom*>::iterator &i);
130
131 //! \name Methods for handling generic data
132 //@{
133 bool HasData(std::string &);
134 bool HasData(const char *);
135 bool HasData(unsigned int type);
136 void DeleteData(unsigned int type);
137 void DeleteData(OBGenericData*);
138 void DeleteData(std::vector<OBGenericData*>&);
139 void SetData(OBGenericData *d)
140 { _vdata.push_back(d); }
141 //! \return the number of OBGenericData items attached to this residue
142 unsigned int DataSize()
143 { return(_vdata.size()); }
144 OBGenericData *GetData(unsigned int type);
145 OBGenericData *GetData(std::string&);
146 OBGenericData *GetData(const char *);
147 std::vector<OBGenericData*> &GetData()
148 { return(_vdata); }
149 std::vector<OBGenericData*>::iterator BeginData()
150 { return(_vdata.begin()); }
151 std::vector<OBGenericData*>::iterator EndData()
152 { return(_vdata.end()); }
153 //@}
154
155 protected: // members
156
157 unsigned int _idx; //!< Residue index (i.e., internal index in an OBMol)
158 char _chain; //!< Chain ID
159 unsigned int _aakey; //!< Amino Acid key ID -- see SetResidueKeys()
160 unsigned int _reskey;//!< Residue key ID -- see SetResidueKeys()
161 unsigned int _resnum;//!< Residue number (i.e., in file)
162 std::string _resname;//!< Residue text name
163
164 std::vector<bool> _hetatm;//!< Is a given atom a HETAM
165 std::vector<std::string> _atomid;//!< Residue atom text IDs
166 std::vector<OBAtom*> _atoms; //!< List of OBAtom in this residue
167 std::vector<unsigned int> _sernum;//!< List of serial numbers
168 std::vector<OBGenericData*> _vdata; //!< Custom data
169 }; // OBResidue
170
171
172 //ATOM Property Macros (flags)
173 //! Atom is in a 4-membered ring
174 #define OB_4RING_ATOM (1<<1)
175 //! Atom is in a 3-membered ring
176 #define OB_3RING_ATOM (1<<2)
177 //! Atom is aromatic
178 #define OB_AROMATIC_ATOM (1<<3)
179 //! Atom is in a ring
180 #define OB_RING_ATOM (1<<4)
181 //! Atom has clockwise SMILES chiral stereochemistry (i.e., "@@")
182 #define OB_CSTEREO_ATOM (1<<5)
183 //! Atom has anticlockwise SMILES chiral stereochemistry (i.e., "@")
184 #define OB_ACSTEREO_ATOM (1<<6)
185 //! Atom is an electron donor
186 #define OB_DONOR_ATOM (1<<7)
187 //! Atom is an electron acceptor
188 #define OB_ACCEPTOR_ATOM (1<<8)
189 //! Atom is chiral
190 #define OB_CHIRAL_ATOM (1<<9)
191 //! Atom has + chiral volume
192 #define OB_POS_CHIRAL_ATOM (1<<10)
193 //! Atom has - chiral volume
194 #define OB_NEG_CHIRAL_ATOM (1<<11)
195 // 12-16 currently unused
196
197 // Class OBAtom
198 // class introduction in atom.cpp
199 class OBAPI OBAtom : public OBNodeBase
200 {
201 protected:
202 char _ele; //!< atomic number (type char to minimize space -- allows for 0..255 elements)
203 char _impval; //!< implicit valence
204 char _type[6]; //!< atomic type
205 short _fcharge; //!< formal charge
206 unsigned short _isotope; //!< isotope (0 = most abundant)
207 short _spinmultiplicity;//!< atomic spin, e.g., 2 for radical 1 or 3 for carbene
208
209 //unsigned short int _idx; //!< index in parent (inherited)
210 unsigned short _cidx; //!< index into coordinate array
211 unsigned short _hyb; //!< hybridization
212 unsigned short _flags; //!< bitwise flags (e.g. aromaticity)
213 double _pcharge; //!< partial charge
214 double **_c; //!< coordinate array in double*
215 vector3 _v; //!< coordinate vector
216 OBResidue *_residue; //!< parent residue (if applicable)
217 //OBMol *_parent; //!< parent molecule (inherited)
218 //vector<OBBond*> _bond; //!< connections (inherited)
219 std::vector<OBGenericData*> _vdata; //!< custom data
220
221 int GetFlag() const { return(_flags); }
222 void SetFlag(int flag) { _flags |= flag; }
223 bool HasFlag(int flag) { return((_flags & flag) ? true : false); }
224
225 public:
226
227 //! Constructor
228 OBAtom();
229 //! Destructor
230 virtual ~OBAtom();
231 //! Assignment
232 OBAtom &operator = (OBAtom &);
233 //! Clear all data
234 void Clear();
235
236 //! \name Methods to set atomic information
237 //@{
238 //! Set atom index (i.e., in an OBMol)
239 void SetIdx(int idx) { _idx = idx; _cidx = (idx-1)*3; }
240 //! Set atom hybridization (i.e., 1 = sp, 2 = sp2, 3 = sp3 ...)
241 void SetHyb(int hyb) { _hyb = hyb; }
242 //! Set atomic number
243 void SetAtomicNum(int atomicnum) { _ele = (char)atomicnum; }
244 //! Set isotope number (actual atomic weight is tabulated automatically, 0 = most abundant)
245 void SetIsotope(unsigned int iso);
246 void SetImplicitValence(int val) { _impval = (char)val; }
247 void IncrementImplicitValence() { _impval++; }
248 void DecrementImplicitValence() { _impval--; }
249 void SetFormalCharge(int fcharge) { _fcharge = fcharge; }
250 void SetSpinMultiplicity(short spin){ _spinmultiplicity = spin; }
251 void SetType(char *type);
252 void SetType(std::string &type);
253 void SetPartialCharge(double pcharge){ _pcharge = pcharge; }
254 void SetVector(vector3 &v);
255 void SetVector(const double x,const double y,const double z);
256 //! Set the position of this atom from a pointer-driven array of coordinates
257 void SetCoordPtr(double **c) { _c = c; _cidx = (GetIdx()-1)*3; }
258 //! Set the position of this atom based on the internal pointer array (i.e. from SetCoordPtr() )
259 void SetVector();
260 void SetResidue(OBResidue *res) { _residue=res; }
261 // void SetParent(OBMol *ptr) { _parent=ptr; } // inherited
262 void SetAromatic() { SetFlag(OB_AROMATIC_ATOM); }
263 void UnsetAromatic() { _flags &= (~(OB_AROMATIC_ATOM)); }
264 //! Mark atom as having SMILES clockwise stereochemistry (i.e., "@@")
265 void SetClockwiseStereo() { SetFlag(OB_CSTEREO_ATOM|OB_CHIRAL_ATOM); }
266 //! Mark atom as having SMILES anticlockwise stereochemistry (i.e., "@")
267 void SetAntiClockwiseStereo() { SetFlag(OB_ACSTEREO_ATOM|OB_CHIRAL_ATOM); }
268 //! Mark an atom as having + chiral volume
269 void SetPositiveStereo() { SetFlag(OB_POS_CHIRAL_ATOM|OB_CHIRAL_ATOM); }
270 //! Mark an atom as having - chiral volume
271 void SetNegativeStereo() { SetFlag(OB_NEG_CHIRAL_ATOM|OB_CHIRAL_ATOM); }
272 //! Clear all stereochemistry information
273 void UnsetStereo()
274 {
275 _flags &= ~(OB_ACSTEREO_ATOM);
276 _flags &= ~(OB_CSTEREO_ATOM);
277 _flags &= ~(OB_POS_CHIRAL_ATOM);
278 _flags &= ~(OB_NEG_CHIRAL_ATOM);
279 _flags &= ~(OB_CHIRAL_ATOM);
280 }
281 //! Mark an atom as belonging to at least one ring
282 void SetInRing() { SetFlag(OB_RING_ATOM); }
283 //! Mark an atom as being chiral with unknown stereochemistry
284 void SetChiral() { SetFlag(OB_CHIRAL_ATOM); }
285 //! Clear the internal coordinate pointer
286 void ClearCoordPtr() { _c = NULL; _cidx=0; }
287 //@}
288
289 //! \name Methods to retrieve atomic information
290 //@{
291 //int GetStereo() const { return((int)_stereo);}
292 int GetFormalCharge() const { return(_fcharge); }
293 unsigned int GetAtomicNum() const { return((unsigned int)_ele); }
294 unsigned short int GetIsotope() const { return(_isotope); }
295 int GetSpinMultiplicity() const { return(_spinmultiplicity); }
296 //! The atomic mass of this atom given by standard IUPAC average molar mass
297 double GetAtomicMass() const;
298 //! The atomic mass of given by the isotope (default of 0 s most abundant isotope)
299 double GetExactMass() const;
300 unsigned int GetIdx() const { return((int)_idx); }
301 unsigned int GetCoordinateIdx() const { return((int)_cidx); }
302 //! \deprecated Use GetCoordinateIdx() instead
303 unsigned int GetCIdx() const { return((int)_cidx); }
304 //! The current number of explicit connections
305 unsigned int GetValence() const
306 {
307 return((_vbond.empty()) ? 0 : _vbond.size());
308 }
309 //! The hybridization of this atom (i.e. 1 for sp, 2 for sp2, 3 for sp3)
310 unsigned int GetHyb() const;
311 //! The implicit valence of this atom type (i.e. maximum number of connections expected)
312 unsigned int GetImplicitValence() const;
313 //! The number of non-hydrogens connected to this atom
314 unsigned int GetHvyValence() const;
315 //! The number of heteroatoms connected to an atom
316 unsigned int GetHeteroValence() const;
317 char *GetType();
318
319 //! The x coordinate
320 double GetX() { return(x()); }
321 //! The y coordinate
322 double GetY() { return(y()); }
323 //! The z coordinate
324 double GetZ() { return(z()); }
325 double x()
326 {
327 if (_c)
328 return((*_c)[_cidx]);
329 else
330 return _v.x();
331 }
332 double y()
333 {
334 if (_c)
335 return((*_c)[_cidx+1]);
336 else
337 return _v.y();
338 }
339 double z()
340 {
341 if (_c)
342 return((*_c)[_cidx+2]);
343 else
344 return _v.z();
345 }
346 //! \return the coordinates as a double*
347 double *GetCoordinate()
348 {
349 if (_c)
350 return(&(*_c)[_cidx]);
351 else
352 return NULL;
353 }
354 //! \return the coordinates as a vector3 object
355 vector3 &GetVector();
356 //! \return the partial charge of this atom, calculating a Gasteiger charge if needed
357 double GetPartialCharge();
358 OBResidue *GetResidue();
359 //OBMol *GetParent() {return((OBMol*)_parent);}
360 //! Create a vector for a new bond from this atom, with length given by the supplied parameter
361 bool GetNewBondVector(vector3 &v,double length);
362 OBBond *GetBond(OBAtom *);
363 OBAtom *GetNextAtom();
364 //@}
365
366 //! \name Iterator methods
367 //@{
368 //! \deprecated Use FOR_BONDS_OF_ATOM and OBAtomBondIter instead
369 std::vector<OBEdgeBase*>::iterator BeginBonds()
370 { return(_vbond.begin()); }
371 //! \deprecated Use FOR_BONDS_OF_ATOM and OBAtomBondIter instead
372 std::vector<OBEdgeBase*>::iterator EndBonds()
373 { return(_vbond.end()); }
374 //! \deprecated Use FOR_BONDS_OF_ATOM and OBAtomBondIter instead
375 OBBond *BeginBond(std::vector<OBEdgeBase*>::iterator &i);
376 //! \deprecated Use FOR_BONDS_OF_ATOM and OBAtomBondIter instead
377 OBBond *NextBond(std::vector<OBEdgeBase*>::iterator &i);
378 //! \deprecated Use FOR_NBORS_OF_ATOM and OBAtomAtomIter instead
379 OBAtom *BeginNbrAtom(std::vector<OBEdgeBase*>::iterator &);
380 //! \deprecated Use FOR_NBORS_OF_ATOM and OBAtomAtomIter instead
381 OBAtom *NextNbrAtom(std::vector<OBEdgeBase*>::iterator &);
382 //@}
383
384 //! \return the distance to the atom defined by OBMol::GetAtom()
385 double GetDistance(int index);
386 //! \return the distance to the supplied OBAtom
387 double GetDistance(OBAtom*);
388 //! \return the angle defined by this atom -> b (vertex) -> c
389 double GetAngle(int b, int c);
390 //! \return the angle defined by this atom -> b (vertex) -> c
391 double GetAngle(OBAtom *b, OBAtom *c);
392
393 //! \name Addition of residue/bond info. for an atom
394 //@{
395 void NewResidue()
396 {
397 if (!_residue)
398 _residue = new OBResidue;
399 }
400 void DeleteResidue()
401 {
402 if (_residue)
403 delete _residue;
404 }
405 void AddBond(OBBond *bond)
406 {
407 _vbond.push_back((OBEdgeBase*)bond);
408 }
409 void InsertBond(std::vector<OBEdgeBase*>::iterator &i, OBBond *bond)
410 {
411 _vbond.insert(i, (OBEdgeBase*)bond);
412 }
413 bool DeleteBond(OBBond*);
414 void ClearBond() {_vbond.clear();}
415 //@}
416
417 //! \name Requests for atomic property information
418 //@{
419 //! The number of oxygen atoms connected that only have one heavy valence
420 unsigned int CountFreeOxygens() const;
421 //! The number of hydrogens needed to fill the implicit valence of this atom
422 unsigned int ImplicitHydrogenCount() const;
423 //! The number of hydrogens explicitly bound to this atom currently
424 unsigned int ExplicitHydrogenCount() const;
425 //! The number of rings that contain this atom
426 unsigned int MemberOfRingCount() const;
427 //! The size of the smallest ring that contains this atom (0 if not in a ring)
428 unsigned int MemberOfRingSize() const;
429 //! The smallest angle of bonds to this atom
430 double SmallestBondAngle();
431 //! The average angle of bonds to this atom
432 double AverageBondAngle();
433 //! The sum of the bond orders of the bonds to the atom (i.e. double bond = 2...)
434 unsigned int BOSum() const;
435 //! The sum of the bond orders of bonds to the atom, considering only KDouble, KTriple bonds
436 unsigned int KBOSum() const;
437 //@}
438
439 //! \name Builder utilities
440 //@{
441 //! If this is a hydrogen atom, transform into a methyl group
442 bool HtoMethyl();
443 //! Change the hybridization of this atom and modify the geometry accordingly
444 bool SetHybAndGeom(int);
445 //@}
446
447 //! \name Property information
448 //@{
449 //! Is there any residue information?
450 bool HasResidue() { return(_residue != NULL); }
451 bool IsHydrogen() { return(GetAtomicNum() == 1); }
452 bool IsCarbon() { return(GetAtomicNum() == 6); }
453 bool IsNitrogen() { return(GetAtomicNum() == 7); }
454 bool IsOxygen() { return(GetAtomicNum() == 8); }
455 bool IsSulfur() { return(GetAtomicNum() == 16);}
456 bool IsPhosphorus() { return(GetAtomicNum() == 15);}
457 bool IsAromatic() const;
458 bool IsInRing() const;
459 bool IsInRingSize(int) const;
460 //! Is this atom an element in the 15th or 16th main groups (i.e., N, O, P, S ...) ?
461 bool IsHeteroatom();
462 //! Is this atom any element except carbon or hydrogen?
463 bool IsNotCorH();
464 //! Is this atom connected to the supplied OBAtom?
465 bool IsConnected(OBAtom*);
466 //! Is this atom related to the supplied OBAtom in a 1,3 bonding pattern?
467 bool IsOneThree(OBAtom*);
468 //! Is this atom related to the supplied OBAtom in a 1,4 bonding pattern?
469 bool IsOneFour(OBAtom*);
470 //! Is this atom an oxygen in a carboxyl (-CO2 or CO2H) group?
471 bool IsCarboxylOxygen();
472 //! Is this atom an oxygen in a phosphate (R-PO3) group?
473 bool IsPhosphateOxygen();
474 //! Is this atom an oxygen in a sulfate (-SO3) group?
475 bool IsSulfateOxygen();
476 //! Is this atom an oxygen in a nitro (-NO2) group?
477 bool IsNitroOxygen();
478 bool IsAmideNitrogen();
479 bool IsPolarHydrogen();
480 bool IsNonPolarHydrogen();
481 bool IsAromaticNOxide();
482 //! Is this atom chiral?
483 bool IsChiral();
484 bool IsAxial();
485 //! Does this atom have SMILES-specified clockwise "@@" stereochemistry?
486 bool IsClockwise() { return(HasFlag(OB_CSTEREO_ATOM)); }
487 //! Does this atom have SMILES-specified anticlockwise "@" stereochemistry?
488 bool IsAntiClockwise() { return(HasFlag(OB_ACSTEREO_ATOM)); }
489 //! Does this atom have a positive chiral volume?
490 bool IsPositiveStereo() { return(HasFlag(OB_POS_CHIRAL_ATOM)); }
491 //! Does this atom have a negative chiral volume?
492 bool IsNegativeStereo() { return(HasFlag(OB_NEG_CHIRAL_ATOM)); }
493 //! Does this atom have SMILES-specified stereochemistry?
494 bool HasChiralitySpecified()
495 { return(HasFlag(OB_CSTEREO_ATOM|OB_ACSTEREO_ATOM)); }
496 //! Does this atom have a specified chiral volume?
497 bool HasChiralVolume()
498 { return(HasFlag(OB_POS_CHIRAL_ATOM|OB_NEG_CHIRAL_ATOM)); }
499 //! Is this atom a hydrogen-bond acceptor (receptor)?
500 bool IsHbondAcceptor();
501 //! Is this atom a hydrogen-bond donor?
502 bool IsHbondDonor();
503 //! Is this a hydrogen atom attached to a hydrogen-bond donor?
504 bool IsHbondDonorH();
505 bool HasAlphaBetaUnsat(bool includePandS=true);
506 bool HasBondOfOrder(unsigned int);
507 int CountBondsOfOrder(unsigned int);
508 bool HasNonSingleBond();
509 bool HasSingleBond() { return(HasBondOfOrder(1)); }
510 bool HasDoubleBond() { return(HasBondOfOrder(2)); }
511 bool HasAromaticBond() { return(HasBondOfOrder(5)); }
512 //! Determines if this atom matches the first atom in a given SMARTS pattern
513 bool MatchesSMARTS(const char *);
514 //@}
515
516 //! \name Methods for handling generic data
517 //@{
518 bool HasData(std::string &);
519 bool HasData(const char *);
520 bool HasData(unsigned int type);
521 void DeleteData(unsigned int type);
522 void DeleteData(OBGenericData*);
523 void DeleteData(std::vector<OBGenericData*>&);
524 void SetData(OBGenericData *d)
525 { _vdata.push_back(d); }
526 //! \return the number of OBGenericData items attached to this atom
527 unsigned int DataSize()
528 { return(_vdata.size()); }
529 OBGenericData *GetData(unsigned int type);
530 OBGenericData *GetData(std::string&);
531 OBGenericData *GetData(const char *);
532 std::vector<OBGenericData*> &GetData() { return(_vdata); }
533 std::vector<OBGenericData*>::iterator BeginData()
534 { return(_vdata.begin()); }
535 std::vector<OBGenericData*>::iterator EndData()
536 { return(_vdata.end()); }
537 //@}
538 }; // class OBAtom
539
540
541 // Class OBBond
542
543 //BOND Property Macros (flags)
544 //! An aromatic bond (regardless of bond order)
545 #define OB_AROMATIC_BOND (1<<1)
546 //! A solid black wedge in 2D representations -- i.e., "up" from the 2D plane
547 #define OB_WEDGE_BOND (1<<2)
548 //! A dashed "hash" bond in 2D representations -- i.e., "down" from the 2D plane
549 #define OB_HASH_BOND (1<<3)
550 //! A bond in a ring
551 #define OB_RING_BOND (1<<4)
552 //! The "upper" bond in a double bond cis/trans isomer (i.e., "/" in SMILES)
553 #define OB_TORUP_BOND (1<<5)
554 //! The "down" bond in a double bond cis/trans isomer (i.e., "\" in SMILES)
555 #define OB_TORDOWN_BOND (1<<6)
556 //! A Kekule single bond
557 #define OB_KSINGLE_BOND (1<<7)
558 //! A Kekule double bond
559 #define OB_KDOUBLE_BOND (1<<8)
560 //! A Kekule triple bond
561 #define OB_KTRIPLE_BOND (1<<9)
562 #define OB_CLOSURE_BOND (1<<10)
563 // 11-16 currently unused
564
565 // class introduction in bond.cpp
566 class OBAPI OBBond : public OBEdgeBase
567 {
568 protected:
569 char _order; //!< Bond order (1, 2, 3, 5=aromatic)
570 unsigned short int _flags; //!< Any flags for this bond
571 //OBAtom *_bgn; //!< Not needed, inherited from OBEdgeBase
572 //OBAtom *_end; //!< Not needed, inherited from OBEdgeBase
573 //OBMol *_parent;//!< Not needed, inherited from OBEdgeBase
574 //unsigned short int _idx; //!< Not needed, inherited from OBEdgeBase
575 std::vector<OBGenericData*> _vdata; //!< Generic data for custom information
576
577 bool HasFlag(int flag) { return((_flags & flag) != 0); }
578 void SetFlag(int flag) { _flags |= flag; }
579
580 public:
581 //! Constructor
582 OBBond();
583 //! Destructor
584 virtual ~OBBond();
585
586 //! \name Bond modification methods
587 //@{
588 void SetIdx(int idx)
589 {
590 _idx = idx;
591 }
592 void SetBO(int order);
593 void SetBegin(OBAtom *begin)
594 {
595 _bgn = begin;
596 }
597 void SetEnd(OBAtom *end)
598 {
599 _end = end;
600 }
601 // void SetParent(OBMol *ptr) {_parent=ptr;} // (inherited)
602 void SetLength(OBAtom*,double);
603 void Set(int,OBAtom*,OBAtom*,int,int);
604 void SetKSingle();
605 void SetKDouble();
606 void SetKTriple();
607 void SetAromatic() { SetFlag(OB_AROMATIC_BOND); }
608 void SetHash() { SetFlag(OB_HASH_BOND); }
609 void SetWedge() { SetFlag(OB_WEDGE_BOND); }
610 void SetUp() { SetFlag(OB_TORUP_BOND); }
611 void SetDown() { SetFlag(OB_TORDOWN_BOND); }
612 void SetInRing() { SetFlag(OB_RING_BOND); }
613 void SetClosure() { SetFlag(OB_CLOSURE_BOND); }
614
615 void UnsetAromatic() { _flags &= (~(OB_AROMATIC_BOND)); }
616 void UnsetKekule()
617 {
618 _flags &= (~(OB_KSINGLE_BOND|OB_KDOUBLE_BOND|OB_KTRIPLE_BOND));
619 }
620 //@}
621
622 //! \name bond data request methods
623 //@{
624 unsigned int GetBO() const { return((int)_order); }
625 unsigned int GetBondOrder() const { return((int)_order); }
626 unsigned int GetFlags() const { return(_flags); }
627 unsigned int GetBeginAtomIdx() const { return(_bgn->GetIdx()); }
628 unsigned int GetEndAtomIdx() const { return(_end->GetIdx()); }
629 OBAtom *GetBeginAtom() { return((OBAtom*)_bgn); }
630 OBAtom *GetEndAtom() { return((OBAtom*)_end); }
631 OBAtom *GetNbrAtom(OBAtom *ptr)
632 {
633 return((ptr != _bgn)? (OBAtom*)_bgn : (OBAtom*)_end);
634 }
635 // OBMol *GetParent() {return(_parent);} // (inherited)
636 double GetEquibLength();
637 double GetLength();
638 int GetNbrAtomIdx(OBAtom *ptr)
639 {
640 return((ptr!=_bgn)?_bgn->GetIdx():_end->GetIdx());
641 }
642 //@}
643
644 //! \name property request methods
645 //@{
646 bool IsAromatic() const;
647 bool IsInRing() const;
648 //! Is the bond a rotatable bond?
649 //! Currently, this function classifies any bond with at least one heavy
650 //! atom, no sp-hybrid atoms (e.g., a triple bond somewhere) not in a ring
651 //! as a potential rotor. No other bond typing is attempted.
652 bool IsRotor();
653 bool IsAmide();
654 bool IsPrimaryAmide();
655 bool IsSecondaryAmide();
656 bool IsEster();
657 bool IsCarbonyl();
658 bool IsSingle();
659 bool IsDouble();
660 bool IsTriple();
661 bool IsKSingle();
662 bool IsKDouble();
663 bool IsKTriple();
664 bool IsClosure();
665 //! \return whether this is the "upper" bond in a double bond cis/trans
666 //! isomer (i.e., "/" in SMILES)
667 bool IsUp() { return(HasFlag(OB_TORUP_BOND)); }
668 //! \return whether this is the "lower" bond in a double bond cis/trans
669 //! isomer (i.e., "\" in SMILES)
670 bool IsDown() { return(HasFlag(OB_TORDOWN_BOND)); }
671 bool IsWedge() { return(HasFlag(OB_WEDGE_BOND)); }
672 bool IsHash() { return(HasFlag(OB_HASH_BOND)); }
673 //! \return whether the geometry around this bond looks unsaturated
674 bool IsDoubleBondGeometry();
675 //@}
676
677 //! \name Methods for handling generic data
678 //@{
679 bool HasData(std::string &);
680 bool HasData(const char *);
681 bool HasData(unsigned int type);
682 void DeleteData(unsigned int type);
683 void DeleteData(OBGenericData*);
684 void DeleteData(std::vector<OBGenericData*>&);
685 void SetData(OBGenericData *d)
686 {
687 _vdata.push_back(d);
688 }
689 //! \return the number of OBGenericData items attached to this bond
690 unsigned int DataSize()
691 {
692 return(_vdata.size());
693 }
694 OBGenericData *GetData(unsigned int type);
695 OBGenericData *GetData(std::string&);
696 OBGenericData *GetData(const char *);
697 std::vector<OBGenericData*> &GetData()
698 {
699 return(_vdata);
700 }
701 std::vector<OBGenericData*>::iterator BeginData()
702 {
703 return(_vdata.begin());
704 }
705 std::vector<OBGenericData*>::iterator EndData()
706 {
707 return(_vdata.end());
708 }
709 //@}
710 }
711 ; // class OBBond
712
713
714 // Class OBMol
715
716 //MOL Property Macros (flags) -- 32+ bits
717 #define OB_SSSR_MOL (1<<1)
718 #define OB_RINGFLAGS_MOL (1<<2)
719 #define OB_AROMATIC_MOL (1<<3)
720 #define OB_ATOMTYPES_MOL (1<<4)
721 #define OB_CHIRALITY_MOL (1<<5)
722 #define OB_PCHARGE_MOL (1<<6)
723 #define OB_HYBRID_MOL (1<<8)
724 #define OB_IMPVAL_MOL (1<<9)
725 #define OB_KEKULE_MOL (1<<10)
726 #define OB_CLOSURE_MOL (1<<11)
727 #define OB_H_ADDED_MOL (1<<12)
728 #define OB_PH_CORRECTED_MOL (1<<13)
729 #define OB_AROM_CORRECTED_MOL (1<<14)
730 #define OB_CHAINS_MOL (1<<15)
731 #define OB_TCHARGE_MOL (1<<16)
732 #define OB_TSPIN_MOL (1<<17)
733 // flags 18-32 unspecified
734 #define OB_CURRENT_CONFORMER -1
735
736 // class introduction in mol.cpp
737 class OBAPI OBMol : public OBGraphBase
738 {
739 protected:
740 int _flags; //!< bitfield of flags
741 bool _autoPartialCharge; //!< Assign partial charges automatically
742 bool _autoFormalCharge; //!< Assign formal charges automatically
743 std::string _title; //!< Molecule title
744 //vector<OBAtom*> _atom; //!< not needed (inherited)
745 //vector<OBBond*> _bond; //!< not needed (inherited)
746 unsigned short int _dimension; //!< Dimensionality of coordinates
747 double _energy; //!< Molecular heat of formation (if applicable)
748 int _totalCharge; //!< Total charge on the molecule
749 unsigned int _totalSpin; //!< Total spin on the molecule (if not specified, assumes lowest possible spin)
750 double *_c; //!< coordinate array
751 std::vector<double*> _vconf; //!< vector of conformers
752 unsigned short int _natoms; //!< Number of atoms
753 unsigned short int _nbonds; //!< Number of bonds
754 std::vector<OBResidue*> _residue; //!< Residue information (if applicable)
755 std::vector<OBInternalCoord*> _internals; //!< Internal Coordinates (if applicable)
756 std::vector<OBGenericData*> _vdata; //!< Custom data -- see OBGenericData class for more
757 unsigned short int _mod; //!< Number of nested calls to BeginModify()
758
759 bool HasFlag(int flag) { return((_flags & flag) ? true : false); }
760 void SetFlag(int flag) { _flags |= flag; }
761
762 //! \name Internal Kekulization routines -- see kekulize.cpp and NewPerceiveKekuleBonds()
763 //@{
764 void start_kekulize(std::vector <OBAtom*> &cycle, std::vector<int> &electron);
765 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);
766 int getorden(OBAtom *atom);
767 void expandcycle(OBAtom *atom, OBBitVec &avisit);
768 //@}
769
770 public:
771
772 //! \name Initialization and data (re)size methods
773 //@{
774 //! Constructor
775 OBMol();
776 //! Copy constructor
777 //! \warning Currently does not copy all associated OBGenericData
778 //! This requires a (minor) API change, and will thus only be fixed in 2.1
779 //! or later releases.
780 OBMol(const OBMol &);
781 //! Destructor
782 virtual ~OBMol();
783 //! Assignment
784 OBMol &operator=(const OBMol &mol);
785 OBMol &operator+=(const OBMol &mol);
786 void ReserveAtoms(int natoms)
787 {
788 if (natoms && _mod)
789 _vatom.reserve(natoms);
790 }
791 virtual OBAtom *CreateAtom(void);
792 virtual OBBond *CreateBond(void);
793 virtual void DestroyAtom(OBNodeBase*);
794 virtual void DestroyBond(OBEdgeBase*);
795 bool AddAtom(OBAtom&);
796 bool AddBond(int,int,int,int flags=0,int insertpos=-1);
797 bool AddBond(OBBond&);
798 bool AddResidue(OBResidue&);
799 bool InsertAtom(OBAtom &);
800 bool DeleteAtom(OBAtom*);
801 bool DeleteBond(OBBond*);
802 bool DeleteResidue(OBResidue*);
803 OBAtom *NewAtom();
804 OBResidue *NewResidue();
805 //@}
806
807 //! \name Molecule modification methods
808 //@{
809 //! Call when making many modifications -- clears conformer/rotomer data.
810 virtual void BeginModify(void);
811 //! Call when done with modificaions -- re-perceive data as needed.
812 virtual void EndModify(bool nukePerceivedData=true);
813 int GetMod()
814 {
815 return(_mod);
816 }
817 void IncrementMod()
818 {
819 _mod++;
820 }
821 void DecrementMod()
822 {
823 _mod--;
824 }
825 //@}
826
827 //! \name Generic data handling methods (via OBGenericData)
828 //@{
829 //! \returns whether the generic attribute/value pair exists
830 bool HasData(std::string &);
831 //! \returns whether the generic attribute/value pair exists
832 bool HasData(const char *);
833 //! \returns whether the generic attribute/value pair exists
834 bool HasData(unsigned int type);
835 void DeleteData(unsigned int type);
836 void DeleteData(OBGenericData*);
837 void DeleteData(std::vector<OBGenericData*>&);
838 void SetData(OBGenericData *d)
839 {
840 _vdata.push_back(d);
841 }
842 //! \return the number of OBGenericData items attached to this molecule.
843 unsigned int DataSize(){ return(_vdata.size()); }
844 OBGenericData *GetData(unsigned int type);
845 OBGenericData *GetData(std::string&);
846 OBGenericData *GetData(const char *);
847 std::vector<OBGenericData*> &GetData() { return(_vdata); }
848 std::vector<OBGenericData*>::iterator BeginData()
849 {
850 return(_vdata.begin());
851 }
852 std::vector<OBGenericData*>::iterator EndData()
853 {
854 return(_vdata.end());
855 }
856 //@}
857
858 //! \name Data retrieval methods
859 //@{
860 int GetFlags() { return(_flags); }
861 //! \return the title of this molecule (often the filename)
862 const char *GetTitle() const { return(_title.c_str()); }
863 //! \return the number of atoms (i.e. OBAtom children)
864 unsigned int NumAtoms() const { return(_natoms); }
865 //! \return the number of bonds (i.e. OBBond children)
866 unsigned int NumBonds() const { return(_nbonds); }
867 //! \return the number of non-hydrogen atoms
868 unsigned int NumHvyAtoms();
869 //! \return the number of residues (i.e. OBResidue substituents)
870 unsigned int NumResidues() const { return(_residue.size()); }
871 //! \return the number of rotatble bonds. See OBBond::IsRotor() for details
872 unsigned int NumRotors();
873
874 OBAtom *GetAtom(int);
875 OBAtom *GetFirstAtom();
876 OBBond *GetBond(int);
877 OBBond *GetBond(int, int);
878 OBBond *GetBond(OBAtom*,OBAtom*);
879 OBResidue *GetResidue(int);
880 std::vector<OBInternalCoord*> GetInternalCoord();
881 //! \return the dihedral angle between the four atoms supplied a1-a2-a3-a4)
882 double GetTorsion(int,int,int,int);
883 //! \return the dihedral angle between the four atoms supplied a1-a2-a3-a4)
884 double GetTorsion(OBAtom*,OBAtom*,OBAtom*,OBAtom*);
885 //! \return the stochoimetric formula (e.g., C4H6O)
886 std::string GetFormula();
887 //! \return the heat of formation for this molecule (in kcal/mol)
888 double GetEnergy() const { return(_energy); }
889 //! \return the standard molar mass given by IUPAC atomic masses (amu)
890 double GetMolWt();
891 //! \return the mass given by isotopes (or most abundant isotope, if not specified)
892 double GetExactMass();
893 //! \return the total charge on this molecule (i.e., 0 = neutral, +1, -1...)
894 int GetTotalCharge();
895 //! \return the total spin on this molecule (i.e., 1 = singlet, 2 = doublet...)
896 unsigned int GetTotalSpinMultiplicity();
897 //! \return the dimensionality of coordinates (i.e., 0 = unknown or no coord, 2=2D, 3=3D)
898 unsigned short int GetDimension() const { return _dimension; }
899 double *GetCoordinates() { return(_c); }
900 //! \return the Smallest Set of Smallest Rings has been run (see OBRing class
901 std::vector<OBRing*> &GetSSSR();
902 //! Get the current flag for whether formal charges are set with pH correction
903 bool AutomaticFormalCharge() { return(_autoFormalCharge); }
904 //! Get the current flag for whether partial charges are auto-determined
905 bool AutomaticPartialCharge() { return(_autoPartialCharge); }
906 //@}
907
908
909 //! \name Data modification methods
910 //@{
911 void SetTitle(const char *title);
912 void SetTitle(std::string &title);
913 //! Set the stochiometric formula for this molecule
914 void SetFormula(std::string molFormula);
915 //! Set the heat of formation for this molecule (in kcal/mol)
916 void SetEnergy(double energy) { _energy = energy; }
917 //! Set the dimension of this molecule (i.e., 0, 1 , 2, 3)
918 void SetDimension(unsigned short int d) { _dimension = d; }
919 void SetTotalCharge(int charge);
920 void SetTotalSpinMultiplicity(unsigned int spin);
921 void SetInternalCoord(std::vector<OBInternalCoord*> int_coord)
922 { _internals = int_coord; }
923 //! Set the flag for determining automatic formal charges with pH (default=true)
924 void SetAutomaticFormalCharge(bool val)
925 { _autoFormalCharge=val; }
926 //! Set the flag for determining partial charges automatically (default=true)
927 void SetAutomaticPartialCharge(bool val)
928 { _autoPartialCharge=val; }
929
930 //! Mark that aromaticity has been perceived for this molecule (see OBAromaticTyper)
931 void SetAromaticPerceived() { SetFlag(OB_AROMATIC_MOL); }
932 //! Mark that Smallest Set of Smallest Rings has been run (see OBRing class)
933 void SetSSSRPerceived() { SetFlag(OB_SSSR_MOL); }
934 //! Mark that rings have been perceived (see OBRing class for details)
935 void SetRingAtomsAndBondsPerceived(){SetFlag(OB_RINGFLAGS_MOL);}
936 //! Mark that atom types have been perceived (see OBAtomTyper for details)
937 void SetAtomTypesPerceived() { SetFlag(OB_ATOMTYPES_MOL); }
938 //! Mark that chains and residues have been perceived (see OBChainsParser)
939 void SetChainsPerceived() { SetFlag(OB_CHAINS_MOL); }
940 //! Mark that chirality has been perceived
941 void SetChiralityPerceived() { SetFlag(OB_CHIRALITY_MOL); }
942 //! Mark that partial charges have been assigned
943 void SetPartialChargesPerceived(){ SetFlag(OB_PCHARGE_MOL); }
944 void SetHybridizationPerceived() { SetFlag(OB_HYBRID_MOL); }
945 void SetImplicitValencePerceived(){ SetFlag(OB_IMPVAL_MOL); }
946 void SetKekulePerceived() { SetFlag(OB_KEKULE_MOL); }
947 void SetClosureBondsPerceived(){ SetFlag(OB_CLOSURE_MOL); }
948 void SetHydrogensAdded() { SetFlag(OB_H_ADDED_MOL); }
949 void SetCorrectedForPH() { SetFlag(OB_PH_CORRECTED_MOL);}
950 void SetAromaticCorrected() { SetFlag(OB_AROM_CORRECTED_MOL);}
951 void SetSpinMultiplicityAssigned(){ SetFlag(OB_TSPIN_MOL); }
952 void SetFlags(int flags) { _flags = flags; }
953
954 void UnsetAromaticPerceived() { _flags &= (~(OB_AROMATIC_MOL)); }
955 void UnsetPartialChargesPerceived(){ _flags &= (~(OB_PCHARGE_MOL));}
956 void UnsetImplicitValencePerceived(){_flags &= (~(OB_IMPVAL_MOL)); }
957 void UnsetFlag(int flag) { _flags &= (~(flag)); }
958
959 //! \name Molecule modification methods
960 //@{
961 // Description in transform.cpp
962 virtual OBBase* DoTransformations(const std::map<std::string,std::string>* pOptions);
963 static const char* ClassDescription();
964 //! Clear all information from a molecule
965 bool Clear();
966 //! Renumber the atoms of this molecule according to the order in the supplied vector
967 void RenumberAtoms(std::vector<OBNodeBase*>&);
968 //! Translate one conformer and rotate by a rotation matrix (which is returned) to the inertial frame-of-reference
969 void ToInertialFrame(int conf, double *rmat);
970 //! Translate all conformers to the inertial frame-of-reference
971 void ToInertialFrame();
972 //! Translates all conformers in the molecule by the supplied vector
973 void Translate(const vector3 &v);
974 //! Translates one conformer in the molecule by the supplied vector
975 void Translate(const vector3 &v, int conf);
976 void Rotate(const double u[3][3]);
977 void Rotate(const double m[9]);
978 void Rotate(const double m[9],int nconf);
979 //! Translate to the center of all coordinates (for this conformer)
980 void Center();
981 //! Transform to standard Kekule bond structure (presumably from an aromatic form)
982 bool Kekulize();
983 bool PerceiveKekuleBonds();
984
985 void NewPerceiveKekuleBonds();
986
987 bool DeleteHydrogen(OBAtom*);
988 bool DeleteHydrogens();
989 bool DeleteHydrogens(OBAtom*);
990 bool DeleteNonPolarHydrogens();
991 bool AddHydrogens(bool polaronly=false,bool correctForPH=true);
992 bool AddHydrogens(OBAtom*);
993 bool AddPolarHydrogens();
994
995 //! Deletes all atoms except for the largest contiguous fragment
996 bool StripSalts();
997 //! Converts the charged form of coordinate bonds, e.g.[N+]([O-])=O to N(=O)=O
998 bool ConvertDativeBonds();
999
1000 bool CorrectForPH();
1001 bool AssignSpinMultiplicity();
1002 vector3 Center(int nconf);
1003 //! Set the torsion defined by these atoms, rotating bonded neighbors
1004 void SetTorsion(OBAtom*,OBAtom*,OBAtom*,OBAtom*,double);
1005 //@}
1006
1007 //! \name Molecule utilities and perception methods
1008 //@{
1009 //! Find Smallest Set of Smallest Rings (see OBRing class for more details)
1010 void FindSSSR();
1011 void FindRingAtomsAndBonds();
1012 void FindChiralCenters();
1013 void FindChildren(std::vector<int> &,int,int);
1014 void FindChildren(std::vector<OBAtom*>&,OBAtom*,OBAtom*);
1015 void FindLargestFragment(OBBitVec &);
1016 //! Sort a list of contig fragments by size from largest to smallest
1017 //! Each vector<int> contains the atom numbers of a contig fragment
1018 void ContigFragList(std::vector<std::vector<int> >&);
1019 //! Aligns atom a on p1 and atom b along p1->p2 vector
1020 void Align(OBAtom*,OBAtom*,vector3&,vector3&);
1021 //! Adds single bonds based on atom proximity
1022 void ConnectTheDots();
1023 //! Attempts to perceive multiple bonds based on geometries
1024 void PerceiveBondOrders();
1025 void FindTorsions();
1026 // documented in mol.cpp: graph-theoretical distance for each atom
1027 bool GetGTDVector(std::vector<int> &);
1028 // documented in mol.cpp: graph-invariant index for each atom
1029 void GetGIVector(std::vector<unsigned int> &);
1030 // documented in mol.cpp: calculate symmetry-unique identifiers
1031 void GetGIDVector(std::vector<unsigned int> &);
1032 //@}
1033
1034 //! \name Methods to check for existence of properties
1035 //@{
1036 //! Are there non-zero coordinates in two dimensions (i.e. X and Y)?
1037 bool Has2D();
1038 //! Are there non-zero coordinates in all three dimensions (i.e. X, Y, Z)?
1039 bool Has3D();
1040 //! Are there any non-zero coordinates?
1041 bool HasNonZeroCoords();
1042 bool HasAromaticPerceived() { return(HasFlag(OB_AROMATIC_MOL)); }
1043 bool HasSSSRPerceived() { return(HasFlag(OB_SSSR_MOL)); }
1044 bool HasRingAtomsAndBondsPerceived(){return(HasFlag(OB_RINGFLAGS_MOL));}
1045 bool HasAtomTypesPerceived() { return(HasFlag(OB_ATOMTYPES_MOL));}
1046 bool HasChiralityPerceived() { return(HasFlag(OB_CHIRALITY_MOL));}
1047 bool HasPartialChargesPerceived() { return(HasFlag(OB_PCHARGE_MOL));}
1048 bool HasHybridizationPerceived() { return(HasFlag(OB_HYBRID_MOL)); }
1049 bool HasImplicitValencePerceived() { return(HasFlag(OB_IMPVAL_MOL));}
1050 bool HasKekulePerceived() { return(HasFlag(OB_KEKULE_MOL)); }
1051 bool HasClosureBondsPerceived() { return(HasFlag(OB_CLOSURE_MOL)); }
1052 bool HasChainsPerceived() { return(HasFlag(OB_CHAINS_MOL)); }
1053 bool HasHydrogensAdded() { return(HasFlag(OB_H_ADDED_MOL)); }
1054 bool HasAromaticCorrected() { return(HasFlag(OB_AROM_CORRECTED_MOL));}
1055 bool IsCorrectedForPH() { return(HasFlag(OB_PH_CORRECTED_MOL)); }
1056 bool HasSpinMultiplicityAssigned() { return(HasFlag(OB_TSPIN_MOL)); }
1057 //! Is this molecule chiral?
1058 bool IsChiral();
1059 //! Are there any atoms in this molecule?
1060 bool Empty() { return(_natoms == 0); }
1061 //@}
1062
1063 //! \name Multiple conformer member functions
1064 //@{
1065 int NumConformers() { return((_vconf.empty())?0:_vconf.size()); }
1066 void SetConformers(std::vector<double*> &v);
1067 void AddConformer(double *f) { _vconf.push_back(f); }
1068 void SetConformer(int i) { _c = _vconf[i]; }
1069 void CopyConformer(double*,int);
1070 void DeleteConformer(int);
1071 double *GetConformer(int i) { return(_vconf[i]); }
1072 double *BeginConformer(std::vector<double*>::iterator&i)
1073 { i = _vconf.begin();
1074 return((i == _vconf.end()) ? NULL:*i); }
1075 double *NextConformer(std::vector<double*>::iterator&i)
1076 { i++;
1077 return((i == _vconf.end()) ? NULL:*i); }
1078 std::vector<double*> &GetConformers() { return(_vconf); }
1079 //@}
1080
1081 //! \name Iterator methods
1082 //@{
1083 //! \deprecated Use FOR_ATOMS_OF_MOL and OBMolAtomIter instead
1084 OBAtom *BeginAtom(std::vector<OBNodeBase*>::iterator &i);
1085 //! \deprecated Use FOR_ATOMS_OF_MOL and OBMolAtomIter instead
1086 OBAtom *NextAtom(std::vector<OBNodeBase*>::iterator &i);
1087 //! \deprecated Use FOR_BONDS_OF_MOL and OBMolBondIter instead
1088 OBBond *BeginBond(std::vector<OBEdgeBase*>::iterator &i);
1089 //! \deprecated Use FOR_BONDS_OF_MOL and OBMolBondIter instead
1090 OBBond *NextBond(std::vector<OBEdgeBase*>::iterator &i);
1091 //! \deprecated Use FOR_RESIDUES_OF_MOL and OBResidueIter instead
1092 OBResidue *BeginResidue(std::vector<OBResidue*>::iterator &i)
1093 {
1094 i = _residue.begin();
1095 return((i == _residue.end()) ? NULL:*i);
1096 }
1097 //! \deprecated Use FOR_RESIDUES_OF_MOL and OBResidueIter instead
1098 OBResidue *NextResidue(std::vector<OBResidue*>::iterator &i)
1099 {
1100 i++;
1101 return((i == _residue.end()) ? NULL:*i);
1102 }
1103 OBInternalCoord *BeginInternalCoord(std::vector<OBInternalCoord*>::iterator &i)
1104 {
1105 i = _internals.begin();
1106 return((i == _internals.end()) ? NULL:*i);
1107 }
1108 OBInternalCoord *NextInternalCoord(std::vector<OBInternalCoord*>::iterator &i)
1109 {
1110 i++;
1111 return((i == _internals.end()) ? NULL:*i);
1112 }
1113 //@}
1114
1115 // Removed with OBConversion framework -- see OBConversion class instead
1116 //! \name Convenience functions for I/O
1117 //@{
1118 // friend std::ostream& operator<< ( std::ostream&, OBMol& ) ;
1119 // friend std::istream& operator>> ( std::istream&, OBMol& ) ;
1120 //@}
1121 };
1122
1123 //! \brief Used to transform from z-matrix to cartesian coordinates.
1124 class OBAPI OBInternalCoord
1125 {
1126 public:
1127 //class members
1128 OBAtom *_a,*_b,*_c;
1129 double _dst,_ang,_tor;
1130 //! Constructor
1131 OBInternalCoord(OBAtom *a=(OBAtom*)NULL,
1132 OBAtom *b=(OBAtom*)NULL,
1133 OBAtom *c=(OBAtom*)NULL)
1134 {
1135 _a = a;
1136 _b = b;
1137 _c = c;
1138 _dst = _ang = _tor = 0.0;
1139 }
1140 };
1141
1142 //function prototypes
1143
1144 OBAPI bool tokenize(std::vector<std::string>&, const char *buf, const char *delimstr=" \t\n");
1145 OBAPI bool tokenize(std::vector<std::string>&, std::string&, const char *delimstr=" \t\n", int limit=-1);
1146 //! remove leading and trailing whitespace from a string
1147 OBAPI void Trim(std::string& txt);
1148 //! \deprecated -- use OBMessageHandler class instead
1149 OBAPI void ThrowError(char *str);
1150 //! \deprecated -- use OBMessageHandler class instead
1151 OBAPI void ThrowError(std::string &str);
1152 OBAPI void CartesianToInternal(std::vector<OBInternalCoord*>&,OBMol&);
1153 OBAPI void InternalToCartesian(std::vector<OBInternalCoord*>&,OBMol&);
1154 OBAPI std::string NewExtension(std::string&,char*);
1155 // Now handled by OBConversion class
1156 // OBAPI bool SetInputType(OBMol&,std::string&);
1157 // OBAPI bool SetOutputType(OBMol&,std::string&);
1158
1159 //global definitions
1160 //! Global OBElementTable for element properties
1161 EXTERN OBElementTable etab;
1162 //! Global OBTypeTable for translating between different atom types
1163 //! (e.g., Sybyl <-> MM2)
1164 EXTERN OBTypeTable ttab;
1165 //! Global OBIsotopeTable for isotope properties
1166 EXTERN OBIsotopeTable isotab;
1167 //! Global OBAromaticTyper for detecting aromatic atoms and bonds
1168 EXTERN OBAromaticTyper aromtyper;
1169 //! Global OBAtomTyper for marking internal valence, hybridization,
1170 //! and atom types (for internal and external use)
1171 EXTERN OBAtomTyper atomtyper;
1172 //! Global OBChainsParser for detecting macromolecular chains and residues
1173 EXTERN OBChainsParser chainsparser;
1174 //! Global OBMessageHandler error handler
1175 EXTERN OBMessageHandler obErrorLog;
1176 //! Global OBResidueData biomolecule residue database
1177 EXTERN OBResidueData resdat;
1178
1179 //Utility Macros
1180
1181 #ifndef BUFF_SIZE
1182 #define BUFF_SIZE 32768
1183 #endif
1184
1185 #ifndef EQ
1186 #define EQ(a,b) (!strcmp((a), (b)))
1187 #endif
1188
1189 #ifndef EQn
1190 #define EQn(a,b,n) (!strncmp((a), (b), (n)))
1191 #endif
1192
1193 #ifndef SQUARE
1194 #define SQUARE(x) ((x)*(x))
1195 #endif
1196
1197 #ifndef IsUnsatType
1198 #define IsUnsatType(x) (EQ(x,"Car") || EQ(x,"C2") || EQ(x,"Sox") || EQ(x,"Sac") || EQ(x,"Pac") || EQ(x,"So2"))
1199 #endif
1200
1201 #ifndef __KCC
1202 extern "C"
1203 {
1204 OBAPI void get_rmat(double*,double*,double*,int);
1205 OBAPI void ob_make_rmat(double mat[3][3],double rmat[9]);
1206 OBAPI void qtrfit (double *r,double *f,int size,double u[3][3]);
1207 OBAPI double superimpose(double*,double*,int);
1208 }
1209 #else
1210 OBAPI void get_rmat(double*,double*,double*,int);
1211 OBAPI void ob_make_rmat(double mat[3][3],double rmat[9]);
1212 OBAPI void qtrfit (double *r,double *f,int size,double u[3][3]);
1213 OBAPI double superimpose(double*,double*,int);
1214 #endif // __KCC
1215
1216 } // end namespace OpenBabel
1217
1218 #endif // OB_MOL_H
1219
1220 //! \file
1221 //! \brief Handle molecules. Declarations of OBMol, OBAtom, OBBond, OBResidue.
1222 //! (the main header for Open Babel)