Gromacs  2020.4
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Functions
#include "gromacs/math/vec.h"
#include "gromacs/topology/idef.h"
#include "gromacs/utility/real.h"
+ Include dependency graph for restcbt.h:
+ This graph shows which files directly or indirectly include this file:

Description

This file contains function declarations necessary for computations of forces due to restricted angle, restricted dihedral and combined bending-torsion potentials.

Author
Nicolae Goga

Functions

void compute_factors_restangles (int type, const t_iparams forceparams[], rvec delta_ante, rvec delta_post, double *prefactor, double *ratio_ante, double *ratio_post, real *v)
 This function computes factors needed for restricted angle potentials. More...
 
void compute_factors_restrdihs (int type, const t_iparams forceparams[], rvec delta_ante, rvec delta_crnt, rvec delta_post, real *factor_phi_ai_ante, real *factor_phi_ai_crnt, real *factor_phi_ai_post, real *factor_phi_aj_ante, real *factor_phi_aj_crnt, real *factor_phi_aj_post, real *factor_phi_ak_ante, real *factor_phi_ak_crnt, real *factor_phi_ak_post, real *factor_phi_al_ante, real *factor_phi_al_crnt, real *factor_phi_al_post, real *prefactor_phi, real *v)
 Compute factors for restricted dihedral potentials. More...
 
void compute_factors_cbtdihs (int type, const t_iparams forceparams[], rvec delta_ante, rvec delta_crnt, rvec delta_post, rvec f_phi_ai, rvec f_phi_aj, rvec f_phi_ak, rvec f_phi_al, rvec f_theta_ante_ai, rvec f_theta_ante_aj, rvec f_theta_ante_ak, rvec f_theta_post_aj, rvec f_theta_post_ak, rvec f_theta_post_al, real *v)
 Compute factors for combined bending-torsion (CBT) potentials. More...
 

Function Documentation

void compute_factors_cbtdihs ( int  type,
const t_iparams  forceparams[],
rvec  delta_ante,
rvec  delta_crnt,
rvec  delta_post,
rvec  f_phi_ai,
rvec  f_phi_aj,
rvec  f_phi_ak,
rvec  f_phi_al,
rvec  f_theta_ante_ai,
rvec  f_theta_ante_aj,
rvec  f_theta_ante_ak,
rvec  f_theta_post_aj,
rvec  f_theta_post_ak,
rvec  f_theta_post_al,
real v 
)

Compute factors for combined bending-torsion (CBT) potentials.

The combined bending-torsion potential goes to zero in a very smooth manner, eliminating the numerical instabilities, when three coarse-grained particles align and the dihedral angle and standard dihedral potentials cannot be calculated. The CBT potential is calculated using the formula:

\[V_{\rm CBT}(\theta_{i-1}, \theta_i, \phi_i) = k_{\phi} \sin^3\theta_{i-1} \sin^3\theta_{i} \sum_{n=0}^4 { a_n \cos^n\phi_i}\]

(see section "Proper dihedrals: Combined bending-torsion potential" from the manual). It contains in its expression not only the dihedral angle $\phi$ but also $\theta_{i-1}$ (denoted as theta_ante below) and $\theta_{i}$ (denoted as theta_post below) — the adjacent bending angles. The derivative of the CBT potential is calculated as:

\[\frac{\partial V_{\rm CBT}(\theta_{i-1},\theta_i,\phi_i)} {\partial \vec r_{l}} = \frac{\partial V_ {\rm CBT}}{\partial \theta_{i-1}} \frac{\partial \theta_{i-1}}{\partial \vec r_{l}} + \frac{\partial V_{\rm CBT}}{\partial \phi_{i }} \frac{\partial \phi_{i }}{\partial \vec r_{l}}\]

where all the derivatives of the angles with respect to Cartesian coordinates are calculated as in Allen & Tildesley (pp. 330-332). Factors f_phi_* come from the derivatives of the torsion angle with respect to the beads ai, aj, ak, al (four) coordinates; f_theta_ante_* come from the derivatives of the bending angle theta_ante (theta_{i-1} in formula above) with respect to the beads ai, aj, ak (three particles) coordinates and f_theta_post_* come from the derivatives of the bending angle theta_post (theta_{i} in formula above) with respect to the beads aj, ak, al (three particles) coordinates.

Parameters
[in]typetype of force parameters
[in]forceparamsarray of parameters for force computations
[in]delta_antedistance between the first two particles
[in]delta_crntdistance between the middle pair of particles
[in]delta_postdistance between the last two particles
[out]f_phi_aiforce for particle ai due to derivative in phi angle
[out]f_phi_ajforce for particle aj due to derivative in phi angle
[out]f_phi_akforce for particle ak due to derivative in phi angle
[out]f_phi_alforce for particle al due to derivative in phi angle
[out]f_theta_ante_aiforce for particle ai due to derivative in theta_ante angle
[out]f_theta_ante_ajforce for particle aj due to derivative in theta_ante angle
[out]f_theta_ante_akforce for particle ak due to derivative in theta_ante angle
[out]f_theta_post_ajforce for particle aj due to derivative in theta_post angle
[out]f_theta_post_akforce for particle ak due to derivative in theta_post angle
[out]f_theta_post_alforce for particle al due to derivative in theta_psot angle
[out]vcontribution to energy (see formula above)
void compute_factors_restangles ( int  type,
const t_iparams  forceparams[],
rvec  delta_ante,
rvec  delta_post,
double *  prefactor,
double *  ratio_ante,
double *  ratio_post,
real v 
)

This function computes factors needed for restricted angle potentials.

The restricted angle potential is used in coarse-grained simulations to avoid singularities when three particles align and the dihedral angle and dihedral potential cannot be calculated. This potential is calculated using the formula:

\[V_{\rm ReB}(\theta_i) = \frac{1}{2} k_{\theta} \frac{(\cos\theta_i - \cos\theta_0)^2}{\sin^2\theta_i}\]

(see section "Restricted Bending Potential" from the manual). The derivative of the restricted angle potential is calculated as:

\[\frac{\partial V_{\rm ReB}(\theta_i)} {\partial \vec r_{k}} = \frac{dV_{\rm ReB}(\theta_i)}{dcos\theta_i} \frac{\partial cos\theta_i}{\partial \vec r_{k}}\]

where all the derivatives of the bending angle with respect to Cartesian coordinates are calculated as in Allen & Tildesley (pp. 330-332)

Parameters
[in]typetype of force parameters
[in]forceparamsarray of parameters for force computations
[in]delta_antedistance between the first two particles
[in]delta_postdistance between the last two particles
[out]prefactorcommon term that comes in front of each force
[out]ratio_anteratio of scalar products of delta_ante with delta_post and delta_ante with delta_ante
[out]ratio_postratio of scalar products of delta_ante with delta_post and delta_post with delta_ante
[out]vcontribution to energy (see formula above)
void compute_factors_restrdihs ( int  type,
const t_iparams  forceparams[],
rvec  delta_ante,
rvec  delta_crnt,
rvec  delta_post,
real factor_phi_ai_ante,
real factor_phi_ai_crnt,
real factor_phi_ai_post,
real factor_phi_aj_ante,
real factor_phi_aj_crnt,
real factor_phi_aj_post,
real factor_phi_ak_ante,
real factor_phi_ak_crnt,
real factor_phi_ak_post,
real factor_phi_al_ante,
real factor_phi_al_crnt,
real factor_phi_al_post,
real prefactor_phi,
real v 
)

Compute factors for restricted dihedral potentials.

The restricted dihedral potential is the equivalent of the restricted bending potential for the dihedral angle. It imposes the dihedral angle to have only one equilibrium value. This potential is calculated using the formula:

\[V_{\rm ReT}(\phi_i) = \frac{1}{2} k_{\phi} \frac{(\cos\phi_i - \cos\phi_0)^2}{\sin^2\phi_i}\]

(see section "Proper dihedrals: Restricted torsion potential" from the manual). The derivative of the restricted dihedral potential is calculated as:

\[\frac{\partial V_{\rm ReT}(\phi_i)} {\partial \vec r_{k}} = \frac{dV_{\rm ReT}(\phi_i)}{dcos\phi_i} \frac{\partial cos\phi_i}{\partial \vec r_{k}}\]

where all the derivatives of the dihedral angle with respect to Cartesian coordinates are calculated as in Allen & Tildesley (pp. 330-332). Factors factor_phi_* are coming from the derivatives of the torsion angle (phi) with respect to the beads ai, aj, ak, al, (four) coordinates and they are multiplied in the force computations with the particle distance stored in parameters delta_ante, delta_crnt, delta_post.

Parameters
[in]typetype of force parameters
[in]forceparamsarray of parameters for force computations
[in]delta_antedistance between the first two particles
[in]delta_crntdistance between the middle pair of particles
[in]delta_postdistance between the last two particles
[out]factor_phi_ai_anteforce factor for particle ai multiplied with delta_ante
[out]factor_phi_ai_crntforce factor for particle ai multiplied with delta_crnt
[out]factor_phi_ai_postforce factor for particle ai multiplied with delta_post
[out]factor_phi_aj_anteforce factor for particle aj multiplied with delta_ante
[out]factor_phi_aj_crntforce factor for particle aj multiplied with delta_crnt
[out]factor_phi_aj_postforce factor for particle aj multiplied with delta_post
[out]factor_phi_ak_anteforce factor for particle ak multiplied with delta_ante
[out]factor_phi_ak_crntforce factor for particle ak multiplied with delta_crnt
[out]factor_phi_ak_postforce factor for particle ak multiplied with delta_post
[out]factor_phi_al_anteforce factor for particle al multiplied with delta_ante
[out]factor_phi_al_crntforce factor for particle al multiplied with delta_crnt
[out]factor_phi_al_postforce factor for particle al multiplied with delta_post
[out]prefactor_phimultiplication constant of the torsion force
[out]vcontribution to energy (see formula above)