OBForceField Class Reference

Base class for molecular mechanics force fields. More...

#include <openbabel/forcefield.h>

List of all members.

Public Member Functions

virtual ~OBForceField ()
virtual std::string GetUnit ()
virtual bool Setup (OBMol &mol)
bool UpdateCoordinates (OBMol &mol)
bool UpdateConformers (OBMol &mol)
void OBFFLog (std::string msg)
void OBFFLog (const char *msg)
void kludge ()
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
bool SetLogFile (std::ostream *pos)
bool SetLogLevel (int level)
int GetLogLevel ()
Methods for structure generation
void SystematicRotorSearch ()
Methods for energy minimization
vector3 LineSearch (OBAtom *atom, vector3 &direction)
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 forcefield validation
virtual bool Validate ()
virtual bool ValidateGradients ()
vector3 ValidateGradientError (vector3 &numgrad, vector3 &anagrad)

Static Public Member Functions

static OBForceFieldFindForceField (const std::string &ID)
static OBForceFieldFindForceField (const char *ID)
Methods for vector analysis (used by OBFFXXXXCalculationYYYY)
static double VectorLengthDerivative (vector3 &a, vector3 &b)
static double VectorAngleDerivative (vector3 &a, vector3 &b, vector3 &c)
static double VectorTorsionDerivative (vector3 &a, vector3 &b, vector3 &c, vector3 &d)

Protected Member Functions

OBFFParameterGetParameter (int a, int b, int c, int d, std::vector< OBFFParameter > &parameter)
OBFFParameterGetParameter (const char *a, const char *b, const char *c, const char *d, std::vector< OBFFParameter > &parameter)
vector3 NumericalDerivative (OBAtom *a, int terms=OBFF_ENERGY)
vector3 NumericalSecondDerivative (OBAtom *a, int terms=OBFF_ENERGY)
virtual vector3 GetGradient (OBAtom *a, int terms=OBFF_ENERGY)
bool IsInSameRing (OBAtom *a, OBAtom *b)

Protected Attributes

OBMol _mol
std::ostream * logos
char logbuf [200]
int loglvl
int current_conformer
double _econv
double _e_n1
int _method
int _cstep
int _nsteps
std::vector< vector3_grad1
std::vector< vector3_dir1


Detailed Description

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, ...). Other classes such as OBFFParameter, 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 can always be used to find available forcefields:

      FOR_EACH(OBForceField, iter) {
      cout << "forcefield ID: " << iter.ID() << endl;
      }

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("Ghemical");

      // 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("Ghemical");

      // 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);


Constructor & Destructor Documentation

virtual ~OBForceField (  )  [inline, virtual]

Destructor.


Member Function Documentation

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]

see GetParameter(int a, int b, int c, int d, std::vector<OBFFParameter> &parameter)

vector3 NumericalDerivative ( OBAtom a,
int  terms = OBFF_ENERGY 
) [protected]

Calculate the potential energy function derivative numerically with repect to the coordinates of atom with index a (this vector is the gradient).

Parameters:
a provides coordinates
terms OBFF_ENERGY, OBFF_EBOND, OBFF_EANGLE, OBFF_ESTRBND, OBFF_ETORSION, OBFF_EOOP, OBFF_EVDW, OBFF_ELECTROSTATIC
Returns:
the negative gradient of atom a

vector3 NumericalSecondDerivative ( OBAtom a,
int  terms = OBFF_ENERGY 
) [protected]

OB 3.0.

virtual vector3 GetGradient ( OBAtom a,
int  terms = OBFF_ENERGY 
) [inline, protected, virtual]

Calculate the potential energy function derivative analyticaly with repect to the coordinates of atom with index a (this vector is the gradient)

If the currently selected forcefield doesn't have analytical gradients, we can still call this function which will return the result of NumericalDerivative()

Parameters:
a provides coordinates
terms OBFF_ENERGY, OBFF_EBOND, OBFF_EANGLE, OBFF_ESTRBND, OBFF_ETORSION, OBFF_EOOP, OBFF_EVDW, OBFF_ELECTROSTATIC
Returns:
the negative gradient of atom a

