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