Open Babel  3.0
Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Static Protected Member Functions | Protected Attributes | Static Protected Attributes | List of all members
OBForceField Class Referenceabstract

#include <openbabel/forcefield.h>

Inheritance diagram for OBForceField:
OBPlugin

Public Types

typedef std::map< const char *, OBPlugin *, CharPtrLessPluginMapType
 
typedef PluginMapType::const_iterator PluginIterator
 

Public Member Functions

virtual OBForceFieldMakeNewInstance ()=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)
 
OBGridDataGetGrid (double step, double padding, const char *type, double pchg)
 
virtual vector3 GetGradient (OBAtom *a, int=OBFF_ENERGY)
 
double * GetGradientPtr ()
 
virtual const char * Description ()
 
virtual bool Display (std::string &txt, const char *param, const char *ID=NULL)
 
virtual OBPluginMakeInstance (const std::vector< std::string > &)
 
virtual void Init ()
 
const char * GetID () const
 
virtual PluginMapTypeGetMap () 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 SetDielectricConstant (double epsilon)
 
double GetDielectricConstant ()
 
void SetUpdateFrequency (int f)
 
int GetUpdateFrequency ()
 
void UpdatePairsSimple ()
 
unsigned int GetNumPairs ()
 
unsigned int GetNumElectrostaticPairs ()
 
unsigned int GetNumVDWPairs ()
 
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, bool sampleRingBonds=false)
 
int SystematicRotorSearchInitialize (unsigned int geomSteps=2500, bool sampleRingBonds=false)
 
bool SystematicRotorSearchNextConformer (unsigned int geomSteps=2500)
 
void RandomRotorSearch (unsigned int conformers, unsigned int geomSteps=2500, bool sampleRingBonds=false)
 
void RandomRotorSearchInitialize (unsigned int conformers, unsigned int geomSteps=2500, bool sampleRingBonds=false)
 
bool RandomRotorSearchNextConformer (unsigned int geomSteps=2500)
 
void WeightedRotorSearch (unsigned int conformers, unsigned int geomSteps, bool sampleRingBonds=false)
 
int FastRotorSearch (bool permute=true)
 
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 OBForceFieldFindForceField (const std::string &ID)
 
static OBForceFieldFindForceField (const char *ID)
 
static OBPluginGetPlugin (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)
 
static void LoadAllPlugins ()
 
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

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)
 
int GetParameterIdx (int a, int b, int c, int d, std::vector< OBFFParameter > &parameter)
 
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 void ClearGradients ()
 
bool IsInSameRing (OBAtom *a, OBAtom *b)
 

Static Protected Member Functions

static PluginMapTypePluginMap ()
 
static PluginMapTypeGetTypeMap (const char *PluginID)
 
static OBPluginBaseFindType (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 _gconv
 
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
 
double _epsilon
 
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 unsigned int _fixAtom = 0
 
static unsigned int _ignoreAtom = 0
 
static int AllPluginsLoaded = 0
 

Methods for constraints

OBFFConstraintsGetConstraints ()
 
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)
 

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, 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/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/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/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->ConjugateGradients(1000);

Minimize a ligand molecule in a binding pocket.

#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.BitIsSet(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->ConjugateGradients(1000);
Examples:
obforcefield_energy.cpp.

Constructor & Destructor Documentation

◆ ~OBForceField()

virtual ~OBForceField ( )
inlinevirtual

Destructor.

Member Function Documentation

◆ GetParameter() [1/2]

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, ...

◆ GetParameter() [2/2]

OBFFParameter * GetParameter ( const char *  a,
const char *  b,
const char *  c,
const char *  d,
std::vector< OBFFParameter > &  parameter 
)
protected

◆ GetParameterIdx()

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

Get index for vector<OBFFParameter> ...

◆ NumericalDerivative()

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
aprovides coordinates
termsOBFF_ENERGY, OBFF_EBOND, OBFF_EANGLE, OBFF_ESTRBND, OBFF_ETORSION, OBFF_EOOP, OBFF_EVDW, OBFF_ELECTROSTATIC
Returns
the negative gradient of atom a

◆ NumericalSecondDerivative()

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

OB 3.0.

◆ SetGradient()

void SetGradient ( double *  grad,
int  idx 
)
inlineprotected

Set the gradient for atom with index idx to grad

◆ AddGradient()

void AddGradient ( double *  grad,
int  idx 
)
inlineprotected

Add grad to the gradient for atom with index idx

◆ ClearGradients()

virtual void ClearGradients ( )
inlineprotectedvirtual

Set all gradients to zero

◆ IsInSameRing()

bool IsInSameRing ( OBAtom a,
OBAtom b 
)
protected

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.]

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

◆ MakeNewInstance()

virtual OBForceField* MakeNewInstance ( )
pure virtual

Clone the current instance. May be desirable in multithreaded environments, Should be deleted after use

◆ TypeID()

const char* TypeID ( )
inlinevirtual
Returns
Plugin type ("forcefields")

Reimplemented from OBPlugin.

◆ FindForceField() [1/2]

static OBForceField* FindForceField ( const std::string &  ID)
inlinestatic
Parameters
IDforcefield id (Ghemical, MMFF94, UFF, ...).
Returns
A pointer to a forcefield (the default if ID is empty), or NULL if not available.

◆ FindForceField() [2/2]

static OBForceField* FindForceField ( const char *  ID)
inlinestatic
Parameters
IDforcefield id (Ghemical, MMFF94, UFF, ...).
Returns
A pointer to a forcefield (the default if ID is empty), or NULL if not available.

◆ SetParameterFile()

void SetParameterFile ( const std::string &  filename)
inline

◆ GetUnit()

virtual std::string GetUnit ( )
inlinevirtual
Returns
The unit (kcal/mol, kJ/mol, ...) in which the energy is expressed as std::string.
Examples:
obforcefield_energy.cpp.

◆ HasAnalyticalGradients()

virtual bool HasAnalyticalGradients ( )
inlinevirtual

◆ Setup() [1/2]

bool Setup ( OBMol mol)

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

Parameters
molThe OBMol object that contains the atoms and bonds.
Returns
True if succesfull.
Examples:
obforcefield_energy.cpp.

Referenced by OBEnergyConformerScore::Score(), OBMinimizingEnergyConformerScore::Score(), and OBMinimizingRMSDConformerScore::Score().

◆ Setup() [2/2]

bool Setup ( OBMol mol,
OBFFConstraints constraints 
)

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

Parameters
molThe OBMol object that contains the atoms and bonds.
constraintsThe OBFFConstraints object that contains the constraints.
Returns
True if succesfull.

◆ ParseParamFile()

virtual bool ParseParamFile ( )
inlinevirtual

Load the parameters (this function is overloaded by the individual forcefields, and is called autoamically from OBForceField::Setup()).

◆ SetTypes()

virtual bool SetTypes ( )
inlinevirtual

Set the atom types (this function is overloaded by the individual forcefields, and is called autoamically from OBForceField::Setup()).

◆ SetFormalCharges()

virtual bool SetFormalCharges ( )
inlinevirtual

Set the formal charges (this function is overloaded by the individual forcefields, and is called autoamically from OBForceField::Setup()).

◆ SetPartialCharges()

virtual bool SetPartialCharges ( )
inlinevirtual

Set the partial charges (this function is overloaded by the individual forcefields, and is called autoamically from OBForceField::Setup()).

◆ SetupCalculations()

virtual bool SetupCalculations ( )
inlinevirtual

Setup the calculations (this function is overloaded by the individual forcefields, and is called autoamically from OBForceField::Setup()).

◆ SetupPointers()

virtual bool SetupPointers ( )
inlinevirtual

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).

◆ IsSetupNeeded()

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.

Returns
True if Setup needs to be called.

◆ GetAtomTypes()

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;
}
...

◆ GetPartialCharges()

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;
}
...

◆ GetCoordinates()

bool GetCoordinates ( OBMol mol)

Get coordinates for current conformer and attach OBConformerData with energies, forces, ... to mol.

Parameters
molThe OBMol object to copy the coordinates to (from OBForceField::_mol).
Returns
True if succesfull.

◆ UpdateCoordinates()

bool UpdateCoordinates ( OBMol mol)
inline
Deprecated:
Use GetCooordinates instead.

◆ GetConformers()

bool GetConformers ( OBMol mol)

Get coordinates for all conformers and attach OBConformerData with energies, forces, ... to mol.

Parameters
molThe OBMol object to copy the coordinates to (from OBForceField::_mol).
Returns
True if succesfull.

◆ UpdateConformers()

bool UpdateConformers ( OBMol mol)
inline
Deprecated:
Use GetConformers instead.

◆ SetCoordinates()

bool SetCoordinates ( OBMol mol)

Set coordinates for current conformer.

Parameters
molthe OBMol object to copy the coordinates from (to OBForceField::_mol).
Returns
true if succesfull.

◆ SetConformers()

bool SetConformers ( OBMol mol)

Set coordinates for all conformers.

Parameters
molThe OBMol object to copy the coordinates from (to OBForceField::_mol).
Returns
True if succesfull.

◆ GetGrid()

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.

Parameters
stepThe grid step size in A..
paddingThe padding for the grid in A.
typeThe force field atom type for the probe.
pchgThe partial charge for the probe atom.
Returns
Pointer to the grid constaining the results.

◆ AddIntraGroup()

void AddIntraGroup ( OBBitVec group)

Enable intra-molecular interactions for group (bonds, angles, strbnd, torsions, oop). This function should be called before Setup().

Parameters
groupOBBitVec with bits set for the indexes of the atoms which make up the group.

◆ AddInterGroup()

void AddInterGroup ( OBBitVec group)

Enable inter-molecular interactions for group (non-bonded: vdw & ele). This function should be called before Setup().

Parameters
groupOBBitVec with bits set for the indexes of the atoms which make up the group.

◆ AddInterGroups()

void AddInterGroups ( OBBitVec group1,
OBBitVec group2 
)

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().

Parameters
group1OBBitVec with bits set for the indexes of the atoms which make up the first group.
group2OBBitVec with bits set for the indexes of the atoms which make up the second group.

◆ ClearGroups()

void ClearGroups ( )

Clear all previously specified groups.

◆ HasGroups()

bool HasGroups ( )
Returns
true if there are groups.

◆ EnableCutOff()

void EnableCutOff ( bool  enable)
inline

Enable or disable Cut-offs. Cut-offs are disabled by default.

Parameters
enableEnable when true, disable when false.

◆ IsCutOffEnabled()

bool IsCutOffEnabled ( )
inline
Returns
True if Cut-off distances are used.

◆ SetVDWCutOff()

void SetVDWCutOff ( double  r)
inline

Set the VDW cut-off distance to r. Note that this does not enable cut-off distances.

Parameters
rThe VDW cut-off distance to be used in A.

◆ GetVDWCutOff()

double GetVDWCutOff ( )
inline

Get the VDW cut-off distance.

Returns
The VDW cut-off distance in A.

◆ SetElectrostaticCutOff()

void SetElectrostaticCutOff ( double  r)
inline

Set the Electrostatic cut-off distance to r. Note that this does not enable cut-off distances.

Parameters
rThe electrostatic cut-off distance to be used in A.

◆ GetElectrostaticCutOff()

double GetElectrostaticCutOff ( )
inline

Get the Electrostatic cut-off distance.

Returns
The electrostatic cut-off distance in A.

◆ SetDielectricConstant()

void SetDielectricConstant ( double  epsilon)
inline

Set the dielectric constant for electrostatic SetupCalculations

Parameters
epsilonThe relative permittivity to use (default = 1.0)

◆ GetDielectricConstant()

double GetDielectricConstant ( )
inline

◆ SetUpdateFrequency()

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.

Parameters
fThe pair list update frequency.

◆ GetUpdateFrequency()

int GetUpdateFrequency ( )
inline

Get the frequency by which non-bonded pairs are updated.

Returns
The pair list update frequency.

◆ UpdatePairsSimple()

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.

Todo:
set the criteria as squared values

◆ GetNumPairs()

unsigned int GetNumPairs ( )

Get the number of non-bonded pairs in _mol.

Returns
The number of atom pairs (ignores cutoff)

◆ GetNumElectrostaticPairs()

unsigned int GetNumElectrostaticPairs ( )

Get the number of enabled electrostatic pairs in _mol.

Returns
The number of pairs currently enabled (within cut-off distance)

◆ GetNumVDWPairs()

unsigned int GetNumVDWPairs ( )

Get the number of enabled VDW pairs in _mol.

Returns
The number of pairs currently enabled (within cut-off distance)

◆ EnableAllPairs()

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.

◆ GetGradient()

virtual vector3 GetGradient ( OBAtom a,
int  = OBFF_ENERGY 
)
inlinevirtual

Get the pointer to the gradients

◆ GetGradientPtr()

double* GetGradientPtr ( )
inline

Get the pointer to the gradients

◆ Energy()

virtual double Energy ( bool  gradients = true)
inlinevirtual
Parameters
gradientsSet 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 individual energy terms
OBFF_LOGLVL_HIGH: energy for individual energy interactions
Examples:
obforcefield_energy.cpp.

Referenced by OBEnergyConformerScore::Score(), OBMinimizingEnergyConformerScore::Score(), and OBMinimizingRMSDConformerScore::Score().

◆ E_Bond()

virtual double E_Bond ( bool  gradients = true)
inlinevirtual
Parameters
gradientsSet 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()

◆ E_Angle()

virtual double E_Angle ( bool  gradients = true)
inlinevirtual
Parameters
gradientsSet 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()

◆ E_StrBnd()

virtual double E_StrBnd ( bool  gradients = true)
inlinevirtual
Parameters
gradientsSet 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()

◆ E_Torsion()

virtual double E_Torsion ( bool  gradients = true)
inlinevirtual
Parameters
gradientsSet to true when the gradients need to be calculated (needs to be done before calling GetGradient()).
Returns
Torsional energy.
Output to log:
see Energy()

◆ E_OOP()

virtual double E_OOP ( bool  gradients = true)
inlinevirtual
Parameters
gradientsSet 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()

◆ E_VDW()

virtual double E_VDW ( bool  gradients = true)
inlinevirtual
Parameters
gradientsSet 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()

◆ E_Electrostatic()

virtual double E_Electrostatic ( bool  gradients = true)
inlinevirtual
Parameters
gradientsSet to true when the gradients need to be calculated (needs to be done before calling GetGradient()).
Returns
Electrostatic energy.
Output to log:
see Energy()

◆ PrintTypes()

void PrintTypes ( )

Print the atom types to the log.

◆ PrintFormalCharges()

void PrintFormalCharges ( )

Print the formal charges to the log (atom.GetPartialCharge(), MMFF94 FC's are not always int).

◆ PrintPartialCharges()

void PrintPartialCharges ( )

Print the partial charges to the log.

◆ PrintVelocities()

void PrintVelocities ( )

Print the velocities to the log.

◆ SetLogFile()

bool SetLogFile ( std::ostream *  pos)

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

Parameters
posStream (when pos is 0, std::cout wil be used).
Returns
True if succesfull.

◆ SetLogLevel()

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:

OBFFLog("this text will NOT be logged...\n");
}
OBFFLog"this text will be logged...\n");
}
OBFFLog("this text will also be logged...\n");
}

◆ GetLogLevel()

int GetLogLevel ( )
inline
Returns
The log level.

◆ OBFFLog() [1/2]

void OBFFLog ( std::string  msg)
inline

Print msg to the logfile.

Parameters
msgThe message to print.

◆ OBFFLog() [2/2]

void OBFFLog ( const char *  msg)
inline

Print msg to the logfile.

Parameters
msgThe message to print.

◆ DistanceGeometry()

void DistanceGeometry ( )

Generate coordinates for the molecule (distance geometry)

Deprecated:
Use OBDistanceGeometry class instead

◆ SystematicRotorSearch()

void SystematicRotorSearch ( unsigned int  geomSteps = 2500,
bool  sampleRingBonds = false 
)

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.

Parameters
geomStepsThe number of steps to take during geometry optimization.
sampleRingBondsWhether to sample ring torsions.
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.

◆ SystematicRotorSearchInitialize()

int SystematicRotorSearchInitialize ( unsigned int  geomSteps = 2500,
bool  sampleRingBonds = false 
)

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.

Parameters
geomStepsThe number of steps to take during geometry optimization.
sampleRingBondsWhether to sample ring torsions.
Returns
The number of conformers.

◆ SystematicRotorSearchNextConformer()

bool SystematicRotorSearchNextConformer ( unsigned int  geomSteps = 2500)

Evaluate the next conformer.

Parameters
geomStepsThe number of steps to take during geometry optimization.
Returns
True if there are more conformers.

◆ RandomRotorSearch()

void RandomRotorSearch ( unsigned int  conformers,
unsigned int  geomSteps = 2500,
bool  sampleRingBonds = false 
)

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.

Parameters
conformersThe number of random conformers to consider during the search.
geomStepsThe number of steps to take during geometry optimization for each conformer.
sampleRingBondsWhether to sample ring torsions.
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.

◆ RandomRotorSearchInitialize()

void RandomRotorSearchInitialize ( unsigned int  conformers,
unsigned int  geomSteps = 2500,
bool  sampleRingBonds = false 
)

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.

Parameters
conformersThe number of random conformers to consider during the search
geomStepsThe number of steps to take during geometry optimization
sampleRingBondsWhether to sample ring torsions.

◆ RandomRotorSearchNextConformer()

bool RandomRotorSearchNextConformer ( unsigned int  geomSteps = 2500)

Evaluate the next conformer.

Parameters
geomStepsThe number of steps to take during geometry optimization.
Returns
True if there are more conformers.

◆ WeightedRotorSearch()

void WeightedRotorSearch ( unsigned int  conformers,
unsigned int  geomSteps,
bool  sampleRingBonds = false 
)

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.

Parameters
conformersThe number of random conformers to consider during the search.
geomStepsThe number of steps to take during geometry optimization for each conformer.
sampleRingBondsWhether to sample ring torsions.
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.

◆ FastRotorSearch()

int FastRotorSearch ( bool  permute = true)

A fast rotor search to find low energy conformations.

Iterate over each of the rotors, and set the torsion angle to that which minimizes the energy (while keeping the rest of the molecule fixed). In general (for molecules with more than one rotatable bond), this procedure will not find the global minimum, but it will at least get rid of any bad clashes, and it do so quickly.

Torsions closer to the center of the molecule will be optimized first as these most likely to generate large clashes.

One possible use of this procedure is to prepare a reasonable 3D structure of a molecule for viewing. Another is to prepare the starting structure for a more systematic rotor search (in which case you should geometry optimize the final structure).

Parameters
permuteWhether or not to permute the order of the 4 most central rotors. Default is true. This does a more thorough search, but takes 4! = 24 times as long.
Since
version 2.4

◆ SetLineSearchType()

void SetLineSearchType ( int  type)
inline

Set the LineSearchType. The default type is LineSearchType::Newton2Num.

Parameters
typeThe LineSearchType to be used in SteepestDescent and ConjugateGradients.

◆ GetLineSearchType()

int GetLineSearchType ( )
inline

Get the LineSearchType.

Returns
The current LineSearchType.

◆ LineSearch() [1/2]

vector3 LineSearch ( OBAtom atom,
vector3 direction 
)

Perform a linesearch starting at atom in direction direction.

Deprecated:
Current code should use LineSearch(double *, double*) instead.

◆ LineSearch() [2/2]

double LineSearch ( double *  currentCoords,
double *  direction 
)

Perform a linesearch for the entire molecule in direction direction. This function is called when using LineSearchType::Simple.

Parameters
currentCoordsStart coordinates.
directionThe search direction.
Returns
alpha, The scale of the step we moved along the direction vector.
Output to log:
OBFF_LOGLVL_NONE: none
OBFF_LOGLVL_LOW: none
OBFF_LOGLVL_MEDIUM: none
OBFF_LOGLVL_HIGH: none

◆ Newton2NumLineSearch()

double Newton2NumLineSearch ( double *  direction)

Perform a linesearch for the entire molecule. This function is called when using LineSearchType::Newton2Num.

Parameters
directionThe search direction.
Returns
alpha, The scale of the step we moved along the direction vector.
Output to log:
OBFF_LOGLVL_NONE: none
OBFF_LOGLVL_LOW: none
OBFF_LOGLVL_MEDIUM: none
OBFF_LOGLVL_HIGH: none

◆ LineSearchTakeStep()

void LineSearchTakeStep ( double *  origCoords,
double *  direction,
double  step 
)

Set the coordinates of the atoms to origCoord + step.

Parameters
origCoordsStart coordinates.
directionThe search direction.
stepThe step to take.

◆ SteepestDescent()

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
stepsThe number of steps.
econvEnergy convergence criteria. (defualt is 1e-6)
methodDeprecated. (see HasAnalyticalGradients())
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

◆ SteepestDescentInitialize()

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
stepsThe number of steps.
econvEnergy convergence criteria. (defualt is 1e-6)
methodDeprecated. (see HasAnalyticalGradients())
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

◆ SteepestDescentTakeNSteps()

bool SteepestDescentTakeNSteps ( int  n)

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

Parameters
nThe 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

◆ ConjugateGradients()

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
stepsThe number of steps.
econvEnergy convergence criteria. (defualt is 1e-6)
methodDeprecated. (see HasAnalyticalGradients())
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

Referenced by OBMinimizingEnergyConformerScore::Score(), and OBMinimizingRMSDConformerScore::Score().

◆ ConjugateGradientsInitialize()

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
stepsThe number of steps.
econvEnergy convergence criteria. (defualt is 1e-6)
methodDeprecated. (see HasAnalyticalGradients())
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

◆ ConjugateGradientsTakeNSteps()

bool ConjugateGradientsTakeNSteps ( int  n)

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

Parameters
nThe 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

◆ GenerateVelocities()

void GenerateVelocities ( )

Generate starting velocities with a Maxwellian distribution.

◆ CorrectVelocities()

void CorrectVelocities ( )

Correct the velocities so that the following is true:

3N
----
0.5 \ m_i * v_i^2 = 0.5 * Ndf * R * 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)
R : gas constant
T : temperature

◆ MolecularDynamicsTakeNSteps()

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, ...)
}
Parameters
nThe number of steps to take.
TAbsolute temperature in Kelvin.
timestepThe time step in picoseconds. (10e-12 s)
methodOBFF_ANALYTICAL_GRADIENTS (default) or OBFF_NUMERICAL_GRADIENTS

◆ GetConstraints()

OBFFConstraints & GetConstraints ( )

Get the current constraints.

Returns
The current constrains stored in the force field.

◆ SetConstraints()

void SetConstraints ( OBFFConstraints constraints)

Set the constraints.

Parameters
constraintsThe new constraints to be used.

◆ SetFixAtom()

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.

Parameters
indexThe index for the atom to fix.

◆ UnsetFixAtom()

void UnsetFixAtom ( )

Undo last SetFixAtom. This function will not remove the fix atom constraint for this atom if set by Setup() or SetConstraints().

◆ SetIgnoreAtom()

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.

Parameters
indexThe index for the atom to ignore.

◆ UnsetIgnoreAtom()

void UnsetIgnoreAtom ( )

Undo last SetIgnoreAtom. This function will not remove the ignore atom constraint for this atom if set by Setup() or SetConstraints().

◆ IgnoreCalculation() [1/3]

bool IgnoreCalculation ( int  a,
int  b 
)
static

internal function

◆ IgnoreCalculation() [2/3]

bool IgnoreCalculation ( int  a,
int  b,
int  c 
)
static

internal function

◆ IgnoreCalculation() [3/3]

bool IgnoreCalculation ( int  a,
int  b,
int  c,
int  d 
)
static

internal function

◆ DetectExplosion()

bool DetectExplosion ( )

(debugging)

◆ ValidateLineSearch()

vector3 ValidateLineSearch ( OBAtom atom,
vector3 direction 
)

(debugging)

◆ ValidateSteepestDescent()

void ValidateSteepestDescent ( int  steps)

(debugging)

◆ ValidateConjugateGradients()

void ValidateConjugateGradients ( int  steps)

(debugging)

◆ Validate()

virtual bool Validate ( )
inlinevirtual

Validate the force field implementation (debugging)

◆ ValidateGradients()

virtual bool ValidateGradients ( )
inlinevirtual

Validate the analytical gradients by comparing them to numerical ones. This function has to be implemented force field specific. (debugging)

◆ ValidateGradientError()

vector3 ValidateGradientError ( vector3 numgrad,
vector3 anagrad 
)

Calculate the error of the analytical gradient (debugging)

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

◆ VectorBondDerivative()

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).

Parameters
pos_aatom a (coordinates)
pos_batom b (coordinates)
force_a- return value for the force on atom a
force_b- return value for the force on atom b
Returns
The distance between a and b (bondlength for bond stretching, separation for vdw, electrostatic)

◆ VectorDistanceDerivative()

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.

◆ VectorLengthDerivative()

double VectorLengthDerivative ( vector3 a,
vector3 b 
)
static

◆ VectorAngleDerivative() [1/2]

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.

Parameters
pos_aatom a (coordinates)
pos_batom b (coordinates)
pos_catom 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
Returns
The angle between a-b-c

◆ VectorAngleDerivative() [2/2]

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

◆ VectorOOPDerivative() [1/2]

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).

Parameters
pos_aatom a (coordinates)
pos_batom b (coordinates)
pos_catom c (coordinates)
pos_datom 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
Returns
The OOP angle for a-b-c-d

◆ VectorOOPDerivative() [2/2]

double VectorOOPDerivative ( vector3 a,
vector3 b,
vector3 c,
vector3 d 
)
static

◆ VectorTorsionDerivative() [1/2]

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).

Parameters
pos_aatom a (coordinates)
pos_batom b (coordinates)
pos_catom c (coordinates)
pos_datom 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
Returns
The tosion angle for a-b-c-d

◆ VectorTorsionDerivative() [2/2]

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

◆ VectorSubtract() [1/2]

static void VectorSubtract ( double *  i,
double *  j,
double *  result 
)
inlinestatic

inline fuction to speed up minimization speed

Parameters
ipointer to i[3]
jpointer to j[3]
resultpointer to result[3], will be set to i - j

◆ VectorSubtract() [2/2]

static void VectorSubtract ( const double *const  i,
const double *const  j,
double *  result 
)
inlinestatic

◆ VectorAdd()

static void VectorAdd ( double *  i,
double *  j,
double *  result 
)
inlinestatic

inline fuction to speed up minimization speed

Parameters
ipointer to i[3]
jpointer to j[3]
resultpointer to result[3], will be set to i + j

◆ VectorDivide()

static void VectorDivide ( double *  i,
double  n,
double *  result 
)
inlinestatic

inline fuction to speed up minimization speed

Parameters
ipointer to i[3]
ndivide x,y,z with n
resultpointer to result[3]

◆ VectorMultiply() [1/2]

static void VectorMultiply ( double *  i,
double  n,
double *  result 
)
inlinestatic

inline fuction to speed up minimization speed

Parameters
ipointer to i[3]
nmultiply x,y,z with n
resultpointer to result[3]

◆ VectorMultiply() [2/2]

static void VectorMultiply ( const double *const  i,
const double  n,
double *  result 
)
inlinestatic

◆ VectorSelfMultiply()

static void VectorSelfMultiply ( double *  i,
double  n 
)
inlinestatic

inline fuction to speed up minimization speed

Parameters
ipointer to i[3], multiply this vector by n and set this vector to the result.
nthe scalar value to be multipled

◆ VectorNormalize()

static void VectorNormalize ( double *  i)
inlinestatic

inline fuction to speed up minimization speed

Parameters
ipointer to i[3] to be normalized

◆ VectorCopy()

static void VectorCopy ( double *  from,
double *  to 
)
inlinestatic

inline fuction to speed up minimization speed

Parameters
frompointer to i[3] to be copied from
topointer to j[3] to be copied to

◆ VectorLength()

static double VectorLength ( double *  i)
inlinestatic

inline fuction to speed up minimization speed

Parameters
ipointer to i[3]
Returns
the vector length

◆ VectorDistance()

static double VectorDistance ( double *  pos_i,
double *  pos_j 
)
inlinestatic

◆ VectorAngle()

double VectorAngle ( double *  i,
double *  j,
double *  k 
)
static

inline fuction to speed up minimization speed

Parameters
ipointer to i[3]
jpointer to j[3]
kpointer to k[3]
Returns
the vector angle ijk (deg)

◆ VectorTorsion()

double VectorTorsion ( double *  i,
double *  j,
double *  k,
double *  l 
)
static

inline fuction to speed up minimization speed

