RDKit
Open-source cheminformatics and machine learning.
Loading...
Searching...
No Matches
UFF/AngleBend.h
Go to the documentation of this file.
1//
2// Copyright (C) 2004-2006 Rational Discovery LLC
3//
4// @@ All Rights Reserved @@
5// This file is part of the RDKit.
6// The contents are covered by the terms of the BSD license
7// which is included in the file license.txt, found at the root
8// of the RDKit source tree.
9//
10#include <RDGeneral/export.h>
11#ifndef __RD_ANGLEBEND_H__
12#define __RD_ANGLEBEND_H__
13
14#include <ForceField/Contrib.h>
15#include <Geometry/point.h>
16
17namespace ForceFields {
18namespace UFF {
19class AtomicParams;
20
21//! The angle-bend term for the Universal Force Field
23 public:
25 //! Constructor
26 /*!
27 The angle is between atom1 - atom2 - atom3
28
29 \param owner pointer to the owning ForceField
30 \param idx1 index of atom1 in the ForceField's positions
31 \param idx2 index of atom2 in the ForceField's positions
32 \param idx3 index of atom3 in the ForceField's positions
33 \param bondOrder12 order of the bond between atoms 1 and 2 (as a double)
34 \param bondOrder23 order of the bond between atoms 2 and 3 (as a double)
35 \param at1Params pointer to the parameters for atom 1
36 \param at2Params pointer to the parameters for atom 2
37 \param at3Params pointer to the parameters for atom 3
38 \param order (optional) the order of the angle term, this is for
39 special cases and should adopt values:
40 - 0: not a special case, use the \c theta0 value from \c at2Params
41 - 2: linear coordination
42 - 3: trigonal planar coordination
43 - 4: square planar or tetrahedral coordination
44
45 */
46 AngleBendContrib(ForceField *owner, unsigned int idx1, unsigned int idx2,
47 unsigned int idx3, double bondOrder12, double bondOrder23,
48 const AtomicParams *at1Params, const AtomicParams *at2Params,
49 const AtomicParams *at3Params, unsigned int order = 0);
50 double getEnergy(double *pos) const override;
51 void getGrad(double *pos, double *grad) const override;
52
53 AngleBendContrib *copy() const override {
54 return new AngleBendContrib(*this);
55 }
56
57 private:
58 int d_at1Idx{-1};
59 int d_at2Idx{-1};
60 int d_at3Idx{-1};
61 unsigned int d_order{0};
62 double d_forceConstant, d_C0, d_C1, d_C2;
63
64 double getEnergyTerm(double cosTheta, double sinThetaSq) const;
65 double getThetaDeriv(double cosTheta, double sinTheta) const;
66};
67
68namespace Utils {
69//! Calculate the force constant for an angle bend
70/*!
71 The angle is between atom1 - atom2 - atom3
72
73 \param theta0 the preferred value of the angle (in radians)
74 \param bondOrder12 order of the bond between atoms 1 and 2 (as a double)
75 \param bondOrder23 order of the bond between atoms 2 and 3 (as a double)
76 \param at1Params pointer to the parameters for atom 1
77 \param at2Params pointer to the parameters for atom 2
78 \param at3Params pointer to the parameters for atom 3
79
80*/
82 double theta0, double bondOrder12, double bondOrder23,
83 const AtomicParams *at1Params, const AtomicParams *at2Params,
84 const AtomicParams *at3Params);
86 double **g, double &dE_dTheta,
87 double &cosTheta,
88 double &sinTheta);
89} // namespace Utils
90} // namespace UFF
91} // namespace ForceFields
92#endif
abstract base class for contributions to ForceFields
Definition Contrib.h:18
A class to store forcefields and handle minimization.
Definition ForceField.h:79
The angle-bend term for the Universal Force Field.
AngleBendContrib(ForceField *owner, unsigned int idx1, unsigned int idx2, unsigned int idx3, double bondOrder12, double bondOrder23, const AtomicParams *at1Params, const AtomicParams *at2Params, const AtomicParams *at3Params, unsigned int order=0)
Constructor.
void getGrad(double *pos, double *grad) const override
calculates our contribution to the gradients of a position
AngleBendContrib * copy() const override
return a copy
double getEnergy(double *pos) const override
returns our contribution to the energy of a position
class to store atomic parameters for the Universal Force Field
Definition UFF/Params.h:73
#define RDKIT_FORCEFIELD_EXPORT
Definition export.h:185
RDKIT_FORCEFIELD_EXPORT void calcAngleBendGrad(RDGeom::Point3D *r, double *dist, double **g, double &dE_dTheta, double &cosTheta, double &sinTheta)
RDKIT_FORCEFIELD_EXPORT double calcAngleForceConstant(double theta0, double bondOrder12, double bondOrder23, const AtomicParams *at1Params, const AtomicParams *at2Params, const AtomicParams *at3Params)
Calculate the force constant for an angle bend.