bool IsInSameRing ( OBAtom a,
OBAtom b 
) [protected]

Check if two atoms are in the same ring

Parameters:
a atom a
b atom b
Returns:
true if atom a and b are in the same ring

static OBForceField* FindForceField ( const std::string &  ID  )  [inline, static]

short description of the force field type.

Parameters:
ID forcefield id (Ghemical, ...)
Returns:
A pointer to a forcefield (the default if ID is empty), or NULL if not available

static OBForceField* FindForceField ( const char *  ID  )  [inline, static]

Parameters:
ID forcefield id (Ghemical, ...)
Returns:
A pointer to a forcefield (the default if ID is empty), or NULL if not available

virtual std::string GetUnit (  )  [inline, virtual]

Returns:
The unit (kcal/mol, kJ/mol, ...) in which the energy is expressed as std::string

virtual bool Setup ( OBMol mol  )  [inline, virtual]

Setup the forcefield for mol (assigns atom types, charges, etc.)

Parameters:
mol the OBMol object that contains the atoms and bonds
Returns:
True if succesfull

bool UpdateCoordinates ( OBMol mol  ) 

Update coordinates for current conformer

Parameters:
mol the OBMol object to copy the coordinates to
Returns:
true if succesfull

bool UpdateConformers ( OBMol mol  ) 

Update coordinates for all conformers

Parameters:
mol the OBMol object to copy the coordinates to
Returns:
true if succesfull

void OBFFLog ( std::string  msg  )  [inline]

Print msg to the logfile

Parameters:
msg the message

void OBFFLog ( const char *  msg  )  [inline]

Print msg to the logfile

Parameters:
msg the message

virtual double Energy ( bool  gradients = true  )  [inline, virtual]

Parameters:
gradients Set to true when the gradients need to be calculated (needs to be done before calling GetGradient())
Returns:
Total energy
Output to log:
OBFF_LOGLVL_NONE: none
OBFF_LOGLVL_LOW: none
OBFF_LOGLVL_MEDIUM: energy for indivudual energy terms
OBFF_LOGLVL_HIGH: energy for individual energy interactions

virtual double E_Bond ( bool  gradients = true  )  [inline, virtual]

Parameters:
gradients Set to true when the gradients need to be calculated (needs to be done before calling GetGradient())
Returns:
Bond stretching energy
Output to log:
see Energy()

virtual double E_Angle ( bool  gradients = true  )  [inline, virtual]

Parameters:
gradients Set to true when the gradients need to be calculated (needs to be done before calling GetGradient())
Returns:
Angle bending energy
Output to log:
see Energy()

virtual double E_StrBnd ( bool  gradients = true  )  [inline, virtual]

Parameters:
gradients Set to true when the gradients need to be calculated (needs to be done before calling GetGradient())
Returns:
Stretch bending energy
Output to log:
see Energy()

virtual double E_Torsion ( bool  gradients = true  )  [inline, virtual]

Parameters:
gradients Set to true when the gradients need to be calculated (needs to be done before calling GetGradient())
Returns:
Torsional energy
Output to log:
see Energy()

virtual double E_OOP ( bool  gradients = true  )  [inline, virtual]

Parameters:
gradients Set to true when the gradients need to be calculated (needs to be done before calling GetGradient())
Returns:
Out-Of-Plane bending energy
Output to log:
see Energy()

virtual double E_VDW ( bool  gradients = true  )  [inline, virtual]

Parameters:
gradients Set to true when the gradients need to be calculated (needs to be done before calling GetGradient())
Returns:
Van der Waals energy
Output to log:
see Energy()

virtual double E_Electrostatic ( bool  gradients = true  )  [inline, virtual]

Parameters:
gradients Set to true when the gradients need to be calculated (needs to be done before calling GetGradient())
Returns:
Electrostatic energy
Output to log:
see Energy()

bool SetLogFile ( std::ostream *  pos  ) 

Set the stream for logging (can also be &cout for logging to screen)

Parameters:
pos stream
Returns:
True if succesfull

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 {
      *logos << "this text will NOT be logged..." << endl
      }
   
      IF_OBFF_LOGLVL_LOW {
      *logos << "this text will be logged..." << endl
      }
  
      IF_OBFF_LOGLVL_MEDIUM {
      *logos << "this text will also be logged..." << endl
      }

int GetLogLevel (  )  [inline]

Returns:
log level

void SystematicRotorSearch (  ) 

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.

Output to log:
This function should only be called with the log level set to OBFF_LOGLVL_NONE or OBFF_LOGLVL_LOW. Otherwise too much information about the energy calculations needed for this function will interfere with the output for this function.

OBFF_LOGLVL_NONE: none
OBFF_LOGLVL_LOW: number of rotatable bonds, energies for the conformers, which one is the lowest, ...
OBFF_LOGLVL_MEDIUM: see note above
OBFF_LOGLVL_HIGH: see note above

vector3 LineSearch ( OBAtom atom,
vector3 direction 
)

Perform a linesearch starting at atom in direction direction

Parameters:
atom start coordinates
direction the search direction
Returns:
vector which starts at atom and stops at the minimum (the length is the ideal stepsize)
Output to log:
OBFF_LOGLVL_NONE: none
OBFF_LOGLVL_LOW: none
OBFF_LOGLVL_MEDIUM: none
OBFF_LOGLVL_HIGH: none

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.

Parameters:
steps the number of steps
econv energy convergence criteria (defualt is 1e-6)
method OBFF_ANALYTICAL_GRADIENTS (default) or OBFF_NUMERICAL_GRADIENTS
Output to log:
This function should only be called with the log level set to OBFF_LOGLVL_NONE or OBFF_LOGLVL_LOW. Otherwise too much information about the energy calculations needed for the minimization will interfere with the list of energies for succesive steps.

OBFF_LOGLVL_NONE: none
OBFF_LOGLVL_LOW: header including number of steps and first step
OBFF_LOGLVL_MEDIUM: see note above
OBFF_LOGLVL_HIGH: see note above

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.

Parameters:
steps the number of steps
econv energy convergence criteria (defualt is 1e-6)
method OBFF_ANALYTICAL_GRADIENTS (default) or OBFF_NUMERICAL_GRADIENTS
Output to log:
This function should only be called with the log level set to OBFF_LOGLVL_NONE or OBFF_LOGLVL_LOW. Otherwise too much information about the energy calculations needed for the minimization will interfere with the list of energies for succesive steps.

OBFF_LOGLVL_NONE: none
OBFF_LOGLVL_LOW: header including number of steps
OBFF_LOGLVL_MEDIUM: see note above
OBFF_LOGLVL_HIGH: see note above

bool SteepestDescentTakeNSteps ( int  n  ) 

Take n steps in a steepestdescent optimalization that was previously initialized with SteepestDescentInitialize().

Parameters:
n the number of steps to take
Returns:
false if convergence or the number of steps given by SteepestDescentInitialize() has been reached
Output to log:
This function should only be called with the log level set to OBFF_LOGLVL_NONE or OBFF_LOGLVL_LOW. Otherwise too much information about the energy calculations needed for the minimization will interfere with the list of energies for succesive steps.

OBFF_LOGLVL_NONE: none
OBFF_LOGLVL_LOW: step number, energy and energy for the previous step
OBFF_LOGLVL_MEDIUM: see note above
OBFF_LOGLVL_HIGH: see note above

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.

Parameters:
steps the number of steps
econv energy convergence criteria (defualt is 1e-6)
method OBFF_ANALYTICAL_GRADIENTS (default) or OBFF_NUMERICAL_GRADIENTS
Output to log:
This function should only be called with the log level set to OBFF_LOGLVL_NONE or OBFF_LOGLVL_LOW. Otherwise too much information about the energy calculations needed for the minimization will interfere with the list of energies for succesive steps.

OBFF_LOGLVL_NONE: none
OBFF_LOGLVL_LOW: information about the progress of the minimization
OBFF_LOGLVL_MEDIUM: see note above
OBFF_LOGLVL_HIGH: see note above

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.

Parameters:
steps the number of steps
econv energy convergence criteria (defualt is 1e-6)
method OBFF_ANALYTICAL_GRADIENTS (default) or OBFF_NUMERICAL_GRADIENTS
Output to log:
This function should only be called with the log level set to OBFF_LOGLVL_NONE or OBFF_LOGLVL_LOW. Otherwise too much information about the energy calculations needed for the minimization will interfere with the list of energies for succesive steps.

OBFF_LOGLVL_NONE: none
OBFF_LOGLVL_LOW: header including number of steps and first step
OBFF_LOGLVL_MEDIUM: see note above
OBFF_LOGLVL_HIGH: see note above

bool ConjugateGradientsTakeNSteps ( int  n  ) 

Take n steps in a conjugate gradient optimalization that was previously initialized with ConjugateGradientsInitialize().

Parameters:
n the number of steps to take
Returns:
false if convergence or the number of steps given by ConjugateGradientsInitialize() has been reached
Output to log:
This function should only be called with the log level set to OBFF_LOGLVL_NONE or OBFF_LOGLVL_LOW. Otherwise too much information about the energy calculations needed for the minimization will interfere with the list of energies for succesive steps.

OBFF_LOGLVL_NONE: none
OBFF_LOGLVL_LOW: step number, energy and energy for the previous step
OBFF_LOGLVL_MEDIUM: see note above
OBFF_LOGLVL_HIGH: see note above

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)

vector3 ValidateGradientError ( vector3 numgrad,
vector3 anagrad 
)

Calculate the error of the analytical gradient (debugging)

Returns:
error = fabs(numgrad - anagrad) / anagrad * 100%

double VectorLengthDerivative ( vector3 a,
vector3 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).

Parameters:
a atom a (coordinates), will be changed to -drab/da
b atom b (coordinates), will be changed to -drab/db
Returns:
The distance between a and b (bondlength for bond stretching, separation for vdw, electrostatic)

double VectorAngleDerivative ( vector3 a,
vector3 b,
vector3 c 
) [static]

Calculate the derivative of a angle a-b-c. The angle is given by dot(ab,cb)/rab*rcb.

Parameters:
a atom a (coordinates), will be changed to -dtheta/da
b atom b (coordinates), will be changed to -dtheta/db
c atom c (coordinates), will be changed to -dtheta/dc
Returns:
The angle between a-b-c

double VectorTorsionDerivative ( vector3 a,
vector3 b,
vector3 c,
vector3 d 
) [static]

Calculate the derivative of a torsion angle a-b-c-d. The torsion angle is given by dot(corss(ab,bc),cross(bc,cd)/rabbc*rbccd.

Parameters:
a atom a (coordinates), will be changed to -dtheta/da
b atom b (coordinates), will be changed to -dtheta/db
c atom c (coordinates), will be changed to -dtheta/dc
d atom d (coordinates), will be changed to -dtheta/dd
Returns:
The tosion angle for a-b-c-d

void kludge (  )  [inline]

Do not use. This function contains rubbish merely to ensure the compiler instantiates some templated functions which are needed for the Windows Python build. TODO Find the proper way of doing this.


Member Data Documentation

OBMol _mol [protected]

Molecule to be evaluated or minimized.

std::ostream* logos [protected]

Output for logfile.

char logbuf[200] [protected]

int loglvl [protected]

Log level for output.

int current_conformer [protected]

used to hold i for current conformer (needed by UpdateConformers)

double _econv [protected]

Used for conjugate gradients and steepest descent(Initialize and TakeNSteps).

double _e_n1 [protected]

int _method [protected]

int _cstep [protected]

int _nsteps [protected]

std::vector<vector3> _grad1 [protected]

std::vector<vector3> _dir1 [protected]


The documentation for this class was generated from the following files: