2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
39 * \brief This file defines functions necessary for mdrun and tools to
40 * compute energies and forces for bonded interactions.
42 * \author Mark Abraham <mark.j.abraham@gmail.com>
44 * \ingroup module_listed-forces
58 #include "gromacs/legacyheaders/disre.h"
59 #include "gromacs/legacyheaders/force.h"
60 #include "gromacs/legacyheaders/macros.h"
61 #include "gromacs/legacyheaders/names.h"
62 #include "gromacs/legacyheaders/ns.h"
63 #include "gromacs/legacyheaders/orires.h"
64 #include "gromacs/legacyheaders/txtdump.h"
65 #include "gromacs/math/units.h"
66 #include "gromacs/math/utilities.h"
67 #include "gromacs/math/vec.h"
68 #include "gromacs/pbcutil/ishift.h"
69 #include "gromacs/pbcutil/mshift.h"
70 #include "gromacs/pbcutil/pbc.h"
71 #include "gromacs/simd/simd.h"
72 #include "gromacs/simd/simd_math.h"
73 #include "gromacs/simd/vector_operations.h"
74 #include "gromacs/utility/fatalerror.h"
75 #include "gromacs/utility/smalloc.h"
77 #include "listed-internal.h"
81 /*! \brief Mysterious CMAP coefficient matrix */
82 const int cmap_coeff_matrix[] = {
83 1, 0, -3, 2, 0, 0, 0, 0, -3, 0, 9, -6, 2, 0, -6, 4,
84 0, 0, 0, 0, 0, 0, 0, 0, 3, 0, -9, 6, -2, 0, 6, -4,
85 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 9, -6, 0, 0, -6, 4,
86 0, 0, 3, -2, 0, 0, 0, 0, 0, 0, -9, 6, 0, 0, 6, -4,
87 0, 0, 0, 0, 1, 0, -3, 2, -2, 0, 6, -4, 1, 0, -3, 2,
88 0, 0, 0, 0, 0, 0, 0, 0, -1, 0, 3, -2, 1, 0, -3, 2,
89 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -3, 2, 0, 0, 3, -2,
90 0, 0, 0, 0, 0, 0, 3, -2, 0, 0, -6, 4, 0, 0, 3, -2,
91 0, 1, -2, 1, 0, 0, 0, 0, 0, -3, 6, -3, 0, 2, -4, 2,
92 0, 0, 0, 0, 0, 0, 0, 0, 0, 3, -6, 3, 0, -2, 4, -2,
93 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -3, 3, 0, 0, 2, -2,
94 0, 0, -1, 1, 0, 0, 0, 0, 0, 0, 3, -3, 0, 0, -2, 2,
95 0, 0, 0, 0, 0, 1, -2, 1, 0, -2, 4, -2, 0, 1, -2, 1,
96 0, 0, 0, 0, 0, 0, 0, 0, 0, -1, 2, -1, 0, 1, -2, 1,
97 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, -1, 0, 0, -1, 1,
98 0, 0, 0, 0, 0, 0, -1, 1, 0, 0, 2, -2, 0, 0, -1, 1
102 /*! \brief Compute dx = xi - xj, modulo PBC if non-NULL
104 * \todo This kind of code appears in many places. Consolidate it */
105 static int pbc_rvec_sub(const t_pbc *pbc, const rvec xi, const rvec xj, rvec dx)
109 return pbc_dx_aiuc(pbc, xi, xj, dx);
113 rvec_sub(xi, xj, dx);
118 #ifdef GMX_SIMD_HAVE_REAL
120 /* SIMD PBC data structure, containing 1/boxdiag and the box vectors */
122 gmx_simd_real_t inv_bzz;
123 gmx_simd_real_t inv_byy;
124 gmx_simd_real_t inv_bxx;
133 /*! \brief Set the SIMD pbc data from a normal t_pbc struct */
134 static void set_pbc_simd(const t_pbc *pbc, pbc_simd_t *pbc_simd)
139 /* Setting inv_bdiag to 0 effectively turns off PBC */
140 clear_rvec(inv_bdiag);
143 for (d = 0; d < pbc->ndim_ePBC; d++)
145 inv_bdiag[d] = 1.0/pbc->box[d][d];
149 pbc_simd->inv_bzz = gmx_simd_set1_r(inv_bdiag[ZZ]);
150 pbc_simd->inv_byy = gmx_simd_set1_r(inv_bdiag[YY]);
151 pbc_simd->inv_bxx = gmx_simd_set1_r(inv_bdiag[XX]);
155 pbc_simd->bzx = gmx_simd_set1_r(pbc->box[ZZ][XX]);
156 pbc_simd->bzy = gmx_simd_set1_r(pbc->box[ZZ][YY]);
157 pbc_simd->bzz = gmx_simd_set1_r(pbc->box[ZZ][ZZ]);
158 pbc_simd->byx = gmx_simd_set1_r(pbc->box[YY][XX]);
159 pbc_simd->byy = gmx_simd_set1_r(pbc->box[YY][YY]);
160 pbc_simd->bxx = gmx_simd_set1_r(pbc->box[XX][XX]);
164 pbc_simd->bzx = gmx_simd_setzero_r();
165 pbc_simd->bzy = gmx_simd_setzero_r();
166 pbc_simd->bzz = gmx_simd_setzero_r();
167 pbc_simd->byx = gmx_simd_setzero_r();
168 pbc_simd->byy = gmx_simd_setzero_r();
169 pbc_simd->bxx = gmx_simd_setzero_r();
173 /*! \brief Correct distance vector *dx,*dy,*dz for PBC using SIMD */
174 static gmx_inline void
175 pbc_dx_simd(gmx_simd_real_t *dx, gmx_simd_real_t *dy, gmx_simd_real_t *dz,
176 const pbc_simd_t *pbc)
180 sh = gmx_simd_round_r(gmx_simd_mul_r(*dz, pbc->inv_bzz));
181 *dx = gmx_simd_fnmadd_r(sh, pbc->bzx, *dx);
182 *dy = gmx_simd_fnmadd_r(sh, pbc->bzy, *dy);
183 *dz = gmx_simd_fnmadd_r(sh, pbc->bzz, *dz);
185 sh = gmx_simd_round_r(gmx_simd_mul_r(*dy, pbc->inv_byy));
186 *dx = gmx_simd_fnmadd_r(sh, pbc->byx, *dx);
187 *dy = gmx_simd_fnmadd_r(sh, pbc->byy, *dy);
189 sh = gmx_simd_round_r(gmx_simd_mul_r(*dx, pbc->inv_bxx));
190 *dx = gmx_simd_fnmadd_r(sh, pbc->bxx, *dx);
193 #endif /* GMX_SIMD_HAVE_REAL */
195 /*! \brief Morse potential bond
197 * By Frank Everdij. Three parameters needed:
199 * b0 = equilibrium distance in nm
200 * be = beta in nm^-1 (actually, it's nu_e*Sqrt(2*pi*pi*mu/D_e))
201 * cb = well depth in kJ/mol
203 * Note: the potential is referenced to be +cb at infinite separation
204 * and zero at the equilibrium distance!
206 real morse_bonds(int nbonds,
207 const t_iatom forceatoms[], const t_iparams forceparams[],
208 const rvec x[], rvec f[], rvec fshift[],
209 const t_pbc *pbc, const t_graph *g,
210 real lambda, real *dvdlambda,
211 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
212 int gmx_unused *global_atom_index)
214 const real one = 1.0;
215 const real two = 2.0;
216 real dr, dr2, temp, omtemp, cbomtemp, fbond, vbond, fij, vtot;
217 real b0, be, cb, b0A, beA, cbA, b0B, beB, cbB, L1;
219 int i, m, ki, type, ai, aj;
223 for (i = 0; (i < nbonds); )
225 type = forceatoms[i++];
226 ai = forceatoms[i++];
227 aj = forceatoms[i++];
229 b0A = forceparams[type].morse.b0A;
230 beA = forceparams[type].morse.betaA;
231 cbA = forceparams[type].morse.cbA;
233 b0B = forceparams[type].morse.b0B;
234 beB = forceparams[type].morse.betaB;
235 cbB = forceparams[type].morse.cbB;
237 L1 = one-lambda; /* 1 */
238 b0 = L1*b0A + lambda*b0B; /* 3 */
239 be = L1*beA + lambda*beB; /* 3 */
240 cb = L1*cbA + lambda*cbB; /* 3 */
242 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
243 dr2 = iprod(dx, dx); /* 5 */
244 dr = dr2*gmx_invsqrt(dr2); /* 10 */
245 temp = exp(-be*(dr-b0)); /* 12 */
249 /* bonds are constrainted. This may _not_ include bond constraints if they are lambda dependent */
250 *dvdlambda += cbB-cbA;
254 omtemp = one-temp; /* 1 */
255 cbomtemp = cb*omtemp; /* 1 */
256 vbond = cbomtemp*omtemp; /* 1 */
257 fbond = -two*be*temp*cbomtemp*gmx_invsqrt(dr2); /* 9 */
258 vtot += vbond; /* 1 */
260 *dvdlambda += (cbB - cbA) * omtemp * omtemp - (2-2*omtemp)*omtemp * cb * ((b0B-b0A)*be - (beB-beA)*(dr-b0)); /* 15 */
264 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
268 for (m = 0; (m < DIM); m++) /* 15 */
273 fshift[ki][m] += fij;
274 fshift[CENTRAL][m] -= fij;
281 real cubic_bonds(int nbonds,
282 const t_iatom forceatoms[], const t_iparams forceparams[],
283 const rvec x[], rvec f[], rvec fshift[],
284 const t_pbc *pbc, const t_graph *g,
285 real gmx_unused lambda, real gmx_unused *dvdlambda,
286 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
287 int gmx_unused *global_atom_index)
289 const real three = 3.0;
290 const real two = 2.0;
292 real dr, dr2, dist, kdist, kdist2, fbond, vbond, fij, vtot;
294 int i, m, ki, type, ai, aj;
298 for (i = 0; (i < nbonds); )
300 type = forceatoms[i++];
301 ai = forceatoms[i++];
302 aj = forceatoms[i++];
304 b0 = forceparams[type].cubic.b0;
305 kb = forceparams[type].cubic.kb;
306 kcub = forceparams[type].cubic.kcub;
308 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
309 dr2 = iprod(dx, dx); /* 5 */
316 dr = dr2*gmx_invsqrt(dr2); /* 10 */
321 vbond = kdist2 + kcub*kdist2*dist;
322 fbond = -(two*kdist + three*kdist2*kcub)/dr;
324 vtot += vbond; /* 21 */
328 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
331 for (m = 0; (m < DIM); m++) /* 15 */
336 fshift[ki][m] += fij;
337 fshift[CENTRAL][m] -= fij;
343 real FENE_bonds(int nbonds,
344 const t_iatom forceatoms[], const t_iparams forceparams[],
345 const rvec x[], rvec f[], rvec fshift[],
346 const t_pbc *pbc, const t_graph *g,
347 real gmx_unused lambda, real gmx_unused *dvdlambda,
348 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
349 int *global_atom_index)
351 const real half = 0.5;
352 const real one = 1.0;
354 real dr2, bm2, omdr2obm2, fbond, vbond, fij, vtot;
356 int i, m, ki, type, ai, aj;
360 for (i = 0; (i < nbonds); )
362 type = forceatoms[i++];
363 ai = forceatoms[i++];
364 aj = forceatoms[i++];
366 bm = forceparams[type].fene.bm;
367 kb = forceparams[type].fene.kb;
369 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
370 dr2 = iprod(dx, dx); /* 5 */
382 "r^2 (%f) >= bm^2 (%f) in FENE bond between atoms %d and %d",
384 glatnr(global_atom_index, ai),
385 glatnr(global_atom_index, aj));
388 omdr2obm2 = one - dr2/bm2;
390 vbond = -half*kb*bm2*log(omdr2obm2);
391 fbond = -kb/omdr2obm2;
393 vtot += vbond; /* 35 */
397 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
400 for (m = 0; (m < DIM); m++) /* 15 */
405 fshift[ki][m] += fij;
406 fshift[CENTRAL][m] -= fij;
412 real harmonic(real kA, real kB, real xA, real xB, real x, real lambda,
415 const real half = 0.5;
416 real L1, kk, x0, dx, dx2;
417 real v, f, dvdlambda;
420 kk = L1*kA+lambda*kB;
421 x0 = L1*xA+lambda*xB;
428 dvdlambda = half*(kB-kA)*dx2 + (xA-xB)*kk*dx;
435 /* That was 19 flops */
439 real bonds(int nbonds,
440 const t_iatom forceatoms[], const t_iparams forceparams[],
441 const rvec x[], rvec f[], rvec fshift[],
442 const t_pbc *pbc, const t_graph *g,
443 real lambda, real *dvdlambda,
444 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
445 int gmx_unused *global_atom_index)
447 int i, m, ki, ai, aj, type;
448 real dr, dr2, fbond, vbond, fij, vtot;
453 for (i = 0; (i < nbonds); )
455 type = forceatoms[i++];
456 ai = forceatoms[i++];
457 aj = forceatoms[i++];
459 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
460 dr2 = iprod(dx, dx); /* 5 */
461 dr = dr2*gmx_invsqrt(dr2); /* 10 */
463 *dvdlambda += harmonic(forceparams[type].harmonic.krA,
464 forceparams[type].harmonic.krB,
465 forceparams[type].harmonic.rA,
466 forceparams[type].harmonic.rB,
467 dr, lambda, &vbond, &fbond); /* 19 */
475 vtot += vbond; /* 1*/
476 fbond *= gmx_invsqrt(dr2); /* 6 */
480 fprintf(debug, "BONDS: dr = %10g vbond = %10g fbond = %10g\n",
486 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
489 for (m = 0; (m < DIM); m++) /* 15 */
494 fshift[ki][m] += fij;
495 fshift[CENTRAL][m] -= fij;
501 real restraint_bonds(int nbonds,
502 const t_iatom forceatoms[], const t_iparams forceparams[],
503 const rvec x[], rvec f[], rvec fshift[],
504 const t_pbc *pbc, const t_graph *g,
505 real lambda, real *dvdlambda,
506 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
507 int gmx_unused *global_atom_index)
509 int i, m, ki, ai, aj, type;
510 real dr, dr2, fbond, vbond, fij, vtot;
512 real low, dlow, up1, dup1, up2, dup2, k, dk;
520 for (i = 0; (i < nbonds); )
522 type = forceatoms[i++];
523 ai = forceatoms[i++];
524 aj = forceatoms[i++];
526 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
527 dr2 = iprod(dx, dx); /* 5 */
528 dr = dr2*gmx_invsqrt(dr2); /* 10 */
530 low = L1*forceparams[type].restraint.lowA + lambda*forceparams[type].restraint.lowB;
531 dlow = -forceparams[type].restraint.lowA + forceparams[type].restraint.lowB;
532 up1 = L1*forceparams[type].restraint.up1A + lambda*forceparams[type].restraint.up1B;
533 dup1 = -forceparams[type].restraint.up1A + forceparams[type].restraint.up1B;
534 up2 = L1*forceparams[type].restraint.up2A + lambda*forceparams[type].restraint.up2B;
535 dup2 = -forceparams[type].restraint.up2A + forceparams[type].restraint.up2B;
536 k = L1*forceparams[type].restraint.kA + lambda*forceparams[type].restraint.kB;
537 dk = -forceparams[type].restraint.kA + forceparams[type].restraint.kB;
546 *dvdlambda += 0.5*dk*drh2 - k*dlow*drh;
559 *dvdlambda += 0.5*dk*drh2 - k*dup1*drh;
564 vbond = k*(up2 - up1)*(0.5*(up2 - up1) + drh);
565 fbond = -k*(up2 - up1);
566 *dvdlambda += dk*(up2 - up1)*(0.5*(up2 - up1) + drh)
567 + k*(dup2 - dup1)*(up2 - up1 + drh)
568 - k*(up2 - up1)*dup2;
576 vtot += vbond; /* 1*/
577 fbond *= gmx_invsqrt(dr2); /* 6 */
581 fprintf(debug, "BONDS: dr = %10g vbond = %10g fbond = %10g\n",
587 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
590 for (m = 0; (m < DIM); m++) /* 15 */
595 fshift[ki][m] += fij;
596 fshift[CENTRAL][m] -= fij;
603 real polarize(int nbonds,
604 const t_iatom forceatoms[], const t_iparams forceparams[],
605 const rvec x[], rvec f[], rvec fshift[],
606 const t_pbc *pbc, const t_graph *g,
607 real lambda, real *dvdlambda,
608 const t_mdatoms *md, t_fcdata gmx_unused *fcd,
609 int gmx_unused *global_atom_index)
611 int i, m, ki, ai, aj, type;
612 real dr, dr2, fbond, vbond, fij, vtot, ksh;
617 for (i = 0; (i < nbonds); )
619 type = forceatoms[i++];
620 ai = forceatoms[i++];
621 aj = forceatoms[i++];
622 ksh = sqr(md->chargeA[aj])*ONE_4PI_EPS0/forceparams[type].polarize.alpha;
625 fprintf(debug, "POL: local ai = %d aj = %d ksh = %.3f\n", ai, aj, ksh);
628 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
629 dr2 = iprod(dx, dx); /* 5 */
630 dr = dr2*gmx_invsqrt(dr2); /* 10 */
632 *dvdlambda += harmonic(ksh, ksh, 0, 0, dr, lambda, &vbond, &fbond); /* 19 */
639 vtot += vbond; /* 1*/
640 fbond *= gmx_invsqrt(dr2); /* 6 */
644 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
647 for (m = 0; (m < DIM); m++) /* 15 */
652 fshift[ki][m] += fij;
653 fshift[CENTRAL][m] -= fij;
659 real anharm_polarize(int nbonds,
660 const t_iatom forceatoms[], const t_iparams forceparams[],
661 const rvec x[], rvec f[], rvec fshift[],
662 const t_pbc *pbc, const t_graph *g,
663 real lambda, real *dvdlambda,
664 const t_mdatoms *md, t_fcdata gmx_unused *fcd,
665 int gmx_unused *global_atom_index)
667 int i, m, ki, ai, aj, type;
668 real dr, dr2, fbond, vbond, fij, vtot, ksh, khyp, drcut, ddr, ddr3;
673 for (i = 0; (i < nbonds); )
675 type = forceatoms[i++];
676 ai = forceatoms[i++];
677 aj = forceatoms[i++];
678 ksh = sqr(md->chargeA[aj])*ONE_4PI_EPS0/forceparams[type].anharm_polarize.alpha; /* 7*/
679 khyp = forceparams[type].anharm_polarize.khyp;
680 drcut = forceparams[type].anharm_polarize.drcut;
683 fprintf(debug, "POL: local ai = %d aj = %d ksh = %.3f\n", ai, aj, ksh);
686 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
687 dr2 = iprod(dx, dx); /* 5 */
688 dr = dr2*gmx_invsqrt(dr2); /* 10 */
690 *dvdlambda += harmonic(ksh, ksh, 0, 0, dr, lambda, &vbond, &fbond); /* 19 */
701 vbond += khyp*ddr*ddr3;
702 fbond -= 4*khyp*ddr3;
704 fbond *= gmx_invsqrt(dr2); /* 6 */
705 vtot += vbond; /* 1*/
709 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
712 for (m = 0; (m < DIM); m++) /* 15 */
717 fshift[ki][m] += fij;
718 fshift[CENTRAL][m] -= fij;
724 real water_pol(int nbonds,
725 const t_iatom forceatoms[], const t_iparams forceparams[],
726 const rvec x[], rvec f[], rvec gmx_unused fshift[],
727 const t_pbc gmx_unused *pbc, const t_graph gmx_unused *g,
728 real gmx_unused lambda, real gmx_unused *dvdlambda,
729 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
730 int gmx_unused *global_atom_index)
732 /* This routine implements anisotropic polarizibility for water, through
733 * a shell connected to a dummy with spring constant that differ in the
734 * three spatial dimensions in the molecular frame.
736 int i, m, aO, aH1, aH2, aD, aS, type, type0, ki;
738 rvec dOH1, dOH2, dHH, dOD, dDS, nW, kk, dx, kdx, proj;
742 real vtot, fij, r_HH, r_OD, r_nW, tx, ty, tz, qS;
747 type0 = forceatoms[0];
749 qS = md->chargeA[aS];
750 kk[XX] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_x;
751 kk[YY] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_y;
752 kk[ZZ] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_z;
753 r_HH = 1.0/forceparams[type0].wpol.rHH;
756 fprintf(debug, "WPOL: qS = %10.5f aS = %5d\n", qS, aS);
757 fprintf(debug, "WPOL: kk = %10.3f %10.3f %10.3f\n",
758 kk[XX], kk[YY], kk[ZZ]);
759 fprintf(debug, "WPOL: rOH = %10.3f rHH = %10.3f rOD = %10.3f\n",
760 forceparams[type0].wpol.rOH,
761 forceparams[type0].wpol.rHH,
762 forceparams[type0].wpol.rOD);
764 for (i = 0; (i < nbonds); i += 6)
766 type = forceatoms[i];
769 gmx_fatal(FARGS, "Sorry, type = %d, type0 = %d, file = %s, line = %d",
770 type, type0, __FILE__, __LINE__);
772 aO = forceatoms[i+1];
773 aH1 = forceatoms[i+2];
774 aH2 = forceatoms[i+3];
775 aD = forceatoms[i+4];
776 aS = forceatoms[i+5];
778 /* Compute vectors describing the water frame */
779 pbc_rvec_sub(pbc, x[aH1], x[aO], dOH1);
780 pbc_rvec_sub(pbc, x[aH2], x[aO], dOH2);
781 pbc_rvec_sub(pbc, x[aH2], x[aH1], dHH);
782 pbc_rvec_sub(pbc, x[aD], x[aO], dOD);
783 ki = pbc_rvec_sub(pbc, x[aS], x[aD], dDS);
784 cprod(dOH1, dOH2, nW);
786 /* Compute inverse length of normal vector
787 * (this one could be precomputed, but I'm too lazy now)
789 r_nW = gmx_invsqrt(iprod(nW, nW));
790 /* This is for precision, but does not make a big difference,
793 r_OD = gmx_invsqrt(iprod(dOD, dOD));
795 /* Normalize the vectors in the water frame */
797 svmul(r_HH, dHH, dHH);
798 svmul(r_OD, dOD, dOD);
800 /* Compute displacement of shell along components of the vector */
801 dx[ZZ] = iprod(dDS, dOD);
802 /* Compute projection on the XY plane: dDS - dx[ZZ]*dOD */
803 for (m = 0; (m < DIM); m++)
805 proj[m] = dDS[m]-dx[ZZ]*dOD[m];
808 /*dx[XX] = iprod(dDS,nW);
809 dx[YY] = iprod(dDS,dHH);*/
810 dx[XX] = iprod(proj, nW);
811 for (m = 0; (m < DIM); m++)
813 proj[m] -= dx[XX]*nW[m];
815 dx[YY] = iprod(proj, dHH);
820 fprintf(debug, "WPOL: dx2=%10g dy2=%10g dz2=%10g sum=%10g dDS^2=%10g\n",
821 sqr(dx[XX]), sqr(dx[YY]), sqr(dx[ZZ]), iprod(dx, dx), iprod(dDS, dDS));
822 fprintf(debug, "WPOL: dHH=(%10g,%10g,%10g)\n", dHH[XX], dHH[YY], dHH[ZZ]);
823 fprintf(debug, "WPOL: dOD=(%10g,%10g,%10g), 1/r_OD = %10g\n",
824 dOD[XX], dOD[YY], dOD[ZZ], 1/r_OD);
825 fprintf(debug, "WPOL: nW =(%10g,%10g,%10g), 1/r_nW = %10g\n",
826 nW[XX], nW[YY], nW[ZZ], 1/r_nW);
827 fprintf(debug, "WPOL: dx =%10g, dy =%10g, dz =%10g\n",
828 dx[XX], dx[YY], dx[ZZ]);
829 fprintf(debug, "WPOL: dDSx=%10g, dDSy=%10g, dDSz=%10g\n",
830 dDS[XX], dDS[YY], dDS[ZZ]);
833 /* Now compute the forces and energy */
834 kdx[XX] = kk[XX]*dx[XX];
835 kdx[YY] = kk[YY]*dx[YY];
836 kdx[ZZ] = kk[ZZ]*dx[ZZ];
837 vtot += iprod(dx, kdx);
841 ivec_sub(SHIFT_IVEC(g, aS), SHIFT_IVEC(g, aD), dt);
845 for (m = 0; (m < DIM); m++)
847 /* This is a tensor operation but written out for speed */
857 fshift[ki][m] += fij;
858 fshift[CENTRAL][m] -= fij;
863 fprintf(debug, "WPOL: vwpol=%g\n", 0.5*iprod(dx, kdx));
864 fprintf(debug, "WPOL: df = (%10g, %10g, %10g)\n", df[XX], df[YY], df[ZZ]);
872 static real do_1_thole(const rvec xi, const rvec xj, rvec fi, rvec fj,
873 const t_pbc *pbc, real qq,
874 rvec fshift[], real afac)
877 real r12sq, r12_1, r12bar, v0, v1, fscal, ebar, fff;
880 t = pbc_rvec_sub(pbc, xi, xj, r12); /* 3 */
882 r12sq = iprod(r12, r12); /* 5 */
883 r12_1 = gmx_invsqrt(r12sq); /* 5 */
884 r12bar = afac/r12_1; /* 5 */
885 v0 = qq*ONE_4PI_EPS0*r12_1; /* 2 */
886 ebar = exp(-r12bar); /* 5 */
887 v1 = (1-(1+0.5*r12bar)*ebar); /* 4 */
888 fscal = ((v0*r12_1)*v1 - v0*0.5*afac*ebar*(r12bar+1))*r12_1; /* 9 */
891 fprintf(debug, "THOLE: v0 = %.3f v1 = %.3f r12= % .3f r12bar = %.3f fscal = %.3f ebar = %.3f\n", v0, v1, 1/r12_1, r12bar, fscal, ebar);
894 for (m = 0; (m < DIM); m++)
900 fshift[CENTRAL][m] -= fff;
903 return v0*v1; /* 1 */
907 real thole_pol(int nbonds,
908 const t_iatom forceatoms[], const t_iparams forceparams[],
909 const rvec x[], rvec f[], rvec fshift[],
910 const t_pbc *pbc, const t_graph gmx_unused *g,
911 real gmx_unused lambda, real gmx_unused *dvdlambda,
912 const t_mdatoms *md, t_fcdata gmx_unused *fcd,
913 int gmx_unused *global_atom_index)
915 /* Interaction between two pairs of particles with opposite charge */
916 int i, type, a1, da1, a2, da2;
917 real q1, q2, qq, a, al1, al2, afac;
919 const real minusOneOnSix = -1.0/6.0;
921 for (i = 0; (i < nbonds); )
923 type = forceatoms[i++];
924 a1 = forceatoms[i++];
925 da1 = forceatoms[i++];
926 a2 = forceatoms[i++];
927 da2 = forceatoms[i++];
928 q1 = md->chargeA[da1];
929 q2 = md->chargeA[da2];
930 a = forceparams[type].thole.a;
931 al1 = forceparams[type].thole.alpha1;
932 al2 = forceparams[type].thole.alpha2;
934 afac = a*pow(al1*al2, minusOneOnSix);
935 V += do_1_thole(x[a1], x[a2], f[a1], f[a2], pbc, qq, fshift, afac);
936 V += do_1_thole(x[da1], x[a2], f[da1], f[a2], pbc, -qq, fshift, afac);
937 V += do_1_thole(x[a1], x[da2], f[a1], f[da2], pbc, -qq, fshift, afac);
938 V += do_1_thole(x[da1], x[da2], f[da1], f[da2], pbc, qq, fshift, afac);
944 real bond_angle(const rvec xi, const rvec xj, const rvec xk, const t_pbc *pbc,
945 rvec r_ij, rvec r_kj, real *costh,
947 /* Return value is the angle between the bonds i-j and j-k */
952 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
953 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
955 *costh = cos_angle(r_ij, r_kj); /* 25 */
956 th = acos(*costh); /* 10 */
961 real angles(int nbonds,
962 const t_iatom forceatoms[], const t_iparams forceparams[],
963 const rvec x[], rvec f[], rvec fshift[],
964 const t_pbc *pbc, const t_graph *g,
965 real lambda, real *dvdlambda,
966 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
967 int gmx_unused *global_atom_index)
969 int i, ai, aj, ak, t1, t2, type;
971 real cos_theta, cos_theta2, theta, dVdt, va, vtot;
972 ivec jt, dt_ij, dt_kj;
975 for (i = 0; i < nbonds; )
977 type = forceatoms[i++];
978 ai = forceatoms[i++];
979 aj = forceatoms[i++];
980 ak = forceatoms[i++];
982 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
983 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
985 *dvdlambda += harmonic(forceparams[type].harmonic.krA,
986 forceparams[type].harmonic.krB,
987 forceparams[type].harmonic.rA*DEG2RAD,
988 forceparams[type].harmonic.rB*DEG2RAD,
989 theta, lambda, &va, &dVdt); /* 21 */
992 cos_theta2 = sqr(cos_theta);
1002 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1003 sth = st*cos_theta; /* 1 */
1007 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1008 theta*RAD2DEG, va, dVdt);
1011 nrij2 = iprod(r_ij, r_ij); /* 5 */
1012 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1014 nrij_1 = gmx_invsqrt(nrij2); /* 10 */
1015 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1017 cik = st*nrij_1*nrkj_1; /* 2 */
1018 cii = sth*nrij_1*nrij_1; /* 2 */
1019 ckk = sth*nrkj_1*nrkj_1; /* 2 */
1021 for (m = 0; m < DIM; m++)
1023 f_i[m] = -(cik*r_kj[m] - cii*r_ij[m]);
1024 f_k[m] = -(cik*r_ij[m] - ckk*r_kj[m]);
1025 f_j[m] = -f_i[m] - f_k[m];
1032 copy_ivec(SHIFT_IVEC(g, aj), jt);
1034 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1035 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1036 t1 = IVEC2IS(dt_ij);
1037 t2 = IVEC2IS(dt_kj);
1039 rvec_inc(fshift[t1], f_i);
1040 rvec_inc(fshift[CENTRAL], f_j);
1041 rvec_inc(fshift[t2], f_k);
1048 #ifdef GMX_SIMD_HAVE_REAL
1050 /* As angles, but using SIMD to calculate many dihedrals at once.
1051 * This routines does not calculate energies and shift forces.
1053 static gmx_inline void
1054 angles_noener_simd(int nbonds,
1055 const t_iatom forceatoms[], const t_iparams forceparams[],
1056 const rvec x[], rvec f[],
1057 const t_pbc *pbc, const t_graph gmx_unused *g,
1058 real gmx_unused lambda,
1059 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1060 int gmx_unused *global_atom_index)
1064 int type, ai[GMX_SIMD_REAL_WIDTH], aj[GMX_SIMD_REAL_WIDTH];
1065 int ak[GMX_SIMD_REAL_WIDTH];
1066 real coeff_array[2*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *coeff;
1067 real dr_array[2*DIM*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *dr;
1068 real f_buf_array[6*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *f_buf;
1069 gmx_simd_real_t k_S, theta0_S;
1070 gmx_simd_real_t rijx_S, rijy_S, rijz_S;
1071 gmx_simd_real_t rkjx_S, rkjy_S, rkjz_S;
1072 gmx_simd_real_t one_S;
1073 gmx_simd_real_t min_one_plus_eps_S;
1074 gmx_simd_real_t rij_rkj_S;
1075 gmx_simd_real_t nrij2_S, nrij_1_S;
1076 gmx_simd_real_t nrkj2_S, nrkj_1_S;
1077 gmx_simd_real_t cos_S, invsin_S;
1078 gmx_simd_real_t theta_S;
1079 gmx_simd_real_t st_S, sth_S;
1080 gmx_simd_real_t cik_S, cii_S, ckk_S;
1081 gmx_simd_real_t f_ix_S, f_iy_S, f_iz_S;
1082 gmx_simd_real_t f_kx_S, f_ky_S, f_kz_S;
1083 pbc_simd_t pbc_simd;
1085 /* Ensure register memory alignment */
1086 coeff = gmx_simd_align_r(coeff_array);
1087 dr = gmx_simd_align_r(dr_array);
1088 f_buf = gmx_simd_align_r(f_buf_array);
1090 set_pbc_simd(pbc, &pbc_simd);
1092 one_S = gmx_simd_set1_r(1.0);
1094 /* The smallest number > -1 */
1095 min_one_plus_eps_S = gmx_simd_set1_r(-1.0 + 2*GMX_REAL_EPS);
1097 /* nbonds is the number of angles times nfa1, here we step GMX_SIMD_REAL_WIDTH angles */
1098 for (i = 0; (i < nbonds); i += GMX_SIMD_REAL_WIDTH*nfa1)
1100 /* Collect atoms for GMX_SIMD_REAL_WIDTH angles.
1101 * iu indexes into forceatoms, we should not let iu go beyond nbonds.
1104 for (s = 0; s < GMX_SIMD_REAL_WIDTH; s++)
1106 type = forceatoms[iu];
1107 ai[s] = forceatoms[iu+1];
1108 aj[s] = forceatoms[iu+2];
1109 ak[s] = forceatoms[iu+3];
1111 coeff[s] = forceparams[type].harmonic.krA;
1112 coeff[GMX_SIMD_REAL_WIDTH+s] = forceparams[type].harmonic.rA*DEG2RAD;
1114 /* If you can't use pbc_dx_simd below for PBC, e.g. because
1115 * you can't round in SIMD, use pbc_rvec_sub here.
1117 /* Store the non PBC corrected distances packed and aligned */
1118 for (m = 0; m < DIM; m++)
1120 dr[s + m *GMX_SIMD_REAL_WIDTH] = x[ai[s]][m] - x[aj[s]][m];
1121 dr[s + (DIM+m)*GMX_SIMD_REAL_WIDTH] = x[ak[s]][m] - x[aj[s]][m];
1124 /* At the end fill the arrays with identical entries */
1125 if (iu + nfa1 < nbonds)
1131 k_S = gmx_simd_load_r(coeff);
1132 theta0_S = gmx_simd_load_r(coeff+GMX_SIMD_REAL_WIDTH);
1134 rijx_S = gmx_simd_load_r(dr + 0*GMX_SIMD_REAL_WIDTH);
1135 rijy_S = gmx_simd_load_r(dr + 1*GMX_SIMD_REAL_WIDTH);
1136 rijz_S = gmx_simd_load_r(dr + 2*GMX_SIMD_REAL_WIDTH);
1137 rkjx_S = gmx_simd_load_r(dr + 3*GMX_SIMD_REAL_WIDTH);
1138 rkjy_S = gmx_simd_load_r(dr + 4*GMX_SIMD_REAL_WIDTH);
1139 rkjz_S = gmx_simd_load_r(dr + 5*GMX_SIMD_REAL_WIDTH);
1141 pbc_dx_simd(&rijx_S, &rijy_S, &rijz_S, &pbc_simd);
1142 pbc_dx_simd(&rkjx_S, &rkjy_S, &rkjz_S, &pbc_simd);
1144 rij_rkj_S = gmx_simd_iprod_r(rijx_S, rijy_S, rijz_S,
1145 rkjx_S, rkjy_S, rkjz_S);
1147 nrij2_S = gmx_simd_norm2_r(rijx_S, rijy_S, rijz_S);
1148 nrkj2_S = gmx_simd_norm2_r(rkjx_S, rkjy_S, rkjz_S);
1150 nrij_1_S = gmx_simd_invsqrt_r(nrij2_S);
1151 nrkj_1_S = gmx_simd_invsqrt_r(nrkj2_S);
1153 cos_S = gmx_simd_mul_r(rij_rkj_S, gmx_simd_mul_r(nrij_1_S, nrkj_1_S));
1155 /* To allow for 180 degrees, we take the max of cos and -1 + 1bit,
1156 * so we can safely get the 1/sin from 1/sqrt(1 - cos^2).
1157 * This also ensures that rounding errors would cause the argument
1158 * of gmx_simd_acos_r to be < -1.
1159 * Note that we do not take precautions for cos(0)=1, so the outer
1160 * atoms in an angle should not be on top of each other.
1162 cos_S = gmx_simd_max_r(cos_S, min_one_plus_eps_S);
1164 theta_S = gmx_simd_acos_r(cos_S);
1166 invsin_S = gmx_simd_invsqrt_r(gmx_simd_sub_r(one_S, gmx_simd_mul_r(cos_S, cos_S)));
1168 st_S = gmx_simd_mul_r(gmx_simd_mul_r(k_S, gmx_simd_sub_r(theta0_S, theta_S)),
1170 sth_S = gmx_simd_mul_r(st_S, cos_S);
1172 cik_S = gmx_simd_mul_r(st_S, gmx_simd_mul_r(nrij_1_S, nrkj_1_S));
1173 cii_S = gmx_simd_mul_r(sth_S, gmx_simd_mul_r(nrij_1_S, nrij_1_S));
1174 ckk_S = gmx_simd_mul_r(sth_S, gmx_simd_mul_r(nrkj_1_S, nrkj_1_S));
1176 f_ix_S = gmx_simd_mul_r(cii_S, rijx_S);
1177 f_ix_S = gmx_simd_fnmadd_r(cik_S, rkjx_S, f_ix_S);
1178 f_iy_S = gmx_simd_mul_r(cii_S, rijy_S);
1179 f_iy_S = gmx_simd_fnmadd_r(cik_S, rkjy_S, f_iy_S);
1180 f_iz_S = gmx_simd_mul_r(cii_S, rijz_S);
1181 f_iz_S = gmx_simd_fnmadd_r(cik_S, rkjz_S, f_iz_S);
1182 f_kx_S = gmx_simd_mul_r(ckk_S, rkjx_S);
1183 f_kx_S = gmx_simd_fnmadd_r(cik_S, rijx_S, f_kx_S);
1184 f_ky_S = gmx_simd_mul_r(ckk_S, rkjy_S);
1185 f_ky_S = gmx_simd_fnmadd_r(cik_S, rijy_S, f_ky_S);
1186 f_kz_S = gmx_simd_mul_r(ckk_S, rkjz_S);
1187 f_kz_S = gmx_simd_fnmadd_r(cik_S, rijz_S, f_kz_S);
1189 gmx_simd_store_r(f_buf + 0*GMX_SIMD_REAL_WIDTH, f_ix_S);
1190 gmx_simd_store_r(f_buf + 1*GMX_SIMD_REAL_WIDTH, f_iy_S);
1191 gmx_simd_store_r(f_buf + 2*GMX_SIMD_REAL_WIDTH, f_iz_S);
1192 gmx_simd_store_r(f_buf + 3*GMX_SIMD_REAL_WIDTH, f_kx_S);
1193 gmx_simd_store_r(f_buf + 4*GMX_SIMD_REAL_WIDTH, f_ky_S);
1194 gmx_simd_store_r(f_buf + 5*GMX_SIMD_REAL_WIDTH, f_kz_S);
1200 for (m = 0; m < DIM; m++)
1202 f[ai[s]][m] += f_buf[s + m*GMX_SIMD_REAL_WIDTH];
1203 f[aj[s]][m] -= f_buf[s + m*GMX_SIMD_REAL_WIDTH] + f_buf[s + (DIM+m)*GMX_SIMD_REAL_WIDTH];
1204 f[ak[s]][m] += f_buf[s + (DIM+m)*GMX_SIMD_REAL_WIDTH];
1209 while (s < GMX_SIMD_REAL_WIDTH && iu < nbonds);
1213 #endif /* GMX_SIMD_HAVE_REAL */
1215 real linear_angles(int nbonds,
1216 const t_iatom forceatoms[], const t_iparams forceparams[],
1217 const rvec x[], rvec f[], rvec fshift[],
1218 const t_pbc *pbc, const t_graph *g,
1219 real lambda, real *dvdlambda,
1220 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1221 int gmx_unused *global_atom_index)
1223 int i, m, ai, aj, ak, t1, t2, type;
1225 real L1, kA, kB, aA, aB, dr, dr2, va, vtot, a, b, klin;
1226 ivec jt, dt_ij, dt_kj;
1227 rvec r_ij, r_kj, r_ik, dx;
1231 for (i = 0; (i < nbonds); )
1233 type = forceatoms[i++];
1234 ai = forceatoms[i++];
1235 aj = forceatoms[i++];
1236 ak = forceatoms[i++];
1238 kA = forceparams[type].linangle.klinA;
1239 kB = forceparams[type].linangle.klinB;
1240 klin = L1*kA + lambda*kB;
1242 aA = forceparams[type].linangle.aA;
1243 aB = forceparams[type].linangle.aB;
1244 a = L1*aA+lambda*aB;
1247 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
1248 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
1249 rvec_sub(r_ij, r_kj, r_ik);
1252 for (m = 0; (m < DIM); m++)
1254 dr = -a * r_ij[m] - b * r_kj[m];
1259 f_j[m] = -(f_i[m]+f_k[m]);
1265 *dvdlambda += 0.5*(kB-kA)*dr2 + klin*(aB-aA)*iprod(dx, r_ik);
1271 copy_ivec(SHIFT_IVEC(g, aj), jt);
1273 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1274 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1275 t1 = IVEC2IS(dt_ij);
1276 t2 = IVEC2IS(dt_kj);
1278 rvec_inc(fshift[t1], f_i);
1279 rvec_inc(fshift[CENTRAL], f_j);
1280 rvec_inc(fshift[t2], f_k);
1285 real urey_bradley(int nbonds,
1286 const t_iatom forceatoms[], const t_iparams forceparams[],
1287 const rvec x[], rvec f[], rvec fshift[],
1288 const t_pbc *pbc, const t_graph *g,
1289 real lambda, real *dvdlambda,
1290 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1291 int gmx_unused *global_atom_index)
1293 int i, m, ai, aj, ak, t1, t2, type, ki;
1294 rvec r_ij, r_kj, r_ik;
1295 real cos_theta, cos_theta2, theta;
1296 real dVdt, va, vtot, dr, dr2, vbond, fbond, fik;
1297 real kthA, th0A, kUBA, r13A, kthB, th0B, kUBB, r13B;
1298 ivec jt, dt_ij, dt_kj, dt_ik;
1301 for (i = 0; (i < nbonds); )
1303 type = forceatoms[i++];
1304 ai = forceatoms[i++];
1305 aj = forceatoms[i++];
1306 ak = forceatoms[i++];
1307 th0A = forceparams[type].u_b.thetaA*DEG2RAD;
1308 kthA = forceparams[type].u_b.kthetaA;
1309 r13A = forceparams[type].u_b.r13A;
1310 kUBA = forceparams[type].u_b.kUBA;
1311 th0B = forceparams[type].u_b.thetaB*DEG2RAD;
1312 kthB = forceparams[type].u_b.kthetaB;
1313 r13B = forceparams[type].u_b.r13B;
1314 kUBB = forceparams[type].u_b.kUBB;
1316 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
1317 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
1319 *dvdlambda += harmonic(kthA, kthB, th0A, th0B, theta, lambda, &va, &dVdt); /* 21 */
1322 ki = pbc_rvec_sub(pbc, x[ai], x[ak], r_ik); /* 3 */
1323 dr2 = iprod(r_ik, r_ik); /* 5 */
1324 dr = dr2*gmx_invsqrt(dr2); /* 10 */
1326 *dvdlambda += harmonic(kUBA, kUBB, r13A, r13B, dr, lambda, &vbond, &fbond); /* 19 */
1328 cos_theta2 = sqr(cos_theta); /* 1 */
1336 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1337 sth = st*cos_theta; /* 1 */
1341 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1342 theta*RAD2DEG, va, dVdt);
1345 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1346 nrij2 = iprod(r_ij, r_ij);
1348 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
1349 cii = sth/nrij2; /* 10 */
1350 ckk = sth/nrkj2; /* 10 */
1352 for (m = 0; (m < DIM); m++) /* 39 */
1354 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
1355 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
1356 f_j[m] = -f_i[m]-f_k[m];
1363 copy_ivec(SHIFT_IVEC(g, aj), jt);
1365 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1366 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1367 t1 = IVEC2IS(dt_ij);
1368 t2 = IVEC2IS(dt_kj);
1370 rvec_inc(fshift[t1], f_i);
1371 rvec_inc(fshift[CENTRAL], f_j);
1372 rvec_inc(fshift[t2], f_k);
1374 /* Time for the bond calculations */
1380 vtot += vbond; /* 1*/
1381 fbond *= gmx_invsqrt(dr2); /* 6 */
1385 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ak), dt_ik);
1386 ki = IVEC2IS(dt_ik);
1388 for (m = 0; (m < DIM); m++) /* 15 */
1390 fik = fbond*r_ik[m];
1393 fshift[ki][m] += fik;
1394 fshift[CENTRAL][m] -= fik;
1400 real quartic_angles(int nbonds,
1401 const t_iatom forceatoms[], const t_iparams forceparams[],
1402 const rvec x[], rvec f[], rvec fshift[],
1403 const t_pbc *pbc, const t_graph *g,
1404 real gmx_unused lambda, real gmx_unused *dvdlambda,
1405 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1406 int gmx_unused *global_atom_index)
1408 int i, j, ai, aj, ak, t1, t2, type;
1410 real cos_theta, cos_theta2, theta, dt, dVdt, va, dtp, c, vtot;
1411 ivec jt, dt_ij, dt_kj;
1414 for (i = 0; (i < nbonds); )
1416 type = forceatoms[i++];
1417 ai = forceatoms[i++];
1418 aj = forceatoms[i++];
1419 ak = forceatoms[i++];
1421 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
1422 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
1424 dt = theta - forceparams[type].qangle.theta*DEG2RAD; /* 2 */
1427 va = forceparams[type].qangle.c[0];
1429 for (j = 1; j <= 4; j++)
1431 c = forceparams[type].qangle.c[j];
1440 cos_theta2 = sqr(cos_theta); /* 1 */
1449 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1450 sth = st*cos_theta; /* 1 */
1454 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1455 theta*RAD2DEG, va, dVdt);
1458 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1459 nrij2 = iprod(r_ij, r_ij);
1461 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
1462 cii = sth/nrij2; /* 10 */
1463 ckk = sth/nrkj2; /* 10 */
1465 for (m = 0; (m < DIM); m++) /* 39 */
1467 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
1468 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
1469 f_j[m] = -f_i[m]-f_k[m];
1476 copy_ivec(SHIFT_IVEC(g, aj), jt);
1478 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1479 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1480 t1 = IVEC2IS(dt_ij);
1481 t2 = IVEC2IS(dt_kj);
1483 rvec_inc(fshift[t1], f_i);
1484 rvec_inc(fshift[CENTRAL], f_j);
1485 rvec_inc(fshift[t2], f_k);
1491 real dih_angle(const rvec xi, const rvec xj, const rvec xk, const rvec xl,
1493 rvec r_ij, rvec r_kj, rvec r_kl, rvec m, rvec n,
1494 real *sign, int *t1, int *t2, int *t3)
1498 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
1499 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
1500 *t3 = pbc_rvec_sub(pbc, xk, xl, r_kl); /* 3 */
1502 cprod(r_ij, r_kj, m); /* 9 */
1503 cprod(r_kj, r_kl, n); /* 9 */
1504 phi = gmx_angle(m, n); /* 49 (assuming 25 for atan2) */
1505 ipr = iprod(r_ij, n); /* 5 */
1506 (*sign) = (ipr < 0.0) ? -1.0 : 1.0;
1507 phi = (*sign)*phi; /* 1 */
1513 #ifdef GMX_SIMD_HAVE_REAL
1515 /* As dih_angle above, but calculates 4 dihedral angles at once using SIMD,
1516 * also calculates the pre-factor required for the dihedral force update.
1517 * Note that bv and buf should be register aligned.
1519 static gmx_inline void
1520 dih_angle_simd(const rvec *x,
1521 const int *ai, const int *aj, const int *ak, const int *al,
1522 const pbc_simd_t *pbc,
1524 gmx_simd_real_t *phi_S,
1525 gmx_simd_real_t *mx_S, gmx_simd_real_t *my_S, gmx_simd_real_t *mz_S,
1526 gmx_simd_real_t *nx_S, gmx_simd_real_t *ny_S, gmx_simd_real_t *nz_S,
1527 gmx_simd_real_t *nrkj_m2_S,
1528 gmx_simd_real_t *nrkj_n2_S,
1533 gmx_simd_real_t rijx_S, rijy_S, rijz_S;
1534 gmx_simd_real_t rkjx_S, rkjy_S, rkjz_S;
1535 gmx_simd_real_t rklx_S, rkly_S, rklz_S;
1536 gmx_simd_real_t cx_S, cy_S, cz_S;
1537 gmx_simd_real_t cn_S;
1538 gmx_simd_real_t s_S;
1539 gmx_simd_real_t ipr_S;
1540 gmx_simd_real_t iprm_S, iprn_S;
1541 gmx_simd_real_t nrkj2_S, nrkj_1_S, nrkj_2_S, nrkj_S;
1542 gmx_simd_real_t toler_S;
1543 gmx_simd_real_t p_S, q_S;
1544 gmx_simd_real_t nrkj2_min_S;
1545 gmx_simd_real_t real_eps_S;
1547 /* Used to avoid division by zero.
1548 * We take into acount that we multiply the result by real_eps_S.
1550 nrkj2_min_S = gmx_simd_set1_r(GMX_REAL_MIN/(2*GMX_REAL_EPS));
1552 /* The value of the last significant bit (GMX_REAL_EPS is half of that) */
1553 real_eps_S = gmx_simd_set1_r(2*GMX_REAL_EPS);
1555 for (s = 0; s < GMX_SIMD_REAL_WIDTH; s++)
1557 /* If you can't use pbc_dx_simd below for PBC, e.g. because
1558 * you can't round in SIMD, use pbc_rvec_sub here.
1560 for (m = 0; m < DIM; m++)
1562 dr[s + (0*DIM + m)*GMX_SIMD_REAL_WIDTH] = x[ai[s]][m] - x[aj[s]][m];
1563 dr[s + (1*DIM + m)*GMX_SIMD_REAL_WIDTH] = x[ak[s]][m] - x[aj[s]][m];
1564 dr[s + (2*DIM + m)*GMX_SIMD_REAL_WIDTH] = x[ak[s]][m] - x[al[s]][m];
1568 rijx_S = gmx_simd_load_r(dr + 0*GMX_SIMD_REAL_WIDTH);
1569 rijy_S = gmx_simd_load_r(dr + 1*GMX_SIMD_REAL_WIDTH);
1570 rijz_S = gmx_simd_load_r(dr + 2*GMX_SIMD_REAL_WIDTH);
1571 rkjx_S = gmx_simd_load_r(dr + 3*GMX_SIMD_REAL_WIDTH);
1572 rkjy_S = gmx_simd_load_r(dr + 4*GMX_SIMD_REAL_WIDTH);
1573 rkjz_S = gmx_simd_load_r(dr + 5*GMX_SIMD_REAL_WIDTH);
1574 rklx_S = gmx_simd_load_r(dr + 6*GMX_SIMD_REAL_WIDTH);
1575 rkly_S = gmx_simd_load_r(dr + 7*GMX_SIMD_REAL_WIDTH);
1576 rklz_S = gmx_simd_load_r(dr + 8*GMX_SIMD_REAL_WIDTH);
1578 pbc_dx_simd(&rijx_S, &rijy_S, &rijz_S, pbc);
1579 pbc_dx_simd(&rkjx_S, &rkjy_S, &rkjz_S, pbc);
1580 pbc_dx_simd(&rklx_S, &rkly_S, &rklz_S, pbc);
1582 gmx_simd_cprod_r(rijx_S, rijy_S, rijz_S,
1583 rkjx_S, rkjy_S, rkjz_S,
1586 gmx_simd_cprod_r(rkjx_S, rkjy_S, rkjz_S,
1587 rklx_S, rkly_S, rklz_S,
1590 gmx_simd_cprod_r(*mx_S, *my_S, *mz_S,
1591 *nx_S, *ny_S, *nz_S,
1592 &cx_S, &cy_S, &cz_S);
1594 cn_S = gmx_simd_sqrt_r(gmx_simd_norm2_r(cx_S, cy_S, cz_S));
1596 s_S = gmx_simd_iprod_r(*mx_S, *my_S, *mz_S, *nx_S, *ny_S, *nz_S);
1598 /* Determine the dihedral angle, the sign might need correction */
1599 *phi_S = gmx_simd_atan2_r(cn_S, s_S);
1601 ipr_S = gmx_simd_iprod_r(rijx_S, rijy_S, rijz_S,
1602 *nx_S, *ny_S, *nz_S);
1604 iprm_S = gmx_simd_norm2_r(*mx_S, *my_S, *mz_S);
1605 iprn_S = gmx_simd_norm2_r(*nx_S, *ny_S, *nz_S);
1607 nrkj2_S = gmx_simd_norm2_r(rkjx_S, rkjy_S, rkjz_S);
1609 /* Avoid division by zero. When zero, the result is multiplied by 0
1610 * anyhow, so the 3 max below do not affect the final result.
1612 nrkj2_S = gmx_simd_max_r(nrkj2_S, nrkj2_min_S);
1613 nrkj_1_S = gmx_simd_invsqrt_r(nrkj2_S);
1614 nrkj_2_S = gmx_simd_mul_r(nrkj_1_S, nrkj_1_S);
1615 nrkj_S = gmx_simd_mul_r(nrkj2_S, nrkj_1_S);
1617 toler_S = gmx_simd_mul_r(nrkj2_S, real_eps_S);
1619 /* Here the plain-C code uses a conditional, but we can't do that in SIMD.
1620 * So we take a max with the tolerance instead. Since we multiply with
1621 * m or n later, the max does not affect the results.
1623 iprm_S = gmx_simd_max_r(iprm_S, toler_S);
1624 iprn_S = gmx_simd_max_r(iprn_S, toler_S);
1625 *nrkj_m2_S = gmx_simd_mul_r(nrkj_S, gmx_simd_inv_r(iprm_S));
1626 *nrkj_n2_S = gmx_simd_mul_r(nrkj_S, gmx_simd_inv_r(iprn_S));
1628 /* Set sign of phi_S with the sign of ipr_S; phi_S is currently positive */
1629 *phi_S = gmx_simd_xor_sign_r(*phi_S, ipr_S);
1630 p_S = gmx_simd_iprod_r(rijx_S, rijy_S, rijz_S,
1631 rkjx_S, rkjy_S, rkjz_S);
1632 p_S = gmx_simd_mul_r(p_S, nrkj_2_S);
1634 q_S = gmx_simd_iprod_r(rklx_S, rkly_S, rklz_S,
1635 rkjx_S, rkjy_S, rkjz_S);
1636 q_S = gmx_simd_mul_r(q_S, nrkj_2_S);
1638 gmx_simd_store_r(p, p_S);
1639 gmx_simd_store_r(q, q_S);
1642 #endif /* GMX_SIMD_HAVE_REAL */
1645 void do_dih_fup(int i, int j, int k, int l, real ddphi,
1646 rvec r_ij, rvec r_kj, rvec r_kl,
1647 rvec m, rvec n, rvec f[], rvec fshift[],
1648 const t_pbc *pbc, const t_graph *g,
1649 const rvec x[], int t1, int t2, int t3)
1652 rvec f_i, f_j, f_k, f_l;
1653 rvec uvec, vvec, svec, dx_jl;
1654 real iprm, iprn, nrkj, nrkj2, nrkj_1, nrkj_2;
1655 real a, b, p, q, toler;
1656 ivec jt, dt_ij, dt_kj, dt_lj;
1658 iprm = iprod(m, m); /* 5 */
1659 iprn = iprod(n, n); /* 5 */
1660 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1661 toler = nrkj2*GMX_REAL_EPS;
1662 if ((iprm > toler) && (iprn > toler))
1664 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1665 nrkj_2 = nrkj_1*nrkj_1; /* 1 */
1666 nrkj = nrkj2*nrkj_1; /* 1 */
1667 a = -ddphi*nrkj/iprm; /* 11 */
1668 svmul(a, m, f_i); /* 3 */
1669 b = ddphi*nrkj/iprn; /* 11 */
1670 svmul(b, n, f_l); /* 3 */
1671 p = iprod(r_ij, r_kj); /* 5 */
1672 p *= nrkj_2; /* 1 */
1673 q = iprod(r_kl, r_kj); /* 5 */
1674 q *= nrkj_2; /* 1 */
1675 svmul(p, f_i, uvec); /* 3 */
1676 svmul(q, f_l, vvec); /* 3 */
1677 rvec_sub(uvec, vvec, svec); /* 3 */
1678 rvec_sub(f_i, svec, f_j); /* 3 */
1679 rvec_add(f_l, svec, f_k); /* 3 */
1680 rvec_inc(f[i], f_i); /* 3 */
1681 rvec_dec(f[j], f_j); /* 3 */
1682 rvec_dec(f[k], f_k); /* 3 */
1683 rvec_inc(f[l], f_l); /* 3 */
1687 copy_ivec(SHIFT_IVEC(g, j), jt);
1688 ivec_sub(SHIFT_IVEC(g, i), jt, dt_ij);
1689 ivec_sub(SHIFT_IVEC(g, k), jt, dt_kj);
1690 ivec_sub(SHIFT_IVEC(g, l), jt, dt_lj);
1691 t1 = IVEC2IS(dt_ij);
1692 t2 = IVEC2IS(dt_kj);
1693 t3 = IVEC2IS(dt_lj);
1697 t3 = pbc_rvec_sub(pbc, x[l], x[j], dx_jl);
1704 rvec_inc(fshift[t1], f_i);
1705 rvec_dec(fshift[CENTRAL], f_j);
1706 rvec_dec(fshift[t2], f_k);
1707 rvec_inc(fshift[t3], f_l);
1712 /* As do_dih_fup above, but without shift forces */
1714 do_dih_fup_noshiftf(int i, int j, int k, int l, real ddphi,
1715 rvec r_ij, rvec r_kj, rvec r_kl,
1716 rvec m, rvec n, rvec f[])
1718 rvec f_i, f_j, f_k, f_l;
1719 rvec uvec, vvec, svec;
1720 real iprm, iprn, nrkj, nrkj2, nrkj_1, nrkj_2;
1721 real a, b, p, q, toler;
1723 iprm = iprod(m, m); /* 5 */
1724 iprn = iprod(n, n); /* 5 */
1725 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1726 toler = nrkj2*GMX_REAL_EPS;
1727 if ((iprm > toler) && (iprn > toler))
1729 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1730 nrkj_2 = nrkj_1*nrkj_1; /* 1 */
1731 nrkj = nrkj2*nrkj_1; /* 1 */
1732 a = -ddphi*nrkj/iprm; /* 11 */
1733 svmul(a, m, f_i); /* 3 */
1734 b = ddphi*nrkj/iprn; /* 11 */
1735 svmul(b, n, f_l); /* 3 */
1736 p = iprod(r_ij, r_kj); /* 5 */
1737 p *= nrkj_2; /* 1 */
1738 q = iprod(r_kl, r_kj); /* 5 */
1739 q *= nrkj_2; /* 1 */
1740 svmul(p, f_i, uvec); /* 3 */
1741 svmul(q, f_l, vvec); /* 3 */
1742 rvec_sub(uvec, vvec, svec); /* 3 */
1743 rvec_sub(f_i, svec, f_j); /* 3 */
1744 rvec_add(f_l, svec, f_k); /* 3 */
1745 rvec_inc(f[i], f_i); /* 3 */
1746 rvec_dec(f[j], f_j); /* 3 */
1747 rvec_dec(f[k], f_k); /* 3 */
1748 rvec_inc(f[l], f_l); /* 3 */
1752 /* As do_dih_fup_noshiftf above, but with pre-calculated pre-factors */
1753 static gmx_inline void
1754 do_dih_fup_noshiftf_precalc(int i, int j, int k, int l,
1756 real f_i_x, real f_i_y, real f_i_z,
1757 real mf_l_x, real mf_l_y, real mf_l_z,
1760 rvec f_i, f_j, f_k, f_l;
1761 rvec uvec, vvec, svec;
1769 svmul(p, f_i, uvec);
1770 svmul(q, f_l, vvec);
1771 rvec_sub(uvec, vvec, svec);
1772 rvec_sub(f_i, svec, f_j);
1773 rvec_add(f_l, svec, f_k);
1774 rvec_inc(f[i], f_i);
1775 rvec_dec(f[j], f_j);
1776 rvec_dec(f[k], f_k);
1777 rvec_inc(f[l], f_l);
1781 real dopdihs(real cpA, real cpB, real phiA, real phiB, int mult,
1782 real phi, real lambda, real *V, real *F)
1784 real v, dvdlambda, mdphi, v1, sdphi, ddphi;
1785 real L1 = 1.0 - lambda;
1786 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1787 real dph0 = (phiB - phiA)*DEG2RAD;
1788 real cp = L1*cpA + lambda*cpB;
1790 mdphi = mult*phi - ph0;
1792 ddphi = -cp*mult*sdphi;
1793 v1 = 1.0 + cos(mdphi);
1796 dvdlambda = (cpB - cpA)*v1 + cp*dph0*sdphi;
1803 /* That was 40 flops */
1807 dopdihs_noener(real cpA, real cpB, real phiA, real phiB, int mult,
1808 real phi, real lambda, real *F)
1810 real mdphi, sdphi, ddphi;
1811 real L1 = 1.0 - lambda;
1812 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1813 real cp = L1*cpA + lambda*cpB;
1815 mdphi = mult*phi - ph0;
1817 ddphi = -cp*mult*sdphi;
1821 /* That was 20 flops */
1825 dopdihs_mdphi(real cpA, real cpB, real phiA, real phiB, int mult,
1826 real phi, real lambda, real *cp, real *mdphi)
1828 real L1 = 1.0 - lambda;
1829 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1831 *cp = L1*cpA + lambda*cpB;
1833 *mdphi = mult*phi - ph0;
1836 static real dopdihs_min(real cpA, real cpB, real phiA, real phiB, int mult,
1837 real phi, real lambda, real *V, real *F)
1838 /* similar to dopdihs, except for a minus sign *
1839 * and a different treatment of mult/phi0 */
1841 real v, dvdlambda, mdphi, v1, sdphi, ddphi;
1842 real L1 = 1.0 - lambda;
1843 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1844 real dph0 = (phiB - phiA)*DEG2RAD;
1845 real cp = L1*cpA + lambda*cpB;
1847 mdphi = mult*(phi-ph0);
1849 ddphi = cp*mult*sdphi;
1850 v1 = 1.0-cos(mdphi);
1853 dvdlambda = (cpB-cpA)*v1 + cp*dph0*sdphi;
1860 /* That was 40 flops */
1863 real pdihs(int nbonds,
1864 const t_iatom forceatoms[], const t_iparams forceparams[],
1865 const rvec x[], rvec f[], rvec fshift[],
1866 const t_pbc *pbc, const t_graph *g,
1867 real lambda, real *dvdlambda,
1868 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1869 int gmx_unused *global_atom_index)
1871 int i, type, ai, aj, ak, al;
1873 rvec r_ij, r_kj, r_kl, m, n;
1874 real phi, sign, ddphi, vpd, vtot;
1878 for (i = 0; (i < nbonds); )
1880 type = forceatoms[i++];
1881 ai = forceatoms[i++];
1882 aj = forceatoms[i++];
1883 ak = forceatoms[i++];
1884 al = forceatoms[i++];
1886 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
1887 &sign, &t1, &t2, &t3); /* 84 */
1888 *dvdlambda += dopdihs(forceparams[type].pdihs.cpA,
1889 forceparams[type].pdihs.cpB,
1890 forceparams[type].pdihs.phiA,
1891 forceparams[type].pdihs.phiB,
1892 forceparams[type].pdihs.mult,
1893 phi, lambda, &vpd, &ddphi);
1896 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
1897 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
1900 fprintf(debug, "pdih: (%d,%d,%d,%d) phi=%g\n",
1901 ai, aj, ak, al, phi);
1908 void make_dp_periodic(real *dp) /* 1 flop? */
1910 /* dp cannot be outside (-pi,pi) */
1915 else if (*dp < -M_PI)
1922 /* As pdihs above, but without calculating energies and shift forces */
1924 pdihs_noener(int nbonds,
1925 const t_iatom forceatoms[], const t_iparams forceparams[],
1926 const rvec x[], rvec f[],
1927 const t_pbc gmx_unused *pbc, const t_graph gmx_unused *g,
1929 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1930 int gmx_unused *global_atom_index)
1932 int i, type, ai, aj, ak, al;
1934 rvec r_ij, r_kj, r_kl, m, n;
1935 real phi, sign, ddphi_tot, ddphi;
1937 for (i = 0; (i < nbonds); )
1939 ai = forceatoms[i+1];
1940 aj = forceatoms[i+2];
1941 ak = forceatoms[i+3];
1942 al = forceatoms[i+4];
1944 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
1945 &sign, &t1, &t2, &t3);
1949 /* Loop over dihedrals working on the same atoms,
1950 * so we avoid recalculating angles and force distributions.
1954 type = forceatoms[i];
1955 dopdihs_noener(forceparams[type].pdihs.cpA,
1956 forceparams[type].pdihs.cpB,
1957 forceparams[type].pdihs.phiA,
1958 forceparams[type].pdihs.phiB,
1959 forceparams[type].pdihs.mult,
1960 phi, lambda, &ddphi);
1965 while (i < nbonds &&
1966 forceatoms[i+1] == ai &&
1967 forceatoms[i+2] == aj &&
1968 forceatoms[i+3] == ak &&
1969 forceatoms[i+4] == al);
1971 do_dih_fup_noshiftf(ai, aj, ak, al, ddphi_tot, r_ij, r_kj, r_kl, m, n, f);
1976 #ifdef GMX_SIMD_HAVE_REAL
1978 /* As pdihs_noner above, but using SIMD to calculate many dihedrals at once */
1980 pdihs_noener_simd(int nbonds,
1981 const t_iatom forceatoms[], const t_iparams forceparams[],
1982 const rvec x[], rvec f[],
1983 const t_pbc *pbc, const t_graph gmx_unused *g,
1984 real gmx_unused lambda,
1985 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1986 int gmx_unused *global_atom_index)
1990 int type, ai[GMX_SIMD_REAL_WIDTH], aj[GMX_SIMD_REAL_WIDTH], ak[GMX_SIMD_REAL_WIDTH], al[GMX_SIMD_REAL_WIDTH];
1991 real dr_array[3*DIM*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *dr;
1992 real buf_array[7*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *buf;
1993 real *cp, *phi0, *mult, *p, *q;
1994 gmx_simd_real_t phi0_S, phi_S;
1995 gmx_simd_real_t mx_S, my_S, mz_S;
1996 gmx_simd_real_t nx_S, ny_S, nz_S;
1997 gmx_simd_real_t nrkj_m2_S, nrkj_n2_S;
1998 gmx_simd_real_t cp_S, mdphi_S, mult_S;
1999 gmx_simd_real_t sin_S, cos_S;
2000 gmx_simd_real_t mddphi_S;
2001 gmx_simd_real_t sf_i_S, msf_l_S;
2002 pbc_simd_t pbc_simd;
2004 /* Ensure SIMD register alignment */
2005 dr = gmx_simd_align_r(dr_array);
2006 buf = gmx_simd_align_r(buf_array);
2008 /* Extract aligned pointer for parameters and variables */
2009 cp = buf + 0*GMX_SIMD_REAL_WIDTH;
2010 phi0 = buf + 1*GMX_SIMD_REAL_WIDTH;
2011 mult = buf + 2*GMX_SIMD_REAL_WIDTH;
2012 p = buf + 3*GMX_SIMD_REAL_WIDTH;
2013 q = buf + 4*GMX_SIMD_REAL_WIDTH;
2015 set_pbc_simd(pbc, &pbc_simd);
2017 /* nbonds is the number of dihedrals times nfa1, here we step GMX_SIMD_REAL_WIDTH dihs */
2018 for (i = 0; (i < nbonds); i += GMX_SIMD_REAL_WIDTH*nfa1)
2020 /* Collect atoms quadruplets for GMX_SIMD_REAL_WIDTH dihedrals.
2021 * iu indexes into forceatoms, we should not let iu go beyond nbonds.
2024 for (s = 0; s < GMX_SIMD_REAL_WIDTH; s++)
2026 type = forceatoms[iu];
2027 ai[s] = forceatoms[iu+1];
2028 aj[s] = forceatoms[iu+2];
2029 ak[s] = forceatoms[iu+3];
2030 al[s] = forceatoms[iu+4];
2032 cp[s] = forceparams[type].pdihs.cpA;
2033 phi0[s] = forceparams[type].pdihs.phiA*DEG2RAD;
2034 mult[s] = forceparams[type].pdihs.mult;
2036 /* At the end fill the arrays with identical entries */
2037 if (iu + nfa1 < nbonds)
2043 /* Caclulate GMX_SIMD_REAL_WIDTH dihedral angles at once */
2044 dih_angle_simd(x, ai, aj, ak, al, &pbc_simd,
2047 &mx_S, &my_S, &mz_S,
2048 &nx_S, &ny_S, &nz_S,
2053 cp_S = gmx_simd_load_r(cp);
2054 phi0_S = gmx_simd_load_r(phi0);
2055 mult_S = gmx_simd_load_r(mult);
2057 mdphi_S = gmx_simd_sub_r(gmx_simd_mul_r(mult_S, phi_S), phi0_S);
2059 /* Calculate GMX_SIMD_REAL_WIDTH sines at once */
2060 gmx_simd_sincos_r(mdphi_S, &sin_S, &cos_S);
2061 mddphi_S = gmx_simd_mul_r(gmx_simd_mul_r(cp_S, mult_S), sin_S);
2062 sf_i_S = gmx_simd_mul_r(mddphi_S, nrkj_m2_S);
2063 msf_l_S = gmx_simd_mul_r(mddphi_S, nrkj_n2_S);
2065 /* After this m?_S will contain f[i] */
2066 mx_S = gmx_simd_mul_r(sf_i_S, mx_S);
2067 my_S = gmx_simd_mul_r(sf_i_S, my_S);
2068 mz_S = gmx_simd_mul_r(sf_i_S, mz_S);
2070 /* After this m?_S will contain -f[l] */
2071 nx_S = gmx_simd_mul_r(msf_l_S, nx_S);
2072 ny_S = gmx_simd_mul_r(msf_l_S, ny_S);
2073 nz_S = gmx_simd_mul_r(msf_l_S, nz_S);
2075 gmx_simd_store_r(dr + 0*GMX_SIMD_REAL_WIDTH, mx_S);
2076 gmx_simd_store_r(dr + 1*GMX_SIMD_REAL_WIDTH, my_S);
2077 gmx_simd_store_r(dr + 2*GMX_SIMD_REAL_WIDTH, mz_S);
2078 gmx_simd_store_r(dr + 3*GMX_SIMD_REAL_WIDTH, nx_S);
2079 gmx_simd_store_r(dr + 4*GMX_SIMD_REAL_WIDTH, ny_S);
2080 gmx_simd_store_r(dr + 5*GMX_SIMD_REAL_WIDTH, nz_S);
2086 do_dih_fup_noshiftf_precalc(ai[s], aj[s], ak[s], al[s],
2088 dr[ XX *GMX_SIMD_REAL_WIDTH+s],
2089 dr[ YY *GMX_SIMD_REAL_WIDTH+s],
2090 dr[ ZZ *GMX_SIMD_REAL_WIDTH+s],
2091 dr[(DIM+XX)*GMX_SIMD_REAL_WIDTH+s],
2092 dr[(DIM+YY)*GMX_SIMD_REAL_WIDTH+s],
2093 dr[(DIM+ZZ)*GMX_SIMD_REAL_WIDTH+s],
2098 while (s < GMX_SIMD_REAL_WIDTH && iu < nbonds);
2102 /* This is mostly a copy of pdihs_noener_simd above, but with using
2103 * the RB potential instead of a harmonic potential.
2104 * This function can replace rbdihs() when no energy and virial are needed.
2107 rbdihs_noener_simd(int nbonds,
2108 const t_iatom forceatoms[], const t_iparams forceparams[],
2109 const rvec x[], rvec f[],
2110 const t_pbc *pbc, const t_graph gmx_unused *g,
2111 real gmx_unused lambda,
2112 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2113 int gmx_unused *global_atom_index)
2117 int type, ai[GMX_SIMD_REAL_WIDTH], aj[GMX_SIMD_REAL_WIDTH], ak[GMX_SIMD_REAL_WIDTH], al[GMX_SIMD_REAL_WIDTH];
2118 real dr_array[3*DIM*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *dr;
2119 real buf_array[(NR_RBDIHS + 4)*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *buf;
2122 gmx_simd_real_t phi_S;
2123 gmx_simd_real_t ddphi_S, cosfac_S;
2124 gmx_simd_real_t mx_S, my_S, mz_S;
2125 gmx_simd_real_t nx_S, ny_S, nz_S;
2126 gmx_simd_real_t nrkj_m2_S, nrkj_n2_S;
2127 gmx_simd_real_t parm_S, c_S;
2128 gmx_simd_real_t sin_S, cos_S;
2129 gmx_simd_real_t sf_i_S, msf_l_S;
2130 pbc_simd_t pbc_simd;
2132 gmx_simd_real_t pi_S = gmx_simd_set1_r(M_PI);
2133 gmx_simd_real_t one_S = gmx_simd_set1_r(1.0);
2135 /* Ensure SIMD register alignment */
2136 dr = gmx_simd_align_r(dr_array);
2137 buf = gmx_simd_align_r(buf_array);
2139 /* Extract aligned pointer for parameters and variables */
2141 p = buf + (NR_RBDIHS + 0)*GMX_SIMD_REAL_WIDTH;
2142 q = buf + (NR_RBDIHS + 1)*GMX_SIMD_REAL_WIDTH;
2144 set_pbc_simd(pbc, &pbc_simd);
2146 /* nbonds is the number of dihedrals times nfa1, here we step GMX_SIMD_REAL_WIDTH dihs */
2147 for (i = 0; (i < nbonds); i += GMX_SIMD_REAL_WIDTH*nfa1)
2149 /* Collect atoms quadruplets for GMX_SIMD_REAL_WIDTH dihedrals.
2150 * iu indexes into forceatoms, we should not let iu go beyond nbonds.
2153 for (s = 0; s < GMX_SIMD_REAL_WIDTH; s++)
2155 type = forceatoms[iu];
2156 ai[s] = forceatoms[iu+1];
2157 aj[s] = forceatoms[iu+2];
2158 ak[s] = forceatoms[iu+3];
2159 al[s] = forceatoms[iu+4];
2161 /* We don't need the first parameter, since that's a constant
2162 * which only affects the energies, not the forces.
2164 for (j = 1; j < NR_RBDIHS; j++)
2166 parm[j*GMX_SIMD_REAL_WIDTH + s] =
2167 forceparams[type].rbdihs.rbcA[j];
2170 /* At the end fill the arrays with identical entries */
2171 if (iu + nfa1 < nbonds)
2177 /* Caclulate GMX_SIMD_REAL_WIDTH dihedral angles at once */
2178 dih_angle_simd(x, ai, aj, ak, al, &pbc_simd,
2181 &mx_S, &my_S, &mz_S,
2182 &nx_S, &ny_S, &nz_S,
2187 /* Change to polymer convention */
2188 phi_S = gmx_simd_sub_r(phi_S, pi_S);
2190 gmx_simd_sincos_r(phi_S, &sin_S, &cos_S);
2192 ddphi_S = gmx_simd_setzero_r();
2195 for (j = 1; j < NR_RBDIHS; j++)
2197 parm_S = gmx_simd_load_r(parm + j*GMX_SIMD_REAL_WIDTH);
2198 ddphi_S = gmx_simd_fmadd_r(gmx_simd_mul_r(c_S, parm_S), cosfac_S, ddphi_S);
2199 cosfac_S = gmx_simd_mul_r(cosfac_S, cos_S);
2200 c_S = gmx_simd_add_r(c_S, one_S);
2203 /* Note that here we do not use the minus sign which is present
2204 * in the normal RB code. This is corrected for through (m)sf below.
2206 ddphi_S = gmx_simd_mul_r(ddphi_S, sin_S);
2208 sf_i_S = gmx_simd_mul_r(ddphi_S, nrkj_m2_S);
2209 msf_l_S = gmx_simd_mul_r(ddphi_S, nrkj_n2_S);
2211 /* After this m?_S will contain f[i] */
2212 mx_S = gmx_simd_mul_r(sf_i_S, mx_S);
2213 my_S = gmx_simd_mul_r(sf_i_S, my_S);
2214 mz_S = gmx_simd_mul_r(sf_i_S, mz_S);
2216 /* After this m?_S will contain -f[l] */
2217 nx_S = gmx_simd_mul_r(msf_l_S, nx_S);
2218 ny_S = gmx_simd_mul_r(msf_l_S, ny_S);
2219 nz_S = gmx_simd_mul_r(msf_l_S, nz_S);
2221 gmx_simd_store_r(dr + 0*GMX_SIMD_REAL_WIDTH, mx_S);
2222 gmx_simd_store_r(dr + 1*GMX_SIMD_REAL_WIDTH, my_S);
2223 gmx_simd_store_r(dr + 2*GMX_SIMD_REAL_WIDTH, mz_S);
2224 gmx_simd_store_r(dr + 3*GMX_SIMD_REAL_WIDTH, nx_S);
2225 gmx_simd_store_r(dr + 4*GMX_SIMD_REAL_WIDTH, ny_S);
2226 gmx_simd_store_r(dr + 5*GMX_SIMD_REAL_WIDTH, nz_S);
2232 do_dih_fup_noshiftf_precalc(ai[s], aj[s], ak[s], al[s],
2234 dr[ XX *GMX_SIMD_REAL_WIDTH+s],
2235 dr[ YY *GMX_SIMD_REAL_WIDTH+s],
2236 dr[ ZZ *GMX_SIMD_REAL_WIDTH+s],
2237 dr[(DIM+XX)*GMX_SIMD_REAL_WIDTH+s],
2238 dr[(DIM+YY)*GMX_SIMD_REAL_WIDTH+s],
2239 dr[(DIM+ZZ)*GMX_SIMD_REAL_WIDTH+s],
2244 while (s < GMX_SIMD_REAL_WIDTH && iu < nbonds);
2248 #endif /* GMX_SIMD_HAVE_REAL */
2251 real idihs(int nbonds,
2252 const t_iatom forceatoms[], const t_iparams forceparams[],
2253 const rvec x[], rvec f[], rvec fshift[],
2254 const t_pbc *pbc, const t_graph *g,
2255 real lambda, real *dvdlambda,
2256 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2257 int gmx_unused *global_atom_index)
2259 int i, type, ai, aj, ak, al;
2261 real phi, phi0, dphi0, ddphi, sign, vtot;
2262 rvec r_ij, r_kj, r_kl, m, n;
2263 real L1, kk, dp, dp2, kA, kB, pA, pB, dvdl_term;
2268 for (i = 0; (i < nbonds); )
2270 type = forceatoms[i++];
2271 ai = forceatoms[i++];
2272 aj = forceatoms[i++];
2273 ak = forceatoms[i++];
2274 al = forceatoms[i++];
2276 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2277 &sign, &t1, &t2, &t3); /* 84 */
2279 /* phi can jump if phi0 is close to Pi/-Pi, which will cause huge
2280 * force changes if we just apply a normal harmonic.
2281 * Instead, we first calculate phi-phi0 and take it modulo (-Pi,Pi).
2282 * This means we will never have the periodicity problem, unless
2283 * the dihedral is Pi away from phiO, which is very unlikely due to
2286 kA = forceparams[type].harmonic.krA;
2287 kB = forceparams[type].harmonic.krB;
2288 pA = forceparams[type].harmonic.rA;
2289 pB = forceparams[type].harmonic.rB;
2291 kk = L1*kA + lambda*kB;
2292 phi0 = (L1*pA + lambda*pB)*DEG2RAD;
2293 dphi0 = (pB - pA)*DEG2RAD;
2297 make_dp_periodic(&dp);
2304 dvdl_term += 0.5*(kB - kA)*dp2 - kk*dphi0*dp;
2306 do_dih_fup(ai, aj, ak, al, -ddphi, r_ij, r_kj, r_kl, m, n,
2307 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2312 fprintf(debug, "idih: (%d,%d,%d,%d) phi=%g\n",
2313 ai, aj, ak, al, phi);
2318 *dvdlambda += dvdl_term;
2323 /*! \brief returns dx, rdist, and dpdl for functions posres() and fbposres()
2325 static void posres_dx(const rvec x, const rvec pos0A, const rvec pos0B,
2326 const rvec comA_sc, const rvec comB_sc,
2328 t_pbc *pbc, int refcoord_scaling, int npbcdim,
2329 rvec dx, rvec rdist, rvec dpdl)
2332 real posA, posB, L1, ref = 0.;
2337 for (m = 0; m < DIM; m++)
2343 switch (refcoord_scaling)
2347 rdist[m] = L1*posA + lambda*posB;
2348 dpdl[m] = posB - posA;
2351 /* Box relative coordinates are stored for dimensions with pbc */
2352 posA *= pbc->box[m][m];
2353 posB *= pbc->box[m][m];
2354 assert(npbcdim <= DIM);
2355 for (d = m+1; d < npbcdim; d++)
2357 posA += pos0A[d]*pbc->box[d][m];
2358 posB += pos0B[d]*pbc->box[d][m];
2360 ref = L1*posA + lambda*posB;
2362 dpdl[m] = posB - posA;
2365 ref = L1*comA_sc[m] + lambda*comB_sc[m];
2366 rdist[m] = L1*posA + lambda*posB;
2367 dpdl[m] = comB_sc[m] - comA_sc[m] + posB - posA;
2370 gmx_fatal(FARGS, "No such scaling method implemented");
2375 ref = L1*posA + lambda*posB;
2377 dpdl[m] = posB - posA;
2380 /* We do pbc_dx with ref+rdist,
2381 * since with only ref we can be up to half a box vector wrong.
2383 pos[m] = ref + rdist[m];
2388 pbc_dx(pbc, x, pos, dx);
2392 rvec_sub(x, pos, dx);
2396 /*! \brief Computes forces and potential for flat-bottom cylindrical restraints.
2397 * Returns the flat-bottom potential. */
2398 static real do_fbposres_cylinder(int fbdim, rvec fm, rvec dx, real rfb, real kk, gmx_bool bInvert)
2401 real dr, dr2, invdr, v, rfb2;
2407 for (d = 0; d < DIM; d++)
2416 ( (dr2 > rfb2 && bInvert == FALSE ) || (dr2 < rfb2 && bInvert == TRUE ) )
2421 v = 0.5*kk*sqr(dr - rfb);
2422 for (d = 0; d < DIM; d++)
2426 fm[d] = -kk*(dr-rfb)*dx[d]*invdr; /* Force pointing to the center */
2434 /*! \brief Adds forces of flat-bottomed positions restraints to f[]
2435 * and fixes vir_diag. Returns the flat-bottomed potential. */
2436 real fbposres(int nbonds,
2437 const t_iatom forceatoms[], const t_iparams forceparams[],
2438 const rvec x[], rvec f[], rvec vir_diag,
2440 int refcoord_scaling, int ePBC, rvec com)
2441 /* compute flat-bottomed positions restraints */
2443 int i, ai, m, d, type, npbcdim = 0, fbdim;
2444 const t_iparams *pr;
2446 real dr, dr2, rfb, rfb2, fact;
2447 rvec com_sc, rdist, dx, dpdl, fm;
2450 npbcdim = ePBC2npbcdim(ePBC);
2452 if (refcoord_scaling == erscCOM)
2455 for (m = 0; m < npbcdim; m++)
2457 assert(npbcdim <= DIM);
2458 for (d = m; d < npbcdim; d++)
2460 com_sc[m] += com[d]*pbc->box[d][m];
2466 for (i = 0; (i < nbonds); )
2468 type = forceatoms[i++];
2469 ai = forceatoms[i++];
2470 pr = &forceparams[type];
2472 /* same calculation as for normal posres, but with identical A and B states, and lambda==0 */
2473 posres_dx(x[ai], forceparams[type].fbposres.pos0, forceparams[type].fbposres.pos0,
2474 com_sc, com_sc, 0.0,
2475 pbc, refcoord_scaling, npbcdim,
2481 kk = pr->fbposres.k;
2482 rfb = pr->fbposres.r;
2485 /* with rfb<0, push particle out of the sphere/cylinder/layer */
2493 switch (pr->fbposres.geom)
2495 case efbposresSPHERE:
2496 /* spherical flat-bottom posres */
2499 ( (dr2 > rfb2 && bInvert == FALSE ) || (dr2 < rfb2 && bInvert == TRUE ) )
2503 v = 0.5*kk*sqr(dr - rfb);
2504 fact = -kk*(dr-rfb)/dr; /* Force pointing to the center pos0 */
2505 svmul(fact, dx, fm);
2508 case efbposresCYLINDERX:
2509 /* cylindrical flat-bottom posres in y-z plane. fm[XX] = 0. */
2511 v = do_fbposres_cylinder(fbdim, fm, dx, rfb, kk, bInvert);
2513 case efbposresCYLINDERY:
2514 /* cylindrical flat-bottom posres in x-z plane. fm[YY] = 0. */
2516 v = do_fbposres_cylinder(fbdim, fm, dx, rfb, kk, bInvert);
2518 case efbposresCYLINDER:
2519 /* equivalent to efbposresCYLINDERZ for backwards compatibility */
2520 case efbposresCYLINDERZ:
2521 /* cylindrical flat-bottom posres in x-y plane. fm[ZZ] = 0. */
2523 v = do_fbposres_cylinder(fbdim, fm, dx, rfb, kk, bInvert);
2525 case efbposresX: /* fbdim=XX */
2526 case efbposresY: /* fbdim=YY */
2527 case efbposresZ: /* fbdim=ZZ */
2528 /* 1D flat-bottom potential */
2529 fbdim = pr->fbposres.geom - efbposresX;
2531 if ( ( dr > rfb && bInvert == FALSE ) || ( 0 < dr && dr < rfb && bInvert == TRUE ) )
2533 v = 0.5*kk*sqr(dr - rfb);
2534 fm[fbdim] = -kk*(dr - rfb);
2536 else if ( (dr < (-rfb) && bInvert == FALSE ) || ( (-rfb) < dr && dr < 0 && bInvert == TRUE ))
2538 v = 0.5*kk*sqr(dr + rfb);
2539 fm[fbdim] = -kk*(dr + rfb);
2546 for (m = 0; (m < DIM); m++)
2549 /* Here we correct for the pbc_dx which included rdist */
2550 vir_diag[m] -= 0.5*(dx[m] + rdist[m])*fm[m];
2558 real posres(int nbonds,
2559 const t_iatom forceatoms[], const t_iparams forceparams[],
2560 const rvec x[], rvec f[], rvec vir_diag,
2562 real lambda, real *dvdlambda,
2563 int refcoord_scaling, int ePBC, rvec comA, rvec comB)
2565 int i, ai, m, d, type, npbcdim = 0;
2566 const t_iparams *pr;
2569 rvec comA_sc, comB_sc, rdist, dpdl, dx;
2570 gmx_bool bForceValid = TRUE;
2572 if ((f == NULL) || (vir_diag == NULL)) /* should both be null together! */
2574 bForceValid = FALSE;
2577 npbcdim = ePBC2npbcdim(ePBC);
2579 if (refcoord_scaling == erscCOM)
2581 clear_rvec(comA_sc);
2582 clear_rvec(comB_sc);
2583 for (m = 0; m < npbcdim; m++)
2585 assert(npbcdim <= DIM);
2586 for (d = m; d < npbcdim; d++)
2588 comA_sc[m] += comA[d]*pbc->box[d][m];
2589 comB_sc[m] += comB[d]*pbc->box[d][m];
2597 for (i = 0; (i < nbonds); )
2599 type = forceatoms[i++];
2600 ai = forceatoms[i++];
2601 pr = &forceparams[type];
2603 /* return dx, rdist, and dpdl */
2604 posres_dx(x[ai], forceparams[type].posres.pos0A, forceparams[type].posres.pos0B,
2605 comA_sc, comB_sc, lambda,
2606 pbc, refcoord_scaling, npbcdim,
2609 for (m = 0; (m < DIM); m++)
2611 kk = L1*pr->posres.fcA[m] + lambda*pr->posres.fcB[m];
2613 vtot += 0.5*kk*dx[m]*dx[m];
2615 0.5*(pr->posres.fcB[m] - pr->posres.fcA[m])*dx[m]*dx[m]
2618 /* Here we correct for the pbc_dx which included rdist */
2622 vir_diag[m] -= 0.5*(dx[m] + rdist[m])*fm;
2630 static real low_angres(int nbonds,
2631 const t_iatom forceatoms[], const t_iparams forceparams[],
2632 const rvec x[], rvec f[], rvec fshift[],
2633 const t_pbc *pbc, const t_graph *g,
2634 real lambda, real *dvdlambda,
2637 int i, m, type, ai, aj, ak, al;
2639 real phi, cos_phi, cos_phi2, vid, vtot, dVdphi;
2640 rvec r_ij, r_kl, f_i, f_k = {0, 0, 0};
2641 real st, sth, nrij2, nrkl2, c, cij, ckl;
2644 t2 = 0; /* avoid warning with gcc-3.3. It is never used uninitialized */
2647 ak = al = 0; /* to avoid warnings */
2648 for (i = 0; i < nbonds; )
2650 type = forceatoms[i++];
2651 ai = forceatoms[i++];
2652 aj = forceatoms[i++];
2653 t1 = pbc_rvec_sub(pbc, x[aj], x[ai], r_ij); /* 3 */
2656 ak = forceatoms[i++];
2657 al = forceatoms[i++];
2658 t2 = pbc_rvec_sub(pbc, x[al], x[ak], r_kl); /* 3 */
2667 cos_phi = cos_angle(r_ij, r_kl); /* 25 */
2668 phi = acos(cos_phi); /* 10 */
2670 *dvdlambda += dopdihs_min(forceparams[type].pdihs.cpA,
2671 forceparams[type].pdihs.cpB,
2672 forceparams[type].pdihs.phiA,
2673 forceparams[type].pdihs.phiB,
2674 forceparams[type].pdihs.mult,
2675 phi, lambda, &vid, &dVdphi); /* 40 */
2679 cos_phi2 = sqr(cos_phi); /* 1 */
2682 st = -dVdphi*gmx_invsqrt(1 - cos_phi2); /* 12 */
2683 sth = st*cos_phi; /* 1 */
2684 nrij2 = iprod(r_ij, r_ij); /* 5 */
2685 nrkl2 = iprod(r_kl, r_kl); /* 5 */
2687 c = st*gmx_invsqrt(nrij2*nrkl2); /* 11 */
2688 cij = sth/nrij2; /* 10 */
2689 ckl = sth/nrkl2; /* 10 */
2691 for (m = 0; m < DIM; m++) /* 18+18 */
2693 f_i[m] = (c*r_kl[m]-cij*r_ij[m]);
2698 f_k[m] = (c*r_ij[m]-ckl*r_kl[m]);
2706 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
2709 rvec_inc(fshift[t1], f_i);
2710 rvec_dec(fshift[CENTRAL], f_i);
2715 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, al), dt);
2718 rvec_inc(fshift[t2], f_k);
2719 rvec_dec(fshift[CENTRAL], f_k);
2724 return vtot; /* 184 / 157 (bZAxis) total */
2727 real angres(int nbonds,
2728 const t_iatom forceatoms[], const t_iparams forceparams[],
2729 const rvec x[], rvec f[], rvec fshift[],
2730 const t_pbc *pbc, const t_graph *g,
2731 real lambda, real *dvdlambda,
2732 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2733 int gmx_unused *global_atom_index)
2735 return low_angres(nbonds, forceatoms, forceparams, x, f, fshift, pbc, g,
2736 lambda, dvdlambda, FALSE);
2739 real angresz(int nbonds,
2740 const t_iatom forceatoms[], const t_iparams forceparams[],
2741 const rvec x[], rvec f[], rvec fshift[],
2742 const t_pbc *pbc, const t_graph *g,
2743 real lambda, real *dvdlambda,
2744 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2745 int gmx_unused *global_atom_index)
2747 return low_angres(nbonds, forceatoms, forceparams, x, f, fshift, pbc, g,
2748 lambda, dvdlambda, TRUE);
2751 real dihres(int nbonds,
2752 const t_iatom forceatoms[], const t_iparams forceparams[],
2753 const rvec x[], rvec f[], rvec fshift[],
2754 const t_pbc *pbc, const t_graph *g,
2755 real lambda, real *dvdlambda,
2756 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2757 int gmx_unused *global_atom_index)
2760 int ai, aj, ak, al, i, k, type, t1, t2, t3;
2761 real phi0A, phi0B, dphiA, dphiB, kfacA, kfacB, phi0, dphi, kfac;
2762 real phi, ddphi, ddp, ddp2, dp, sign, d2r, L1;
2763 rvec r_ij, r_kj, r_kl, m, n;
2770 for (i = 0; (i < nbonds); )
2772 type = forceatoms[i++];
2773 ai = forceatoms[i++];
2774 aj = forceatoms[i++];
2775 ak = forceatoms[i++];
2776 al = forceatoms[i++];
2778 phi0A = forceparams[type].dihres.phiA*d2r;
2779 dphiA = forceparams[type].dihres.dphiA*d2r;
2780 kfacA = forceparams[type].dihres.kfacA;
2782 phi0B = forceparams[type].dihres.phiB*d2r;
2783 dphiB = forceparams[type].dihres.dphiB*d2r;
2784 kfacB = forceparams[type].dihres.kfacB;
2786 phi0 = L1*phi0A + lambda*phi0B;
2787 dphi = L1*dphiA + lambda*dphiB;
2788 kfac = L1*kfacA + lambda*kfacB;
2790 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2791 &sign, &t1, &t2, &t3);
2796 fprintf(debug, "dihres[%d]: %d %d %d %d : phi=%f, dphi=%f, kfac=%f\n",
2797 k++, ai, aj, ak, al, phi0, dphi, kfac);
2799 /* phi can jump if phi0 is close to Pi/-Pi, which will cause huge
2800 * force changes if we just apply a normal harmonic.
2801 * Instead, we first calculate phi-phi0 and take it modulo (-Pi,Pi).
2802 * This means we will never have the periodicity problem, unless
2803 * the dihedral is Pi away from phiO, which is very unlikely due to
2807 make_dp_periodic(&dp);
2813 else if (dp < -dphi)
2825 vtot += 0.5*kfac*ddp2;
2828 *dvdlambda += 0.5*(kfacB - kfacA)*ddp2;
2829 /* lambda dependence from changing restraint distances */
2832 *dvdlambda -= kfac*ddp*((dphiB - dphiA)+(phi0B - phi0A));
2836 *dvdlambda += kfac*ddp*((dphiB - dphiA)-(phi0B - phi0A));
2838 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
2839 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2846 real unimplemented(int gmx_unused nbonds,
2847 const t_iatom gmx_unused forceatoms[], const t_iparams gmx_unused forceparams[],
2848 const rvec gmx_unused x[], rvec gmx_unused f[], rvec gmx_unused fshift[],
2849 const t_pbc gmx_unused *pbc, const t_graph gmx_unused *g,
2850 real gmx_unused lambda, real gmx_unused *dvdlambda,
2851 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2852 int gmx_unused *global_atom_index)
2854 gmx_impl("*** you are using a not implemented function");
2856 return 0.0; /* To make the compiler happy */
2859 real restrangles(int nbonds,
2860 const t_iatom forceatoms[], const t_iparams forceparams[],
2861 const rvec x[], rvec f[], rvec fshift[],
2862 const t_pbc *pbc, const t_graph *g,
2863 real gmx_unused lambda, real gmx_unused *dvdlambda,
2864 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2865 int gmx_unused *global_atom_index)
2867 int i, d, ai, aj, ak, type, m;
2870 ivec jt, dt_ij, dt_kj;
2872 real prefactor, ratio_ante, ratio_post;
2873 rvec delta_ante, delta_post, vec_temp;
2876 for (i = 0; (i < nbonds); )
2878 type = forceatoms[i++];
2879 ai = forceatoms[i++];
2880 aj = forceatoms[i++];
2881 ak = forceatoms[i++];
2883 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], vec_temp);
2884 pbc_rvec_sub(pbc, x[aj], x[ai], delta_ante);
2885 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], delta_post);
2888 /* This function computes factors needed for restricted angle potential.
2889 * The restricted angle potential is used in coarse-grained simulations to avoid singularities
2890 * when three particles align and the dihedral angle and dihedral potential
2891 * cannot be calculated. This potential is calculated using the formula:
2892 real restrangles(int nbonds,
2893 const t_iatom forceatoms[],const t_iparams forceparams[],
2894 const rvec x[],rvec f[],rvec fshift[],
2895 const t_pbc *pbc,const t_graph *g,
2896 real gmx_unused lambda,real gmx_unused *dvdlambda,
2897 const t_mdatoms gmx_unused *md,t_fcdata gmx_unused *fcd,
2898 int gmx_unused *global_atom_index)
2900 int i, d, ai, aj, ak, type, m;
2904 ivec jt,dt_ij,dt_kj;
2906 real prefactor, ratio_ante, ratio_post;
2907 rvec delta_ante, delta_post, vec_temp;
2910 for(i=0; (i<nbonds); )
2912 type = forceatoms[i++];
2913 ai = forceatoms[i++];
2914 aj = forceatoms[i++];
2915 ak = forceatoms[i++];
2917 * \f[V_{\rm ReB}(\theta_i) = \frac{1}{2} k_{\theta} \frac{(\cos\theta_i - \cos\theta_0)^2}
2918 * {\sin^2\theta_i}\f] ({eq:ReB} and ref \cite{MonicaGoga2013} from the manual).
2919 * For more explanations see comments file "restcbt.h". */
2921 compute_factors_restangles(type, forceparams, delta_ante, delta_post,
2922 &prefactor, &ratio_ante, &ratio_post, &v);
2924 /* Forces are computed per component */
2925 for (d = 0; d < DIM; d++)
2927 f_i[d] = prefactor * (ratio_ante * delta_ante[d] - delta_post[d]);
2928 f_j[d] = prefactor * ((ratio_post + 1.0) * delta_post[d] - (ratio_ante + 1.0) * delta_ante[d]);
2929 f_k[d] = prefactor * (delta_ante[d] - ratio_post * delta_post[d]);
2932 /* Computation of potential energy */
2938 for (m = 0; (m < DIM); m++)
2947 copy_ivec(SHIFT_IVEC(g, aj), jt);
2948 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
2949 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
2950 t1 = IVEC2IS(dt_ij);
2951 t2 = IVEC2IS(dt_kj);
2954 rvec_inc(fshift[t1], f_i);
2955 rvec_inc(fshift[CENTRAL], f_j);
2956 rvec_inc(fshift[t2], f_k);
2962 real restrdihs(int nbonds,
2963 const t_iatom forceatoms[], const t_iparams forceparams[],
2964 const rvec x[], rvec f[], rvec fshift[],
2965 const t_pbc *pbc, const t_graph *g,
2966 real gmx_unused lambda, real gmx_unused *dvlambda,
2967 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2968 int gmx_unused *global_atom_index)
2970 int i, d, type, ai, aj, ak, al;
2971 rvec f_i, f_j, f_k, f_l;
2973 ivec jt, dt_ij, dt_kj, dt_lj;
2976 rvec delta_ante, delta_crnt, delta_post, vec_temp;
2977 real factor_phi_ai_ante, factor_phi_ai_crnt, factor_phi_ai_post;
2978 real factor_phi_aj_ante, factor_phi_aj_crnt, factor_phi_aj_post;
2979 real factor_phi_ak_ante, factor_phi_ak_crnt, factor_phi_ak_post;
2980 real factor_phi_al_ante, factor_phi_al_crnt, factor_phi_al_post;
2985 for (i = 0; (i < nbonds); )
2987 type = forceatoms[i++];
2988 ai = forceatoms[i++];
2989 aj = forceatoms[i++];
2990 ak = forceatoms[i++];
2991 al = forceatoms[i++];
2993 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], vec_temp);
2994 pbc_rvec_sub(pbc, x[aj], x[ai], delta_ante);
2995 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], delta_crnt);
2996 pbc_rvec_sub(pbc, x[ak], x[al], vec_temp);
2997 pbc_rvec_sub(pbc, x[al], x[ak], delta_post);
2999 /* This function computes factors needed for restricted angle potential.
3000 * The restricted angle potential is used in coarse-grained simulations to avoid singularities
3001 * when three particles align and the dihedral angle and dihedral potential cannot be calculated.
3002 * This potential is calculated using the formula:
3003 * \f[V_{\rm ReB}(\theta_i) = \frac{1}{2} k_{\theta}
3004 * \frac{(\cos\theta_i - \cos\theta_0)^2}{\sin^2\theta_i}\f]
3005 * ({eq:ReB} and ref \cite{MonicaGoga2013} from the manual).
3006 * For more explanations see comments file "restcbt.h" */
3008 compute_factors_restrdihs(type, forceparams,
3009 delta_ante, delta_crnt, delta_post,
3010 &factor_phi_ai_ante, &factor_phi_ai_crnt, &factor_phi_ai_post,
3011 &factor_phi_aj_ante, &factor_phi_aj_crnt, &factor_phi_aj_post,
3012 &factor_phi_ak_ante, &factor_phi_ak_crnt, &factor_phi_ak_post,
3013 &factor_phi_al_ante, &factor_phi_al_crnt, &factor_phi_al_post,
3014 &prefactor_phi, &v);
3017 /* Computation of forces per component */
3018 for (d = 0; d < DIM; d++)
3020 f_i[d] = prefactor_phi * (factor_phi_ai_ante * delta_ante[d] + factor_phi_ai_crnt * delta_crnt[d] + factor_phi_ai_post * delta_post[d]);
3021 f_j[d] = prefactor_phi * (factor_phi_aj_ante * delta_ante[d] + factor_phi_aj_crnt * delta_crnt[d] + factor_phi_aj_post * delta_post[d]);
3022 f_k[d] = prefactor_phi * (factor_phi_ak_ante * delta_ante[d] + factor_phi_ak_crnt * delta_crnt[d] + factor_phi_ak_post * delta_post[d]);
3023 f_l[d] = prefactor_phi * (factor_phi_al_ante * delta_ante[d] + factor_phi_al_crnt * delta_crnt[d] + factor_phi_al_post * delta_post[d]);
3025 /* Computation of the energy */
3031 /* Updating the forces */
3033 rvec_inc(f[ai], f_i);
3034 rvec_inc(f[aj], f_j);
3035 rvec_inc(f[ak], f_k);
3036 rvec_inc(f[al], f_l);
3039 /* Updating the fshift forces for the pressure coupling */
3042 copy_ivec(SHIFT_IVEC(g, aj), jt);
3043 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3044 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3045 ivec_sub(SHIFT_IVEC(g, al), jt, dt_lj);
3046 t1 = IVEC2IS(dt_ij);
3047 t2 = IVEC2IS(dt_kj);
3048 t3 = IVEC2IS(dt_lj);
3052 t3 = pbc_rvec_sub(pbc, x[al], x[aj], dx_jl);
3059 rvec_inc(fshift[t1], f_i);
3060 rvec_inc(fshift[CENTRAL], f_j);
3061 rvec_inc(fshift[t2], f_k);
3062 rvec_inc(fshift[t3], f_l);
3070 real cbtdihs(int nbonds,
3071 const t_iatom forceatoms[], const t_iparams forceparams[],
3072 const rvec x[], rvec f[], rvec fshift[],
3073 const t_pbc *pbc, const t_graph *g,
3074 real gmx_unused lambda, real gmx_unused *dvdlambda,
3075 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3076 int gmx_unused *global_atom_index)
3078 int type, ai, aj, ak, al, i, d;
3082 rvec f_i, f_j, f_k, f_l;
3083 ivec jt, dt_ij, dt_kj, dt_lj;
3085 rvec delta_ante, delta_crnt, delta_post;
3086 rvec f_phi_ai, f_phi_aj, f_phi_ak, f_phi_al;
3087 rvec f_theta_ante_ai, f_theta_ante_aj, f_theta_ante_ak;
3088 rvec f_theta_post_aj, f_theta_post_ak, f_theta_post_al;
3094 for (i = 0; (i < nbonds); )
3096 type = forceatoms[i++];
3097 ai = forceatoms[i++];
3098 aj = forceatoms[i++];
3099 ak = forceatoms[i++];
3100 al = forceatoms[i++];
3103 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], vec_temp);
3104 pbc_rvec_sub(pbc, x[aj], x[ai], delta_ante);
3105 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], vec_temp);
3106 pbc_rvec_sub(pbc, x[ak], x[aj], delta_crnt);
3107 pbc_rvec_sub(pbc, x[ak], x[al], vec_temp);
3108 pbc_rvec_sub(pbc, x[al], x[ak], delta_post);
3110 /* \brief Compute factors for CBT potential
3111 * The combined bending-torsion potential goes to zero in a very smooth manner, eliminating the numerical
3112 * instabilities, when three coarse-grained particles align and the dihedral angle and standard
3113 * dihedral potentials cannot be calculated. The CBT potential is calculated using the formula:
3114 * \f[V_{\rm CBT}(\theta_{i-1}, \theta_i, \phi_i) = k_{\phi} \sin^3\theta_{i-1} \sin^3\theta_{i}
3115 * \sum_{n=0}^4 { a_n \cos^n\phi_i}\f] ({eq:CBT} and ref \cite{MonicaGoga2013} from the manual).
3116 * It contains in its expression not only the dihedral angle \f$\phi\f$
3117 * but also \f[\theta_{i-1}\f] (theta_ante bellow) and \f[\theta_{i}\f] (theta_post bellow)
3118 * --- the adjacent bending angles.
3119 * For more explanations see comments file "restcbt.h". */
3121 compute_factors_cbtdihs(type, forceparams, delta_ante, delta_crnt, delta_post,
3122 f_phi_ai, f_phi_aj, f_phi_ak, f_phi_al,
3123 f_theta_ante_ai, f_theta_ante_aj, f_theta_ante_ak,
3124 f_theta_post_aj, f_theta_post_ak, f_theta_post_al,
3128 /* Acumulate the resuts per beads */
3129 for (d = 0; d < DIM; d++)
3131 f_i[d] = f_phi_ai[d] + f_theta_ante_ai[d];
3132 f_j[d] = f_phi_aj[d] + f_theta_ante_aj[d] + f_theta_post_aj[d];
3133 f_k[d] = f_phi_ak[d] + f_theta_ante_ak[d] + f_theta_post_ak[d];
3134 f_l[d] = f_phi_al[d] + f_theta_post_al[d];
3137 /* Compute the potential energy */
3142 /* Updating the forces */
3143 rvec_inc(f[ai], f_i);
3144 rvec_inc(f[aj], f_j);
3145 rvec_inc(f[ak], f_k);
3146 rvec_inc(f[al], f_l);
3149 /* Updating the fshift forces for the pressure coupling */
3152 copy_ivec(SHIFT_IVEC(g, aj), jt);
3153 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3154 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3155 ivec_sub(SHIFT_IVEC(g, al), jt, dt_lj);
3156 t1 = IVEC2IS(dt_ij);
3157 t2 = IVEC2IS(dt_kj);
3158 t3 = IVEC2IS(dt_lj);
3162 t3 = pbc_rvec_sub(pbc, x[al], x[aj], dx_jl);
3169 rvec_inc(fshift[t1], f_i);
3170 rvec_inc(fshift[CENTRAL], f_j);
3171 rvec_inc(fshift[t2], f_k);
3172 rvec_inc(fshift[t3], f_l);
3178 real rbdihs(int nbonds,
3179 const t_iatom forceatoms[], const t_iparams forceparams[],
3180 const rvec x[], rvec f[], rvec fshift[],
3181 const t_pbc *pbc, const t_graph *g,
3182 real lambda, real *dvdlambda,
3183 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3184 int gmx_unused *global_atom_index)
3186 const real c0 = 0.0, c1 = 1.0, c2 = 2.0, c3 = 3.0, c4 = 4.0, c5 = 5.0;
3187 int type, ai, aj, ak, al, i, j;
3189 rvec r_ij, r_kj, r_kl, m, n;
3190 real parmA[NR_RBDIHS];
3191 real parmB[NR_RBDIHS];
3192 real parm[NR_RBDIHS];
3193 real cos_phi, phi, rbp, rbpBA;
3194 real v, sign, ddphi, sin_phi;
3196 real L1 = 1.0-lambda;
3200 for (i = 0; (i < nbonds); )
3202 type = forceatoms[i++];
3203 ai = forceatoms[i++];
3204 aj = forceatoms[i++];
3205 ak = forceatoms[i++];
3206 al = forceatoms[i++];
3208 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
3209 &sign, &t1, &t2, &t3); /* 84 */
3211 /* Change to polymer convention */
3218 phi -= M_PI; /* 1 */
3222 /* Beware of accuracy loss, cannot use 1-sqrt(cos^2) ! */
3225 for (j = 0; (j < NR_RBDIHS); j++)
3227 parmA[j] = forceparams[type].rbdihs.rbcA[j];
3228 parmB[j] = forceparams[type].rbdihs.rbcB[j];
3229 parm[j] = L1*parmA[j]+lambda*parmB[j];
3231 /* Calculate cosine powers */
3232 /* Calculate the energy */
3233 /* Calculate the derivative */
3236 dvdl_term += (parmB[0]-parmA[0]);
3241 rbpBA = parmB[1]-parmA[1];
3242 ddphi += rbp*cosfac;
3245 dvdl_term += cosfac*rbpBA;
3247 rbpBA = parmB[2]-parmA[2];
3248 ddphi += c2*rbp*cosfac;
3251 dvdl_term += cosfac*rbpBA;
3253 rbpBA = parmB[3]-parmA[3];
3254 ddphi += c3*rbp*cosfac;
3257 dvdl_term += cosfac*rbpBA;
3259 rbpBA = parmB[4]-parmA[4];
3260 ddphi += c4*rbp*cosfac;
3263 dvdl_term += cosfac*rbpBA;
3265 rbpBA = parmB[5]-parmA[5];
3266 ddphi += c5*rbp*cosfac;
3269 dvdl_term += cosfac*rbpBA;
3271 ddphi = -ddphi*sin_phi; /* 11 */
3273 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
3274 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
3277 *dvdlambda += dvdl_term;
3284 /*! \brief Mysterious undocumented function */
3286 cmap_setup_grid_index(int ip, int grid_spacing, int *ipm1, int *ipp1, int *ipp2)
3292 ip = ip + grid_spacing - 1;
3294 else if (ip > grid_spacing)
3296 ip = ip - grid_spacing - 1;
3305 im1 = grid_spacing - 1;
3307 else if (ip == grid_spacing-2)
3311 else if (ip == grid_spacing-1)
3325 /*! \brief Compute CMAP dihedral energies and forces */
3327 cmap_dihs(int nbonds,
3328 const t_iatom forceatoms[], const t_iparams forceparams[],
3329 const gmx_cmap_t *cmap_grid,
3330 const rvec x[], rvec f[], rvec fshift[],
3331 const t_pbc *pbc, const t_graph *g,
3332 real gmx_unused lambda, real gmx_unused *dvdlambda,
3333 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3334 int gmx_unused *global_atom_index)
3336 int i, j, k, n, idx;
3337 int ai, aj, ak, al, am;
3338 int a1i, a1j, a1k, a1l, a2i, a2j, a2k, a2l;
3340 int t11, t21, t31, t12, t22, t32;
3341 int iphi1, ip1m1, ip1p1, ip1p2;
3342 int iphi2, ip2m1, ip2p1, ip2p2;
3344 int pos1, pos2, pos3, pos4;
3346 real ty[4], ty1[4], ty2[4], ty12[4], tc[16], tx[16];
3347 real phi1, cos_phi1, sin_phi1, sign1, xphi1;
3348 real phi2, cos_phi2, sin_phi2, sign2, xphi2;
3349 real dx, xx, tt, tu, e, df1, df2, vtot;
3350 real ra21, rb21, rg21, rg1, rgr1, ra2r1, rb2r1, rabr1;
3351 real ra22, rb22, rg22, rg2, rgr2, ra2r2, rb2r2, rabr2;
3352 real fg1, hg1, fga1, hgb1, gaa1, gbb1;
3353 real fg2, hg2, fga2, hgb2, gaa2, gbb2;
3356 rvec r1_ij, r1_kj, r1_kl, m1, n1;
3357 rvec r2_ij, r2_kj, r2_kl, m2, n2;
3358 rvec f1_i, f1_j, f1_k, f1_l;
3359 rvec f2_i, f2_j, f2_k, f2_l;
3360 rvec a1, b1, a2, b2;
3361 rvec f1, g1, h1, f2, g2, h2;
3362 rvec dtf1, dtg1, dth1, dtf2, dtg2, dth2;
3363 ivec jt1, dt1_ij, dt1_kj, dt1_lj;
3364 ivec jt2, dt2_ij, dt2_kj, dt2_lj;
3368 int loop_index[4][4] = {
3375 /* Total CMAP energy */
3378 for (n = 0; n < nbonds; )
3380 /* Five atoms are involved in the two torsions */
3381 type = forceatoms[n++];
3382 ai = forceatoms[n++];
3383 aj = forceatoms[n++];
3384 ak = forceatoms[n++];
3385 al = forceatoms[n++];
3386 am = forceatoms[n++];
3388 /* Which CMAP type is this */
3389 cmapA = forceparams[type].cmap.cmapA;
3390 cmapd = cmap_grid->cmapdata[cmapA].cmap;
3398 phi1 = dih_angle(x[a1i], x[a1j], x[a1k], x[a1l], pbc, r1_ij, r1_kj, r1_kl, m1, n1,
3399 &sign1, &t11, &t21, &t31); /* 84 */
3401 cos_phi1 = cos(phi1);
3403 a1[0] = r1_ij[1]*r1_kj[2]-r1_ij[2]*r1_kj[1];
3404 a1[1] = r1_ij[2]*r1_kj[0]-r1_ij[0]*r1_kj[2];
3405 a1[2] = r1_ij[0]*r1_kj[1]-r1_ij[1]*r1_kj[0]; /* 9 */
3407 b1[0] = r1_kl[1]*r1_kj[2]-r1_kl[2]*r1_kj[1];
3408 b1[1] = r1_kl[2]*r1_kj[0]-r1_kl[0]*r1_kj[2];
3409 b1[2] = r1_kl[0]*r1_kj[1]-r1_kl[1]*r1_kj[0]; /* 9 */
3411 pbc_rvec_sub(pbc, x[a1l], x[a1k], h1);
3413 ra21 = iprod(a1, a1); /* 5 */
3414 rb21 = iprod(b1, b1); /* 5 */
3415 rg21 = iprod(r1_kj, r1_kj); /* 5 */
3421 rabr1 = sqrt(ra2r1*rb2r1);
3423 sin_phi1 = rg1 * rabr1 * iprod(a1, h1) * (-1);
3425 if (cos_phi1 < -0.5 || cos_phi1 > 0.5)
3427 phi1 = asin(sin_phi1);
3437 phi1 = -M_PI - phi1;
3443 phi1 = acos(cos_phi1);
3451 xphi1 = phi1 + M_PI; /* 1 */
3453 /* Second torsion */
3459 phi2 = dih_angle(x[a2i], x[a2j], x[a2k], x[a2l], pbc, r2_ij, r2_kj, r2_kl, m2, n2,
3460 &sign2, &t12, &t22, &t32); /* 84 */
3462 cos_phi2 = cos(phi2);
3464 a2[0] = r2_ij[1]*r2_kj[2]-r2_ij[2]*r2_kj[1];
3465 a2[1] = r2_ij[2]*r2_kj[0]-r2_ij[0]*r2_kj[2];
3466 a2[2] = r2_ij[0]*r2_kj[1]-r2_ij[1]*r2_kj[0]; /* 9 */
3468 b2[0] = r2_kl[1]*r2_kj[2]-r2_kl[2]*r2_kj[1];
3469 b2[1] = r2_kl[2]*r2_kj[0]-r2_kl[0]*r2_kj[2];
3470 b2[2] = r2_kl[0]*r2_kj[1]-r2_kl[1]*r2_kj[0]; /* 9 */
3472 pbc_rvec_sub(pbc, x[a2l], x[a2k], h2);
3474 ra22 = iprod(a2, a2); /* 5 */
3475 rb22 = iprod(b2, b2); /* 5 */
3476 rg22 = iprod(r2_kj, r2_kj); /* 5 */
3482 rabr2 = sqrt(ra2r2*rb2r2);
3484 sin_phi2 = rg2 * rabr2 * iprod(a2, h2) * (-1);
3486 if (cos_phi2 < -0.5 || cos_phi2 > 0.5)
3488 phi2 = asin(sin_phi2);
3498 phi2 = -M_PI - phi2;
3504 phi2 = acos(cos_phi2);
3512 xphi2 = phi2 + M_PI; /* 1 */
3514 /* Range mangling */
3517 xphi1 = xphi1 + 2*M_PI;
3519 else if (xphi1 >= 2*M_PI)
3521 xphi1 = xphi1 - 2*M_PI;
3526 xphi2 = xphi2 + 2*M_PI;
3528 else if (xphi2 >= 2*M_PI)
3530 xphi2 = xphi2 - 2*M_PI;
3533 /* Number of grid points */
3534 dx = 2*M_PI / cmap_grid->grid_spacing;
3536 /* Where on the grid are we */
3537 iphi1 = static_cast<int>(xphi1/dx);
3538 iphi2 = static_cast<int>(xphi2/dx);
3540 iphi1 = cmap_setup_grid_index(iphi1, cmap_grid->grid_spacing, &ip1m1, &ip1p1, &ip1p2);
3541 iphi2 = cmap_setup_grid_index(iphi2, cmap_grid->grid_spacing, &ip2m1, &ip2p1, &ip2p2);
3543 pos1 = iphi1*cmap_grid->grid_spacing+iphi2;
3544 pos2 = ip1p1*cmap_grid->grid_spacing+iphi2;
3545 pos3 = ip1p1*cmap_grid->grid_spacing+ip2p1;
3546 pos4 = iphi1*cmap_grid->grid_spacing+ip2p1;
3548 ty[0] = cmapd[pos1*4];
3549 ty[1] = cmapd[pos2*4];
3550 ty[2] = cmapd[pos3*4];
3551 ty[3] = cmapd[pos4*4];
3553 ty1[0] = cmapd[pos1*4+1];
3554 ty1[1] = cmapd[pos2*4+1];
3555 ty1[2] = cmapd[pos3*4+1];
3556 ty1[3] = cmapd[pos4*4+1];
3558 ty2[0] = cmapd[pos1*4+2];
3559 ty2[1] = cmapd[pos2*4+2];
3560 ty2[2] = cmapd[pos3*4+2];
3561 ty2[3] = cmapd[pos4*4+2];
3563 ty12[0] = cmapd[pos1*4+3];
3564 ty12[1] = cmapd[pos2*4+3];
3565 ty12[2] = cmapd[pos3*4+3];
3566 ty12[3] = cmapd[pos4*4+3];
3568 /* Switch to degrees */
3569 dx = 360.0 / cmap_grid->grid_spacing;
3570 xphi1 = xphi1 * RAD2DEG;
3571 xphi2 = xphi2 * RAD2DEG;
3573 for (i = 0; i < 4; i++) /* 16 */
3576 tx[i+4] = ty1[i]*dx;
3577 tx[i+8] = ty2[i]*dx;
3578 tx[i+12] = ty12[i]*dx*dx;
3582 for (i = 0; i < 4; i++) /* 1056 */
3584 for (j = 0; j < 4; j++)
3587 for (k = 0; k < 16; k++)
3589 xx = xx + cmap_coeff_matrix[k*16+idx]*tx[k];
3597 tt = (xphi1-iphi1*dx)/dx;
3598 tu = (xphi2-iphi2*dx)/dx;
3604 for (i = 3; i >= 0; i--)
3606 l1 = loop_index[i][3];
3607 l2 = loop_index[i][2];
3608 l3 = loop_index[i][1];
3610 e = tt * e + ((tc[i*4+3]*tu+tc[i*4+2])*tu + tc[i*4+1])*tu+tc[i*4];
3611 df1 = tu * df1 + (3.0*tc[l1]*tt+2.0*tc[l2])*tt+tc[l3];
3612 df2 = tt * df2 + (3.0*tc[i*4+3]*tu+2.0*tc[i*4+2])*tu+tc[i*4+1];
3622 /* Do forces - first torsion */
3623 fg1 = iprod(r1_ij, r1_kj);
3624 hg1 = iprod(r1_kl, r1_kj);
3625 fga1 = fg1*ra2r1*rgr1;
3626 hgb1 = hg1*rb2r1*rgr1;
3630 for (i = 0; i < DIM; i++)
3632 dtf1[i] = gaa1 * a1[i];
3633 dtg1[i] = fga1 * a1[i] - hgb1 * b1[i];
3634 dth1[i] = gbb1 * b1[i];
3636 f1[i] = df1 * dtf1[i];
3637 g1[i] = df1 * dtg1[i];
3638 h1[i] = df1 * dth1[i];
3641 f1_j[i] = -f1[i] - g1[i];
3642 f1_k[i] = h1[i] + g1[i];
3645 f[a1i][i] = f[a1i][i] + f1_i[i];
3646 f[a1j][i] = f[a1j][i] + f1_j[i]; /* - f1[i] - g1[i] */
3647 f[a1k][i] = f[a1k][i] + f1_k[i]; /* h1[i] + g1[i] */
3648 f[a1l][i] = f[a1l][i] + f1_l[i]; /* h1[i] */
3651 /* Do forces - second torsion */
3652 fg2 = iprod(r2_ij, r2_kj);
3653 hg2 = iprod(r2_kl, r2_kj);
3654 fga2 = fg2*ra2r2*rgr2;
3655 hgb2 = hg2*rb2r2*rgr2;
3659 for (i = 0; i < DIM; i++)
3661 dtf2[i] = gaa2 * a2[i];
3662 dtg2[i] = fga2 * a2[i] - hgb2 * b2[i];
3663 dth2[i] = gbb2 * b2[i];
3665 f2[i] = df2 * dtf2[i];
3666 g2[i] = df2 * dtg2[i];
3667 h2[i] = df2 * dth2[i];
3670 f2_j[i] = -f2[i] - g2[i];
3671 f2_k[i] = h2[i] + g2[i];
3674 f[a2i][i] = f[a2i][i] + f2_i[i]; /* f2[i] */
3675 f[a2j][i] = f[a2j][i] + f2_j[i]; /* - f2[i] - g2[i] */
3676 f[a2k][i] = f[a2k][i] + f2_k[i]; /* h2[i] + g2[i] */
3677 f[a2l][i] = f[a2l][i] + f2_l[i]; /* - h2[i] */
3683 copy_ivec(SHIFT_IVEC(g, a1j), jt1);
3684 ivec_sub(SHIFT_IVEC(g, a1i), jt1, dt1_ij);
3685 ivec_sub(SHIFT_IVEC(g, a1k), jt1, dt1_kj);
3686 ivec_sub(SHIFT_IVEC(g, a1l), jt1, dt1_lj);
3687 t11 = IVEC2IS(dt1_ij);
3688 t21 = IVEC2IS(dt1_kj);
3689 t31 = IVEC2IS(dt1_lj);
3691 copy_ivec(SHIFT_IVEC(g, a2j), jt2);
3692 ivec_sub(SHIFT_IVEC(g, a2i), jt2, dt2_ij);
3693 ivec_sub(SHIFT_IVEC(g, a2k), jt2, dt2_kj);
3694 ivec_sub(SHIFT_IVEC(g, a2l), jt2, dt2_lj);
3695 t12 = IVEC2IS(dt2_ij);
3696 t22 = IVEC2IS(dt2_kj);
3697 t32 = IVEC2IS(dt2_lj);
3701 t31 = pbc_rvec_sub(pbc, x[a1l], x[a1j], h1);
3702 t32 = pbc_rvec_sub(pbc, x[a2l], x[a2j], h2);
3710 rvec_inc(fshift[t11], f1_i);
3711 rvec_inc(fshift[CENTRAL], f1_j);
3712 rvec_inc(fshift[t21], f1_k);
3713 rvec_inc(fshift[t31], f1_l);
3715 rvec_inc(fshift[t21], f2_i);
3716 rvec_inc(fshift[CENTRAL], f2_j);
3717 rvec_inc(fshift[t22], f2_k);
3718 rvec_inc(fshift[t32], f2_l);
3725 /***********************************************************
3727 * G R O M O S 9 6 F U N C T I O N S
3729 ***********************************************************/
3730 real g96harmonic(real kA, real kB, real xA, real xB, real x, real lambda,
3733 const real half = 0.5;
3734 real L1, kk, x0, dx, dx2;
3735 real v, f, dvdlambda;
3738 kk = L1*kA+lambda*kB;
3739 x0 = L1*xA+lambda*xB;
3746 dvdlambda = half*(kB-kA)*dx2 + (xA-xB)*kk*dx;
3753 /* That was 21 flops */
3756 real g96bonds(int nbonds,
3757 const t_iatom forceatoms[], const t_iparams forceparams[],
3758 const rvec x[], rvec f[], rvec fshift[],
3759 const t_pbc *pbc, const t_graph *g,
3760 real lambda, real *dvdlambda,
3761 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3762 int gmx_unused *global_atom_index)
3764 int i, m, ki, ai, aj, type;
3765 real dr2, fbond, vbond, fij, vtot;
3770 for (i = 0; (i < nbonds); )
3772 type = forceatoms[i++];
3773 ai = forceatoms[i++];
3774 aj = forceatoms[i++];
3776 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
3777 dr2 = iprod(dx, dx); /* 5 */
3779 *dvdlambda += g96harmonic(forceparams[type].harmonic.krA,
3780 forceparams[type].harmonic.krB,
3781 forceparams[type].harmonic.rA,
3782 forceparams[type].harmonic.rB,
3783 dr2, lambda, &vbond, &fbond);
3785 vtot += 0.5*vbond; /* 1*/
3789 fprintf(debug, "G96-BONDS: dr = %10g vbond = %10g fbond = %10g\n",
3790 sqrt(dr2), vbond, fbond);
3796 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
3799 for (m = 0; (m < DIM); m++) /* 15 */
3804 fshift[ki][m] += fij;
3805 fshift[CENTRAL][m] -= fij;
3811 real g96bond_angle(const rvec xi, const rvec xj, const rvec xk, const t_pbc *pbc,
3812 rvec r_ij, rvec r_kj,
3814 /* Return value is the angle between the bonds i-j and j-k */
3818 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
3819 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
3821 costh = cos_angle(r_ij, r_kj); /* 25 */
3826 real g96angles(int nbonds,
3827 const t_iatom forceatoms[], const t_iparams forceparams[],
3828 const rvec x[], rvec f[], rvec fshift[],
3829 const t_pbc *pbc, const t_graph *g,
3830 real lambda, real *dvdlambda,
3831 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3832 int gmx_unused *global_atom_index)
3834 int i, ai, aj, ak, type, m, t1, t2;
3836 real cos_theta, dVdt, va, vtot;
3837 real rij_1, rij_2, rkj_1, rkj_2, rijrkj_1;
3839 ivec jt, dt_ij, dt_kj;
3842 for (i = 0; (i < nbonds); )
3844 type = forceatoms[i++];
3845 ai = forceatoms[i++];
3846 aj = forceatoms[i++];
3847 ak = forceatoms[i++];
3849 cos_theta = g96bond_angle(x[ai], x[aj], x[ak], pbc, r_ij, r_kj, &t1, &t2);
3851 *dvdlambda += g96harmonic(forceparams[type].harmonic.krA,
3852 forceparams[type].harmonic.krB,
3853 forceparams[type].harmonic.rA,
3854 forceparams[type].harmonic.rB,
3855 cos_theta, lambda, &va, &dVdt);
3858 rij_1 = gmx_invsqrt(iprod(r_ij, r_ij));
3859 rkj_1 = gmx_invsqrt(iprod(r_kj, r_kj));
3860 rij_2 = rij_1*rij_1;
3861 rkj_2 = rkj_1*rkj_1;
3862 rijrkj_1 = rij_1*rkj_1; /* 23 */
3867 fprintf(debug, "G96ANGLES: costheta = %10g vth = %10g dV/dct = %10g\n",
3868 cos_theta, va, dVdt);
3871 for (m = 0; (m < DIM); m++) /* 42 */
3873 f_i[m] = dVdt*(r_kj[m]*rijrkj_1 - r_ij[m]*rij_2*cos_theta);
3874 f_k[m] = dVdt*(r_ij[m]*rijrkj_1 - r_kj[m]*rkj_2*cos_theta);
3875 f_j[m] = -f_i[m]-f_k[m];
3883 copy_ivec(SHIFT_IVEC(g, aj), jt);
3885 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3886 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3887 t1 = IVEC2IS(dt_ij);
3888 t2 = IVEC2IS(dt_kj);
3890 rvec_inc(fshift[t1], f_i);
3891 rvec_inc(fshift[CENTRAL], f_j);
3892 rvec_inc(fshift[t2], f_k); /* 9 */
3898 real cross_bond_bond(int nbonds,
3899 const t_iatom forceatoms[], const t_iparams forceparams[],
3900 const rvec x[], rvec f[], rvec fshift[],
3901 const t_pbc *pbc, const t_graph *g,
3902 real gmx_unused lambda, real gmx_unused *dvdlambda,
3903 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3904 int gmx_unused *global_atom_index)
3906 /* Potential from Lawrence and Skimmer, Chem. Phys. Lett. 372 (2003)
3909 int i, ai, aj, ak, type, m, t1, t2;
3911 real vtot, vrr, s1, s2, r1, r2, r1e, r2e, krr;
3913 ivec jt, dt_ij, dt_kj;
3916 for (i = 0; (i < nbonds); )
3918 type = forceatoms[i++];
3919 ai = forceatoms[i++];
3920 aj = forceatoms[i++];
3921 ak = forceatoms[i++];
3922 r1e = forceparams[type].cross_bb.r1e;
3923 r2e = forceparams[type].cross_bb.r2e;
3924 krr = forceparams[type].cross_bb.krr;
3926 /* Compute distance vectors ... */
3927 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
3928 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
3930 /* ... and their lengths */
3934 /* Deviations from ideality */
3938 /* Energy (can be negative!) */
3943 svmul(-krr*s2/r1, r_ij, f_i);
3944 svmul(-krr*s1/r2, r_kj, f_k);
3946 for (m = 0; (m < DIM); m++) /* 12 */
3948 f_j[m] = -f_i[m] - f_k[m];
3957 copy_ivec(SHIFT_IVEC(g, aj), jt);
3959 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3960 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3961 t1 = IVEC2IS(dt_ij);
3962 t2 = IVEC2IS(dt_kj);
3964 rvec_inc(fshift[t1], f_i);
3965 rvec_inc(fshift[CENTRAL], f_j);
3966 rvec_inc(fshift[t2], f_k); /* 9 */
3972 real cross_bond_angle(int nbonds,
3973 const t_iatom forceatoms[], const t_iparams forceparams[],
3974 const rvec x[], rvec f[], rvec fshift[],
3975 const t_pbc *pbc, const t_graph *g,
3976 real gmx_unused lambda, real gmx_unused *dvdlambda,
3977 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3978 int gmx_unused *global_atom_index)
3980 /* Potential from Lawrence and Skimmer, Chem. Phys. Lett. 372 (2003)
3983 int i, ai, aj, ak, type, m, t1, t2;
3984 rvec r_ij, r_kj, r_ik;
3985 real vtot, vrt, s1, s2, s3, r1, r2, r3, r1e, r2e, r3e, krt, k1, k2, k3;
3987 ivec jt, dt_ij, dt_kj;
3990 for (i = 0; (i < nbonds); )
3992 type = forceatoms[i++];
3993 ai = forceatoms[i++];
3994 aj = forceatoms[i++];
3995 ak = forceatoms[i++];
3996 r1e = forceparams[type].cross_ba.r1e;
3997 r2e = forceparams[type].cross_ba.r2e;
3998 r3e = forceparams[type].cross_ba.r3e;
3999 krt = forceparams[type].cross_ba.krt;
4001 /* Compute distance vectors ... */
4002 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
4003 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
4004 pbc_rvec_sub(pbc, x[ai], x[ak], r_ik);
4006 /* ... and their lengths */
4011 /* Deviations from ideality */
4016 /* Energy (can be negative!) */
4017 vrt = krt*s3*(s1+s2);
4023 k3 = -krt*(s1+s2)/r3;
4024 for (m = 0; (m < DIM); m++)
4026 f_i[m] = k1*r_ij[m] + k3*r_ik[m];
4027 f_k[m] = k2*r_kj[m] - k3*r_ik[m];
4028 f_j[m] = -f_i[m] - f_k[m];
4031 for (m = 0; (m < DIM); m++) /* 12 */
4041 copy_ivec(SHIFT_IVEC(g, aj), jt);
4043 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
4044 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
4045 t1 = IVEC2IS(dt_ij);
4046 t2 = IVEC2IS(dt_kj);
4048 rvec_inc(fshift[t1], f_i);
4049 rvec_inc(fshift[CENTRAL], f_j);
4050 rvec_inc(fshift[t2], f_k); /* 9 */
4056 static real bonded_tab(const char *type, int table_nr,
4057 const bondedtable_t *table, real kA, real kB, real r,
4058 real lambda, real *V, real *F)
4060 real k, tabscale, *VFtab, rt, eps, eps2, Yt, Ft, Geps, Heps2, Fp, VV, FF;
4064 k = (1.0 - lambda)*kA + lambda*kB;
4066 tabscale = table->scale;
4067 VFtab = table->data;
4070 n0 = static_cast<int>(rt);
4073 gmx_fatal(FARGS, "A tabulated %s interaction table number %d is out of the table range: r %f, between table indices %d and %d, table length %d",
4074 type, table_nr, r, n0, n0+1, table->n);
4081 Geps = VFtab[nnn+2]*eps;
4082 Heps2 = VFtab[nnn+3]*eps2;
4083 Fp = Ft + Geps + Heps2;
4085 FF = Fp + Geps + 2.0*Heps2;
4087 *F = -k*FF*tabscale;
4089 dvdlambda = (kB - kA)*VV;
4093 /* That was 22 flops */
4096 real tab_bonds(int nbonds,
4097 const t_iatom forceatoms[], const t_iparams forceparams[],
4098 const rvec x[], rvec f[], rvec fshift[],
4099 const t_pbc *pbc, const t_graph *g,
4100 real lambda, real *dvdlambda,
4101 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
4102 int gmx_unused *global_atom_index)
4104 int i, m, ki, ai, aj, type, table;
4105 real dr, dr2, fbond, vbond, fij, vtot;
4110 for (i = 0; (i < nbonds); )
4112 type = forceatoms[i++];
4113 ai = forceatoms[i++];
4114 aj = forceatoms[i++];
4116 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
4117 dr2 = iprod(dx, dx); /* 5 */
4118 dr = dr2*gmx_invsqrt(dr2); /* 10 */
4120 table = forceparams[type].tab.table;
4122 *dvdlambda += bonded_tab("bond", table,
4123 &fcd->bondtab[table],
4124 forceparams[type].tab.kA,
4125 forceparams[type].tab.kB,
4126 dr, lambda, &vbond, &fbond); /* 22 */
4134 vtot += vbond; /* 1*/
4135 fbond *= gmx_invsqrt(dr2); /* 6 */
4139 fprintf(debug, "TABBONDS: dr = %10g vbond = %10g fbond = %10g\n",
4145 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
4148 for (m = 0; (m < DIM); m++) /* 15 */
4153 fshift[ki][m] += fij;
4154 fshift[CENTRAL][m] -= fij;
4160 real tab_angles(int nbonds,
4161 const t_iatom forceatoms[], const t_iparams forceparams[],
4162 const rvec x[], rvec f[], rvec fshift[],
4163 const t_pbc *pbc, const t_graph *g,
4164 real lambda, real *dvdlambda,
4165 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
4166 int gmx_unused *global_atom_index)
4168 int i, ai, aj, ak, t1, t2, type, table;
4170 real cos_theta, cos_theta2, theta, dVdt, va, vtot;
4171 ivec jt, dt_ij, dt_kj;
4174 for (i = 0; (i < nbonds); )
4176 type = forceatoms[i++];
4177 ai = forceatoms[i++];
4178 aj = forceatoms[i++];
4179 ak = forceatoms[i++];
4181 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
4182 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
4184 table = forceparams[type].tab.table;
4186 *dvdlambda += bonded_tab("angle", table,
4187 &fcd->angletab[table],
4188 forceparams[type].tab.kA,
4189 forceparams[type].tab.kB,
4190 theta, lambda, &va, &dVdt); /* 22 */
4193 cos_theta2 = sqr(cos_theta); /* 1 */
4202 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
4203 sth = st*cos_theta; /* 1 */
4207 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
4208 theta*RAD2DEG, va, dVdt);
4211 nrkj2 = iprod(r_kj, r_kj); /* 5 */
4212 nrij2 = iprod(r_ij, r_ij);
4214 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
4215 cii = sth/nrij2; /* 10 */
4216 ckk = sth/nrkj2; /* 10 */
4218 for (m = 0; (m < DIM); m++) /* 39 */
4220 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
4221 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
4222 f_j[m] = -f_i[m]-f_k[m];
4229 copy_ivec(SHIFT_IVEC(g, aj), jt);
4231 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
4232 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
4233 t1 = IVEC2IS(dt_ij);
4234 t2 = IVEC2IS(dt_kj);
4236 rvec_inc(fshift[t1], f_i);
4237 rvec_inc(fshift[CENTRAL], f_j);
4238 rvec_inc(fshift[t2], f_k);
4244 real tab_dihs(int nbonds,
4245 const t_iatom forceatoms[], const t_iparams forceparams[],
4246 const rvec x[], rvec f[], rvec fshift[],
4247 const t_pbc *pbc, const t_graph *g,
4248 real lambda, real *dvdlambda,
4249 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
4250 int gmx_unused *global_atom_index)
4252 int i, type, ai, aj, ak, al, table;
4254 rvec r_ij, r_kj, r_kl, m, n;
4255 real phi, sign, ddphi, vpd, vtot;
4258 for (i = 0; (i < nbonds); )
4260 type = forceatoms[i++];
4261 ai = forceatoms[i++];
4262 aj = forceatoms[i++];
4263 ak = forceatoms[i++];
4264 al = forceatoms[i++];
4266 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
4267 &sign, &t1, &t2, &t3); /* 84 */
4269 table = forceparams[type].tab.table;
4271 /* Hopefully phi+M_PI never results in values < 0 */
4272 *dvdlambda += bonded_tab("dihedral", table,
4273 &fcd->dihtab[table],
4274 forceparams[type].tab.kA,
4275 forceparams[type].tab.kB,
4276 phi+M_PI, lambda, &vpd, &ddphi);
4279 do_dih_fup(ai, aj, ak, al, -ddphi, r_ij, r_kj, r_kl, m, n,
4280 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
4283 fprintf(debug, "pdih: (%d,%d,%d,%d) phi=%g\n",
4284 ai, aj, ak, al, phi);
4293 /*! \brief Return true if ftype is an explicit pair-listed LJ or
4294 * COULOMB interaction type: bonded LJ (usually 1-4), or special
4295 * listed non-bonded for FEP. */
4297 isPairInteraction(int ftype)
4299 return ((ftype) >= F_LJ14 && (ftype) <= F_LJC_PAIRS_NB);
4303 ftype_is_bonded_potential(int ftype)
4306 (interaction_function[ftype].flags & IF_BOND) &&
4307 !(ftype == F_CONNBONDS || ftype == F_POSRES || ftype == F_FBPOSRES) &&
4308 (ftype < F_GB12 || ftype > F_GB14);
4311 /*! \brief Zero thread-local force-output buffers */
4312 static void zero_thread_forces(f_thread_t *f_t, int n,
4313 int nblock, int blocksize)
4315 int b, a0, a1, a, i, j;
4317 if (n > f_t->f_nalloc)
4319 f_t->f_nalloc = over_alloc_large(n);
4320 srenew(f_t->f, f_t->f_nalloc);
4323 if (f_t->red_mask != 0)
4325 for (b = 0; b < nblock; b++)
4327 if (f_t->red_mask && (1U<<b))
4330 a1 = std::min((b+1)*blocksize, n);
4331 for (a = a0; a < a1; a++)
4333 clear_rvec(f_t->f[a]);
4338 for (i = 0; i < SHIFTS; i++)
4340 clear_rvec(f_t->fshift[i]);
4342 for (i = 0; i < F_NRE; i++)
4346 for (i = 0; i < egNR; i++)
4348 for (j = 0; j < f_t->grpp.nener; j++)
4350 f_t->grpp.ener[i][j] = 0;
4353 for (i = 0; i < efptNR; i++)
4359 /*! \brief The max thread number is arbitrary, we used a fixed number
4360 * to avoid memory management. Using more than 16 threads is probably
4361 * never useful performance wise. */
4362 #define MAX_BONDED_THREADS 256
4364 /*! \brief Reduce thread-local force buffers */
4365 static void reduce_thread_force_buffer(int n, rvec *f,
4366 int nthreads, f_thread_t *f_t,
4367 int nblock, int block_size)
4371 if (nthreads > MAX_BONDED_THREADS)
4373 gmx_fatal(FARGS, "Can not reduce bonded forces on more than %d threads",
4374 MAX_BONDED_THREADS);
4377 /* This reduction can run on any number of threads,
4378 * independently of nthreads.
4380 #pragma omp parallel for num_threads(nthreads) schedule(static)
4381 for (b = 0; b < nblock; b++)
4383 rvec *fp[MAX_BONDED_THREADS];
4387 /* Determine which threads contribute to this block */
4389 for (ft = 1; ft < nthreads; ft++)
4391 if (f_t[ft].red_mask & (1U<<b))
4393 fp[nfb++] = f_t[ft].f;
4398 /* Reduce force buffers for threads that contribute */
4400 a1 = (b+1)*block_size;
4401 a1 = std::min(a1, n);
4402 for (a = a0; a < a1; a++)
4404 for (fb = 0; fb < nfb; fb++)
4406 rvec_inc(f[a], fp[fb][a]);
4413 /*! \brief Reduce thread-local forces */
4414 static void reduce_thread_forces(int n, rvec *f, rvec *fshift,
4415 real *ener, gmx_grppairener_t *grpp, real *dvdl,
4416 int nthreads, f_thread_t *f_t,
4417 int nblock, int block_size,
4418 gmx_bool bCalcEnerVir,
4423 /* Reduce the bonded force buffer */
4424 reduce_thread_force_buffer(n, f, nthreads, f_t, nblock, block_size);
4427 /* When necessary, reduce energy and virial using one thread only */
4432 for (i = 0; i < SHIFTS; i++)
4434 for (t = 1; t < nthreads; t++)
4436 rvec_inc(fshift[i], f_t[t].fshift[i]);
4439 for (i = 0; i < F_NRE; i++)
4441 for (t = 1; t < nthreads; t++)
4443 ener[i] += f_t[t].ener[i];
4446 for (i = 0; i < egNR; i++)
4448 for (j = 0; j < f_t[1].grpp.nener; j++)
4450 for (t = 1; t < nthreads; t++)
4453 grpp->ener[i][j] += f_t[t].grpp.ener[i][j];
4459 for (i = 0; i < efptNR; i++)
4462 for (t = 1; t < nthreads; t++)
4464 dvdl[i] += f_t[t].dvdl[i];
4471 /*! \brief Calculate one element of the list of bonded interactions
4473 static real calc_one_bond(int thread,
4474 int ftype, const t_idef *idef,
4475 const rvec x[], rvec f[], rvec fshift[],
4477 const t_pbc *pbc, const t_graph *g,
4478 gmx_grppairener_t *grpp,
4480 real *lambda, real *dvdl,
4481 const t_mdatoms *md, t_fcdata *fcd,
4482 gmx_bool bCalcEnerVir,
4483 int *global_atom_index)
4485 int nat1, nbonds, efptFTYPE;
4490 if (IS_RESTRAINT_TYPE(ftype))
4492 efptFTYPE = efptRESTRAINT;
4496 efptFTYPE = efptBONDED;
4499 nat1 = interaction_function[ftype].nratoms + 1;
4500 nbonds = idef->il[ftype].nr/nat1;
4501 iatoms = idef->il[ftype].iatoms;
4503 nb0 = idef->il_thread_division[ftype*(idef->nthreads+1)+thread];
4504 nbn = idef->il_thread_division[ftype*(idef->nthreads+1)+thread+1] - nb0;
4506 if (!isPairInteraction(ftype))
4508 if (ftype == F_CMAP)
4510 v = cmap_dihs(nbn, iatoms+nb0,
4511 idef->iparams, &idef->cmap_grid,
4513 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
4514 md, fcd, global_atom_index);
4516 #ifdef GMX_SIMD_HAVE_REAL
4517 else if (ftype == F_ANGLES &&
4518 !bCalcEnerVir && fr->efep == efepNO)
4520 /* No energies, shift forces, dvdl */
4521 angles_noener_simd(nbn, idef->il[ftype].iatoms+nb0,
4524 pbc, g, lambda[efptFTYPE], md, fcd,
4529 else if (ftype == F_PDIHS &&
4530 !bCalcEnerVir && fr->efep == efepNO)
4532 /* No energies, shift forces, dvdl */
4533 #ifdef GMX_SIMD_HAVE_REAL
4538 (nbn, idef->il[ftype].iatoms+nb0,
4541 pbc, g, lambda[efptFTYPE], md, fcd,
4545 #ifdef GMX_SIMD_HAVE_REAL
4546 else if (ftype == F_RBDIHS &&
4547 !bCalcEnerVir && fr->efep == efepNO)
4549 /* No energies, shift forces, dvdl */
4550 rbdihs_noener_simd(nbn, idef->il[ftype].iatoms+nb0,
4553 pbc, g, lambda[efptFTYPE], md, fcd,
4560 v = interaction_function[ftype].ifunc(nbn, iatoms+nb0,
4563 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
4564 md, fcd, global_atom_index);
4569 v = do_pairs(ftype, nbn, iatoms+nb0, idef->iparams, x, f, fshift,
4570 pbc, g, lambda, dvdl, md, fr, grpp, global_atom_index);
4575 inc_nrnb(nrnb, interaction_function[ftype].nrnb_ind, nbonds);
4581 void calc_bonds(const gmx_multisim_t *ms,
4583 const rvec x[], history_t *hist,
4584 rvec f[], t_forcerec *fr,
4585 const struct t_pbc *pbc, const struct t_graph *g,
4586 gmx_enerdata_t *enerd, t_nrnb *nrnb,
4588 const t_mdatoms *md,
4589 t_fcdata *fcd, int *global_atom_index,
4592 gmx_bool bCalcEnerVir;
4594 real dvdl[efptNR]; /* The dummy array is to have a place to store the dhdl at other values
4595 of lambda, which will be thrown away in the end*/
4596 const t_pbc *pbc_null;
4599 assert(fr->nthreads == idef->nthreads);
4601 bCalcEnerVir = (force_flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY));
4603 for (i = 0; i < efptNR; i++)
4619 p_graph(debug, "Bondage is fun", g);
4623 /* Do pre force calculation stuff which might require communication */
4624 if (idef->il[F_ORIRES].nr)
4626 enerd->term[F_ORIRESDEV] =
4627 calc_orires_dev(ms, idef->il[F_ORIRES].nr,
4628 idef->il[F_ORIRES].iatoms,
4629 idef->iparams, md, x,
4630 pbc_null, fcd, hist);
4632 if (idef->il[F_DISRES].nr)
4634 calc_disres_R_6(idef->il[F_DISRES].nr,
4635 idef->il[F_DISRES].iatoms,
4636 idef->iparams, x, pbc_null,
4639 if (fcd->disres.nsystems > 1)
4641 gmx_sum_sim(2*fcd->disres.nres, fcd->disres.Rt_6, ms);
4646 #pragma omp parallel for num_threads(fr->nthreads) schedule(static)
4647 for (thread = 0; thread < fr->nthreads; thread++)
4654 gmx_grppairener_t *grpp;
4659 fshift = fr->fshift;
4661 grpp = &enerd->grpp;
4666 zero_thread_forces(&fr->f_t[thread], fr->natoms_force,
4667 fr->red_nblock, 1<<fr->red_ashift);
4669 ft = fr->f_t[thread].f;
4670 fshift = fr->f_t[thread].fshift;
4671 epot = fr->f_t[thread].ener;
4672 grpp = &fr->f_t[thread].grpp;
4673 dvdlt = fr->f_t[thread].dvdl;
4675 /* Loop over all bonded force types to calculate the bonded forces */
4676 for (ftype = 0; (ftype < F_NRE); ftype++)
4678 if (idef->il[ftype].nr > 0 && ftype_is_bonded_potential(ftype))
4680 v = calc_one_bond(thread, ftype, idef, x,
4681 ft, fshift, fr, pbc_null, g, grpp,
4682 nrnb, lambda, dvdlt,
4683 md, fcd, bCalcEnerVir,
4689 if (fr->nthreads > 1)
4691 reduce_thread_forces(fr->natoms_force, f, fr->fshift,
4692 enerd->term, &enerd->grpp, dvdl,
4693 fr->nthreads, fr->f_t,
4694 fr->red_nblock, 1<<fr->red_ashift,
4696 force_flags & GMX_FORCE_DHDL);
4698 if (force_flags & GMX_FORCE_DHDL)
4700 for (i = 0; i < efptNR; i++)
4702 enerd->dvdl_nonlin[i] += dvdl[i];
4706 /* Copy the sum of violations for the distance restraints from fcd */
4709 enerd->term[F_DISRESVIOL] = fcd->disres.sumviol;
4714 void calc_bonds_lambda(const t_idef *idef,
4717 const struct t_pbc *pbc, const struct t_graph *g,
4718 gmx_grppairener_t *grpp, real *epot, t_nrnb *nrnb,
4720 const t_mdatoms *md,
4722 int *global_atom_index)
4724 int ftype, nr_nonperturbed, nr;
4726 real dvdl_dum[efptNR];
4728 const t_pbc *pbc_null;
4740 /* Copy the whole idef, so we can modify the contents locally */
4742 idef_fe.nthreads = 1;
4743 snew(idef_fe.il_thread_division, F_NRE*(idef_fe.nthreads+1));
4745 /* We already have the forces, so we use temp buffers here */
4746 snew(f, fr->natoms_force);
4747 snew(fshift, SHIFTS);
4749 /* Loop over all bonded force types to calculate the bonded energies */
4750 for (ftype = 0; (ftype < F_NRE); ftype++)
4752 if (ftype_is_bonded_potential(ftype))
4754 /* Set the work range of thread 0 to the perturbed bondeds only */
4755 nr_nonperturbed = idef->il[ftype].nr_nonperturbed;
4756 nr = idef->il[ftype].nr;
4757 idef_fe.il_thread_division[ftype*2+0] = nr_nonperturbed;
4758 idef_fe.il_thread_division[ftype*2+1] = nr;
4760 /* This is only to get the flop count correct */
4761 idef_fe.il[ftype].nr = nr - nr_nonperturbed;
4763 if (nr - nr_nonperturbed > 0)
4765 v = calc_one_bond(0, ftype, &idef_fe,
4766 x, f, fshift, fr, pbc_null, g,
4767 grpp, nrnb, lambda, dvdl_dum,
4778 sfree(idef_fe.il_thread_division);