forcefield.h
Go to the documentation of this file.
00001 /**********************************************************************
00002 forcefield.h - Handle OBForceField class.
00003 
00004 Copyright (C) 2006-2007 by Tim Vandermeersch <tim.vandermeersch@gmail.com>
00005 
00006 This file is part of the Open Babel project.
00007 For more information, see <http://openbabel.org/>
00008 
00009 This program is free software; you can redistribute it and/or modify
00010 it under the terms of the GNU General Public License as published by
00011 the Free Software Foundation version 2 of the License.
00012 
00013 This program is distributed in the hope that it will be useful,
00014 but WITHOUT ANY WARRANTY; without even the implied warranty of
00015 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00016 GNU General Public License for more details.
00017 ***********************************************************************/
00018 
00019 #ifndef OB_FORCEFIELD_H
00020 #define OB_FORCEFIELD_H
00021 
00022 #include <vector>
00023 #include <string>
00024 #include <map>
00025 
00026 #include <list>
00027 #include <set>
00028 #include <openbabel/babelconfig.h>
00029 #include <openbabel/base.h>
00030 #include <openbabel/mol.h>
00031 #include <openbabel/plugin.h>
00032 #include <openbabel/grid.h>
00033 #include <openbabel/griddata.h>
00034 #include <float.h>
00035 
00036 namespace OpenBabel
00037 {
00038   // log levels
00039 #define OBFF_LOGLVL_NONE        0   //!< no output
00040 #define OBFF_LOGLVL_LOW         1   //!< SteepestDescent progress... (no output from Energy())
00041 #define OBFF_LOGLVL_MEDIUM      2   //!< individual energy terms
00042 #define OBFF_LOGLVL_HIGH        3   //!< individual calculations and parameters
00043 
00044   // terms
00045 #define OBFF_ENERGY             (1 << 0)   //!< all terms
00046 #define OBFF_EBOND              (1 << 1)   //!< bond term
00047 #define OBFF_EANGLE             (1 << 2)   //!< angle term
00048 #define OBFF_ESTRBND            (1 << 3)   //!< strbnd term
00049 #define OBFF_ETORSION           (1 << 4)   //!< torsion term
00050 #define OBFF_EOOP               (1 << 5)   //!< oop term
00051 #define OBFF_EVDW               (1 << 6)   //!< vdw term
00052 #define OBFF_EELECTROSTATIC     (1 << 7)   //!< electrostatic term
00053 
00054   // constraint types
00055 #define OBFF_CONST_IGNORE       (1 << 0)   //!< ignore the atom while setting up calculations
00056 #define OBFF_CONST_ATOM         (1 << 1)   //!< fix the atom position
00057 #define OBFF_CONST_ATOM_X       (1 << 2)   //!< fix the x coordinate of the atom position
00058 #define OBFF_CONST_ATOM_Y       (1 << 3)   //!< fix the y coordinate of the atom position
00059 #define OBFF_CONST_ATOM_Z       (1 << 4)   //!< fix the z coordinate of the atom position
00060 #define OBFF_CONST_DISTANCE     (1 << 5)   //!< constrain distance length
00061 #define OBFF_CONST_ANGLE        (1 << 6)   //!< constrain angle
00062 #define OBFF_CONST_TORSION      (1 << 7)   //!< constrain torsion
00063 #define OBFF_CONST_CHIRAL       (1 << 8)   //!< constrain chiral volume
00064 
00065   // mode arguments for SteepestDescent, ConjugateGradients, ...
00066 #define OBFF_NUMERICAL_GRADIENT         (1 << 0)  //!< use numerical gradients
00067 #define OBFF_ANALYTICAL_GRADIENT        (1 << 1)  //!< use analytical gradients
00068 
00069 #define KCAL_TO_KJ      4.1868
00070 
00071   // inline if statements for logging.
00072 #define IF_OBFF_LOGLVL_LOW    if(_loglvl >= OBFF_LOGLVL_LOW)
00073 #define IF_OBFF_LOGLVL_MEDIUM if(_loglvl >= OBFF_LOGLVL_MEDIUM)
00074 #define IF_OBFF_LOGLVL_HIGH   if(_loglvl >= OBFF_LOGLVL_HIGH)
00075 
00077   struct LineSearchType
00078   {
00079     enum {
00080       Simple, Newton2Num
00081     };
00082   };
00083   /*
00084   struct ConstraintType
00085   {
00086     enum {
00087       Ignore, Atom, AtomX, AtomY, AtomZ, Distance, Angle, Torsion, Chiral
00088     };
00089   };
00090   */
00091 
00094   class OBFPRT OBFFParameter {
00095   public:
00097     int         a, b, c, d;
00099     std::string _a, _b, _c, _d;
00100 
00102     std::vector<int>    _ipar;
00104     std::vector<double> _dpar;
00105 
00107     OBFFParameter& operator=(const OBFFParameter &ai)
00108       {
00109         if (this != &ai) {
00110           a = ai.a;
00111           b = ai.b;
00112           c = ai.c;
00113           d = ai.d;
00114           _a = ai._a;
00115           _b = ai._b;
00116           _c = ai._c;
00117           _d = ai._d;
00118           _ipar = ai._ipar;
00119           _dpar = ai._dpar;
00120         }
00121 
00122         return *this;
00123       }
00124 
00126     void clear ()
00127     {
00128       a = b = c = d = 0;
00129       _ipar.clear();
00130       _dpar.clear();
00131     }
00132   }; // class OBFFParameter
00133 
00134   // specific class introductions in forcefieldYYYY.cpp (for YYYY calculations)
00135 
00138   class OBFPRT OBFFCalculation2
00139   {
00140   public:
00142     double energy;
00144     OBAtom *a, *b;
00146     int idx_a, idx_b;
00148     double *pos_a, *pos_b;
00150     double force_a[3], force_b[3];
00152     virtual ~OBFFCalculation2()
00153     {
00154     }
00157     virtual void SetupPointers()
00158     {
00159       if (!a || !b) return;
00160       pos_a = a->GetCoordinate();
00161       idx_a = a->GetIdx();
00162       pos_b = b->GetCoordinate();
00163       idx_b = b->GetIdx();
00164     }
00165   };
00166 
00169   class OBFPRT OBFFCalculation3: public OBFFCalculation2
00170   {
00171   public:
00173     OBAtom *c;
00175     int idx_c;
00177     double *pos_c;
00179     double force_c[3];
00181     virtual ~OBFFCalculation3()
00182     {
00183     }
00186     virtual void SetupPointers()
00187     {
00188       if (!a || !b || !c) return;
00189       pos_a = a->GetCoordinate();
00190       idx_a = a->GetIdx();
00191       pos_b = b->GetCoordinate();
00192       idx_b = b->GetIdx();
00193       pos_c = c->GetCoordinate();
00194       idx_c = c->GetIdx();
00195     }
00196   };
00197 
00200   class OBFPRT OBFFCalculation4: public OBFFCalculation3
00201   {
00202   public:
00204     OBAtom *d;
00206     int idx_d;
00208     double *pos_d;
00210     double force_d[3];
00212     virtual ~OBFFCalculation4()
00213     {
00214     }
00217     void SetupPointers()
00218     {
00219       if (!a || !b || !c || !d) return;
00220       pos_a = a->GetCoordinate();
00221       idx_a = a->GetIdx();
00222       pos_b = b->GetCoordinate();
00223       idx_b = b->GetIdx();
00224       pos_c = c->GetCoordinate();
00225       idx_c = c->GetIdx();
00226       pos_d = d->GetCoordinate();
00227       idx_d = d->GetIdx();
00228     }
00229   };
00230 
00234   class OBFPRT OBFFConstraint
00235   {
00236   public:
00238     double factor, constraint_value;
00239     double rab0, rbc0;
00241     int type, ia, ib, ic, id;
00243     OBAtom *a, *b, *c, *d;
00245     vector3 grada, gradb, gradc, gradd;
00246 
00248     OBFFConstraint()
00249       {
00250         a = b = c = d = NULL;
00251         ia = ib = ic = id = 0;
00252         constraint_value = 0.0;
00253         factor = 0.0;
00254       }
00256     ~OBFFConstraint()
00257       {
00258       }
00259 
00260     vector3 GetGradient(int a)
00261     {
00262       if (a == ia)
00263         return grada;
00264       else if (a == ib)
00265         return gradb;
00266       else if (a == ic)
00267         return gradc;
00268       else if (a == id)
00269         return gradd;
00270       else
00271         return  VZero;
00272     }
00273   };
00274 
00278   class OBFPRT OBFFConstraints
00279   {
00280   public:
00282     OBFFConstraints();
00284     ~OBFFConstraints()
00285       {
00286         _constraints.clear();
00287         _ignored.Clear();
00288         _fixed.Clear();
00289         _Xfixed.Clear();
00290         _Yfixed.Clear();
00291         _Zfixed.Clear();
00292       }
00294     void Clear();
00296     double GetConstraintEnergy();
00298     vector3 GetGradient(int a);
00300     OBFFConstraints& operator=(const OBFFConstraints &ai)
00301       {
00302         if (this != &ai) {
00303           _constraints = ai._constraints;
00304           _ignored = ai._ignored;
00305           _fixed = ai._fixed;
00306           _Xfixed = ai._Xfixed;
00307           _Yfixed = ai._Yfixed;
00308           _Zfixed = ai._Zfixed;
00309         }
00310         return *this;
00311       }
00312 
00316     void Setup(OBMol &mol);
00317 
00319     // Set Constraints                                                     //
00322 
00323 
00324     void SetFactor(double factor);
00326     void AddIgnore(int a);
00328     void AddAtomConstraint(int a);
00330     void AddAtomXConstraint(int a);
00332     void AddAtomYConstraint(int a);
00334     void AddAtomZConstraint(int a);
00336     void AddDistanceConstraint(int a, int b, double length);
00338     void AddAngleConstraint(int a, int b, int c, double angle);
00340     void AddTorsionConstraint(int a, int b, int c, int d, double torsion);
00343     void DeleteConstraint(int index);
00345 
00346     // Get Constraints                                                     //
00349 
00350 
00351     double GetFactor();
00353     int Size() const;
00363     int GetConstraintType(int index) const;
00367     double GetConstraintValue(int index) const;
00370     int GetConstraintAtomA(int index) const;
00373     int GetConstraintAtomB(int index) const;
00376     int GetConstraintAtomC(int index) const;
00379     int GetConstraintAtomD(int index) const;
00382     bool IsIgnored(int a);
00385     bool IsFixed(int a);
00388     bool IsXFixed(int a);
00391     bool IsYFixed(int a);
00394     bool IsZFixed(int a);
00398     OBBitVec GetIgnoredBitVec() { return _ignored; }
00401     OBBitVec GetFixedBitVec() { return _fixed; }
00403 
00404   private:
00405     std::vector<OBFFConstraint> _constraints;
00406     OBBitVec    _ignored;
00407     OBBitVec    _fixed;
00408     OBBitVec    _Xfixed;
00409     OBBitVec    _Yfixed;
00410     OBBitVec    _Zfixed;
00411     double _factor;
00412   };
00413 
00414   // Class OBForceField
00415   // class introduction in forcefield.cpp
00416   class OBFPRT OBForceField : public OBPlugin
00417   {
00418     MAKE_PLUGIN(OBForceField)
00419 
00420     protected:
00421 
00463     OBFFParameter* GetParameter(int a, int b, int c, int d, std::vector<OBFFParameter> &parameter);
00465     OBFFParameter* GetParameter(const char* a, const char* b, const char* c, const char* d,
00466         std::vector<OBFFParameter> &parameter);
00468     int GetParameterIdx(int a, int b, int c, int d, std::vector<OBFFParameter> &parameter);
00469 
00478     vector3 NumericalDerivative(OBAtom *a, int terms = OBFF_ENERGY);
00480     vector3 NumericalSecondDerivative(OBAtom *a, int terms = OBFF_ENERGY);
00481 
00482     /*
00483      *   NEW gradients functions
00484      */
00485 
00488     void SetGradient(double *grad, int idx)
00489     {
00490       const int coordIdx = (idx - 1) * 3;
00491       for (unsigned int i = 0; i < 3; ++i) {
00492         _gradientPtr[coordIdx + i] = grad[i];
00493       }
00494     }
00495 
00498     void AddGradient(double *grad, int idx)
00499     {
00500       const int coordIdx = (idx - 1) * 3;
00501       for (unsigned int i = 0; i < 3; ++i) {
00502         _gradientPtr[coordIdx + i] += grad[i];
00503       }
00504     }
00505 
00508     virtual vector3 GetGradient(OBAtom *a, int /*terms*/ = OBFF_ENERGY)
00509     {
00510       const int coordIdx = (a->GetIdx() - 1) * 3;
00511       return _gradientPtr + coordIdx;
00512     }
00513 
00516     double* GetGradientPtr()
00517     {
00518       return _gradientPtr;
00519     }
00520 
00523     virtual void ClearGradients()
00524     {
00525       // We cannot use memset because IEEE floating point representations
00526       // are not guaranteed by C/C++ standard, but this loop can be
00527       // unrolled or vectorized by compilers
00528       for (unsigned int i = 0; i < _ncoords; ++i)
00529         _gradientPtr[i] = 0.0;
00530       //      memset(_gradientPtr, '\0', sizeof(double)*_ncoords);
00531     }
00532 
00540     bool IsInSameRing(OBAtom* a, OBAtom* b);
00541 
00542     // general variables
00543     OBMol       _mol; 
00544     bool        _init; 
00545     std::string _parFile; 
00546     bool        _validSetup; 
00547     double      *_gradientPtr; 
00548     // logging variables
00549     std::ostream* _logos; 
00550     char        _logbuf[BUFF_SIZE+1]; 
00551     int         _loglvl; 
00552     int         _origLogLevel;
00553     // conformer genereation (rotor search) variables
00554     int         _current_conformer; 
00555     std::vector<double> _energies; 
00556     // minimization variables
00557     double      _econv, _e_n1; 
00558     int         _cstep, _nsteps; 
00559     double      *_grad1; 
00560     unsigned int _ncoords; 
00561     int         _linesearch; 
00562     // molecular dynamics variables
00563     double      _timestep; 
00564     double      _temp; 
00565     double      *_velocityPtr; 
00566     // contraint varibles
00567     static OBFFConstraints _constraints; 
00568     static int _fixAtom; 
00569     static int _ignoreAtom; 
00570     // cut-off variables
00571     bool        _cutoff; 
00572     double      _rvdw; 
00573     double      _rele; 
00574     OBBitVec    _vdwpairs; 
00575     OBBitVec    _elepairs; 
00576     int         _pairfreq; 
00577     // group variables
00578     std::vector<OBBitVec> _intraGroup; 
00579     std::vector<OBBitVec> _interGroup; 
00580     std::vector<std::pair<OBBitVec, OBBitVec> > _interGroups; 
00581 
00582   public:
00586     virtual OBForceField* MakeNewInstance()=0;
00587 
00589     virtual ~OBForceField()
00590     {
00591       if (_grad1 != NULL) {
00592         delete [] _grad1;
00593         _grad1 = NULL;
00594       }
00595       if (_gradientPtr != NULL) {
00596         delete [] _gradientPtr;
00597         _gradientPtr = NULL;
00598       }
00599     }
00600 
00602     const char* TypeID()
00603     {
00604       return "forcefields";
00605     }
00606 
00610     static OBForceField* FindForceField(const std::string& ID)
00611     {
00612       return FindType(ID.c_str());
00613     }
00617     static OBForceField* FindForceField(const char *ID)
00618     {
00619       return FindType(ID);
00620     }
00621     /*
00622      *
00623      */
00624     void SetParameterFile(const std::string &filename)
00625     {
00626       _parFile = filename;
00627       _init = false;
00628     }
00631     virtual std::string GetUnit() { return std::string("au"); }
00632     /* Does this force field have analytical gradients defined for all
00633      * calculation components (bonds, angles, non-bonded, etc.)
00634      * If this is true, code should default to using OBFF_ANALYTICAL_GRADIENT
00635      * for SteepestDescent() or ConjugateGradients().
00636      * \return True if all analytical gradients are implemented.
00637      */
00638     virtual bool HasAnalyticalGradients() { return false; }
00643     bool Setup(OBMol &mol);
00649     bool Setup(OBMol &mol, OBFFConstraints &constraints);
00653     // move to protected in future version
00654     virtual bool ParseParamFile() { return false; }
00658     // move to protected in future version
00659     virtual bool SetTypes() { return false; }
00663     // move to protected in future version
00664     virtual bool SetFormalCharges() { return false; }
00668     // move to protected in future version
00669     virtual bool SetPartialCharges() { return false; }
00673     // move to protected in future version
00674     virtual bool SetupCalculations() { return false; }
00679     // move to protected in future version
00680     virtual bool SetupPointers() { return false; }
00686     bool IsSetupNeeded(OBMol &mol);
00702     bool GetAtomTypes(OBMol &mol);
00718     bool GetPartialCharges(OBMol &mol);
00719 
00720 
00721 
00726     bool GetCoordinates(OBMol &mol);
00728     bool UpdateCoordinates(OBMol &mol) {return GetCoordinates(mol); }
00733     bool GetConformers(OBMol &mol);
00735     bool UpdateConformers(OBMol &mol) { return GetConformers(mol); }
00740     bool SetCoordinates(OBMol &mol);
00745     bool SetConformers(OBMol &mol);
00755     OBGridData *GetGrid(double step, double padding, const char *type, double pchg);
00756 
00758     // Interacting groups                                                  //
00760 
00762 
00763 
00767     void AddIntraGroup(OBBitVec &group);
00772     void AddInterGroup(OBBitVec &group);
00780     void AddInterGroups(OBBitVec &group1, OBBitVec &group2);
00783     void ClearGroups();
00786     bool HasGroups();
00788 
00790     // Cut-off                                                             //
00792 
00794 
00795 
00798     void EnableCutOff(bool enable)
00799     {
00800       _cutoff = enable;
00801     }
00804     bool IsCutOffEnabled()
00805     {
00806       return _cutoff;
00807     }
00811     void SetVDWCutOff(double r)
00812     {
00813       _rvdw = r;
00814     }
00818     double GetVDWCutOff()
00819     {
00820       return _rvdw;
00821     }
00826     void SetElectrostaticCutOff(double r)
00827     {
00828       _rele = r;
00829     }
00833     double GetElectrostaticCutOff()
00834     {
00835       return _rele;
00836     }
00842     void SetUpdateFrequency(int f)
00843     {
00844       _pairfreq = f;
00845     }
00849     int GetUpdateFrequency()
00850     {
00851       return _pairfreq;
00852     }
00857     void UpdatePairsSimple();
00858 
00859     //void UpdatePairsGroup(); TODO
00860 
00864     unsigned int GetNumPairs();
00868     void EnableAllPairs()
00869     {
00870       // TODO: OBBitVec doesn't seem to be allocating it's memory correctly
00871       //_vdwpairs.SetRangeOn(0, _numpairs-1);
00872       //_elepairs.SetRangeOn(0, _numpairs-1);
00873     }
00875 
00877     // Energy Evaluation                                                   //
00879 
00881 
00882 
00891     virtual double Energy(bool UNUSED(gradients) = true) { return 0.0f; }
00898     virtual double E_Bond(bool UNUSED(gradients) = true) { return 0.0f; }
00905     virtual double E_Angle(bool UNUSED(gradients) = true) { return 0.0f; }
00912     virtual double E_StrBnd(bool UNUSED(gradients) = true) { return 0.0f; }
00919     virtual double E_Torsion(bool UNUSED(gradients) = true) { return 0.0f; }
00926     virtual double E_OOP(bool UNUSED(gradients) = true) { return 0.0f; }
00933     virtual double E_VDW(bool UNUSED(gradients) = true) { return 0.0f; }
00940     virtual double E_Electrostatic(bool UNUSED(gradients) = true) { return 0.0f; }
00942 
00944     // Logging                                                             //
00946 
00948 
00949 
00951     void PrintTypes();
00955     void PrintFormalCharges();
00958     void PrintPartialCharges();
00961     void PrintVelocities();
00966     bool SetLogFile(std::ostream *pos);
00991     bool SetLogLevel(int level);
00994     int GetLogLevel() { return _loglvl; }
00998     void OBFFLog(std::string msg)
00999     {
01000       if (!_logos)
01001         return;
01002 
01003       *_logos << msg;
01004     }
01008     void OBFFLog(const char *msg)
01009     {
01010       if (!_logos)
01011         return;
01012 
01013       *_logos << msg;
01014     }
01016 
01018     // Structure Generation                                                //
01020 
01022 
01023 
01024     void DistanceGeometry();
01044     void SystematicRotorSearch(unsigned int geomSteps = 2500);
01062     int SystematicRotorSearchInitialize(unsigned int geomSteps = 2500);
01067     bool SystematicRotorSearchNextConformer(unsigned int geomSteps = 2500);
01088     void RandomRotorSearch(unsigned int conformers, unsigned int geomSteps = 2500);
01106     void RandomRotorSearchInitialize(unsigned int conformers, unsigned int geomSteps = 2500);
01111     bool RandomRotorSearchNextConformer(unsigned int geomSteps = 2500);
01133     void WeightedRotorSearch(unsigned int conformers, unsigned int geomSteps);
01134 
01136     // Energy Minimization                                                 //
01138 
01140 
01141 
01144     void SetLineSearchType(int type)
01145     {
01146       _linesearch = type;
01147     }
01151     int GetLineSearchType()
01152     {
01153       return _linesearch;
01154     }
01158     vector3 LineSearch(OBAtom *atom, vector3 &direction);
01172     double LineSearch(double *currentCoords, double *direction);
01185     double Newton2NumLineSearch(double *direction);
01191     void   LineSearchTakeStep(double *origCoords, double *direction, double step);
01207     void SteepestDescent(int steps, double econv = 1e-6f, int method = OBFF_ANALYTICAL_GRADIENT);
01234     void SteepestDescentInitialize(int steps = 1000, double econv = 1e-6f, int method = OBFF_ANALYTICAL_GRADIENT);
01249     bool SteepestDescentTakeNSteps(int n);
01265     void ConjugateGradients(int steps, double econv = 1e-6f, int method = OBFF_ANALYTICAL_GRADIENT);
01293     void ConjugateGradientsInitialize(int steps = 1000, double econv = 1e-6f, int method = OBFF_ANALYTICAL_GRADIENT);
01309     bool ConjugateGradientsTakeNSteps(int n);
01311 
01313     // Molecular Dynamics                                                  //
01315 
01317 
01318 
01320     void GenerateVelocities();
01340     void CorrectVelocities();
01356     void MolecularDynamicsTakeNSteps(int n, double T, double timestep = 0.001, int method = OBFF_ANALYTICAL_GRADIENT);
01358 
01360     // Constraints                                                         //
01362 
01364 
01365 
01368     OBFFConstraints& GetConstraints();
01372     void SetConstraints(OBFFConstraints& constraints);
01380     void SetFixAtom(int index);
01384     void UnsetFixAtom();
01393     void SetIgnoreAtom(int index);
01397     void UnsetIgnoreAtom();
01398 
01400     static bool IgnoreCalculation(int a, int b);
01402     static bool IgnoreCalculation(int a, int b, int c);
01404     static bool IgnoreCalculation(int a, int b, int c, int d);
01406 
01407 
01409     // Validation                                                          //
01411 
01413 
01414 
01415     bool DetectExplosion();
01417     vector3 ValidateLineSearch(OBAtom *atom, vector3 &direction);
01419     void ValidateSteepestDescent(int steps);
01421     void ValidateConjugateGradients(int steps);
01423     virtual bool Validate() { return false; }
01428     virtual bool ValidateGradients() { return false; }
01433     vector3 ValidateGradientError(vector3 &numgrad, vector3 &anagrad);
01435 
01437     // Vector Analysis                                                     //
01439 
01441 
01442 
01450     static double VectorBondDerivative(double *pos_a, double *pos_b,
01451                                        double *force_a, double *force_b);
01455     static double VectorDistanceDerivative(const double* const pos_i, const double* const pos_j,
01456                                            double *force_i, double *force_j);
01458     static double VectorLengthDerivative(vector3 &a, vector3 &b);
01459 
01470     static double VectorAngleDerivative(double *pos_a, double *pos_b, double *pos_c,
01471                                         double *force_a, double *force_b, double *force_c);
01473     static double VectorAngleDerivative(vector3 &a, vector3 &b, vector3 &c);
01486     static double VectorOOPDerivative(double *pos_a, double *pos_b, double *pos_c, double *pos_d,
01487                                       double *force_a, double *force_b, double *force_c, double *force_d);
01489     static double VectorOOPDerivative(vector3 &a, vector3 &b, vector3 &c, vector3 &d);
01501     static double VectorTorsionDerivative(double *pos_a, double *pos_b, double *pos_c, double *pos_d,
01502                                           double *force_a, double *force_b, double *force_c, double *force_d);
01504     static double VectorTorsionDerivative(vector3 &a, vector3 &b, vector3 &c, vector3 &d);
01505 
01511     static void VectorSubtract(double *i, double *j, double *result)
01512     {
01513       for (unsigned int c = 0; c < 3; ++c)
01514         result[c] = i[c] - j[c];
01515     }
01516 
01517     static void VectorSubtract(const double* const i, const double* const j, double *result)
01518     {
01519       for (unsigned int c = 0; c < 3; ++c)
01520         result[c] = i[c] - j[c];
01521     }
01522 
01528     static void VectorAdd(double *i, double *j, double *result)
01529     {
01530       for (unsigned int c = 0; c < 3; ++c)
01531         result[c] = i[c] + j[c];
01532     }
01533 
01539     static void VectorDivide(double *i, double n, double *result)
01540     {
01541       for (unsigned int c = 0; c < 3; ++c)
01542         result[c] = i[c] / n;
01543     }
01544 
01550     static void VectorMultiply(double *i, double n, double *result)
01551     {
01552       for (unsigned int c = 0; c < 3; ++c)
01553         result[c] = i[c] * n;
01554     }
01555 
01556     static void VectorMultiply(const double* const i, const double n, double *result)
01557     {
01558       for (unsigned int c = 0; c < 3; ++c)
01559         result[c] = i[c] * n;
01560     }
01561 
01566     static void VectorSelfMultiply(double *i, double n)
01567     {
01568       for (unsigned int c = 0; c < 3; ++c)
01569         i[c] *= n;
01570     }
01571 
01575     static void VectorNormalize(double *i)
01576     {
01577       double length = VectorLength(i);
01578       for (unsigned int c = 0; c < 3; ++c)
01579         i[c] /= length;
01580     }
01581 
01586     static void VectorCopy(double *from, double *to)
01587     {
01588       for (unsigned int c = 0; c < 3; ++c)
01589         to[c] = from[c];
01590     }
01591 
01596     static double VectorLength(double *i)
01597     {
01598       return sqrt( i[0]*i[0] + i[1]*i[1] + i[2]*i[2] );
01599     }
01600 
01601     static double VectorDistance(double *pos_i, double *pos_j)
01602     {
01603       double ij[3];
01604       VectorSubtract(pos_i, pos_j, ij);
01605       const double rij = VectorLength(ij);
01606       return rij;
01607     }
01608 
01615     static double VectorAngle(double *i, double *j, double *k);
01616 
01624     static double VectorTorsion(double *i, double *j, double *k, double *l);
01625 
01633     static double VectorOOP(double *i, double *j, double *k, double *l);
01634 
01638     static void VectorClear(double *i)
01639     {
01640       for (unsigned int c = 0; c < 3; ++c)
01641         i[c] = 0.0;
01642     }
01643 
01649     static double VectorDot(double *i, double *j)
01650     {
01651       double result = 0.0;
01652       // Written as a loop for vectorization
01653       // Loop will be unrolled by compiler otherwise
01654       for (unsigned int c = 0; c < 3; ++c)
01655         result += i[c]*j[c];
01656       return result;
01657     }
01658 
01664     static void VectorCross(double *i, double *j, double *result)
01665     {
01666       result[0] =   i[1]*j[2] - i[2]*j[1];
01667       result[1] = - i[0]*j[2] + i[2]*j[0];
01668       result[2] =   i[0]*j[1] - i[1]*j[0];
01669     }
01670 
01671     static void PrintVector(double *i)
01672     {
01673       std::cout << "<" << i[0] << ", " << i[1] << ", " << i[2] << ">" << std::endl;
01674     }
01676 
01677   }; // class OBForceField
01678 
01679 }// namespace OpenBabel
01680 
01681 #endif   // OB_FORCEFIELD_H
01682 
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines