2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
37 #include "gromacs/math/units.h"
38 #include "gromacs/math/utilities.h"
39 #include "gromacs/math/vec.h"
40 #include "gromacs/topology/idef.h"
42 /* This function computes factors needed for restricted angle potential.
43 * For explanations on formula used see file "restcbt.h" */
45 void compute_factors_restangles(int type, const t_iparams forceparams[],
46 rvec delta_ante, rvec delta_post,
47 real *prefactor, real *ratio_ante, real *ratio_post, real *v)
49 real theta_equil, k_bending;
50 real cosine_theta_equil;
51 real c_ante, c_cros, c_post;
53 real delta_cosine, cosine_theta;
55 real term_theta_theta_equil;
57 k_bending = forceparams[type].harmonic.krA;
58 theta_equil = forceparams[type].harmonic.rA*DEG2RAD;
59 theta_equil = M_PI - theta_equil;
60 cosine_theta_equil = cos(theta_equil);
62 c_ante = iprod(delta_ante, delta_ante);
63 c_cros = iprod(delta_ante, delta_post);
64 c_post = iprod(delta_post, delta_post);
66 norm = gmx_invsqrt(c_ante * c_post);
67 cosine_theta = c_cros * norm;
68 sine_theta_sq = 1 - cosine_theta * cosine_theta;
70 *ratio_ante = c_cros / c_ante;
71 *ratio_post = c_cros / c_post;
73 delta_cosine = cosine_theta - cosine_theta_equil;
74 term_theta_theta_equil = 1 - cosine_theta * cosine_theta_equil;
75 *prefactor = -(k_bending) * delta_cosine * norm * term_theta_theta_equil / (sine_theta_sq * sine_theta_sq);
77 *v = k_bending * 0.5 * delta_cosine * delta_cosine / sine_theta_sq;
82 /* Compute factors for restricted dihedral potential
83 * For explanations on formula used see file "restcbt.h" */
84 void compute_factors_restrdihs(int type, const t_iparams forceparams[],
85 rvec delta_ante, rvec delta_crnt, rvec delta_post,
86 real *factor_phi_ai_ante, real *factor_phi_ai_crnt, real *factor_phi_ai_post,
87 real *factor_phi_aj_ante, real *factor_phi_aj_crnt, real *factor_phi_aj_post,
88 real *factor_phi_ak_ante, real *factor_phi_ak_crnt, real *factor_phi_ak_post,
89 real *factor_phi_al_ante, real *factor_phi_al_crnt, real *factor_phi_al_post,
90 real *prefactor_phi, real *v)
93 real phi0, sine_phi0, cosine_phi0;
95 real c_self_ante, c_self_crnt, c_self_post;
96 real c_cros_ante, c_cros_acrs, c_cros_post;
97 real c_prod, d_post, d_ante;
98 real sine_phi_sq, cosine_phi;
99 real delta_cosine, term_phi_phi0;
100 real ratio_phi_ante, ratio_phi_post;
101 real cos_phi, norm_phi;
103 /* Read parameters phi0 and k_torsion */
104 phi0 = forceparams[type].pdihs.phiA * DEG2RAD;
105 cosine_phi0 = cos(phi0);
106 sine_phi0 = sin(phi0);
107 k_torsion = forceparams[type].pdihs.cpA;
109 /* Computation of the cosine of the dihedral angle. The scalar ("dot") product method
110 * is used. c_*_* cummulate the scalar products of the differences of particles
111 * positions while c_prod, d_ante and d_post are differences of products of scalar
112 * terms that are parts of the derivatives of forces */
113 c_self_ante = iprod(delta_ante, delta_ante);
114 c_self_crnt = iprod(delta_crnt, delta_crnt);
115 c_self_post = iprod(delta_post, delta_post);
116 c_cros_ante = iprod(delta_ante, delta_crnt);
117 c_cros_acrs = iprod(delta_ante, delta_post);
118 c_cros_post = iprod(delta_crnt, delta_post);
119 c_prod = c_cros_ante * c_cros_post - c_self_crnt * c_cros_acrs;
120 d_ante = c_self_ante * c_self_crnt - c_cros_ante * c_cros_ante;
121 d_post = c_self_post * c_self_crnt - c_cros_post * c_cros_post;
123 /* When three consecutive beads align, we obtain values close to zero.
124 * Here we avoid small values to prevent round-off errors. */
125 if (d_ante < GMX_REAL_EPS)
127 d_ante = GMX_REAL_EPS;
129 if (d_post < GMX_REAL_EPS)
131 d_post = GMX_REAL_EPS;
134 /* Computes the square of the sinus of phi in sine_phi_sq */
135 norm_phi = gmx_invsqrt(d_ante * d_post);
136 cosine_phi = c_prod * norm_phi;
137 sine_phi_sq = 1.0 - cosine_phi * cosine_phi;
139 /* It is possible that cosine_phi is slightly bigger than 1.0 due to round-off errors. */
140 if (sine_phi_sq < 0.0)
145 /* Computation of the differences of cosines (delta_cosine) and a term (term_phi_phi0)
146 * that is part of the common prefactor_phi */
148 delta_cosine = cosine_phi - cosine_phi0;
149 term_phi_phi0 = 1 - cosine_phi * cosine_phi0;
152 /* Computation of ratios */
153 ratio_phi_ante = c_prod / d_ante;
154 ratio_phi_post = c_prod / d_post;
156 /* Computation of the prefactor - common term for all forces */
157 *prefactor_phi = -(k_torsion) * delta_cosine * norm_phi * term_phi_phi0 / (sine_phi_sq * sine_phi_sq);
159 /* Computation of force factors. Factors factor_phi_* are coming from the
160 * derivatives of the torsion angle (phi) with respect to the beads ai, aj, al, ak,
161 * (four) coordinates and they are multiplied in the force computations with the
162 * differences of the particles positions stored in parameters delta_ante,
163 * delta_crnt, delta_post. For formulas see file "restcbt.h" */
165 *factor_phi_ai_ante = ratio_phi_ante * c_self_crnt;
166 *factor_phi_ai_crnt = -c_cros_post - ratio_phi_ante * c_cros_ante;
167 *factor_phi_ai_post = c_self_crnt;
168 *factor_phi_aj_ante = -c_cros_post - ratio_phi_ante * (c_self_crnt + c_cros_ante);
169 *factor_phi_aj_crnt = c_cros_post + c_cros_acrs * 2.0 + ratio_phi_ante * (c_self_ante + c_cros_ante) + ratio_phi_post * c_self_post;
170 *factor_phi_aj_post = -(c_cros_ante + c_self_crnt) - ratio_phi_post * c_cros_post;
171 *factor_phi_ak_ante = c_cros_post + c_self_crnt + ratio_phi_ante * c_cros_ante;
172 *factor_phi_ak_crnt = -(c_cros_ante + c_cros_acrs * 2.0)- ratio_phi_ante * c_self_ante - ratio_phi_post * (c_self_post + c_cros_post);
173 *factor_phi_ak_post = c_cros_ante + ratio_phi_post * (c_self_crnt + c_cros_post);
174 *factor_phi_al_ante = -c_self_crnt;
175 *factor_phi_al_crnt = c_cros_ante + ratio_phi_post * c_cros_post;
176 *factor_phi_al_post = -ratio_phi_post * c_self_crnt;
178 /* Contribution to energy - see formula in file "restcbt.h"*/
179 *v = k_torsion * 0.5 * delta_cosine * delta_cosine / sine_phi_sq;
184 /* Compute factors for CBT potential
185 * For explanations on formula used see file "restcbt.h" */
187 void compute_factors_cbtdihs(int type, const t_iparams forceparams[],
188 rvec delta_ante, rvec delta_crnt, rvec delta_post,
189 rvec f_phi_ai, rvec f_phi_aj, rvec f_phi_ak, rvec f_phi_al,
190 rvec f_theta_ante_ai, rvec f_theta_ante_aj, rvec f_theta_ante_ak,
191 rvec f_theta_post_aj, rvec f_theta_post_ak, rvec f_theta_post_al,
195 real torsion_coef[NR_CBTDIHS];
196 real c_self_ante, c_self_crnt, c_self_post;
197 real c_cros_ante, c_cros_acrs, c_cros_post;
198 real c_prod, d_ante, d_post;
199 real norm_phi, norm_theta_ante, norm_theta_post;
200 real cosine_phi, cosine_theta_ante, cosine_theta_post;
201 real sine_theta_ante_sq, sine_theta_post_sq;
202 real sine_theta_ante, sine_theta_post;
204 real ratio_phi_ante, ratio_phi_post;
206 real factor_phi_ai_ante, factor_phi_ai_crnt, factor_phi_ai_post;
207 real factor_phi_aj_ante, factor_phi_aj_crnt, factor_phi_aj_post;
208 real factor_phi_ak_ante, factor_phi_ak_crnt, factor_phi_ak_post;
209 real factor_phi_al_ante, factor_phi_al_crnt, factor_phi_al_post;
210 real prefactor_theta_ante, ratio_theta_ante_ante, ratio_theta_ante_crnt;
211 real prefactor_theta_post, ratio_theta_post_crnt, ratio_theta_post_post;
213 /* The formula for combined bending-torsion potential (see file "restcbt.h") contains
214 * in its expression not only the dihedral angle \f[\phi\f] but also \f[\theta_{i-1}\f]
215 * (theta_ante bellow) and \f[\theta_{i}\f] (theta_post bellow)--- the adjacent bending
216 * angles. The forces for the particles ai, aj, ak, al have components coming from the
217 * derivatives of the potential with respect to all three angles.
218 * This function is organised in 4 parts
219 * PART 1 - Computes force factors common to all the derivatives for the four particles
220 * PART 2 - Computes the force components due to the derivatives of dihedral angle Phi
221 * PART 3 - Computes the force components due to the derivatives of bending angle Theta_Ante
222 * PART 4 - Computes the force components due to the derivatives of bending angle Theta_Post
223 * Bellow we will respct thuis structure */
226 /* PART 1 - COMPUTES FORCE FACTORS COMMON TO ALL DERIVATIVES FOR THE FOUR PARTICLES */
229 for (j = 0; (j < NR_CBTDIHS); j++)
231 torsion_coef[j] = forceparams[type].cbtdihs.cbtcA[j];
234 /* Computation of the cosine of the dihedral angle. The scalar ("dot") product method
235 * is used. c_*_* cummulate the scalar products of the differences of particles
236 * positions while c_prod, d_ante and d_post are differences of products of scalar
237 * terms that are parts of the derivatives of forces */
239 c_self_ante = iprod(delta_ante, delta_ante);
240 c_self_crnt = iprod(delta_crnt, delta_crnt);
241 c_self_post = iprod(delta_post, delta_post);
242 c_cros_ante = iprod(delta_ante, delta_crnt);
243 c_cros_acrs = iprod(delta_ante, delta_post);
244 c_cros_post = iprod(delta_crnt, delta_post);
245 c_prod = c_cros_ante * c_cros_post - c_self_crnt * c_cros_acrs;
246 d_ante = c_self_ante * c_self_crnt - c_cros_ante * c_cros_ante;
247 d_post = c_self_post * c_self_crnt - c_cros_post * c_cros_post;
249 /* When three consecutive beads align, we obtain values close to zero.
250 Here we avoid small values to prevent round-off errors. */
251 if (d_ante < GMX_REAL_EPS)
253 d_ante = GMX_REAL_EPS;
255 if (d_post < GMX_REAL_EPS)
257 d_post = GMX_REAL_EPS;
260 /* Computations of cosines */
261 norm_phi = gmx_invsqrt(d_ante * d_post);
262 norm_theta_ante = gmx_invsqrt(c_self_ante * c_self_crnt);
263 norm_theta_post = gmx_invsqrt(c_self_crnt * c_self_post);
264 cosine_phi = c_prod * norm_phi;
265 cosine_theta_ante = c_cros_ante * norm_theta_ante;
266 cosine_theta_post = c_cros_post * norm_theta_post;
267 sine_theta_ante_sq = 1 - cosine_theta_ante * cosine_theta_ante;
268 sine_theta_post_sq = 1 - cosine_theta_post * cosine_theta_post;
270 /* It is possible that cosine_theta is slightly bigger than 1.0 due to round-off errors. */
271 if (sine_theta_ante_sq < 0.0)
273 sine_theta_ante_sq = 0.0;
275 if (sine_theta_post_sq < 0.0)
277 sine_theta_post_sq = 0.0;
280 sine_theta_ante = sqrt(sine_theta_ante_sq);
281 sine_theta_post = sqrt(sine_theta_post_sq);
283 /* PART 2 - COMPUTES FORCE COMPONENTS DUE TO DERIVATIVES TO DIHEDRAL ANGLE PHI */
285 /* Computation of ratios */
286 ratio_phi_ante = c_prod / d_ante;
287 ratio_phi_post = c_prod / d_post;
289 /* Computation of the prefactor */
290 /* Computing 2nd power */
293 prefactor_phi = -torsion_coef[0] * norm_phi * (torsion_coef[2] + torsion_coef[3] * 2.0 * cosine_phi + torsion_coef[4] * 3.0 * (r1 * r1) + 4*torsion_coef[5]*r1*r1*r1) *
294 sine_theta_ante_sq * sine_theta_ante * sine_theta_post_sq * sine_theta_post;
296 /* Computation of factors (important for gaining speed). Factors factor_phi_* are coming from the
297 * derivatives of the torsion angle (phi) with respect to the beads ai, aj, al, ak,
298 * (four) coordinates and they are multiplied in the force computations with the
299 * differences of the particles positions stored in parameters delta_ante,
300 * delta_crnt, delta_post. For formulas see file "restcbt.h" */
302 factor_phi_ai_ante = ratio_phi_ante * c_self_crnt;
303 factor_phi_ai_crnt = -c_cros_post - ratio_phi_ante * c_cros_ante;
304 factor_phi_ai_post = c_self_crnt;
305 factor_phi_aj_ante = -c_cros_post - ratio_phi_ante * (c_self_crnt + c_cros_ante);
306 factor_phi_aj_crnt = c_cros_post + c_cros_acrs * 2.0 + ratio_phi_ante * (c_self_ante + c_cros_ante) + ratio_phi_post * c_self_post;
307 factor_phi_aj_post = -(c_cros_ante + c_self_crnt) - ratio_phi_post * c_cros_post;
308 factor_phi_ak_ante = c_cros_post + c_self_crnt + ratio_phi_ante * c_cros_ante;
309 factor_phi_ak_crnt = -(c_cros_ante + c_cros_acrs * 2.0) - ratio_phi_ante * c_self_ante - ratio_phi_post * (c_self_post + c_cros_post);
310 factor_phi_ak_post = c_cros_ante + ratio_phi_post * (c_self_crnt + c_cros_post);
311 factor_phi_al_ante = -c_self_crnt;
312 factor_phi_al_crnt = c_cros_ante + ratio_phi_post * c_cros_post;
313 factor_phi_al_post = -ratio_phi_post * c_self_crnt;
315 /* Computation of forces due to the derivatives of dihedral angle phi*/
316 for (d = 0; d < DIM; d++)
318 f_phi_ai[d] = prefactor_phi * (factor_phi_ai_ante * delta_ante[d] + factor_phi_ai_crnt * delta_crnt[d] + factor_phi_ai_post * delta_post[d]);
319 f_phi_aj[d] = prefactor_phi * (factor_phi_aj_ante * delta_ante[d] + factor_phi_aj_crnt * delta_crnt[d] + factor_phi_aj_post * delta_post[d]);
320 f_phi_ak[d] = prefactor_phi * (factor_phi_ak_ante * delta_ante[d] + factor_phi_ak_crnt * delta_crnt[d] + factor_phi_ak_post * delta_post[d]);
321 f_phi_al[d] = prefactor_phi * (factor_phi_al_ante * delta_ante[d] + factor_phi_al_crnt * delta_crnt[d] + factor_phi_al_post * delta_post[d]);
324 /* PART 3 - COMPUTES THE FORCE COMPONENTS DUE TO THE DERIVATIVES OF BENDING ANGLE THETHA_ANTHE */
325 /* Computation of ratios */
326 ratio_theta_ante_ante = c_cros_ante / c_self_ante;
327 ratio_theta_ante_crnt = c_cros_ante / c_self_crnt;
329 /* Computation of the prefactor */
330 /* Computing 2nd power */
332 /* Computing 3rd power */
335 prefactor_theta_ante = -torsion_coef[0] * norm_theta_ante * ( torsion_coef[1] + torsion_coef[2] * cosine_phi + torsion_coef[3] * (r1 * r1) +
336 torsion_coef[4] * (r2 * (r2 * r2))+ torsion_coef[5] * (r2 * (r2 * (r2 * r2)))) * (-3.0) * cosine_theta_ante * sine_theta_ante * sine_theta_post_sq * sine_theta_post;
339 /* Computation of forces due to the derivatives of bending angle theta_ante */
340 for (d = 0; d < DIM; d++)
342 f_theta_ante_ai[d] = prefactor_theta_ante * (ratio_theta_ante_ante * delta_ante[d] - delta_crnt[d]);
343 f_theta_ante_aj[d] = prefactor_theta_ante * ((ratio_theta_ante_crnt + 1.0) * delta_crnt[d] - (ratio_theta_ante_ante + 1.0) * delta_ante[d]);
344 f_theta_ante_ak[d] = prefactor_theta_ante * (delta_ante[d] - ratio_theta_ante_crnt * delta_crnt[d]);
347 /* PART 4 - COMPUTES THE FORCE COMPONENTS DUE TO THE DERIVATIVES OF THE BENDING ANGLE THETA_POST */
349 /* Computation of ratios */
350 ratio_theta_post_crnt = c_cros_post / c_self_crnt;
351 ratio_theta_post_post = c_cros_post / c_self_post;
353 /* Computation of the prefactor */
354 /* Computing 2nd power */
356 /* Computing 3rd power */
359 prefactor_theta_post = -torsion_coef[0] * norm_theta_post * (torsion_coef[1] + torsion_coef[2] * cosine_phi + torsion_coef[3] * (r1 * r1) +
360 torsion_coef[4] * (r2 * (r2 * r2)) + torsion_coef[5] * (r2 * (r2 * (r2 * r2)))) * sine_theta_ante_sq * sine_theta_ante * (-3.0) * cosine_theta_post * sine_theta_post;
363 /* Computation of forces due to the derivatives of bending angle Theta_Post */
364 for (d = 0; d < DIM; d++)
366 f_theta_post_aj[d] = prefactor_theta_post * (ratio_theta_post_crnt * delta_crnt[d] - delta_post[d]);
367 f_theta_post_ak[d] = prefactor_theta_post * ((ratio_theta_post_post + 1.0) * delta_post[d] - (ratio_theta_post_crnt + 1.0) * delta_crnt[d]);
368 f_theta_post_al[d] = prefactor_theta_post * (delta_crnt[d] - ratio_theta_post_post * delta_post[d]);
373 /* Contribution to energy - for formula see file "restcbt.h" */
374 *v = torsion_coef[0] * (torsion_coef[1] + torsion_coef[2] * cosine_phi + torsion_coef[3] * (r1 * r1) +
375 torsion_coef[4] * (r2 * (r2 * r2)) + torsion_coef[5] * (r2 * (r2 * (r2 * r2)))) * sine_theta_ante_sq *
376 sine_theta_ante * sine_theta_post_sq * sine_theta_post;