Gromacs
2020.4
|
This file contains function declarations necessary for computations of forces due to restricted angle, restricted dihedral and combined bending-torsion potentials.
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... | |
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:
(see section "Proper dihedrals: Combined bending-torsion potential" from the manual). It contains in its expression not only the dihedral angle but also (denoted as theta_ante below) and (denoted as theta_post below) — the adjacent bending angles. The derivative of the CBT potential is calculated as:
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.
[in] | type | type of force parameters |
[in] | forceparams | array of parameters for force computations |
[in] | delta_ante | distance between the first two particles |
[in] | delta_crnt | distance between the middle pair of particles |
[in] | delta_post | distance between the last two particles |
[out] | f_phi_ai | force for particle ai due to derivative in phi angle |
[out] | f_phi_aj | force for particle aj due to derivative in phi angle |
[out] | f_phi_ak | force for particle ak due to derivative in phi angle |
[out] | f_phi_al | force for particle al due to derivative in phi angle |
[out] | f_theta_ante_ai | force for particle ai due to derivative in theta_ante angle |
[out] | f_theta_ante_aj | force for particle aj due to derivative in theta_ante angle |
[out] | f_theta_ante_ak | force for particle ak due to derivative in theta_ante angle |
[out] | f_theta_post_aj | force for particle aj due to derivative in theta_post angle |
[out] | f_theta_post_ak | force for particle ak due to derivative in theta_post angle |
[out] | f_theta_post_al | force for particle al due to derivative in theta_psot angle |
[out] | v | contribution 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:
(see section "Restricted Bending Potential" from the manual). The derivative of the restricted angle potential is calculated as:
where all the derivatives of the bending angle with respect to Cartesian coordinates are calculated as in Allen & Tildesley (pp. 330-332)
[in] | type | type of force parameters |
[in] | forceparams | array of parameters for force computations |
[in] | delta_ante | distance between the first two particles |
[in] | delta_post | distance between the last two particles |
[out] | prefactor | common term that comes in front of each force |
[out] | ratio_ante | ratio of scalar products of delta_ante with delta_post and delta_ante with delta_ante |
[out] | ratio_post | ratio of scalar products of delta_ante with delta_post and delta_post with delta_ante |
[out] | v | contribution 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:
(see section "Proper dihedrals: Restricted torsion potential" from the manual). The derivative of the restricted dihedral potential is calculated as:
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.
[in] | type | type of force parameters |
[in] | forceparams | array of parameters for force computations |
[in] | delta_ante | distance between the first two particles |
[in] | delta_crnt | distance between the middle pair of particles |
[in] | delta_post | distance between the last two particles |
[out] | factor_phi_ai_ante | force factor for particle ai multiplied with delta_ante |
[out] | factor_phi_ai_crnt | force factor for particle ai multiplied with delta_crnt |
[out] | factor_phi_ai_post | force factor for particle ai multiplied with delta_post |
[out] | factor_phi_aj_ante | force factor for particle aj multiplied with delta_ante |
[out] | factor_phi_aj_crnt | force factor for particle aj multiplied with delta_crnt |
[out] | factor_phi_aj_post | force factor for particle aj multiplied with delta_post |
[out] | factor_phi_ak_ante | force factor for particle ak multiplied with delta_ante |
[out] | factor_phi_ak_crnt | force factor for particle ak multiplied with delta_crnt |
[out] | factor_phi_ak_post | force factor for particle ak multiplied with delta_post |
[out] | factor_phi_al_ante | force factor for particle al multiplied with delta_ante |
[out] | factor_phi_al_crnt | force factor for particle al multiplied with delta_crnt |
[out] | factor_phi_al_post | force factor for particle al multiplied with delta_post |
[out] | prefactor_phi | multiplication constant of the torsion force |
[out] | v | contribution to energy (see formula above) |