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