Base class for molecular mechanics force fields. More...
#include <openbabel/forcefield.h>
Public Types | |
typedef std::map< const char *, OBPlugin *, CharPtrLess > | PluginMapType |
typedef PluginMapType::const_iterator | PluginIterator |
Public Member Functions | |
virtual OBForceField * | MakeNewInstance ()=0 |
virtual | ~OBForceField () |
const char * | TypeID () |
void | SetParameterFile (const std::string &filename) |
virtual std::string | GetUnit () |
virtual bool | HasAnalyticalGradients () |
bool | Setup (OBMol &mol) |
bool | Setup (OBMol &mol, OBFFConstraints &constraints) |
virtual bool | ParseParamFile () |
virtual bool | SetTypes () |
virtual bool | SetFormalCharges () |
virtual bool | SetPartialCharges () |
virtual bool | SetupCalculations () |
virtual bool | SetupPointers () |
bool | IsSetupNeeded (OBMol &mol) |
bool | GetAtomTypes (OBMol &mol) |
bool | GetPartialCharges (OBMol &mol) |
bool | GetCoordinates (OBMol &mol) |
bool | UpdateCoordinates (OBMol &mol) |
bool | GetConformers (OBMol &mol) |
bool | UpdateConformers (OBMol &mol) |
bool | SetCoordinates (OBMol &mol) |
bool | SetConformers (OBMol &mol) |
OBGridData * | GetGrid (double step, double padding, const char *type, double pchg) |
virtual const char * | Description () |
virtual bool | Display (std::string &txt, const char *param, const char *ID=NULL) |
virtual OBPlugin * | MakeInstance (const std::vector< std::string > &) |
virtual void | Init () |
const char * | GetID () const |
virtual PluginMapType & | GetMap () const =0 |
Methods for specifying interaction groups | |
void | AddIntraGroup (OBBitVec &group) |
void | AddInterGroup (OBBitVec &group) |
void | AddInterGroups (OBBitVec &group1, OBBitVec &group2) |
void | ClearGroups () |
bool | HasGroups () |
Methods for Cut-off distances | |
void | EnableCutOff (bool enable) |
bool | IsCutOffEnabled () |
void | SetVDWCutOff (double r) |
double | GetVDWCutOff () |
void | SetElectrostaticCutOff (double r) |
double | GetElectrostaticCutOff () |
void | SetUpdateFrequency (int f) |
int | GetUpdateFrequency () |
void | UpdatePairsSimple () |
unsigned int | GetNumPairs () |
void | EnableAllPairs () |
Methods for energy evaluation | |
virtual double | Energy (bool gradients=true) |
virtual double | E_Bond (bool gradients=true) |
virtual double | E_Angle (bool gradients=true) |
virtual double | E_StrBnd (bool gradients=true) |
virtual double | E_Torsion (bool gradients=true) |
virtual double | E_OOP (bool gradients=true) |
virtual double | E_VDW (bool gradients=true) |
virtual double | E_Electrostatic (bool gradients=true) |
Methods for logging | |
void | PrintTypes () |
void | PrintFormalCharges () |
void | PrintPartialCharges () |
void | PrintVelocities () |
bool | SetLogFile (std::ostream *pos) |
bool | SetLogLevel (int level) |
int | GetLogLevel () |
void | OBFFLog (std::string msg) |
void | OBFFLog (const char *msg) |
Methods for structure generation | |
void | DistanceGeometry () |
void | SystematicRotorSearch (unsigned int geomSteps=2500) |
int | SystematicRotorSearchInitialize (unsigned int geomSteps=2500) |
bool | SystematicRotorSearchNextConformer (unsigned int geomSteps=2500) |
void | RandomRotorSearch (unsigned int conformers, unsigned int geomSteps=2500) |
void | RandomRotorSearchInitialize (unsigned int conformers, unsigned int geomSteps=2500) |
bool | RandomRotorSearchNextConformer (unsigned int geomSteps=2500) |
void | WeightedRotorSearch (unsigned int conformers, unsigned int geomSteps) |
Methods for energy minimization | |
void | SetLineSearchType (int type) |
int | GetLineSearchType () |
vector3 | LineSearch (OBAtom *atom, vector3 &direction) |
double | LineSearch (double *currentCoords, double *direction) |
double | Newton2NumLineSearch (double *direction) |
void | LineSearchTakeStep (double *origCoords, double *direction, double step) |
void | SteepestDescent (int steps, double econv=1e-6f, int method=OBFF_ANALYTICAL_GRADIENT) |
void | SteepestDescentInitialize (int steps=1000, double econv=1e-6f, int method=OBFF_ANALYTICAL_GRADIENT) |
bool | SteepestDescentTakeNSteps (int n) |
void | ConjugateGradients (int steps, double econv=1e-6f, int method=OBFF_ANALYTICAL_GRADIENT) |
void | ConjugateGradientsInitialize (int steps=1000, double econv=1e-6f, int method=OBFF_ANALYTICAL_GRADIENT) |
bool | ConjugateGradientsTakeNSteps (int n) |
Methods for molecular dynamics | |
void | GenerateVelocities () |
void | CorrectVelocities () |
void | MolecularDynamicsTakeNSteps (int n, double T, double timestep=0.001, int method=OBFF_ANALYTICAL_GRADIENT) |
Methods for forcefield validation | |
bool | DetectExplosion () |
vector3 | ValidateLineSearch (OBAtom *atom, vector3 &direction) |
void | ValidateSteepestDescent (int steps) |
void | ValidateConjugateGradients (int steps) |
virtual bool | Validate () |
virtual bool | ValidateGradients () |
vector3 | ValidateGradientError (vector3 &numgrad, vector3 &anagrad) |
Static Public Member Functions | |
static OBForceField * | FindForceField (const std::string &ID) |
static OBForceField * | FindForceField (const char *ID) |
static OBPlugin * | GetPlugin (const char *Type, const char *ID) |
static bool | ListAsVector (const char *PluginID, const char *param, std::vector< std::string > &vlist) |
static void | List (const char *PluginID, const char *param=NULL, std::ostream *os=&std::cout) |
static std::string | ListAsString (const char *PluginID, const char *param=NULL) |
static std::string | FirstLine (const char *txt) |
static PluginIterator | Begin (const char *PluginID) |
static PluginIterator | End (const char *PluginID) |
Methods for vector analysis (used by OBFFXXXXCalculationYYYY) | |
static double | VectorBondDerivative (double *pos_a, double *pos_b, double *force_a, double *force_b) |
static double | VectorDistanceDerivative (const double *const pos_i, const double *const pos_j, double *force_i, double *force_j) |
static double | VectorLengthDerivative (vector3 &a, vector3 &b) |
static double | VectorAngleDerivative (double *pos_a, double *pos_b, double *pos_c, double *force_a, double *force_b, double *force_c) |
static double | VectorAngleDerivative (vector3 &a, vector3 &b, vector3 &c) |
static double | VectorOOPDerivative (double *pos_a, double *pos_b, double *pos_c, double *pos_d, double *force_a, double *force_b, double *force_c, double *force_d) |
static double | VectorOOPDerivative (vector3 &a, vector3 &b, vector3 &c, vector3 &d) |
static double | VectorTorsionDerivative (double *pos_a, double *pos_b, double *pos_c, double *pos_d, double *force_a, double *force_b, double *force_c, double *force_d) |
static double | VectorTorsionDerivative (vector3 &a, vector3 &b, vector3 &c, vector3 &d) |
static void | VectorSubtract (double *i, double *j, double *result) |
static void | VectorSubtract (const double *const i, const double *const j, double *result) |
static void | VectorAdd (double *i, double *j, double *result) |
static void | VectorDivide (double *i, double n, double *result) |
static void | VectorMultiply (double *i, double n, double *result) |
static void | VectorMultiply (const double *const i, const double n, double *result) |
static void | VectorSelfMultiply (double *i, double n) |
static void | VectorNormalize (double *i) |
static void | VectorCopy (double *from, double *to) |
static double | VectorLength (double *i) |
static double | VectorDistance (double *pos_i, double *pos_j) |
static double | VectorAngle (double *i, double *j, double *k) |
static double | VectorTorsion (double *i, double *j, double *k, double *l) |
static double | VectorOOP (double *i, double *j, double *k, double *l) |
static void | VectorClear (double *i) |
static double | VectorDot (double *i, double *j) |
static void | VectorCross (double *i, double *j, double *result) |
static void | PrintVector (double *i) |
Protected Member Functions | |
OBFFParameter * | GetParameter (int a, int b, int c, int d, std::vector< OBFFParameter > ¶meter) |
OBFFParameter * | GetParameter (const char *a, const char *b, const char *c, const char *d, std::vector< OBFFParameter > ¶meter) |
int | GetParameterIdx (int a, int b, int c, int d, std::vector< OBFFParameter > ¶meter) |
vector3 | NumericalDerivative (OBAtom *a, int terms=OBFF_ENERGY) |
vector3 | NumericalSecondDerivative (OBAtom *a, int terms=OBFF_ENERGY) |
void | SetGradient (double *grad, int idx) |
void | AddGradient (double *grad, int idx) |
virtual vector3 | GetGradient (OBAtom *a, int=OBFF_ENERGY) |
double * | GetGradientPtr () |
virtual void | ClearGradients () |
bool | IsInSameRing (OBAtom *a, OBAtom *b) |
Static Protected Member Functions | |
static PluginMapType & | PluginMap () |
static PluginMapType & | GetTypeMap (const char *PluginID) |
static OBPlugin * | BaseFindType (PluginMapType &Map, const char *ID) |
Protected Attributes | |
OBMol | _mol |
bool | _init |
std::string | _parFile |
bool | _validSetup |
double * | _gradientPtr |
std::ostream * | _logos |
char | _logbuf [BUFF_SIZE+1] |
int | _loglvl |
int | _origLogLevel |
int | _current_conformer |
std::vector< double > | _energies |
double | _econv |
double | _e_n1 |
int | _cstep |
int | _nsteps |
double * | _grad1 |
unsigned int | _ncoords |
int | _linesearch |
double | _timestep |
double | _temp |
double * | _velocityPtr |
bool | _cutoff |
double | _rvdw |
double | _rele |
OBBitVec | _vdwpairs |
OBBitVec | _elepairs |
int | _pairfreq |
std::vector< OBBitVec > | _intraGroup |
std::vector< OBBitVec > | _interGroup |
std::vector< std::pair < OBBitVec, OBBitVec > > | _interGroups |
const char * | _id |
Static Protected Attributes | |
static OBFFConstraints | _constraints = OBFFConstraints() |
static int | _fixAtom = 0 |
static int | _ignoreAtom = 0 |
Methods for constraints | |
OBFFConstraints & | GetConstraints () |
void | SetConstraints (OBFFConstraints &constraints) |
void | SetFixAtom (int index) |
void | UnsetFixAtom () |
void | SetIgnoreAtom (int index) |
void | UnsetIgnoreAtom () |
static bool | IgnoreCalculation (int a, int b) |
static bool | IgnoreCalculation (int a, int b, int c) |
static bool | IgnoreCalculation (int a, int b, int c, int d) |
Base class for molecular mechanics force fields.
The OBForceField class is the base class for molecular mechanics in Open Babel. Classes derived from the OBForceField implement specific force fields (Ghemical, MMFF94, UFF, ...).Other classes such as OBFFParameter, OBFFConstraint, OBFFCalculation and its derived classes are only for internal use. As a user interested in using the available force fields in Open Babel, you don't need these classes. The rest of this short introduction is aimed at these users. For information on how to implement additional force fields, see the wiki pages or post your questions to the openbabel-devel mailing list.
Before we can start using a force field, we must first select it and set it up. This is illustrated in the first example below. The Setup procedure assigns atom types, charges and parameters. There are several reasons why this may fail, a log message will be written to the logfile before Setup() returns false.
The force field classes use their own logging functions. You can set the logfile using SetLogFile() and set the log level using SetLogLevel(). If needed you can also write to the logfile using OBFFLog(). There are four log levels: BFF_LOGLVL_NONE, OBFF_LOGLVL_LOW, OBFF_LOGLVL_MEDIUM, OBFF_LOGLVL_HIGH. See the API documentation to know what kind of output each function writes to the logfile for the different log levels.
Below are two examples which explain the basics.
This piece of code will output a list of available forcefields to cout:
OBPlugin::List("forcefields");
Calculate the energy for the structure in mol using the Ghemical forcefield.
#include <openbabel/forcefield.h> #include <openbabel/mol.h> // See OBConversion class to fill the mol object. OBMol mol; // Select the forcefield, this returns a pointer that we // will later use to access the forcefield functions. OBForceField* pFF = OBForceField::FindForceField("MMFF94"); // Make sure we have a valid pointer if (!pFF) // exit... // Set the logfile (can also be &cout or &cerr) pFF->SetLogFile(&cerr); // Set the log level. See indivual functions to know // what kind of output each function produces for the // different log levels. pFF->SetLogLevel(OBFF_LOGLVL_HIGH); // We need to setup the forcefield before we can use it. Setup() // returns false if it failes to find the atom types, parameters, ... if (!pFF->Setup(mol)) { cerr << "ERROR: could not setup force field." << endl; } // Calculate the energy. The output will be written to the // logfile specified by SetLogFile() pFF->Energy();
Minimize the structure in mol using conjugate gradients.
#include <openbabel/forcefield.h> #include <openbabel/mol.h> OBMol mol; OBForceField* pFF = OBForceField::FindForceField("MMFF94"); // Make sure we have a valid pointer if (!pFF) // exit... pFF->SetLogFile(&cerr); pFF->SetLogLevel(OBFF_LOGLVL_LOW); if (!pFF->Setup(mol)) { cerr << "ERROR: could not setup force field." << endl; } // Perform the actual minimization, maximum 1000 steps pFF->ConjugateGradients(1000);
Minimize the structure in mol using steepest descent and fix the position of atom with index 1.
#include <openbabel/forcefield.h> #include <openbabel/mol.h> OBMol mol; OBForceField* pFF = OBForceField::FindForceField("MMFF94"); // Make sure we have a valid pointer if (!pFF) // exit... pFF->SetLogFile(&cerr); pFF->SetLogLevel(OBFF_LOGLVL_LOW); // Set the constraints OBFFConstraints constraints; constraints.AddAtomConstraint(1); // We pass the constraints as argument for Setup() if (!pFF->Setup(mol, constraints)) { cerr << "ERROR: could not setup force field." << endl; } // Perform the actual minimization, maximum 1000 steps pFF->SteepestDescent(1000);
Minimize a ligand molecule in a binding pocket.
#include <openbabel/forcefield.h> #include <openbabel/mol.h> OBMol mol; // // Read the pocket + ligand (initial guess for position) into mol... // OBBitVec pocket; // set the bits with atoms indexes for the pocket to 1... OBBitVec ligand; // set the bits with atoms indexes for the ligand to 1... OBForceField* pFF = OBForceField::FindForceField("MMFF94"); // Make sure we have a valid pointer if (!pFF) // exit... pFF->SetLogFile(&cerr); pFF->SetLogLevel(OBFF_LOGLVL_LOW); // Fix the binding pocket atoms OBFFConstraints constraints; FOR_ATOMS_OF_MOL (a, mol) { if (pocket.BitIsOn(a->GetIdx()) constraints.AddAtomConstraint(a->GetIdx()); } // Specify the interacting groups. The pocket atoms are fixed, so there // is no need to calculate intra- and inter-molecular interactions for // the binding pocket. pFF->AddIntraGroup(ligand); // bonded interactions in the ligand pFF->AddInterGroup(ligand); // non-bonded between ligand-ligand atoms pFF->AddInterGroups(ligand, pocket); // non-bonded between ligand and pocket atoms // We pass the constraints as argument for Setup() if (!pFF->Setup(mol, constraints)) { cerr << "ERROR: could not setup force field." << endl; } // Perform the actual minimization, maximum 1000 steps pFF->SteepestDescent(1000);
typedef std::map<const char*, OBPlugin*, CharPtrLess> PluginMapType [inherited] |
virtual ~OBForceField | ( | ) | [inline, virtual] |
Destructor.
OBFFParameter * GetParameter | ( | int | a, |
int | b, | ||
int | c, | ||
int | d, | ||
std::vector< OBFFParameter > & | parameter | ||
) | [protected] |
Get the correct OBFFParameter from a OBFFParameter vector.
vector<OBFFParameter> parameters;
this vector is filled with entries (as OBFFParameter) from a parameter file. This happens in the Setup() function.
GetParameter(a, 0, 0, 0, parameters);
returns the first OBFFParameter from vector<OBFFParameter> parameters where: pa = a (pa = parameter.a)
use: vdw parameters, ...
GetParameter(a, b, 0, 0, parameters);
returns the first OBFFParameter from vector<OBFFParameter> parameters where: pa = a & pb = b (ab) or: pa = b & pb = a (ba)
use: bond parameters, vdw parameters (pairs), ...
GetParameter(a, b, c, 0, parameters);
returns the first OBFFParameter from vector<OBFFParameter> parameters where: pa = a & pb = b & pc = c (abc) or: pa = c & pb = b & pc = a (cba)
use: angle parameters, ...
GetParameter(a, b, c, d, parameters);
returns the first OBFFParameter from vector<OBFFParameter> parameters where: pa = a & pb = b & pc = c & pd = d (abcd) or: pa = d & pb = b & pc = c & pd = a (dbca) or: pa = a & pb = c & pc = b & pd = d (acbd) or: pa = d & pb = c & pc = b & pd = a (dcba)
use: torsion parameters, ...
OBFFParameter * GetParameter | ( | const char * | a, |
const char * | b, | ||
const char * | c, | ||
const char * | d, | ||
std::vector< OBFFParameter > & | parameter | ||
) | [protected] |
int GetParameterIdx | ( | int | a, |
int | b, | ||
int | c, | ||
int | d, | ||
std::vector< OBFFParameter > & | parameter | ||
) | [protected] |
Get index for vector<OBFFParameter> ...
Calculate the potential energy function derivative numerically with repect to the coordinates of atom with index a (this vector is the gradient)
a | provides coordinates |
terms | OBFF_ENERGY, OBFF_EBOND, OBFF_EANGLE, OBFF_ESTRBND, OBFF_ETORSION, OBFF_EOOP, OBFF_EVDW, OBFF_ELECTROSTATIC |
void SetGradient | ( | double * | grad, |
int | idx | ||
) | [inline, protected] |
Set the gradient for atom with index idx to grad
void AddGradient | ( | double * | grad, |
int | idx | ||
) | [inline, protected] |
Add grad to the gradient for atom with index idx
Get the pointer to the gradients
double* GetGradientPtr | ( | ) | [inline, protected] |
Get the pointer to the gradients
virtual void ClearGradients | ( | ) | [inline, protected, virtual] |
Set all gradients to zero
Check if two atoms are in the same ring. [NOTE: this function uses SSSR, this means that not all rings are found for bridged rings. This causes some problems with the MMFF94 validation.]
a | atom a |
b | atom b |
virtual OBForceField* MakeNewInstance | ( | ) | [pure virtual] |
Clone the current instance. May be desirable in multithreaded environments, Should be deleted after use
const char* TypeID | ( | ) | [inline, virtual] |
Reimplemented from OBPlugin.
static OBForceField* FindForceField | ( | const std::string & | ID ) | [inline, static] |
ID | forcefield id (Ghemical, MMFF94, UFF, ...). |
static OBForceField* FindForceField | ( | const char * | ID ) | [inline, static] |
ID | forcefield id (Ghemical, MMFF94, UFF, ...). |
void SetParameterFile | ( | const std::string & | filename ) | [inline] |
virtual std::string GetUnit | ( | ) | [inline, virtual] |
virtual bool HasAnalyticalGradients | ( | ) | [inline, virtual] |
bool Setup | ( | OBMol & | mol ) |
Setup the forcefield for mol (assigns atom types, charges, etc.). Keep current constraints.
mol | The OBMol object that contains the atoms and bonds. |
Referenced by OBEnergyConformerScore::Score().
bool Setup | ( | OBMol & | mol, |
OBFFConstraints & | constraints | ||
) |
Setup the forcefield for mol (assigns atom types, charges, etc.). Use new constraints.
mol | The OBMol object that contains the atoms and bonds. |
constraints | The OBFFConstraints object that contains the constraints. |
virtual bool ParseParamFile | ( | ) | [inline, virtual] |
Load the parameters (this function is overloaded by the individual forcefields, and is called autoamically from OBForceField::Setup()).
virtual bool SetTypes | ( | ) | [inline, virtual] |
Set the atom types (this function is overloaded by the individual forcefields, and is called autoamically from OBForceField::Setup()).
virtual bool SetFormalCharges | ( | ) | [inline, virtual] |
Set the formal charges (this function is overloaded by the individual forcefields, and is called autoamically from OBForceField::Setup()).
virtual bool SetPartialCharges | ( | ) | [inline, virtual] |
Set the partial charges (this function is overloaded by the individual forcefields, and is called autoamically from OBForceField::Setup()).
virtual bool SetupCalculations | ( | ) | [inline, virtual] |
Setup the calculations (this function is overloaded by the individual forcefields, and is called autoamically from OBForceField::Setup()).
virtual bool SetupPointers | ( | ) | [inline, virtual] |
Setup the pointers to the atom positions in the OBFFCalculation objects. This method will iterate over all the calculations and call SetupPointers for each one. (This function should be implemented by the individual force field implementations).
bool IsSetupNeeded | ( | OBMol & | mol ) |
Compare the internal forcefield OBMol object to mol. If the two have the same number of atoms and bonds, and all atomic numbers are the same, this function returns false, and no call to Setup is needed.
bool GetAtomTypes | ( | OBMol & | mol ) |
Get the force atom types. The atom types will be added to the atoms of mol as OBPairData. The attribute will be "FFAtomType".
... pFF->Setup(&mol); pFF->GetAtomTypes(&mol); FOR_ATOMS_OF_MOL (atom, mol) { OBPairData *type = (OBPairData*) atom->GetData("FFAtomType"); if (type) cout << "atom " << atom->GetIdx() << " : " << type->GetValue() << endl; } ...
bool GetPartialCharges | ( | OBMol & | mol ) |
Get the force field formal charges. The formal charges will be added to the atoms of mol as OBPairData. The attribute will be "FFPartialCharge".
... pFF->Setup(&mol); pFF->GetPartialCharges(&mol); FOR_ATOMS_OF_MOL (atom, mol) { OBPairData *chg = (OBPairData*) atom->GetData("FFPartialCharge"); if (chg) cout << "atom " << atom->GetIdx() << " : " << chg->GetValue() << endl; } ...
bool GetCoordinates | ( | OBMol & | mol ) |
Get coordinates for current conformer and attach OBConformerData with energies, forces, ... to mol.
mol | The OBMol object to copy the coordinates to (from OBForceField::_mol). |
bool UpdateCoordinates | ( | OBMol & | mol ) | [inline] |
bool GetConformers | ( | OBMol & | mol ) |
Get coordinates for all conformers and attach OBConformerData with energies, forces, ... to mol.
mol | The OBMol object to copy the coordinates to (from OBForceField::_mol). |
bool UpdateConformers | ( | OBMol & | mol ) | [inline] |
bool SetCoordinates | ( | OBMol & | mol ) |
Set coordinates for current conformer.
mol | the OBMol object to copy the coordinates from (to OBForceField::_mol). |
bool SetConformers | ( | OBMol & | mol ) |
Set coordinates for all conformers.
mol | The OBMol object to copy the coordinates from (to OBForceField::_mol). |
OBGridData * GetGrid | ( | double | step, |
double | padding, | ||
const char * | type, | ||
double | pchg | ||
) |
Create a grid with spacing step
and padding
. Place a probe atom of type probe at every grid point, calculate the energy and store it in the grid. These grids can then be used to create isosurfaces to identify locations where the probe atom has favourable interactions with the molecule.
step | The grid step size in A.. |
padding | The padding for the grid in A. |
type | The force field atom type for the probe. |
pchg | The partial charge for the probe atom. |
void AddIntraGroup | ( | OBBitVec & | group ) |
void AddInterGroup | ( | OBBitVec & | group ) |
Enable inter-molecular interactions between group1 and group2 (non-bonded: vdw & ele). Note that this function doesn't enable bonded interactions in either group. Non-bonded interactions in the groups itself are also not enabled. This function should be called before Setup().
void ClearGroups | ( | ) |
Clear all previously specified groups.
bool HasGroups | ( | ) |
void EnableCutOff | ( | bool | enable ) | [inline] |
Enable or disable Cut-offs. Cut-offs are disabled by default.
enable | Enable when true, disable when false. |
bool IsCutOffEnabled | ( | ) | [inline] |
void SetVDWCutOff | ( | double | r ) | [inline] |
Set the VDW cut-off distance to r. Note that this does not enable cut-off distances.
r | The VDW cut-off distance to be used in A. |
double GetVDWCutOff | ( | ) | [inline] |
Get the VDW cut-off distance.
void SetElectrostaticCutOff | ( | double | r ) | [inline] |
Set the Electrostatic cut-off distance to r. Note that this does not enable cut-off distances.
r | The electrostatic cut-off distance to be used in A. |
double GetElectrostaticCutOff | ( | ) | [inline] |
Get the Electrostatic cut-off distance.
void SetUpdateFrequency | ( | int | f ) | [inline] |
Set the frequency by which non-bonded pairs are updated. Values from 10 to 20 are recommended. Too low will decrease performance, too high will cause non-bonded interactions within cut-off not to be calculated.
f | The pair list update frequency. |
int GetUpdateFrequency | ( | ) | [inline] |
Get the frequency by which non-bonded pairs are updated.
void UpdatePairsSimple | ( | ) |
Set the bits in _vdwpairs and _elepairs to 1 for interactions that are within cut-off distance. This function is called in minimizing algorithms such as SteepestDescent and ConjugateGradients.
unsigned int GetNumPairs | ( | ) |
Get the number of non-bonded pairs in _mol.
void EnableAllPairs | ( | ) | [inline] |
Set bits in range 0..._numpairs-1 to 1. Using this means there will be no cut-off. (not-working: see code for more information.
virtual double Energy | ( | bool | gradients = true ) |
[inline, virtual] |
gradients | Set to true when the gradients need to be calculated (needs to be done before calling GetGradient()). |
Referenced by OBEnergyConformerScore::Score().
virtual double E_Bond | ( | bool | gradients = true ) |
[inline, virtual] |
gradients | Set to true when the gradients need to be calculated (needs to be done before calling GetGradient()). |
virtual double E_Angle | ( | bool | gradients = true ) |
[inline, virtual] |
gradients | Set to true when the gradients need to be calculated (needs to be done before calling GetGradient()). |
virtual double E_StrBnd | ( | bool | gradients = true ) |
[inline, virtual] |
gradients | Set to true when the gradients need to be calculated (needs to be done before calling GetGradient()). |
virtual double E_Torsion | ( | bool | gradients = true ) |
[inline, virtual] |
gradients | Set to true when the gradients need to be calculated (needs to be done before calling GetGradient()). |
virtual double E_OOP | ( | bool | gradients = true ) |
[inline, virtual] |
gradients | Set to true when the gradients need to be calculated (needs to be done before calling GetGradient()). |
virtual double E_VDW | ( | bool | gradients = true ) |
[inline, virtual] |
gradients | Set to true when the gradients need to be calculated (needs to be done before calling GetGradient()). |
virtual double E_Electrostatic | ( | bool | gradients = true ) |
[inline, virtual] |
gradients | Set to true when the gradients need to be calculated (needs to be done before calling GetGradient()). |
void PrintTypes | ( | ) |
Print the atom types to the log.
void PrintFormalCharges | ( | ) |
Print the formal charges to the log (atom.GetPartialCharge(), MMFF94 FC's are not always int).
void PrintPartialCharges | ( | ) |
Print the partial charges to the log.
void PrintVelocities | ( | ) |
Print the velocities to the log.
bool SetLogFile | ( | std::ostream * | pos ) |
Set the stream for logging (can also be &cout for logging to screen).
pos | Stream (when pos is 0, std::cout wil be used). |
bool SetLogLevel | ( | int | level ) |
Set the log level (OBFF_LOGLVL_NONE, OBFF_LOGLVL_LOW, OBFF_LOGLVL_MEDIUM, OBFF_LOGLVL_HIGH). Inline if statements for logging are available:
#define IF_OBFF_LOGLVL_LOW if(_loglvl >= OBFF_LOGLVL_LOW) #define IF_OBFF_LOGLVL_MEDIUM if(_loglvl >= OBFF_LOGLVL_MEDIUM) #define IF_OBFF_LOGLVL_HIGH if(_loglvl >= OBFF_LOGLVL_HIGH)
example:
SetLogLevel(OBFF_LOGLVL_MEDIUM); IF_OBFF_LOGLVL_HIGH { OBFFLog("this text will NOT be logged...\n"); } IF_OBFF_LOGLVL_LOW { OBFFLog"this text will be logged...\n"); } IF_OBFF_LOGLVL_MEDIUM { OBFFLog("this text will also be logged...\n"); }
int GetLogLevel | ( | ) | [inline] |
void OBFFLog | ( | std::string | msg ) | [inline] |
Print msg to the logfile.
msg | The message to print. |
void OBFFLog | ( | const char * | msg ) | [inline] |
Print msg to the logfile.
msg | The message to print. |
void DistanceGeometry | ( | ) |
Generate coordinates for the molecule (distance geometry). (OB 3.0).
void SystematicRotorSearch | ( | unsigned int | geomSteps = 2500 ) |
Generate conformers for the molecule (systematicaly rotating torsions).
The initial starting structure here is important, this structure should be minimized for the best results. SystematicRotorSearch works by rotating around the rotatable bond in a molecule (see OBRotamerList class). This rotating generates multiple conformers. The energy for all these conformers is then evaluated and the lowest energy conformer is selected.
geomSteps | The number of steps to take during geometry optimization. |
int SystematicRotorSearchInitialize | ( | unsigned int | geomSteps = 2500 ) |
Generate conformers for the molecule by systematicaly rotating torsions. To be used in combination with SystematicRotorSearchNexConformer().
example:
// pFF is a pointer to a OBForceField class pFF->SystematicRotorSearchInitialize(300); while (pFF->SystematicRotorSearchNextConformer(300)) { // do some updating in your program (show last generated conformer, ...) }
If you don't need any updating in your program, SystematicRotorSearch() is recommended.
geomSteps | The number of steps to take during geometry optimization. |
bool SystematicRotorSearchNextConformer | ( | unsigned int | geomSteps = 2500 ) |
Evaluate the next conformer.
geomSteps | The number of steps to take during geometry optimization. |
void RandomRotorSearch | ( | unsigned int | conformers, |
unsigned int | geomSteps = 2500 |
||
) |
Generate conformers for the molecule (randomly rotating torsions).
The initial starting structure here is important, this structure should be minimized for the best results. RandomRotorSearch works by randomly rotating around the rotatable bonds in a molecule (see OBRotamerList class). This rotating generates multiple conformers. The energy for all these conformers is then evaluated and the lowest energy conformer is selected.
conformers | The number of random conformers to consider during the search. |
geomSteps | The number of steps to take during geometry optimization for each conformer. |
void RandomRotorSearchInitialize | ( | unsigned int | conformers, |
unsigned int | geomSteps = 2500 |
||
) |
Generate conformers for the molecule by randomly rotating torsions. To be used in combination with RandomRotorSearchNexConformer().
example:
// pFF is a pointer to a OBForceField class pFF->RandomRotorSearchInitialize(300); while (pFF->RandomRotorSearchNextConformer(300)) { // do some updating in your program (show last generated conformer, ...) }
If you don't need any updating in your program, RandomRotorSearch() is recommended.
conformers | The number of random conformers to consider during the search |
geomSteps | The number of steps to take during geometry optimization |
bool RandomRotorSearchNextConformer | ( | unsigned int | geomSteps = 2500 ) |
Evaluate the next conformer.
geomSteps | The number of steps to take during geometry optimization. |
void WeightedRotorSearch | ( | unsigned int | conformers, |
unsigned int | geomSteps | ||
) |
Generate conformers for the molecule (randomly rotating torsions).
The initial starting structure here is important, this structure should be minimized for the best results. WeightedRotorSearch works by randomly rotating around the rotatable bonds in a molecule (see OBRotamerList class). Unlike RandomRotorSearch() the random choice of torsions is reweighted based on the energy of the generated conformer. Over time, the generated conformers for each step should become increasingly better. The lowest energy conformer is selected.
conformers | The number of random conformers to consider during the search. |
geomSteps | The number of steps to take during geometry optimization for each conformer. |
void SetLineSearchType | ( | int | type ) | [inline] |
Set the LineSearchType. The default type is LineSearchType::Simple.
type | The LineSearchType to be used in SteepestDescent and ConjugateGradients. |
int GetLineSearchType | ( | ) | [inline] |
Get the LineSearchType.
Perform a linesearch starting at atom in direction direction.
double LineSearch | ( | double * | currentCoords, |
double * | direction | ||
) |
Perform a linesearch for the entire molecule in direction direction
. This function is called when using LineSearchType::Simple.
currentCoords | Start coordinates. |
direction | The search direction. |
double Newton2NumLineSearch | ( | double * | direction ) |
Perform a linesearch for the entire molecule. This function is called when using LineSearchType::Newton2Num.
direction | The search direction. |
void LineSearchTakeStep | ( | double * | origCoords, |
double * | direction, | ||
double | step | ||
) |
Set the coordinates of the atoms to origCoord + step.
origCoords | Start coordinates. |
direction | The search direction. |
step | The step to take. |
void SteepestDescent | ( | int | steps, |
double | econv = 1e-6f , |
||
int | method = OBFF_ANALYTICAL_GRADIENT |
||
) |
Perform steepest descent optimalization for steps steps or until convergence criteria is reached.
steps | The number of steps. |
econv | Energy convergence criteria. (defualt is 1e-6) |
method | Deprecated. (see HasAnalyticalGradients()) |
void SteepestDescentInitialize | ( | int | steps = 1000 , |
double | econv = 1e-6f , |
||
int | method = OBFF_ANALYTICAL_GRADIENT |
||
) |
Initialize steepest descent optimalization, to be used in combination with SteepestDescentTakeNSteps().
example:
// pFF is a pointer to a OBForceField class pFF->SteepestDescentInitialize(100, 1e-5f); while (pFF->SteepestDescentTakeNSteps(5)) { // do some updating in your program (redraw structure, ...) }
If you don't need any updating in your program, SteepestDescent() is recommended.
steps | The number of steps. |
econv | Energy convergence criteria. (defualt is 1e-6) |
method | Deprecated. (see HasAnalyticalGradients()) |
bool SteepestDescentTakeNSteps | ( | int | n ) |
Take n steps in a steepestdescent optimalization that was previously initialized with SteepestDescentInitialize().
n | The number of steps to take. |
void ConjugateGradients | ( | int | steps, |
double | econv = 1e-6f , |
||
int | method = OBFF_ANALYTICAL_GRADIENT |
||
) |
Perform conjugate gradient optimalization for steps steps or until convergence criteria is reached.
steps | The number of steps. |
econv | Energy convergence criteria. (defualt is 1e-6) |
method | Deprecated. (see HasAnalyticalGradients()) |
void ConjugateGradientsInitialize | ( | int | steps = 1000 , |
double | econv = 1e-6f , |
||
int | method = OBFF_ANALYTICAL_GRADIENT |
||
) |
Initialize conjugate gradient optimalization and take the first step, to be used in combination with ConjugateGradientsTakeNSteps().
example:
// pFF is a pointer to a OBForceField class pFF->ConjugateGradientsInitialize(100, 1e-5f); while (pFF->ConjugateGradientsTakeNSteps(5)) { // do some updating in your program (redraw structure, ...) }
If you don't need any updating in your program, ConjugateGradients() is recommended.
steps | The number of steps. |
econv | Energy convergence criteria. (defualt is 1e-6) |
method | Deprecated. (see HasAnalyticalGradients()) |
bool ConjugateGradientsTakeNSteps | ( | int | n ) |
Take n steps in a conjugate gradient optimalization that was previously initialized with ConjugateGradientsInitialize().
n | The number of steps to take. |
void GenerateVelocities | ( | ) |
Generate starting velocities with a Maxwellian distribution.
void CorrectVelocities | ( | ) |
Correct the velocities so that the following is true:
3N ---- 0.5 \ m_i * v_i^2 = 0.5 * Ndf * kB * T = E_kin / ---- i=1 E_kin : kinetic energy m_i : mass of atom i v_i : velocity of atom i Ndf : number of degrees of freedom (3 * number of atoms) kB : Boltzmann's constant T : temperature
void MolecularDynamicsTakeNSteps | ( | int | n, |
double | T, | ||
double | timestep = 0.001 , |
||
int | method = OBFF_ANALYTICAL_GRADIENT |
||
) |
Take n steps at temperature T. If no velocities are set, they will be generated.
example:
// pFF is a pointer to a OBForceField class while (pFF->MolecularDynamicsTakeNSteps(5, 300)) { // do some updating in your program (redraw structure, ...) }
n | The number of steps to take. |
T | Absolute temperature in Kelvin. |
timestep | The time step in picoseconds. (10e-12 s) |
method | OBFF_ANALYTICAL_GRADIENTS (default) or OBFF_NUMERICAL_GRADIENTS |
OBFFConstraints & GetConstraints | ( | ) |
Get the current constraints.
void SetConstraints | ( | OBFFConstraints & | constraints ) |
Set the constraints.
constraints | The new constraints to be used. |
void SetFixAtom | ( | int | index ) |
Fix the atom position until UnsetFixAtom() is called. This function can be used in programs that allow the user to interact with a molecule that is being minimized without having to check if the atom is already fixed in the constraints set by Setup() or SetConstraints(). Using this makes sure the selected atom follows the mouse cursur.
index | The index for the atom to fix. |
void UnsetFixAtom | ( | ) |
Undo last SetFixAtom. This function will not remove the fix atom constraint for this atom if set by Setup() or SetConstraints().
void SetIgnoreAtom | ( | int | index ) |
Ignore the atom until UnsetIgnoreAtom() is called. This function can be used in programs that allow the user to interact with a molecule that is being minimized without having to check if the atom is already ignored in the constraints set by Setup() or SetConstraints(). Using this makes sure, in drawing mode, you can close rings without your newly created puching the other atoms away.
index | The index for the atom to ignore. |
void UnsetIgnoreAtom | ( | ) |
Undo last SetIgnoreAtom. This function will not remove the ignore atom constraint for this atom if set by Setup() or SetConstraints().
bool IgnoreCalculation | ( | int | a, |
int | b | ||
) | [static] |
internal function
bool IgnoreCalculation | ( | int | a, |
int | b, | ||
int | c | ||
) | [static] |
internal function
bool IgnoreCalculation | ( | int | a, |
int | b, | ||
int | c, | ||
int | d | ||
) | [static] |
internal function
bool DetectExplosion | ( | ) |
(debugging)
void ValidateSteepestDescent | ( | int | steps ) |
(debugging)
void ValidateConjugateGradients | ( | int | steps ) |
(debugging)
virtual bool Validate | ( | ) | [inline, virtual] |
Validate the force field implementation (debugging)
virtual bool ValidateGradients | ( | ) | [inline, virtual] |
Validate the analytical gradients by comparing them to numerical ones. This function has to be implemented force field specific. (debugging)
Calculate the error of the analytical gradient (debugging)
double VectorBondDerivative | ( | double * | pos_a, |
double * | pos_b, | ||
double * | force_a, | ||
double * | force_b | ||
) | [static] |
Calculate the derivative of a vector length. The vector is given by a - b, the length of this vector rab = sqrt(ab.x^2 + ab.y^2 + ab.z^2).
pos_a | atom a (coordinates) |
pos_b | atom b (coordinates) |
force_a | - return value for the force on atom a |
force_b | - return value for the force on atom b |
double VectorDistanceDerivative | ( | const double *const | pos_i, |
const double *const | pos_j, | ||
double * | force_i, | ||
double * | force_j | ||
) | [static] |
To be used for VDW or Electrostatic interactions. This is faster than VectorBondDerivative, but does no error checking.
double VectorAngleDerivative | ( | double * | pos_a, |
double * | pos_b, | ||
double * | pos_c, | ||
double * | force_a, | ||
double * | force_b, | ||
double * | force_c | ||
) | [static] |
Calculate the derivative of a angle a-b-c. The angle is given by dot(ab,cb)/rab*rcb. Used for harmonic (cubic) angle potentials.
pos_a | atom a (coordinates) |
pos_b | atom b (coordinates) |
pos_c | atom c (coordinates) |
force_a | - return value for the force on atom a |
force_b | - return value for the force on atom b |
force_c | - return value for the force on atom c |
double VectorOOPDerivative | ( | double * | pos_a, |
double * | pos_b, | ||
double * | pos_c, | ||
double * | pos_d, | ||
double * | force_a, | ||
double * | force_b, | ||
double * | force_c, | ||
double * | force_d | ||
) | [static] |
Calculate the derivative of a OOP angle a-b-c-d. b is the central atom, and a-b-c is the plane. The OOP angle is given by 90° - arccos(dot(corss(ab,cb),db)/rabbc*rdb).
pos_a | atom a (coordinates) |
pos_b | atom b (coordinates) |
pos_c | atom c (coordinates) |
pos_d | atom d (coordinates) |
force_a | - return value for the force on atom a |
force_b | - return value for the force on atom b |
force_c | - return value for the force on atom c |
force_d | - return value for the force on atom d |
double VectorTorsionDerivative | ( | double * | pos_a, |
double * | pos_b, | ||
double * | pos_c, | ||
double * | pos_d, | ||
double * | force_a, | ||
double * | force_b, | ||
double * | force_c, | ||
double * | force_d | ||
) | [static] |
Calculate the derivative of a torsion angle a-b-c-d. The torsion angle is given by arccos(dot(corss(ab,bc),cross(bc,cd))/rabbc*rbccd).
pos_a | atom a (coordinates) |
pos_b | atom b (coordinates) |
pos_c | atom c (coordinates) |
pos_d | atom d (coordinates) |
force_a | - return value for the force on atom a |
force_b | - return value for the force on atom b |
force_c | - return value for the force on atom c |
force_d | - return value for the force on atom d |
static void VectorSubtract | ( | double * | i, |
double * | j, | ||
double * | result | ||
) | [inline, static] |
inline fuction to speed up minimization speed
i | pointer to i[3] |
j | pointer to j[3] |
result | pointer to result[3], will be set to i - j |
static void VectorSubtract | ( | const double *const | i, |
const double *const | j, | ||
double * | result | ||
) | [inline, static] |
static void VectorAdd | ( | double * | i, |
double * | j, | ||
double * | result | ||
) | [inline, static] |
inline fuction to speed up minimization speed
i | pointer to i[3] |
j | pointer to j[3] |
result | pointer to result[3], will be set to i + j |
static void VectorDivide | ( | double * | i, |
double | n, | ||
double * | result | ||
) | [inline, static] |
inline fuction to speed up minimization speed
i | pointer to i[3] |
n | divide x,y,z with n |
result | pointer to result[3] |
static void VectorMultiply | ( | double * | i, |
double | n, | ||
double * | result | ||
) | [inline, static] |
inline fuction to speed up minimization speed
i | pointer to i[3] |
n | multiply x,y,z with n |
result | pointer to result[3] |
static void VectorMultiply | ( | const double *const | i, |
const double | n, | ||
double * | result | ||
) | [inline, static] |
static void VectorSelfMultiply | ( | double * | i, |
double | n | ||
) | [inline, static] |
inline fuction to speed up minimization speed
i | pointer to i[3], multiply this vector by n and set this vector to the result. |
n | the scalar value to be multipled |
static void VectorNormalize | ( | double * | i ) | [inline, static] |
inline fuction to speed up minimization speed
i | pointer to i[3] to be normalized |
static void VectorCopy | ( | double * | from, |
double * | to | ||
) | [inline, static] |
inline fuction to speed up minimization speed
from | pointer to i[3] to be copied from |
to | pointer to j[3] to be copied to |
static double VectorLength | ( | double * | i ) | [inline, static] |
inline fuction to speed up minimization speed
i | pointer to i[3] |
static double VectorDistance | ( | double * | pos_i, |
double * | pos_j | ||
) | [inline, static] |
double VectorAngle | ( | double * | i, |
double * | j, | ||
double * | k | ||
) | [static] |
inline fuction to speed up minimization speed
i | pointer to i[3] |
j | pointer to j[3] |
k | pointer to k[3] |
double VectorTorsion | ( | double * | i, |
double * | j, | ||
double * | k, | ||
double * | l | ||
) | [static] |
inline fuction to speed up minimization speed
i | pointer to i[3] |
j | pointer to j[3] |
k | pointer to k[3] |
l | pointer to l[3] |
double VectorOOP | ( | double * | i, |
double * | j, | ||
double * | k, | ||
double * | l | ||
) | [static] |
inline fuction to speed up minimization speed
i | pointer to i[3] |
j | pointer to j[3] |
k | pointer to k[3] |
l | pointer to l[3] |
static void VectorClear | ( | double * | i ) | [inline, static] |
inline fuction to speed up minimization speed
i | pointer to i[3], will set x,y,z to 0,0,0 |
static double VectorDot | ( | double * | i, |
double * | j | ||
) | [inline, static] |
inline fuction to speed up minimization speed
i | pointer to i[3] |
j | pointer to j[3] |
static void VectorCross | ( | double * | i, |
double * | j, | ||
double * | result | ||
) | [inline, static] |
inline fuction to speed up minimization speed
i | pointer to i[3] |
j | pointer to j[3] |
result | the dot product (as a return value double[3]) |
static void PrintVector | ( | double * | i ) | [inline, static] |
bool _init [protected] |
Used to make sure we only parse the parameter file once, when needed.
std::string _parFile [protected] |
bool _validSetup [protected] |
< parameter file name
was the last call to Setup succesfull
double* _gradientPtr [protected] |
pointer to the gradients (used by AddGradient(), minimization functions, ...)
std::ostream* _logos [protected] |
Output for logfile.
char _logbuf[BUFF_SIZE+1] [protected] |
Temporary buffer for logfile output.
int _loglvl [protected] |
Log level for output.
int _origLogLevel [protected] |
int _current_conformer [protected] |
used to hold i for current conformer (needed by UpdateConformers)
std::vector<double> _energies [protected] |
used to hold the energies for all conformers
double _econv [protected] |
double _e_n1 [protected] |
Used for conjugate gradients and steepest descent(Initialize and TakeNSteps)
int _cstep [protected] |
int _nsteps [protected] |
Used for conjugate gradients and steepest descent(Initialize and TakeNSteps)
double* _grad1 [protected] |
Used for conjugate gradients and steepest descent(Initialize and TakeNSteps)
unsigned int _ncoords [protected] |
Number of coordinates for conjugate gradients.
int _linesearch [protected] |
LineSearch type.
double _timestep [protected] |
Molecular dynamics time step in picoseconds.
double _temp [protected] |
Molecular dynamics temperature in Kelvin.
double* _velocityPtr [protected] |
pointer to the velocities
OBFFConstraints _constraints = OBFFConstraints() [static, protected] |
Constraints.
int _fixAtom = 0 [static, protected] |
SetFixAtom()/UnsetFixAtom()
int _ignoreAtom = 0 [static, protected] |
SetIgnoreAtom()/UnsetIgnoreAtom()
bool _cutoff [protected] |
true = cut-off enabled
double _rvdw [protected] |
VDW cut-off distance.
double _rele [protected] |
Electrostatic cut-off distance.
int _pairfreq [protected] |
The frequence to update non-bonded pairs.
std::vector<OBBitVec> _intraGroup [protected] |
groups for which intra-molecular interactions should be calculated
std::vector<OBBitVec> _interGroup [protected] |
groups for which intra-molecular interactions should be calculated
std::vector<std::pair<OBBitVec, OBBitVec> > _interGroups [protected] |
groups for which intra-molecular interactions should be calculated