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