Parameters
ipointer to i[3]
jpointer to j[3]
kpointer to k[3]
lpointer to l[3]
Returns
the vector torson ijkl (deg)

◆ VectorOOP()

double VectorOOP ( double *  i,
double *  j,
double *  k,
double *  l 
)
static

inline fuction to speed up minimization speed

Parameters
ipointer to i[3]
jpointer to j[3]
kpointer to k[3]
lpointer to l[3]
Returns
the vector torson ijkl (deg)

◆ VectorClear()

static void VectorClear ( double *  i)
inlinestatic

inline fuction to speed up minimization speed

Parameters
ipointer to i[3], will set x,y,z to 0,0,0

◆ VectorDot()

static double VectorDot ( double *  i,
double *  j 
)
inlinestatic

inline fuction to speed up minimization speed

Parameters
ipointer to i[3]
jpointer to j[3]
Returns
the dot product

◆ VectorCross()

static void VectorCross ( double *  i,
double *  j,
double *  result 
)
inlinestatic

inline fuction to speed up minimization speed

Parameters
ipointer to i[3]
jpointer to j[3]
resultthe dot product (as a return value double[3])

◆ PrintVector()

static void PrintVector ( double *  i)
inlinestatic

Member Data Documentation

◆ _mol

OBMol _mol
protected

Molecule to be evaluated or minimized.

◆ _init

bool _init
protected

Used to make sure we only parse the parameter file once, when needed.

◆ _parFile

std::string _parFile
protected

◆ _validSetup

bool _validSetup
protected

< parameter file name

was the last call to Setup succesfull

◆ _gradientPtr

double* _gradientPtr
protected

pointer to the gradients (used by AddGradient(), minimization functions, ...)

◆ _logos

std::ostream* _logos
protected

Output for logfile.

◆ _logbuf

char _logbuf[BUFF_SIZE+1]
protected

Temporary buffer for logfile output.

◆ _loglvl

int _loglvl
protected

Log level for output.

◆ _origLogLevel

int _origLogLevel
protected

◆ _current_conformer

int _current_conformer
protected

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

◆ _energies

std::vector<double> _energies
protected

used to hold the energies for all conformers

◆ _econv

double _econv
protected

◆ _gconv

double _gconv
protected

◆ _e_n1

double _e_n1
protected

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

◆ _cstep

int _cstep
protected

◆ _nsteps

int _nsteps
protected

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

◆ _grad1

double* _grad1
protected

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

◆ _ncoords

unsigned int _ncoords
protected

Number of coordinates for conjugate gradients.

◆ _linesearch

int _linesearch
protected

LineSearch type.

◆ _timestep

double _timestep
protected

Molecular dynamics time step in picoseconds.

◆ _temp

double _temp
protected

Molecular dynamics temperature in Kelvin.

◆ _velocityPtr

double* _velocityPtr
protected

pointer to the velocities

◆ _constraints

OBFFConstraints _constraints = OBFFConstraints()
staticprotected

Constraints.

◆ _fixAtom

unsigned int _fixAtom = 0
staticprotected

SetFixAtom()/UnsetFixAtom()

◆ _ignoreAtom

unsigned int _ignoreAtom = 0
staticprotected

SetIgnoreAtom()/UnsetIgnoreAtom()

◆ _cutoff

bool _cutoff
protected

true = cut-off enabled

◆ _rvdw

double _rvdw
protected

VDW cut-off distance.

◆ _rele

double _rele
protected

Electrostatic cut-off distance.

◆ _epsilon

double _epsilon
protected

Dielectric constant for electrostatics.

◆ _vdwpairs

OBBitVec _vdwpairs
protected

VDW pairs that should be calculated.

◆ _elepairs

OBBitVec _elepairs
protected

Electrostatic pairs that should be calculated.

◆ _pairfreq

int _pairfreq
protected

The frequence to update non-bonded pairs.

◆ _intraGroup

std::vector<OBBitVec> _intraGroup
protected

groups for which intra-molecular interactions should be calculated

◆ _interGroup

std::vector<OBBitVec> _interGroup
protected

groups for which intra-molecular interactions should be calculated

◆ _interGroups

std::vector<std::pair<OBBitVec, OBBitVec> > _interGroups
protected

groups for which intra-molecular interactions should be calculated


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