f264029e43622540483b83781b1a89dbcc529eb2
[alexxy/gromacs.git] / src / gromacs / listed-forces / restcbt.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2014,2015,2017, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \internal \file
36  *
37  * \brief
38  * This file contains function definitions necessary
39  * for computations of forces due to restricted angle, restricted dihedral and
40  * combined bending-torsion potentials.
41  *
42  * \author Nicolae Goga
43  *
44  * \ingroup module_listed-forces
45  */
46 #include "gmxpre.h"
47
48 #include "restcbt.h"
49
50 #include <cmath>
51
52 #include "gromacs/math/functions.h"
53 #include "gromacs/math/units.h"
54 #include "gromacs/math/utilities.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/topology/idef.h"
57
58 /* This function computes factors needed for restricted angle potential.
59  * For explanations on formula used see file "restcbt.h" */
60
61 void compute_factors_restangles(int type, const t_iparams forceparams[],
62                                 rvec delta_ante,  rvec delta_post,
63                                 real *prefactor, real *ratio_ante, real *ratio_post, real *v)
64 {
65     real theta_equil, k_bending;
66     real cosine_theta_equil;
67     real c_ante, c_cros, c_post;
68     real norm;
69     real delta_cosine, cosine_theta;
70     real sine_theta_sq;
71     real term_theta_theta_equil;
72
73     k_bending          = forceparams[type].harmonic.krA;
74     theta_equil        =  forceparams[type].harmonic.rA*DEG2RAD;
75     theta_equil        = M_PI - theta_equil;
76     cosine_theta_equil = cos(theta_equil);
77
78     c_ante = iprod(delta_ante, delta_ante);
79     c_cros = iprod(delta_ante, delta_post);
80     c_post = iprod(delta_post, delta_post);
81
82     norm          = gmx::invsqrt(c_ante * c_post);
83     cosine_theta  = c_cros * norm;
84     sine_theta_sq = 1 - cosine_theta * cosine_theta;
85
86     *ratio_ante = c_cros / c_ante;
87     *ratio_post = c_cros / c_post;
88
89     delta_cosine           = cosine_theta - cosine_theta_equil;
90     term_theta_theta_equil = 1 - cosine_theta * cosine_theta_equil;
91     *prefactor             = -(k_bending) * delta_cosine * norm * term_theta_theta_equil / (sine_theta_sq * sine_theta_sq);
92
93     *v = k_bending * 0.5 * delta_cosine * delta_cosine / sine_theta_sq;
94
95 }
96
97
98 /* Compute factors for restricted dihedral potential
99  * For explanations on formula used see file "restcbt.h" */
100 void compute_factors_restrdihs(int type,  const t_iparams forceparams[],
101                                rvec delta_ante, rvec delta_crnt, rvec delta_post,
102                                real *factor_phi_ai_ante, real *factor_phi_ai_crnt, real *factor_phi_ai_post,
103                                real *factor_phi_aj_ante, real *factor_phi_aj_crnt, real *factor_phi_aj_post,
104                                real *factor_phi_ak_ante, real *factor_phi_ak_crnt, real *factor_phi_ak_post,
105                                real *factor_phi_al_ante, real *factor_phi_al_crnt, real *factor_phi_al_post,
106                                real *prefactor_phi, real *v)
107 {
108
109     real phi0, cosine_phi0;
110     real k_torsion;
111     real c_self_ante, c_self_crnt, c_self_post;
112     real c_cros_ante, c_cros_acrs, c_cros_post;
113     real c_prod, d_post, d_ante;
114     real sine_phi_sq, cosine_phi;
115     real delta_cosine, term_phi_phi0;
116     real ratio_phi_ante, ratio_phi_post;
117     real norm_phi;
118
119     /* Read parameters phi0 and k_torsion */
120     phi0        = forceparams[type].pdihs.phiA * DEG2RAD;
121     cosine_phi0 = cos(phi0);
122     k_torsion   = forceparams[type].pdihs.cpA;
123
124     /* Computation of the cosine of the dihedral angle. The scalar ("dot") product  method
125      * is used. c_*_* cummulate the scalar products of the differences of particles
126      * positions while c_prod, d_ante and d_post are differences of products of scalar
127      * terms that are parts of the derivatives of forces */
128     c_self_ante = iprod(delta_ante, delta_ante);
129     c_self_crnt = iprod(delta_crnt, delta_crnt);
130     c_self_post = iprod(delta_post, delta_post);
131     c_cros_ante = iprod(delta_ante, delta_crnt);
132     c_cros_acrs = iprod(delta_ante, delta_post);
133     c_cros_post = iprod(delta_crnt, delta_post);
134     c_prod      = c_cros_ante * c_cros_post - c_self_crnt * c_cros_acrs;
135     d_ante      = c_self_ante * c_self_crnt - c_cros_ante * c_cros_ante;
136     d_post      = c_self_post * c_self_crnt - c_cros_post * c_cros_post;
137
138     /*  When three consecutive beads align, we obtain values close to zero.
139      *  Here we avoid small values to prevent round-off errors. */
140     if (d_ante < GMX_REAL_EPS)
141     {
142         d_ante = GMX_REAL_EPS;
143     }
144     if (d_post < GMX_REAL_EPS)
145     {
146         d_post = GMX_REAL_EPS;
147     }
148
149     /* Computes the square of the sinus of phi  in sine_phi_sq */
150     norm_phi    = gmx::invsqrt(d_ante * d_post);
151     cosine_phi  = c_prod * norm_phi;
152     sine_phi_sq = 1.0 - cosine_phi * cosine_phi;
153
154     /*  It is possible that cosine_phi is slightly bigger than 1.0 due to round-off errors. */
155     if (sine_phi_sq < 0.0)
156     {
157         sine_phi_sq = 0.0;
158     }
159
160     /* Computation of the differences of cosines (delta_cosine) and a term (term_phi_phi0)
161      * that is part of the common prefactor_phi */
162
163     delta_cosine  = cosine_phi - cosine_phi0;
164     term_phi_phi0 = 1 - cosine_phi * cosine_phi0;
165
166
167     /*      Computation of ratios */
168     ratio_phi_ante = c_prod / d_ante;
169     ratio_phi_post = c_prod / d_post;
170
171     /*      Computation of the prefactor - common term for all forces */
172     *prefactor_phi = -(k_torsion) * delta_cosine * norm_phi * term_phi_phi0 / (sine_phi_sq * sine_phi_sq);
173
174     /* Computation of force factors.  Factors factor_phi_*  are coming from the
175      * derivatives of the torsion angle (phi) with respect to the beads ai, aj, al, ak,
176      * (four) coordinates and they are multiplied in the force computations with the
177      * differences of the particles positions stored in parameters delta_ante,
178      * delta_crnt, delta_post. For formulas see file "restcbt.h" */
179
180     *factor_phi_ai_ante = ratio_phi_ante * c_self_crnt;
181     *factor_phi_ai_crnt = -c_cros_post - ratio_phi_ante * c_cros_ante;
182     *factor_phi_ai_post = c_self_crnt;
183     *factor_phi_aj_ante = -c_cros_post - ratio_phi_ante * (c_self_crnt + c_cros_ante);
184     *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;
185     *factor_phi_aj_post = -(c_cros_ante + c_self_crnt) - ratio_phi_post * c_cros_post;
186     *factor_phi_ak_ante = c_cros_post + c_self_crnt + ratio_phi_ante * c_cros_ante;
187     *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);
188     *factor_phi_ak_post = c_cros_ante + ratio_phi_post * (c_self_crnt + c_cros_post);
189     *factor_phi_al_ante = -c_self_crnt;
190     *factor_phi_al_crnt = c_cros_ante + ratio_phi_post * c_cros_post;
191     *factor_phi_al_post = -ratio_phi_post * c_self_crnt;
192
193     /* Contribution to energy  - see formula in file "restcbt.h"*/
194     *v = k_torsion * 0.5 * delta_cosine * delta_cosine / sine_phi_sq;
195 }
196
197
198
199 /* Compute factors for CBT potential
200  * For explanations on formula used see file "restcbt.h" */
201
202 void compute_factors_cbtdihs(int type,  const t_iparams forceparams[],
203                              rvec delta_ante, rvec delta_crnt, rvec delta_post,
204                              rvec f_phi_ai, rvec f_phi_aj, rvec f_phi_ak, rvec f_phi_al,
205                              rvec f_theta_ante_ai, rvec f_theta_ante_aj, rvec f_theta_ante_ak,
206                              rvec f_theta_post_aj, rvec f_theta_post_ak, rvec f_theta_post_al,
207                              real * v)
208 {
209     int  j, d;
210     real torsion_coef[NR_CBTDIHS];
211     real c_self_ante, c_self_crnt, c_self_post;
212     real c_cros_ante, c_cros_acrs, c_cros_post;
213     real c_prod, d_ante,  d_post;
214     real norm_phi, norm_theta_ante, norm_theta_post;
215     real cosine_phi, cosine_theta_ante, cosine_theta_post;
216     real sine_theta_ante_sq, sine_theta_post_sq;
217     real sine_theta_ante, sine_theta_post;
218     real prefactor_phi;
219     real ratio_phi_ante, ratio_phi_post;
220     real r1, r2;
221     real factor_phi_ai_ante, factor_phi_ai_crnt, factor_phi_ai_post;
222     real factor_phi_aj_ante, factor_phi_aj_crnt, factor_phi_aj_post;
223     real factor_phi_ak_ante, factor_phi_ak_crnt, factor_phi_ak_post;
224     real factor_phi_al_ante, factor_phi_al_crnt, factor_phi_al_post;
225     real prefactor_theta_ante, ratio_theta_ante_ante, ratio_theta_ante_crnt;
226     real prefactor_theta_post, ratio_theta_post_crnt, ratio_theta_post_post;
227
228     /* The formula for combined bending-torsion potential (see file "restcbt.h") contains
229      * in its expression not only the dihedral angle \f[\phi\f] but also \f[\theta_{i-1}\f]
230      * (theta_ante bellow) and \f[\theta_{i}\f] (theta_post bellow)--- the adjacent bending
231      * angles. The forces for the particles ai, aj, ak, al have components coming from the
232      * derivatives of the potential with respect to all three angles.
233      * This function is organised in 4 parts
234      * PART 1 - Computes force factors common to all the derivatives for the four particles
235      * PART 2 - Computes the force components due to the derivatives of dihedral angle Phi
236      * PART 3 - Computes the force components due to the derivatives of bending angle Theta_Ante
237      * PART 4 - Computes the force components due to the derivatives of bending angle Theta_Post
238      * Bellow we will respct thuis structure */
239
240
241     /* PART 1 - COMPUTES FORCE FACTORS COMMON TO ALL DERIVATIVES FOR THE FOUR PARTICLES */
242
243
244     for (j = 0; (j < NR_CBTDIHS); j++)
245     {
246         torsion_coef[j]  = forceparams[type].cbtdihs.cbtcA[j];
247     }
248
249     /* Computation of the cosine of the dihedral angle. The scalar ("dot") product  method
250      * is used. c_*_* cummulate the scalar products of the differences of particles
251      * positions while c_prod, d_ante and d_post are differences of products of scalar
252      * terms that are parts of the derivatives of forces */
253
254     c_self_ante = iprod(delta_ante, delta_ante);
255     c_self_crnt = iprod(delta_crnt, delta_crnt);
256     c_self_post = iprod(delta_post, delta_post);
257     c_cros_ante = iprod(delta_ante, delta_crnt);
258     c_cros_acrs = iprod(delta_ante, delta_post);
259     c_cros_post = iprod(delta_crnt, delta_post);
260     c_prod      = c_cros_ante * c_cros_post - c_self_crnt * c_cros_acrs;
261     d_ante      = c_self_ante * c_self_crnt - c_cros_ante * c_cros_ante;
262     d_post      = c_self_post * c_self_crnt - c_cros_post * c_cros_post;
263
264     /*  When three consecutive beads align, we obtain values close to zero.
265        Here we avoid small values to prevent round-off errors. */
266     if (d_ante < GMX_REAL_EPS)
267     {
268         d_ante = GMX_REAL_EPS;
269     }
270     if (d_post < GMX_REAL_EPS)
271     {
272         d_post = GMX_REAL_EPS;
273     }
274
275     /* Computations of cosines */
276     norm_phi           = gmx::invsqrt(d_ante * d_post);
277     norm_theta_ante    = gmx::invsqrt(c_self_ante * c_self_crnt);
278     norm_theta_post    = gmx::invsqrt(c_self_crnt * c_self_post);
279     cosine_phi         = c_prod * norm_phi;
280     cosine_theta_ante  = c_cros_ante * norm_theta_ante;
281     cosine_theta_post  = c_cros_post * norm_theta_post;
282     sine_theta_ante_sq = 1 - cosine_theta_ante * cosine_theta_ante;
283     sine_theta_post_sq = 1 - cosine_theta_post * cosine_theta_post;
284
285     /*  It is possible that cosine_theta is slightly bigger than 1.0 due to round-off errors. */
286     if (sine_theta_ante_sq < 0.0)
287     {
288         sine_theta_ante_sq = 0.0;
289     }
290     if (sine_theta_post_sq < 0.0)
291     {
292         sine_theta_post_sq = 0.0;
293     }
294
295     sine_theta_ante = std::sqrt(sine_theta_ante_sq);
296     sine_theta_post = std::sqrt(sine_theta_post_sq);
297
298     /* PART 2 - COMPUTES FORCE COMPONENTS DUE TO DERIVATIVES TO DIHEDRAL ANGLE PHI */
299
300     /*      Computation of ratios */
301     ratio_phi_ante = c_prod / d_ante;
302     ratio_phi_post = c_prod / d_post;
303
304     /*       Computation of the prefactor */
305     /*      Computing 2nd power */
306     r1 = cosine_phi;
307
308     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) *
309         sine_theta_ante_sq * sine_theta_ante * sine_theta_post_sq * sine_theta_post;
310
311     /* Computation of factors (important for gaining speed). Factors factor_phi_*  are coming from the
312      * derivatives of the torsion angle (phi) with respect to the beads ai, aj, al, ak,
313      * (four) coordinates and they are multiplied in the force computations with the
314      * differences of the particles positions stored in parameters delta_ante,
315      * delta_crnt, delta_post. For formulas see file "restcbt.h" */
316
317     factor_phi_ai_ante = ratio_phi_ante * c_self_crnt;
318     factor_phi_ai_crnt = -c_cros_post - ratio_phi_ante * c_cros_ante;
319     factor_phi_ai_post = c_self_crnt;
320     factor_phi_aj_ante = -c_cros_post - ratio_phi_ante * (c_self_crnt + c_cros_ante);
321     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;
322     factor_phi_aj_post = -(c_cros_ante + c_self_crnt) - ratio_phi_post * c_cros_post;
323     factor_phi_ak_ante = c_cros_post + c_self_crnt + ratio_phi_ante * c_cros_ante;
324     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);
325     factor_phi_ak_post = c_cros_ante + ratio_phi_post * (c_self_crnt + c_cros_post);
326     factor_phi_al_ante = -c_self_crnt;
327     factor_phi_al_crnt = c_cros_ante + ratio_phi_post * c_cros_post;
328     factor_phi_al_post = -ratio_phi_post * c_self_crnt;
329
330     /* Computation of forces due to the derivatives of dihedral angle phi*/
331     for (d = 0; d < DIM; d++)
332     {
333         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]);
334         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]);
335         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]);
336         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]);
337     }
338
339     /* PART 3 - COMPUTES THE FORCE COMPONENTS DUE TO THE DERIVATIVES OF BENDING ANGLE THETHA_ANTHE */
340     /*      Computation of ratios */
341     ratio_theta_ante_ante = c_cros_ante / c_self_ante;
342     ratio_theta_ante_crnt = c_cros_ante / c_self_crnt;
343
344     /*      Computation of the prefactor */
345     /*      Computing 2nd power */
346     r1 = cosine_phi;
347     /*      Computing 3rd power */
348     r2 = cosine_phi;
349
350     prefactor_theta_ante = -torsion_coef[0] * norm_theta_ante * ( torsion_coef[1] + torsion_coef[2] * cosine_phi + torsion_coef[3] * (r1 * r1) +
351                                                                   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;
352
353
354     /*      Computation of forces due to the derivatives of bending angle theta_ante */
355     for (d = 0; d < DIM; d++)
356     {
357         f_theta_ante_ai[d] = prefactor_theta_ante * (ratio_theta_ante_ante * delta_ante[d] - delta_crnt[d]);
358         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]);
359         f_theta_ante_ak[d] = prefactor_theta_ante * (delta_ante[d] - ratio_theta_ante_crnt * delta_crnt[d]);
360     }
361
362     /* PART 4 - COMPUTES THE FORCE COMPONENTS DUE TO THE DERIVATIVES OF THE BENDING ANGLE THETA_POST */
363
364     /*      Computation of ratios */
365     ratio_theta_post_crnt = c_cros_post / c_self_crnt;
366     ratio_theta_post_post = c_cros_post / c_self_post;
367
368     /*     Computation of the prefactor */
369     /*      Computing 2nd power */
370     r1 = cosine_phi;
371     /*      Computing 3rd power */
372     r2 = cosine_phi;
373
374     prefactor_theta_post = -torsion_coef[0] * norm_theta_post * (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 * sine_theta_ante * (-3.0) * cosine_theta_post * sine_theta_post;
376
377
378     /*      Computation of forces due to the derivatives of bending angle Theta_Post */
379     for (d = 0; d < DIM; d++)
380     {
381         f_theta_post_aj[d] = prefactor_theta_post * (ratio_theta_post_crnt * delta_crnt[d] - delta_post[d]);
382         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]);
383         f_theta_post_al[d] = prefactor_theta_post * (delta_crnt[d] - ratio_theta_post_post * delta_post[d]);
384     }
385     r1 = cosine_phi;
386     r2 = cosine_phi;
387
388     /* Contribution to energy - for formula see file "restcbt.h" */
389     *v = torsion_coef[0] * (torsion_coef[1] + torsion_coef[2] * cosine_phi + torsion_coef[3] * (r1 * r1) +
390                             torsion_coef[4] * (r2 * (r2 * r2)) + torsion_coef[5] * (r2 * (r2 * (r2 * r2)))) * sine_theta_ante_sq *
391         sine_theta_ante * sine_theta_post_sq * sine_theta_post;
392
393
394 }