Pseudo-improper-dihedral custom potential¶
-
namespace mdk
-
class LambdaPeak
- #include <PseudoImproperDihedral.hpp>
An object representing either of the three lambda functions. It is theoretically a Gaussian function, but we implement it here in terms of computationally faster functions (either cosine or an “algebraic version” thereof, with a cutoff at the far ends).
Public Functions
-
bool supp(double psi) const
Check whether a value $psi$ is in the support of the function.
- Parameters
psi – Value to test.
- Returns
true if $lambda(psi) != 0$, false otherwise.
-
void eval(double psi, double &L, double &dL_dpsi) const
Compute the value of the function.
- Parameters
psi – Argument, i.e. the improper dihedral angle
L – Value of the lambda function at
psi
,dL_dpsi – Derivative of L wrt
psi
Public Members
-
double alpha
The “width” of the lambda function.
-
double psi0
The “center” of the lambda function.
-
bool cosineVersion
Whether we should use the cosine version of the function.
-
bool supp(double psi) const
-
class PseudoImproperDihedral : public mdk::NonlocalForce
- #include <PseudoImproperDihedral.hpp>
The “pseudo-improper-dihedral” custom potential. A more computationally expensive but more accurate and Hamiltonian alternative to the quasi-adiabatic potential.
Public Functions
-
virtual void bind(Simulation &simulation) override
Bind the object to the simulation. Here we also fetch from it the types and chains.
- Parameters
simulation – Simulation to bind to.
-
virtual void asyncPart(Dynamics &dynamics) override
Asynchronous part of the force computation.
- Parameters
dynamics – Dynamics object to add potential energy and forces to.
-
virtual void vlUpdateHook() override
An action to be performed when the Verlet list is updated. We essentially copy the residues but exclude the pairs within bond distance of less than 4, or at the ends of their chains.
Public Members
-
LambdaPeak bb_pos
Lambda function corresponding to the bb+ part of the potential.
-
LennardJones bb_pos_lj
Lennard-Jones potential corresponding to the bb+ part of the potential.
-
LambdaPeak bb_neg
Lambda function corresponding to the bb- part of the potential.
-
LennardJones bb_neg_lj
Lennard-Jones potential corresponding to the bb- part of the potential.
-
LambdaPeak ss
Lambda function corresponding to the ss part of the potential.
-
SidechainLJ ss_ljs[AminoAcid::N][AminoAcid::N]
An array of sidechain L-J potentials for each pair of types of amino acids interacting.
Protected Functions
Private Functions
-
void deriveAngles(vl::PairInfo const &pair, double psi[2], Vector dpsi_dr[2][6]) const
Derive the angle data (specifically the angles and derivatives) pertaining to the pair.
- Parameters
pair – Pair of residues.
psi – An array, where psi[0] is the improper dihedral angle between i_1-1, i_1, i_1+1 and i_2, and psi[1] is the improper dihedral angle between i_2-1, i_2, i_2+1 and i_1.
dpsi_dr – An array, where dpsi_dr[0][0..2] are the derivatives d psi_1/d r_{i_1-1 .. i_1+1} and dspi_dr[0][3..6] are the derivatives d psi_1/d r_{i_2-1 .. i_2+1}; dpsi_dr[1][0..3] is defined as for 0, but with psi_2 instead of psi_1 (order of variables is not changed).
-
virtual void bind(Simulation &simulation) override
-
class LambdaPeak