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.
42 #include "gromacs/math/units.h"
43 #include "gromacs/math/vec.h"
44 #include "gromacs/math/utilities.h"
53 #include "nonbonded.h"
56 #include "gromacs/pbcutil/ishift.h"
57 #include "gromacs/pbcutil/mshift.h"
58 #include "gromacs/pbcutil/pbc.h"
59 #include "gromacs/simd/simd.h"
60 #include "gromacs/simd/simd_math.h"
61 #include "gromacs/simd/vector_operations.h"
62 #include "gromacs/utility/fatalerror.h"
63 #include "gromacs/utility/smalloc.h"
65 /* Find a better place for this? */
66 const int cmap_coeff_matrix[] = {
67 1, 0, -3, 2, 0, 0, 0, 0, -3, 0, 9, -6, 2, 0, -6, 4,
68 0, 0, 0, 0, 0, 0, 0, 0, 3, 0, -9, 6, -2, 0, 6, -4,
69 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 9, -6, 0, 0, -6, 4,
70 0, 0, 3, -2, 0, 0, 0, 0, 0, 0, -9, 6, 0, 0, 6, -4,
71 0, 0, 0, 0, 1, 0, -3, 2, -2, 0, 6, -4, 1, 0, -3, 2,
72 0, 0, 0, 0, 0, 0, 0, 0, -1, 0, 3, -2, 1, 0, -3, 2,
73 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -3, 2, 0, 0, 3, -2,
74 0, 0, 0, 0, 0, 0, 3, -2, 0, 0, -6, 4, 0, 0, 3, -2,
75 0, 1, -2, 1, 0, 0, 0, 0, 0, -3, 6, -3, 0, 2, -4, 2,
76 0, 0, 0, 0, 0, 0, 0, 0, 0, 3, -6, 3, 0, -2, 4, -2,
77 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -3, 3, 0, 0, 2, -2,
78 0, 0, -1, 1, 0, 0, 0, 0, 0, 0, 3, -3, 0, 0, -2, 2,
79 0, 0, 0, 0, 0, 1, -2, 1, 0, -2, 4, -2, 0, 1, -2, 1,
80 0, 0, 0, 0, 0, 0, 0, 0, 0, -1, 2, -1, 0, 1, -2, 1,
81 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, -1, 0, 0, -1, 1,
82 0, 0, 0, 0, 0, 0, -1, 1, 0, 0, 2, -2, 0, 0, -1, 1
87 int glatnr(int *global_atom_index, int i)
91 if (global_atom_index == NULL)
97 atnr = global_atom_index[i] + 1;
103 static int pbc_rvec_sub(const t_pbc *pbc, const rvec xi, const rvec xj, rvec dx)
107 return pbc_dx_aiuc(pbc, xi, xj, dx);
111 rvec_sub(xi, xj, dx);
116 #ifdef GMX_SIMD_HAVE_REAL
118 /* SIMD PBC data structure, containing 1/boxdiag and the box vectors */
120 gmx_simd_real_t inv_bzz;
121 gmx_simd_real_t inv_byy;
122 gmx_simd_real_t inv_bxx;
131 /* Set the SIMD pbc data from a normal t_pbc struct */
132 static void set_pbc_simd(const t_pbc *pbc, pbc_simd_t *pbc_simd)
137 /* Setting inv_bdiag to 0 effectively turns off PBC */
138 clear_rvec(inv_bdiag);
141 for (d = 0; d < pbc->ndim_ePBC; d++)
143 inv_bdiag[d] = 1.0/pbc->box[d][d];
147 pbc_simd->inv_bzz = gmx_simd_set1_r(inv_bdiag[ZZ]);
148 pbc_simd->inv_byy = gmx_simd_set1_r(inv_bdiag[YY]);
149 pbc_simd->inv_bxx = gmx_simd_set1_r(inv_bdiag[XX]);
153 pbc_simd->bzx = gmx_simd_set1_r(pbc->box[ZZ][XX]);
154 pbc_simd->bzy = gmx_simd_set1_r(pbc->box[ZZ][YY]);
155 pbc_simd->bzz = gmx_simd_set1_r(pbc->box[ZZ][ZZ]);
156 pbc_simd->byx = gmx_simd_set1_r(pbc->box[YY][XX]);
157 pbc_simd->byy = gmx_simd_set1_r(pbc->box[YY][YY]);
158 pbc_simd->bxx = gmx_simd_set1_r(pbc->box[XX][XX]);
162 pbc_simd->bzx = gmx_simd_setzero_r();
163 pbc_simd->bzy = gmx_simd_setzero_r();
164 pbc_simd->bzz = gmx_simd_setzero_r();
165 pbc_simd->byx = gmx_simd_setzero_r();
166 pbc_simd->byy = gmx_simd_setzero_r();
167 pbc_simd->bxx = gmx_simd_setzero_r();
171 /* Correct distance vector *dx,*dy,*dz for PBC using SIMD */
172 static gmx_inline void
173 pbc_dx_simd(gmx_simd_real_t *dx, gmx_simd_real_t *dy, gmx_simd_real_t *dz,
174 const pbc_simd_t *pbc)
178 sh = gmx_simd_round_r(gmx_simd_mul_r(*dz, pbc->inv_bzz));
179 *dx = gmx_simd_fnmadd_r(sh, pbc->bzx, *dx);
180 *dy = gmx_simd_fnmadd_r(sh, pbc->bzy, *dy);
181 *dz = gmx_simd_fnmadd_r(sh, pbc->bzz, *dz);
183 sh = gmx_simd_round_r(gmx_simd_mul_r(*dy, pbc->inv_byy));
184 *dx = gmx_simd_fnmadd_r(sh, pbc->byx, *dx);
185 *dy = gmx_simd_fnmadd_r(sh, pbc->byy, *dy);
187 sh = gmx_simd_round_r(gmx_simd_mul_r(*dx, pbc->inv_bxx));
188 *dx = gmx_simd_fnmadd_r(sh, pbc->bxx, *dx);
191 #endif /* GMX_SIMD_HAVE_REAL */
194 * Morse potential bond by Frank Everdij
196 * Three parameters needed:
198 * b0 = equilibrium distance in nm
199 * be = beta in nm^-1 (actually, it's nu_e*Sqrt(2*pi*pi*mu/D_e))
200 * cb = well depth in kJ/mol
202 * Note: the potential is referenced to be +cb at infinite separation
203 * 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;
280 real cubic_bonds(int nbonds,
281 const t_iatom forceatoms[], const t_iparams forceparams[],
282 const rvec x[], rvec f[], rvec fshift[],
283 const t_pbc *pbc, const t_graph *g,
284 real gmx_unused lambda, real gmx_unused *dvdlambda,
285 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
286 int gmx_unused *global_atom_index)
288 const real three = 3.0;
289 const real two = 2.0;
291 real dr, dr2, dist, kdist, kdist2, fbond, vbond, fij, vtot;
293 int i, m, ki, type, ai, aj;
297 for (i = 0; (i < nbonds); )
299 type = forceatoms[i++];
300 ai = forceatoms[i++];
301 aj = forceatoms[i++];
303 b0 = forceparams[type].cubic.b0;
304 kb = forceparams[type].cubic.kb;
305 kcub = forceparams[type].cubic.kcub;
307 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
308 dr2 = iprod(dx, dx); /* 5 */
315 dr = dr2*gmx_invsqrt(dr2); /* 10 */
320 vbond = kdist2 + kcub*kdist2*dist;
321 fbond = -(two*kdist + three*kdist2*kcub)/dr;
323 vtot += vbond; /* 21 */
327 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
330 for (m = 0; (m < DIM); m++) /* 15 */
335 fshift[ki][m] += fij;
336 fshift[CENTRAL][m] -= fij;
342 real FENE_bonds(int nbonds,
343 const t_iatom forceatoms[], const t_iparams forceparams[],
344 const rvec x[], rvec f[], rvec fshift[],
345 const t_pbc *pbc, const t_graph *g,
346 real gmx_unused lambda, real gmx_unused *dvdlambda,
347 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
348 int *global_atom_index)
350 const real half = 0.5;
351 const real one = 1.0;
353 real dr, dr2, bm2, omdr2obm2, fbond, vbond, fij, vtot;
355 int i, m, ki, type, ai, aj;
359 for (i = 0; (i < nbonds); )
361 type = forceatoms[i++];
362 ai = forceatoms[i++];
363 aj = forceatoms[i++];
365 bm = forceparams[type].fene.bm;
366 kb = forceparams[type].fene.kb;
368 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
369 dr2 = iprod(dx, dx); /* 5 */
381 "r^2 (%f) >= bm^2 (%f) in FENE bond between atoms %d and %d",
383 glatnr(global_atom_index, ai),
384 glatnr(global_atom_index, aj));
387 omdr2obm2 = one - dr2/bm2;
389 vbond = -half*kb*bm2*log(omdr2obm2);
390 fbond = -kb/omdr2obm2;
392 vtot += vbond; /* 35 */
396 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
399 for (m = 0; (m < DIM); m++) /* 15 */
404 fshift[ki][m] += fij;
405 fshift[CENTRAL][m] -= fij;
411 real harmonic(real kA, real kB, real xA, real xB, real x, real lambda,
414 const real half = 0.5;
415 real L1, kk, x0, dx, dx2;
416 real v, f, dvdlambda;
419 kk = L1*kA+lambda*kB;
420 x0 = L1*xA+lambda*xB;
427 dvdlambda = half*(kB-kA)*dx2 + (xA-xB)*kk*dx;
434 /* That was 19 flops */
438 real bonds(int nbonds,
439 const t_iatom forceatoms[], const t_iparams forceparams[],
440 const rvec x[], rvec f[], rvec fshift[],
441 const t_pbc *pbc, const t_graph *g,
442 real lambda, real *dvdlambda,
443 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
444 int gmx_unused *global_atom_index)
446 int i, m, ki, ai, aj, type;
447 real dr, dr2, fbond, vbond, fij, vtot;
452 for (i = 0; (i < nbonds); )
454 type = forceatoms[i++];
455 ai = forceatoms[i++];
456 aj = forceatoms[i++];
458 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
459 dr2 = iprod(dx, dx); /* 5 */
460 dr = dr2*gmx_invsqrt(dr2); /* 10 */
462 *dvdlambda += harmonic(forceparams[type].harmonic.krA,
463 forceparams[type].harmonic.krB,
464 forceparams[type].harmonic.rA,
465 forceparams[type].harmonic.rB,
466 dr, lambda, &vbond, &fbond); /* 19 */
474 vtot += vbond; /* 1*/
475 fbond *= gmx_invsqrt(dr2); /* 6 */
479 fprintf(debug, "BONDS: dr = %10g vbond = %10g fbond = %10g\n",
485 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
488 for (m = 0; (m < DIM); m++) /* 15 */
493 fshift[ki][m] += fij;
494 fshift[CENTRAL][m] -= fij;
500 real restraint_bonds(int nbonds,
501 const t_iatom forceatoms[], const t_iparams forceparams[],
502 const rvec x[], rvec f[], rvec fshift[],
503 const t_pbc *pbc, const t_graph *g,
504 real lambda, real *dvdlambda,
505 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
506 int gmx_unused *global_atom_index)
508 int i, m, ki, ai, aj, type;
509 real dr, dr2, fbond, vbond, fij, vtot;
511 real low, dlow, up1, dup1, up2, dup2, k, dk;
519 for (i = 0; (i < nbonds); )
521 type = forceatoms[i++];
522 ai = forceatoms[i++];
523 aj = forceatoms[i++];
525 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
526 dr2 = iprod(dx, dx); /* 5 */
527 dr = dr2*gmx_invsqrt(dr2); /* 10 */
529 low = L1*forceparams[type].restraint.lowA + lambda*forceparams[type].restraint.lowB;
530 dlow = -forceparams[type].restraint.lowA + forceparams[type].restraint.lowB;
531 up1 = L1*forceparams[type].restraint.up1A + lambda*forceparams[type].restraint.up1B;
532 dup1 = -forceparams[type].restraint.up1A + forceparams[type].restraint.up1B;
533 up2 = L1*forceparams[type].restraint.up2A + lambda*forceparams[type].restraint.up2B;
534 dup2 = -forceparams[type].restraint.up2A + forceparams[type].restraint.up2B;
535 k = L1*forceparams[type].restraint.kA + lambda*forceparams[type].restraint.kB;
536 dk = -forceparams[type].restraint.kA + forceparams[type].restraint.kB;
545 *dvdlambda += 0.5*dk*drh2 - k*dlow*drh;
558 *dvdlambda += 0.5*dk*drh2 - k*dup1*drh;
563 vbond = k*(up2 - up1)*(0.5*(up2 - up1) + drh);
564 fbond = -k*(up2 - up1);
565 *dvdlambda += dk*(up2 - up1)*(0.5*(up2 - up1) + drh)
566 + k*(dup2 - dup1)*(up2 - up1 + drh)
567 - k*(up2 - up1)*dup2;
575 vtot += vbond; /* 1*/
576 fbond *= gmx_invsqrt(dr2); /* 6 */
580 fprintf(debug, "BONDS: dr = %10g vbond = %10g fbond = %10g\n",
586 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
589 for (m = 0; (m < DIM); m++) /* 15 */
594 fshift[ki][m] += fij;
595 fshift[CENTRAL][m] -= fij;
602 real polarize(int nbonds,
603 const t_iatom forceatoms[], const t_iparams forceparams[],
604 const rvec x[], rvec f[], rvec fshift[],
605 const t_pbc *pbc, const t_graph *g,
606 real lambda, real *dvdlambda,
607 const t_mdatoms *md, t_fcdata gmx_unused *fcd,
608 int gmx_unused *global_atom_index)
610 int i, m, ki, ai, aj, type;
611 real dr, dr2, fbond, vbond, fij, vtot, ksh;
616 for (i = 0; (i < nbonds); )
618 type = forceatoms[i++];
619 ai = forceatoms[i++];
620 aj = forceatoms[i++];
621 ksh = sqr(md->chargeA[aj])*ONE_4PI_EPS0/forceparams[type].polarize.alpha;
624 fprintf(debug, "POL: local ai = %d aj = %d ksh = %.3f\n", ai, aj, ksh);
627 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
628 dr2 = iprod(dx, dx); /* 5 */
629 dr = dr2*gmx_invsqrt(dr2); /* 10 */
631 *dvdlambda += harmonic(ksh, ksh, 0, 0, dr, lambda, &vbond, &fbond); /* 19 */
638 vtot += vbond; /* 1*/
639 fbond *= gmx_invsqrt(dr2); /* 6 */
643 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
646 for (m = 0; (m < DIM); m++) /* 15 */
651 fshift[ki][m] += fij;
652 fshift[CENTRAL][m] -= fij;
658 real anharm_polarize(int nbonds,
659 const t_iatom forceatoms[], const t_iparams forceparams[],
660 const rvec x[], rvec f[], rvec fshift[],
661 const t_pbc *pbc, const t_graph *g,
662 real lambda, real *dvdlambda,
663 const t_mdatoms *md, t_fcdata gmx_unused *fcd,
664 int gmx_unused *global_atom_index)
666 int i, m, ki, ai, aj, type;
667 real dr, dr2, fbond, vbond, fij, vtot, ksh, khyp, drcut, ddr, ddr3;
672 for (i = 0; (i < nbonds); )
674 type = forceatoms[i++];
675 ai = forceatoms[i++];
676 aj = forceatoms[i++];
677 ksh = sqr(md->chargeA[aj])*ONE_4PI_EPS0/forceparams[type].anharm_polarize.alpha; /* 7*/
678 khyp = forceparams[type].anharm_polarize.khyp;
679 drcut = forceparams[type].anharm_polarize.drcut;
682 fprintf(debug, "POL: local ai = %d aj = %d ksh = %.3f\n", ai, aj, ksh);
685 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
686 dr2 = iprod(dx, dx); /* 5 */
687 dr = dr2*gmx_invsqrt(dr2); /* 10 */
689 *dvdlambda += harmonic(ksh, ksh, 0, 0, dr, lambda, &vbond, &fbond); /* 19 */
700 vbond += khyp*ddr*ddr3;
701 fbond -= 4*khyp*ddr3;
703 fbond *= gmx_invsqrt(dr2); /* 6 */
704 vtot += vbond; /* 1*/
708 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
711 for (m = 0; (m < DIM); m++) /* 15 */
716 fshift[ki][m] += fij;
717 fshift[CENTRAL][m] -= fij;
723 real water_pol(int nbonds,
724 const t_iatom forceatoms[], const t_iparams forceparams[],
725 const rvec x[], rvec f[], rvec gmx_unused fshift[],
726 const t_pbc gmx_unused *pbc, const t_graph gmx_unused *g,
727 real gmx_unused lambda, real gmx_unused *dvdlambda,
728 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
729 int gmx_unused *global_atom_index)
731 /* This routine implements anisotropic polarizibility for water, through
732 * a shell connected to a dummy with spring constant that differ in the
733 * three spatial dimensions in the molecular frame.
735 int i, m, aO, aH1, aH2, aD, aS, type, type0, ki;
737 rvec dOH1, dOH2, dHH, dOD, dDS, nW, kk, dx, kdx, proj;
741 real vtot, fij, r_HH, r_OD, r_nW, tx, ty, tz, qS;
746 type0 = forceatoms[0];
748 qS = md->chargeA[aS];
749 kk[XX] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_x;
750 kk[YY] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_y;
751 kk[ZZ] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_z;
752 r_HH = 1.0/forceparams[type0].wpol.rHH;
753 r_OD = 1.0/forceparams[type0].wpol.rOD;
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, r12n, 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;
920 for (i = 0; (i < nbonds); )
922 type = forceatoms[i++];
923 a1 = forceatoms[i++];
924 da1 = forceatoms[i++];
925 a2 = forceatoms[i++];
926 da2 = forceatoms[i++];
927 q1 = md->chargeA[da1];
928 q2 = md->chargeA[da2];
929 a = forceparams[type].thole.a;
930 al1 = forceparams[type].thole.alpha1;
931 al2 = forceparams[type].thole.alpha2;
933 afac = a*pow(al1*al2, -1.0/6.0);
934 V += do_1_thole(x[a1], x[a2], f[a1], f[a2], pbc, qq, fshift, afac);
935 V += do_1_thole(x[da1], x[a2], f[da1], f[a2], pbc, -qq, fshift, afac);
936 V += do_1_thole(x[a1], x[da2], f[a1], f[da2], pbc, -qq, fshift, afac);
937 V += do_1_thole(x[da1], x[da2], f[da1], f[da2], pbc, qq, fshift, afac);
943 real bond_angle(const rvec xi, const rvec xj, const rvec xk, const t_pbc *pbc,
944 rvec r_ij, rvec r_kj, real *costh,
946 /* Return value is the angle between the bonds i-j and j-k */
951 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
952 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
954 *costh = cos_angle(r_ij, r_kj); /* 25 */
955 th = acos(*costh); /* 10 */
960 real angles(int nbonds,
961 const t_iatom forceatoms[], const t_iparams forceparams[],
962 const rvec x[], rvec f[], rvec fshift[],
963 const t_pbc *pbc, const t_graph *g,
964 real lambda, real *dvdlambda,
965 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
966 int gmx_unused *global_atom_index)
968 int i, ai, aj, ak, t1, t2, type;
970 real cos_theta, cos_theta2, theta, dVdt, va, vtot;
971 ivec jt, dt_ij, dt_kj;
974 for (i = 0; i < nbonds; )
976 type = forceatoms[i++];
977 ai = forceatoms[i++];
978 aj = forceatoms[i++];
979 ak = forceatoms[i++];
981 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
982 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
984 *dvdlambda += harmonic(forceparams[type].harmonic.krA,
985 forceparams[type].harmonic.krB,
986 forceparams[type].harmonic.rA*DEG2RAD,
987 forceparams[type].harmonic.rB*DEG2RAD,
988 theta, lambda, &va, &dVdt); /* 21 */
991 cos_theta2 = sqr(cos_theta);
1001 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1002 sth = st*cos_theta; /* 1 */
1006 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1007 theta*RAD2DEG, va, dVdt);
1010 nrij2 = iprod(r_ij, r_ij); /* 5 */
1011 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1013 nrij_1 = gmx_invsqrt(nrij2); /* 10 */
1014 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1016 cik = st*nrij_1*nrkj_1; /* 2 */
1017 cii = sth*nrij_1*nrij_1; /* 2 */
1018 ckk = sth*nrkj_1*nrkj_1; /* 2 */
1020 for (m = 0; m < DIM; m++)
1022 f_i[m] = -(cik*r_kj[m] - cii*r_ij[m]);
1023 f_k[m] = -(cik*r_ij[m] - ckk*r_kj[m]);
1024 f_j[m] = -f_i[m] - f_k[m];
1031 copy_ivec(SHIFT_IVEC(g, aj), jt);
1033 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1034 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1035 t1 = IVEC2IS(dt_ij);
1036 t2 = IVEC2IS(dt_kj);
1038 rvec_inc(fshift[t1], f_i);
1039 rvec_inc(fshift[CENTRAL], f_j);
1040 rvec_inc(fshift[t2], f_k);
1047 #ifdef GMX_SIMD_HAVE_REAL
1049 /* As angles, but using SIMD to calculate many dihedrals at once.
1050 * This routines does not calculate energies and shift forces.
1052 static gmx_inline void
1053 angles_noener_simd(int nbonds,
1054 const t_iatom forceatoms[], const t_iparams forceparams[],
1055 const rvec x[], rvec f[],
1056 const t_pbc *pbc, const t_graph gmx_unused *g,
1057 real gmx_unused lambda,
1058 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1059 int gmx_unused *global_atom_index)
1063 int type, ai[GMX_SIMD_REAL_WIDTH], aj[GMX_SIMD_REAL_WIDTH];
1064 int ak[GMX_SIMD_REAL_WIDTH];
1065 real coeff_array[2*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *coeff;
1066 real dr_array[2*DIM*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *dr;
1067 real f_buf_array[6*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *f_buf;
1068 gmx_simd_real_t k_S, theta0_S;
1069 gmx_simd_real_t rijx_S, rijy_S, rijz_S;
1070 gmx_simd_real_t rkjx_S, rkjy_S, rkjz_S;
1071 gmx_simd_real_t one_S;
1072 gmx_simd_real_t min_one_plus_eps_S;
1073 gmx_simd_real_t rij_rkj_S;
1074 gmx_simd_real_t nrij2_S, nrij_1_S;
1075 gmx_simd_real_t nrkj2_S, nrkj_1_S;
1076 gmx_simd_real_t cos_S, invsin_S;
1077 gmx_simd_real_t theta_S;
1078 gmx_simd_real_t st_S, sth_S;
1079 gmx_simd_real_t cik_S, cii_S, ckk_S;
1080 gmx_simd_real_t f_ix_S, f_iy_S, f_iz_S;
1081 gmx_simd_real_t f_kx_S, f_ky_S, f_kz_S;
1082 pbc_simd_t pbc_simd;
1084 /* Ensure register memory alignment */
1085 coeff = gmx_simd_align_r(coeff_array);
1086 dr = gmx_simd_align_r(dr_array);
1087 f_buf = gmx_simd_align_r(f_buf_array);
1089 set_pbc_simd(pbc, &pbc_simd);
1091 one_S = gmx_simd_set1_r(1.0);
1093 /* The smallest number > -1 */
1094 min_one_plus_eps_S = gmx_simd_set1_r(-1.0 + 2*GMX_REAL_EPS);
1096 /* nbonds is the number of angles times nfa1, here we step GMX_SIMD_REAL_WIDTH angles */
1097 for (i = 0; (i < nbonds); i += GMX_SIMD_REAL_WIDTH*nfa1)
1099 /* Collect atoms for GMX_SIMD_REAL_WIDTH angles.
1100 * iu indexes into forceatoms, we should not let iu go beyond nbonds.
1103 for (s = 0; s < GMX_SIMD_REAL_WIDTH; s++)
1105 type = forceatoms[iu];
1106 ai[s] = forceatoms[iu+1];
1107 aj[s] = forceatoms[iu+2];
1108 ak[s] = forceatoms[iu+3];
1110 coeff[s] = forceparams[type].harmonic.krA;
1111 coeff[GMX_SIMD_REAL_WIDTH+s] = forceparams[type].harmonic.rA*DEG2RAD;
1113 /* If you can't use pbc_dx_simd below for PBC, e.g. because
1114 * you can't round in SIMD, use pbc_rvec_sub here.
1116 /* Store the non PBC corrected distances packed and aligned */
1117 for (m = 0; m < DIM; m++)
1119 dr[s + m *GMX_SIMD_REAL_WIDTH] = x[ai[s]][m] - x[aj[s]][m];
1120 dr[s + (DIM+m)*GMX_SIMD_REAL_WIDTH] = x[ak[s]][m] - x[aj[s]][m];
1123 /* At the end fill the arrays with identical entries */
1124 if (iu + nfa1 < nbonds)
1130 k_S = gmx_simd_load_r(coeff);
1131 theta0_S = gmx_simd_load_r(coeff+GMX_SIMD_REAL_WIDTH);
1133 rijx_S = gmx_simd_load_r(dr + 0*GMX_SIMD_REAL_WIDTH);
1134 rijy_S = gmx_simd_load_r(dr + 1*GMX_SIMD_REAL_WIDTH);
1135 rijz_S = gmx_simd_load_r(dr + 2*GMX_SIMD_REAL_WIDTH);
1136 rkjx_S = gmx_simd_load_r(dr + 3*GMX_SIMD_REAL_WIDTH);
1137 rkjy_S = gmx_simd_load_r(dr + 4*GMX_SIMD_REAL_WIDTH);
1138 rkjz_S = gmx_simd_load_r(dr + 5*GMX_SIMD_REAL_WIDTH);
1140 pbc_dx_simd(&rijx_S, &rijy_S, &rijz_S, &pbc_simd);
1141 pbc_dx_simd(&rkjx_S, &rkjy_S, &rkjz_S, &pbc_simd);
1143 rij_rkj_S = gmx_simd_iprod_r(rijx_S, rijy_S, rijz_S,
1144 rkjx_S, rkjy_S, rkjz_S);
1146 nrij2_S = gmx_simd_norm2_r(rijx_S, rijy_S, rijz_S);
1147 nrkj2_S = gmx_simd_norm2_r(rkjx_S, rkjy_S, rkjz_S);
1149 nrij_1_S = gmx_simd_invsqrt_r(nrij2_S);
1150 nrkj_1_S = gmx_simd_invsqrt_r(nrkj2_S);
1152 cos_S = gmx_simd_mul_r(rij_rkj_S, gmx_simd_mul_r(nrij_1_S, nrkj_1_S));
1154 /* To allow for 180 degrees, we take the max of cos and -1 + 1bit,
1155 * so we can safely get the 1/sin from 1/sqrt(1 - cos^2).
1156 * This also ensures that rounding errors would cause the argument
1157 * of gmx_simd_acos_r to be < -1.
1158 * Note that we do not take precautions for cos(0)=1, so the outer
1159 * atoms in an angle should not be on top of each other.
1161 cos_S = gmx_simd_max_r(cos_S, min_one_plus_eps_S);
1163 theta_S = gmx_simd_acos_r(cos_S);
1165 invsin_S = gmx_simd_invsqrt_r(gmx_simd_sub_r(one_S, gmx_simd_mul_r(cos_S, cos_S)));
1167 st_S = gmx_simd_mul_r(gmx_simd_mul_r(k_S, gmx_simd_sub_r(theta0_S, theta_S)),
1169 sth_S = gmx_simd_mul_r(st_S, cos_S);
1171 cik_S = gmx_simd_mul_r(st_S, gmx_simd_mul_r(nrij_1_S, nrkj_1_S));
1172 cii_S = gmx_simd_mul_r(sth_S, gmx_simd_mul_r(nrij_1_S, nrij_1_S));
1173 ckk_S = gmx_simd_mul_r(sth_S, gmx_simd_mul_r(nrkj_1_S, nrkj_1_S));
1175 f_ix_S = gmx_simd_mul_r(cii_S, rijx_S);
1176 f_ix_S = gmx_simd_fnmadd_r(cik_S, rkjx_S, f_ix_S);
1177 f_iy_S = gmx_simd_mul_r(cii_S, rijy_S);
1178 f_iy_S = gmx_simd_fnmadd_r(cik_S, rkjy_S, f_iy_S);
1179 f_iz_S = gmx_simd_mul_r(cii_S, rijz_S);
1180 f_iz_S = gmx_simd_fnmadd_r(cik_S, rkjz_S, f_iz_S);
1181 f_kx_S = gmx_simd_mul_r(ckk_S, rkjx_S);
1182 f_kx_S = gmx_simd_fnmadd_r(cik_S, rijx_S, f_kx_S);
1183 f_ky_S = gmx_simd_mul_r(ckk_S, rkjy_S);
1184 f_ky_S = gmx_simd_fnmadd_r(cik_S, rijy_S, f_ky_S);
1185 f_kz_S = gmx_simd_mul_r(ckk_S, rkjz_S);
1186 f_kz_S = gmx_simd_fnmadd_r(cik_S, rijz_S, f_kz_S);
1188 gmx_simd_store_r(f_buf + 0*GMX_SIMD_REAL_WIDTH, f_ix_S);
1189 gmx_simd_store_r(f_buf + 1*GMX_SIMD_REAL_WIDTH, f_iy_S);
1190 gmx_simd_store_r(f_buf + 2*GMX_SIMD_REAL_WIDTH, f_iz_S);
1191 gmx_simd_store_r(f_buf + 3*GMX_SIMD_REAL_WIDTH, f_kx_S);
1192 gmx_simd_store_r(f_buf + 4*GMX_SIMD_REAL_WIDTH, f_ky_S);
1193 gmx_simd_store_r(f_buf + 5*GMX_SIMD_REAL_WIDTH, f_kz_S);
1199 for (m = 0; m < DIM; m++)
1201 f[ai[s]][m] += f_buf[s + m*GMX_SIMD_REAL_WIDTH];
1202 f[aj[s]][m] -= f_buf[s + m*GMX_SIMD_REAL_WIDTH] + f_buf[s + (DIM+m)*GMX_SIMD_REAL_WIDTH];
1203 f[ak[s]][m] += f_buf[s + (DIM+m)*GMX_SIMD_REAL_WIDTH];
1208 while (s < GMX_SIMD_REAL_WIDTH && iu < nbonds);
1212 #endif /* GMX_SIMD_HAVE_REAL */
1214 real linear_angles(int nbonds,
1215 const t_iatom forceatoms[], const t_iparams forceparams[],
1216 const rvec x[], rvec f[], rvec fshift[],
1217 const t_pbc *pbc, const t_graph *g,
1218 real lambda, real *dvdlambda,
1219 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1220 int gmx_unused *global_atom_index)
1222 int i, m, ai, aj, ak, t1, t2, type;
1224 real L1, kA, kB, aA, aB, dr, dr2, va, vtot, a, b, klin;
1225 ivec jt, dt_ij, dt_kj;
1226 rvec r_ij, r_kj, r_ik, dx;
1230 for (i = 0; (i < nbonds); )
1232 type = forceatoms[i++];
1233 ai = forceatoms[i++];
1234 aj = forceatoms[i++];
1235 ak = forceatoms[i++];
1237 kA = forceparams[type].linangle.klinA;
1238 kB = forceparams[type].linangle.klinB;
1239 klin = L1*kA + lambda*kB;
1241 aA = forceparams[type].linangle.aA;
1242 aB = forceparams[type].linangle.aB;
1243 a = L1*aA+lambda*aB;
1246 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
1247 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
1248 rvec_sub(r_ij, r_kj, r_ik);
1251 for (m = 0; (m < DIM); m++)
1253 dr = -a * r_ij[m] - b * r_kj[m];
1258 f_j[m] = -(f_i[m]+f_k[m]);
1264 *dvdlambda += 0.5*(kB-kA)*dr2 + klin*(aB-aA)*iprod(dx, r_ik);
1270 copy_ivec(SHIFT_IVEC(g, aj), jt);
1272 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1273 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1274 t1 = IVEC2IS(dt_ij);
1275 t2 = IVEC2IS(dt_kj);
1277 rvec_inc(fshift[t1], f_i);
1278 rvec_inc(fshift[CENTRAL], f_j);
1279 rvec_inc(fshift[t2], f_k);
1284 real urey_bradley(int nbonds,
1285 const t_iatom forceatoms[], const t_iparams forceparams[],
1286 const rvec x[], rvec f[], rvec fshift[],
1287 const t_pbc *pbc, const t_graph *g,
1288 real lambda, real *dvdlambda,
1289 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1290 int gmx_unused *global_atom_index)
1292 int i, m, ai, aj, ak, t1, t2, type, ki;
1293 rvec r_ij, r_kj, r_ik;
1294 real cos_theta, cos_theta2, theta;
1295 real dVdt, va, vtot, dr, dr2, vbond, fbond, fik;
1296 real kthA, th0A, kUBA, r13A, kthB, th0B, kUBB, r13B;
1297 ivec jt, dt_ij, dt_kj, dt_ik;
1300 for (i = 0; (i < nbonds); )
1302 type = forceatoms[i++];
1303 ai = forceatoms[i++];
1304 aj = forceatoms[i++];
1305 ak = forceatoms[i++];
1306 th0A = forceparams[type].u_b.thetaA*DEG2RAD;
1307 kthA = forceparams[type].u_b.kthetaA;
1308 r13A = forceparams[type].u_b.r13A;
1309 kUBA = forceparams[type].u_b.kUBA;
1310 th0B = forceparams[type].u_b.thetaB*DEG2RAD;
1311 kthB = forceparams[type].u_b.kthetaB;
1312 r13B = forceparams[type].u_b.r13B;
1313 kUBB = forceparams[type].u_b.kUBB;
1315 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
1316 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
1318 *dvdlambda += harmonic(kthA, kthB, th0A, th0B, theta, lambda, &va, &dVdt); /* 21 */
1321 ki = pbc_rvec_sub(pbc, x[ai], x[ak], r_ik); /* 3 */
1322 dr2 = iprod(r_ik, r_ik); /* 5 */
1323 dr = dr2*gmx_invsqrt(dr2); /* 10 */
1325 *dvdlambda += harmonic(kUBA, kUBB, r13A, r13B, dr, lambda, &vbond, &fbond); /* 19 */
1327 cos_theta2 = sqr(cos_theta); /* 1 */
1335 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1336 sth = st*cos_theta; /* 1 */
1340 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1341 theta*RAD2DEG, va, dVdt);
1344 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1345 nrij2 = iprod(r_ij, r_ij);
1347 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
1348 cii = sth/nrij2; /* 10 */
1349 ckk = sth/nrkj2; /* 10 */
1351 for (m = 0; (m < DIM); m++) /* 39 */
1353 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
1354 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
1355 f_j[m] = -f_i[m]-f_k[m];
1362 copy_ivec(SHIFT_IVEC(g, aj), jt);
1364 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1365 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1366 t1 = IVEC2IS(dt_ij);
1367 t2 = IVEC2IS(dt_kj);
1369 rvec_inc(fshift[t1], f_i);
1370 rvec_inc(fshift[CENTRAL], f_j);
1371 rvec_inc(fshift[t2], f_k);
1373 /* Time for the bond calculations */
1379 vtot += vbond; /* 1*/
1380 fbond *= gmx_invsqrt(dr2); /* 6 */
1384 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ak), dt_ik);
1385 ki = IVEC2IS(dt_ik);
1387 for (m = 0; (m < DIM); m++) /* 15 */
1389 fik = fbond*r_ik[m];
1392 fshift[ki][m] += fik;
1393 fshift[CENTRAL][m] -= fik;
1399 real quartic_angles(int nbonds,
1400 const t_iatom forceatoms[], const t_iparams forceparams[],
1401 const rvec x[], rvec f[], rvec fshift[],
1402 const t_pbc *pbc, const t_graph *g,
1403 real gmx_unused lambda, real gmx_unused *dvdlambda,
1404 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1405 int gmx_unused *global_atom_index)
1407 int i, j, ai, aj, ak, t1, t2, type;
1409 real cos_theta, cos_theta2, theta, dt, dVdt, va, dtp, c, vtot;
1410 ivec jt, dt_ij, dt_kj;
1413 for (i = 0; (i < nbonds); )
1415 type = forceatoms[i++];
1416 ai = forceatoms[i++];
1417 aj = forceatoms[i++];
1418 ak = forceatoms[i++];
1420 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
1421 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
1423 dt = theta - forceparams[type].qangle.theta*DEG2RAD; /* 2 */
1426 va = forceparams[type].qangle.c[0];
1428 for (j = 1; j <= 4; j++)
1430 c = forceparams[type].qangle.c[j];
1439 cos_theta2 = sqr(cos_theta); /* 1 */
1448 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1449 sth = st*cos_theta; /* 1 */
1453 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1454 theta*RAD2DEG, va, dVdt);
1457 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1458 nrij2 = iprod(r_ij, r_ij);
1460 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
1461 cii = sth/nrij2; /* 10 */
1462 ckk = sth/nrkj2; /* 10 */
1464 for (m = 0; (m < DIM); m++) /* 39 */
1466 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
1467 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
1468 f_j[m] = -f_i[m]-f_k[m];
1475 copy_ivec(SHIFT_IVEC(g, aj), jt);
1477 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1478 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1479 t1 = IVEC2IS(dt_ij);
1480 t2 = IVEC2IS(dt_kj);
1482 rvec_inc(fshift[t1], f_i);
1483 rvec_inc(fshift[CENTRAL], f_j);
1484 rvec_inc(fshift[t2], f_k);
1490 real dih_angle(const rvec xi, const rvec xj, const rvec xk, const rvec xl,
1492 rvec r_ij, rvec r_kj, rvec r_kl, rvec m, rvec n,
1493 real *sign, int *t1, int *t2, int *t3)
1497 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
1498 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
1499 *t3 = pbc_rvec_sub(pbc, xk, xl, r_kl); /* 3 */
1501 cprod(r_ij, r_kj, m); /* 9 */
1502 cprod(r_kj, r_kl, n); /* 9 */
1503 phi = gmx_angle(m, n); /* 49 (assuming 25 for atan2) */
1504 ipr = iprod(r_ij, n); /* 5 */
1505 (*sign) = (ipr < 0.0) ? -1.0 : 1.0;
1506 phi = (*sign)*phi; /* 1 */
1512 #ifdef GMX_SIMD_HAVE_REAL
1514 /* As dih_angle above, but calculates 4 dihedral angles at once using SIMD,
1515 * also calculates the pre-factor required for the dihedral force update.
1516 * Note that bv and buf should be register aligned.
1518 static gmx_inline void
1519 dih_angle_simd(const rvec *x,
1520 const int *ai, const int *aj, const int *ak, const int *al,
1521 const pbc_simd_t *pbc,
1523 gmx_simd_real_t *phi_S,
1524 gmx_simd_real_t *mx_S, gmx_simd_real_t *my_S, gmx_simd_real_t *mz_S,
1525 gmx_simd_real_t *nx_S, gmx_simd_real_t *ny_S, gmx_simd_real_t *nz_S,
1526 gmx_simd_real_t *nrkj_m2_S,
1527 gmx_simd_real_t *nrkj_n2_S,
1532 gmx_simd_real_t rijx_S, rijy_S, rijz_S;
1533 gmx_simd_real_t rkjx_S, rkjy_S, rkjz_S;
1534 gmx_simd_real_t rklx_S, rkly_S, rklz_S;
1535 gmx_simd_real_t cx_S, cy_S, cz_S;
1536 gmx_simd_real_t cn_S;
1537 gmx_simd_real_t s_S;
1538 gmx_simd_real_t ipr_S;
1539 gmx_simd_real_t iprm_S, iprn_S;
1540 gmx_simd_real_t nrkj2_S, nrkj_1_S, nrkj_2_S, nrkj_S;
1541 gmx_simd_real_t toler_S;
1542 gmx_simd_real_t p_S, q_S;
1543 gmx_simd_real_t nrkj2_min_S;
1544 gmx_simd_real_t real_eps_S;
1546 /* Used to avoid division by zero.
1547 * We take into acount that we multiply the result by real_eps_S.
1549 nrkj2_min_S = gmx_simd_set1_r(GMX_REAL_MIN/(2*GMX_REAL_EPS));
1551 /* The value of the last significant bit (GMX_REAL_EPS is half of that) */
1552 real_eps_S = gmx_simd_set1_r(2*GMX_REAL_EPS);
1554 for (s = 0; s < GMX_SIMD_REAL_WIDTH; s++)
1556 /* If you can't use pbc_dx_simd below for PBC, e.g. because
1557 * you can't round in SIMD, use pbc_rvec_sub here.
1559 for (m = 0; m < DIM; m++)
1561 dr[s + (0*DIM + m)*GMX_SIMD_REAL_WIDTH] = x[ai[s]][m] - x[aj[s]][m];
1562 dr[s + (1*DIM + m)*GMX_SIMD_REAL_WIDTH] = x[ak[s]][m] - x[aj[s]][m];
1563 dr[s + (2*DIM + m)*GMX_SIMD_REAL_WIDTH] = x[ak[s]][m] - x[al[s]][m];
1567 rijx_S = gmx_simd_load_r(dr + 0*GMX_SIMD_REAL_WIDTH);
1568 rijy_S = gmx_simd_load_r(dr + 1*GMX_SIMD_REAL_WIDTH);
1569 rijz_S = gmx_simd_load_r(dr + 2*GMX_SIMD_REAL_WIDTH);
1570 rkjx_S = gmx_simd_load_r(dr + 3*GMX_SIMD_REAL_WIDTH);
1571 rkjy_S = gmx_simd_load_r(dr + 4*GMX_SIMD_REAL_WIDTH);
1572 rkjz_S = gmx_simd_load_r(dr + 5*GMX_SIMD_REAL_WIDTH);
1573 rklx_S = gmx_simd_load_r(dr + 6*GMX_SIMD_REAL_WIDTH);
1574 rkly_S = gmx_simd_load_r(dr + 7*GMX_SIMD_REAL_WIDTH);
1575 rklz_S = gmx_simd_load_r(dr + 8*GMX_SIMD_REAL_WIDTH);
1577 pbc_dx_simd(&rijx_S, &rijy_S, &rijz_S, pbc);
1578 pbc_dx_simd(&rkjx_S, &rkjy_S, &rkjz_S, pbc);
1579 pbc_dx_simd(&rklx_S, &rkly_S, &rklz_S, pbc);
1581 gmx_simd_cprod_r(rijx_S, rijy_S, rijz_S,
1582 rkjx_S, rkjy_S, rkjz_S,
1585 gmx_simd_cprod_r(rkjx_S, rkjy_S, rkjz_S,
1586 rklx_S, rkly_S, rklz_S,
1589 gmx_simd_cprod_r(*mx_S, *my_S, *mz_S,
1590 *nx_S, *ny_S, *nz_S,
1591 &cx_S, &cy_S, &cz_S);
1593 cn_S = gmx_simd_sqrt_r(gmx_simd_norm2_r(cx_S, cy_S, cz_S));
1595 s_S = gmx_simd_iprod_r(*mx_S, *my_S, *mz_S, *nx_S, *ny_S, *nz_S);
1597 /* Determine the dihedral angle, the sign might need correction */
1598 *phi_S = gmx_simd_atan2_r(cn_S, s_S);
1600 ipr_S = gmx_simd_iprod_r(rijx_S, rijy_S, rijz_S,
1601 *nx_S, *ny_S, *nz_S);
1603 iprm_S = gmx_simd_norm2_r(*mx_S, *my_S, *mz_S);
1604 iprn_S = gmx_simd_norm2_r(*nx_S, *ny_S, *nz_S);
1606 nrkj2_S = gmx_simd_norm2_r(rkjx_S, rkjy_S, rkjz_S);
1608 /* Avoid division by zero. When zero, the result is multiplied by 0
1609 * anyhow, so the 3 max below do not affect the final result.
1611 nrkj2_S = gmx_simd_max_r(nrkj2_S, nrkj2_min_S);
1612 nrkj_1_S = gmx_simd_invsqrt_r(nrkj2_S);
1613 nrkj_2_S = gmx_simd_mul_r(nrkj_1_S, nrkj_1_S);
1614 nrkj_S = gmx_simd_mul_r(nrkj2_S, nrkj_1_S);
1616 toler_S = gmx_simd_mul_r(nrkj2_S, real_eps_S);
1618 /* Here the plain-C code uses a conditional, but we can't do that in SIMD.
1619 * So we take a max with the tolerance instead. Since we multiply with
1620 * m or n later, the max does not affect the results.
1622 iprm_S = gmx_simd_max_r(iprm_S, toler_S);
1623 iprn_S = gmx_simd_max_r(iprn_S, toler_S);
1624 *nrkj_m2_S = gmx_simd_mul_r(nrkj_S, gmx_simd_inv_r(iprm_S));
1625 *nrkj_n2_S = gmx_simd_mul_r(nrkj_S, gmx_simd_inv_r(iprn_S));
1627 /* Set sign of phi_S with the sign of ipr_S; phi_S is currently positive */
1628 *phi_S = gmx_simd_xor_sign_r(*phi_S, ipr_S);
1629 p_S = gmx_simd_iprod_r(rijx_S, rijy_S, rijz_S,
1630 rkjx_S, rkjy_S, rkjz_S);
1631 p_S = gmx_simd_mul_r(p_S, nrkj_2_S);
1633 q_S = gmx_simd_iprod_r(rklx_S, rkly_S, rklz_S,
1634 rkjx_S, rkjy_S, rkjz_S);
1635 q_S = gmx_simd_mul_r(q_S, nrkj_2_S);
1637 gmx_simd_store_r(p, p_S);
1638 gmx_simd_store_r(q, q_S);
1641 #endif /* GMX_SIMD_HAVE_REAL */
1644 void do_dih_fup(int i, int j, int k, int l, real ddphi,
1645 rvec r_ij, rvec r_kj, rvec r_kl,
1646 rvec m, rvec n, rvec f[], rvec fshift[],
1647 const t_pbc *pbc, const t_graph *g,
1648 const rvec x[], int t1, int t2, int t3)
1651 rvec f_i, f_j, f_k, f_l;
1652 rvec uvec, vvec, svec, dx_jl;
1653 real iprm, iprn, nrkj, nrkj2, nrkj_1, nrkj_2;
1654 real a, b, p, q, toler;
1655 ivec jt, dt_ij, dt_kj, dt_lj;
1657 iprm = iprod(m, m); /* 5 */
1658 iprn = iprod(n, n); /* 5 */
1659 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1660 toler = nrkj2*GMX_REAL_EPS;
1661 if ((iprm > toler) && (iprn > toler))
1663 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1664 nrkj_2 = nrkj_1*nrkj_1; /* 1 */
1665 nrkj = nrkj2*nrkj_1; /* 1 */
1666 a = -ddphi*nrkj/iprm; /* 11 */
1667 svmul(a, m, f_i); /* 3 */
1668 b = ddphi*nrkj/iprn; /* 11 */
1669 svmul(b, n, f_l); /* 3 */
1670 p = iprod(r_ij, r_kj); /* 5 */
1671 p *= nrkj_2; /* 1 */
1672 q = iprod(r_kl, r_kj); /* 5 */
1673 q *= nrkj_2; /* 1 */
1674 svmul(p, f_i, uvec); /* 3 */
1675 svmul(q, f_l, vvec); /* 3 */
1676 rvec_sub(uvec, vvec, svec); /* 3 */
1677 rvec_sub(f_i, svec, f_j); /* 3 */
1678 rvec_add(f_l, svec, f_k); /* 3 */
1679 rvec_inc(f[i], f_i); /* 3 */
1680 rvec_dec(f[j], f_j); /* 3 */
1681 rvec_dec(f[k], f_k); /* 3 */
1682 rvec_inc(f[l], f_l); /* 3 */
1686 copy_ivec(SHIFT_IVEC(g, j), jt);
1687 ivec_sub(SHIFT_IVEC(g, i), jt, dt_ij);
1688 ivec_sub(SHIFT_IVEC(g, k), jt, dt_kj);
1689 ivec_sub(SHIFT_IVEC(g, l), jt, dt_lj);
1690 t1 = IVEC2IS(dt_ij);
1691 t2 = IVEC2IS(dt_kj);
1692 t3 = IVEC2IS(dt_lj);
1696 t3 = pbc_rvec_sub(pbc, x[l], x[j], dx_jl);
1703 rvec_inc(fshift[t1], f_i);
1704 rvec_dec(fshift[CENTRAL], f_j);
1705 rvec_dec(fshift[t2], f_k);
1706 rvec_inc(fshift[t3], f_l);
1711 /* As do_dih_fup above, but without shift forces */
1713 do_dih_fup_noshiftf(int i, int j, int k, int l, real ddphi,
1714 rvec r_ij, rvec r_kj, rvec r_kl,
1715 rvec m, rvec n, rvec f[])
1717 rvec f_i, f_j, f_k, f_l;
1718 rvec uvec, vvec, svec, dx_jl;
1719 real iprm, iprn, nrkj, nrkj2, nrkj_1, nrkj_2;
1720 real a, b, p, q, toler;
1721 ivec jt, dt_ij, dt_kj, dt_lj;
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 int t1[GMX_SIMD_REAL_WIDTH], t2[GMX_SIMD_REAL_WIDTH], t3[GMX_SIMD_REAL_WIDTH];
1993 real dr_array[3*DIM*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *dr;
1994 real buf_array[7*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *buf;
1995 real *cp, *phi0, *mult, *phi, *p, *q, *sf_i, *msf_l;
1996 gmx_simd_real_t phi0_S, phi_S;
1997 gmx_simd_real_t mx_S, my_S, mz_S;
1998 gmx_simd_real_t nx_S, ny_S, nz_S;
1999 gmx_simd_real_t nrkj_m2_S, nrkj_n2_S;
2000 gmx_simd_real_t cp_S, mdphi_S, mult_S;
2001 gmx_simd_real_t sin_S, cos_S;
2002 gmx_simd_real_t mddphi_S;
2003 gmx_simd_real_t sf_i_S, msf_l_S;
2004 pbc_simd_t pbc_simd;
2006 /* Ensure SIMD register alignment */
2007 dr = gmx_simd_align_r(dr_array);
2008 buf = gmx_simd_align_r(buf_array);
2010 /* Extract aligned pointer for parameters and variables */
2011 cp = buf + 0*GMX_SIMD_REAL_WIDTH;
2012 phi0 = buf + 1*GMX_SIMD_REAL_WIDTH;
2013 mult = buf + 2*GMX_SIMD_REAL_WIDTH;
2014 p = buf + 3*GMX_SIMD_REAL_WIDTH;
2015 q = buf + 4*GMX_SIMD_REAL_WIDTH;
2016 sf_i = buf + 5*GMX_SIMD_REAL_WIDTH;
2017 msf_l = buf + 6*GMX_SIMD_REAL_WIDTH;
2019 set_pbc_simd(pbc, &pbc_simd);
2021 /* nbonds is the number of dihedrals times nfa1, here we step GMX_SIMD_REAL_WIDTH dihs */
2022 for (i = 0; (i < nbonds); i += GMX_SIMD_REAL_WIDTH*nfa1)
2024 /* Collect atoms quadruplets for GMX_SIMD_REAL_WIDTH dihedrals.
2025 * iu indexes into forceatoms, we should not let iu go beyond nbonds.
2028 for (s = 0; s < GMX_SIMD_REAL_WIDTH; s++)
2030 type = forceatoms[iu];
2031 ai[s] = forceatoms[iu+1];
2032 aj[s] = forceatoms[iu+2];
2033 ak[s] = forceatoms[iu+3];
2034 al[s] = forceatoms[iu+4];
2036 cp[s] = forceparams[type].pdihs.cpA;
2037 phi0[s] = forceparams[type].pdihs.phiA*DEG2RAD;
2038 mult[s] = forceparams[type].pdihs.mult;
2040 /* At the end fill the arrays with identical entries */
2041 if (iu + nfa1 < nbonds)
2047 /* Caclulate GMX_SIMD_REAL_WIDTH dihedral angles at once */
2048 dih_angle_simd(x, ai, aj, ak, al, &pbc_simd,
2051 &mx_S, &my_S, &mz_S,
2052 &nx_S, &ny_S, &nz_S,
2057 cp_S = gmx_simd_load_r(cp);
2058 phi0_S = gmx_simd_load_r(phi0);
2059 mult_S = gmx_simd_load_r(mult);
2061 mdphi_S = gmx_simd_sub_r(gmx_simd_mul_r(mult_S, phi_S), phi0_S);
2063 /* Calculate GMX_SIMD_REAL_WIDTH sines at once */
2064 gmx_simd_sincos_r(mdphi_S, &sin_S, &cos_S);
2065 mddphi_S = gmx_simd_mul_r(gmx_simd_mul_r(cp_S, mult_S), sin_S);
2066 sf_i_S = gmx_simd_mul_r(mddphi_S, nrkj_m2_S);
2067 msf_l_S = gmx_simd_mul_r(mddphi_S, nrkj_n2_S);
2069 /* After this m?_S will contain f[i] */
2070 mx_S = gmx_simd_mul_r(sf_i_S, mx_S);
2071 my_S = gmx_simd_mul_r(sf_i_S, my_S);
2072 mz_S = gmx_simd_mul_r(sf_i_S, mz_S);
2074 /* After this m?_S will contain -f[l] */
2075 nx_S = gmx_simd_mul_r(msf_l_S, nx_S);
2076 ny_S = gmx_simd_mul_r(msf_l_S, ny_S);
2077 nz_S = gmx_simd_mul_r(msf_l_S, nz_S);
2079 gmx_simd_store_r(dr + 0*GMX_SIMD_REAL_WIDTH, mx_S);
2080 gmx_simd_store_r(dr + 1*GMX_SIMD_REAL_WIDTH, my_S);
2081 gmx_simd_store_r(dr + 2*GMX_SIMD_REAL_WIDTH, mz_S);
2082 gmx_simd_store_r(dr + 3*GMX_SIMD_REAL_WIDTH, nx_S);
2083 gmx_simd_store_r(dr + 4*GMX_SIMD_REAL_WIDTH, ny_S);
2084 gmx_simd_store_r(dr + 5*GMX_SIMD_REAL_WIDTH, nz_S);
2090 do_dih_fup_noshiftf_precalc(ai[s], aj[s], ak[s], al[s],
2092 dr[ XX *GMX_SIMD_REAL_WIDTH+s],
2093 dr[ YY *GMX_SIMD_REAL_WIDTH+s],
2094 dr[ ZZ *GMX_SIMD_REAL_WIDTH+s],
2095 dr[(DIM+XX)*GMX_SIMD_REAL_WIDTH+s],
2096 dr[(DIM+YY)*GMX_SIMD_REAL_WIDTH+s],
2097 dr[(DIM+ZZ)*GMX_SIMD_REAL_WIDTH+s],
2102 while (s < GMX_SIMD_REAL_WIDTH && iu < nbonds);
2106 #endif /* GMX_SIMD_HAVE_REAL */
2109 real idihs(int nbonds,
2110 const t_iatom forceatoms[], const t_iparams forceparams[],
2111 const rvec x[], rvec f[], rvec fshift[],
2112 const t_pbc *pbc, const t_graph *g,
2113 real lambda, real *dvdlambda,
2114 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2115 int gmx_unused *global_atom_index)
2117 int i, type, ai, aj, ak, al;
2119 real phi, phi0, dphi0, ddphi, sign, vtot;
2120 rvec r_ij, r_kj, r_kl, m, n;
2121 real L1, kk, dp, dp2, kA, kB, pA, pB, dvdl_term;
2126 for (i = 0; (i < nbonds); )
2128 type = forceatoms[i++];
2129 ai = forceatoms[i++];
2130 aj = forceatoms[i++];
2131 ak = forceatoms[i++];
2132 al = forceatoms[i++];
2134 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2135 &sign, &t1, &t2, &t3); /* 84 */
2137 /* phi can jump if phi0 is close to Pi/-Pi, which will cause huge
2138 * force changes if we just apply a normal harmonic.
2139 * Instead, we first calculate phi-phi0 and take it modulo (-Pi,Pi).
2140 * This means we will never have the periodicity problem, unless
2141 * the dihedral is Pi away from phiO, which is very unlikely due to
2144 kA = forceparams[type].harmonic.krA;
2145 kB = forceparams[type].harmonic.krB;
2146 pA = forceparams[type].harmonic.rA;
2147 pB = forceparams[type].harmonic.rB;
2149 kk = L1*kA + lambda*kB;
2150 phi0 = (L1*pA + lambda*pB)*DEG2RAD;
2151 dphi0 = (pB - pA)*DEG2RAD;
2155 make_dp_periodic(&dp);
2162 dvdl_term += 0.5*(kB - kA)*dp2 - kk*dphi0*dp;
2164 do_dih_fup(ai, aj, ak, al, (real)(-ddphi), r_ij, r_kj, r_kl, m, n,
2165 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2170 fprintf(debug, "idih: (%d,%d,%d,%d) phi=%g\n",
2171 ai, aj, ak, al, phi);
2176 *dvdlambda += dvdl_term;
2181 /*! \brief returns dx, rdist, and dpdl for functions posres() and fbposres()
2183 static void posres_dx(const rvec x, const rvec pos0A, const rvec pos0B,
2184 const rvec comA_sc, const rvec comB_sc,
2186 t_pbc *pbc, int refcoord_scaling, int npbcdim,
2187 rvec dx, rvec rdist, rvec dpdl)
2190 real posA, posB, L1, ref = 0.;
2195 for (m = 0; m < DIM; m++)
2201 switch (refcoord_scaling)
2205 rdist[m] = L1*posA + lambda*posB;
2206 dpdl[m] = posB - posA;
2209 /* Box relative coordinates are stored for dimensions with pbc */
2210 posA *= pbc->box[m][m];
2211 posB *= pbc->box[m][m];
2212 assert(npbcdim <= DIM);
2213 for (d = m+1; d < npbcdim; d++)
2215 posA += pos0A[d]*pbc->box[d][m];
2216 posB += pos0B[d]*pbc->box[d][m];
2218 ref = L1*posA + lambda*posB;
2220 dpdl[m] = posB - posA;
2223 ref = L1*comA_sc[m] + lambda*comB_sc[m];
2224 rdist[m] = L1*posA + lambda*posB;
2225 dpdl[m] = comB_sc[m] - comA_sc[m] + posB - posA;
2228 gmx_fatal(FARGS, "No such scaling method implemented");
2233 ref = L1*posA + lambda*posB;
2235 dpdl[m] = posB - posA;
2238 /* We do pbc_dx with ref+rdist,
2239 * since with only ref we can be up to half a box vector wrong.
2241 pos[m] = ref + rdist[m];
2246 pbc_dx(pbc, x, pos, dx);
2250 rvec_sub(x, pos, dx);
2254 /*! \brief Adds forces of flat-bottomed positions restraints to f[]
2255 * and fixes vir_diag. Returns the flat-bottomed potential. */
2256 real fbposres(int nbonds,
2257 const t_iatom forceatoms[], const t_iparams forceparams[],
2258 const rvec x[], rvec f[], rvec vir_diag,
2260 int refcoord_scaling, int ePBC, rvec com)
2261 /* compute flat-bottomed positions restraints */
2263 int i, ai, m, d, type, npbcdim = 0, fbdim;
2264 const t_iparams *pr;
2266 real ref = 0, dr, dr2, rpot, rfb, rfb2, fact, invdr;
2267 rvec com_sc, rdist, pos, dx, dpdl, fm;
2270 npbcdim = ePBC2npbcdim(ePBC);
2272 if (refcoord_scaling == erscCOM)
2275 for (m = 0; m < npbcdim; m++)
2277 assert(npbcdim <= DIM);
2278 for (d = m; d < npbcdim; d++)
2280 com_sc[m] += com[d]*pbc->box[d][m];
2286 for (i = 0; (i < nbonds); )
2288 type = forceatoms[i++];
2289 ai = forceatoms[i++];
2290 pr = &forceparams[type];
2292 /* same calculation as for normal posres, but with identical A and B states, and lambda==0 */
2293 posres_dx(x[ai], forceparams[type].fbposres.pos0, forceparams[type].fbposres.pos0,
2294 com_sc, com_sc, 0.0,
2295 pbc, refcoord_scaling, npbcdim,
2301 kk = pr->fbposres.k;
2302 rfb = pr->fbposres.r;
2305 /* with rfb<0, push particle out of the sphere/cylinder/layer */
2313 switch (pr->fbposres.geom)
2315 case efbposresSPHERE:
2316 /* spherical flat-bottom posres */
2319 ( (dr2 > rfb2 && bInvert == FALSE ) || (dr2 < rfb2 && bInvert == TRUE ) )
2323 v = 0.5*kk*sqr(dr - rfb);
2324 fact = -kk*(dr-rfb)/dr; /* Force pointing to the center pos0 */
2325 svmul(fact, dx, fm);
2328 case efbposresCYLINDER:
2329 /* cylidrical flat-bottom posres in x-y plane. fm[ZZ] = 0. */
2330 dr2 = sqr(dx[XX])+sqr(dx[YY]);
2332 ( (dr2 > rfb2 && bInvert == FALSE ) || (dr2 < rfb2 && bInvert == TRUE ) )
2337 v = 0.5*kk*sqr(dr - rfb);
2338 fm[XX] = -kk*(dr-rfb)*dx[XX]*invdr; /* Force pointing to the center */
2339 fm[YY] = -kk*(dr-rfb)*dx[YY]*invdr;
2342 case efbposresX: /* fbdim=XX */
2343 case efbposresY: /* fbdim=YY */
2344 case efbposresZ: /* fbdim=ZZ */
2345 /* 1D flat-bottom potential */
2346 fbdim = pr->fbposres.geom - efbposresX;
2348 if ( ( dr > rfb && bInvert == FALSE ) || ( 0 < dr && dr < rfb && bInvert == TRUE ) )
2350 v = 0.5*kk*sqr(dr - rfb);
2351 fm[fbdim] = -kk*(dr - rfb);
2353 else if ( (dr < (-rfb) && bInvert == FALSE ) || ( (-rfb) < dr && dr < 0 && bInvert == TRUE ))
2355 v = 0.5*kk*sqr(dr + rfb);
2356 fm[fbdim] = -kk*(dr + rfb);
2363 for (m = 0; (m < DIM); m++)
2366 /* Here we correct for the pbc_dx which included rdist */
2367 vir_diag[m] -= 0.5*(dx[m] + rdist[m])*fm[m];
2375 real posres(int nbonds,
2376 const t_iatom forceatoms[], const t_iparams forceparams[],
2377 const rvec x[], rvec f[], rvec vir_diag,
2379 real lambda, real *dvdlambda,
2380 int refcoord_scaling, int ePBC, rvec comA, rvec comB)
2382 int i, ai, m, d, type, ki, npbcdim = 0;
2383 const t_iparams *pr;
2386 real posA, posB, ref = 0;
2387 rvec comA_sc, comB_sc, rdist, dpdl, pos, dx;
2388 gmx_bool bForceValid = TRUE;
2390 if ((f == NULL) || (vir_diag == NULL)) /* should both be null together! */
2392 bForceValid = FALSE;
2395 npbcdim = ePBC2npbcdim(ePBC);
2397 if (refcoord_scaling == erscCOM)
2399 clear_rvec(comA_sc);
2400 clear_rvec(comB_sc);
2401 for (m = 0; m < npbcdim; m++)
2403 assert(npbcdim <= DIM);
2404 for (d = m; d < npbcdim; d++)
2406 comA_sc[m] += comA[d]*pbc->box[d][m];
2407 comB_sc[m] += comB[d]*pbc->box[d][m];
2415 for (i = 0; (i < nbonds); )
2417 type = forceatoms[i++];
2418 ai = forceatoms[i++];
2419 pr = &forceparams[type];
2421 /* return dx, rdist, and dpdl */
2422 posres_dx(x[ai], forceparams[type].posres.pos0A, forceparams[type].posres.pos0B,
2423 comA_sc, comB_sc, lambda,
2424 pbc, refcoord_scaling, npbcdim,
2427 for (m = 0; (m < DIM); m++)
2429 kk = L1*pr->posres.fcA[m] + lambda*pr->posres.fcB[m];
2431 vtot += 0.5*kk*dx[m]*dx[m];
2433 0.5*(pr->posres.fcB[m] - pr->posres.fcA[m])*dx[m]*dx[m]
2436 /* Here we correct for the pbc_dx which included rdist */
2440 vir_diag[m] -= 0.5*(dx[m] + rdist[m])*fm;
2448 static real low_angres(int nbonds,
2449 const t_iatom forceatoms[], const t_iparams forceparams[],
2450 const rvec x[], rvec f[], rvec fshift[],
2451 const t_pbc *pbc, const t_graph *g,
2452 real lambda, real *dvdlambda,
2455 int i, m, type, ai, aj, ak, al;
2457 real phi, cos_phi, cos_phi2, vid, vtot, dVdphi;
2458 rvec r_ij, r_kl, f_i, f_k = {0, 0, 0};
2459 real st, sth, nrij2, nrkl2, c, cij, ckl;
2462 t2 = 0; /* avoid warning with gcc-3.3. It is never used uninitialized */
2465 ak = al = 0; /* to avoid warnings */
2466 for (i = 0; i < nbonds; )
2468 type = forceatoms[i++];
2469 ai = forceatoms[i++];
2470 aj = forceatoms[i++];
2471 t1 = pbc_rvec_sub(pbc, x[aj], x[ai], r_ij); /* 3 */
2474 ak = forceatoms[i++];
2475 al = forceatoms[i++];
2476 t2 = pbc_rvec_sub(pbc, x[al], x[ak], r_kl); /* 3 */
2485 cos_phi = cos_angle(r_ij, r_kl); /* 25 */
2486 phi = acos(cos_phi); /* 10 */
2488 *dvdlambda += dopdihs_min(forceparams[type].pdihs.cpA,
2489 forceparams[type].pdihs.cpB,
2490 forceparams[type].pdihs.phiA,
2491 forceparams[type].pdihs.phiB,
2492 forceparams[type].pdihs.mult,
2493 phi, lambda, &vid, &dVdphi); /* 40 */
2497 cos_phi2 = sqr(cos_phi); /* 1 */
2500 st = -dVdphi*gmx_invsqrt(1 - cos_phi2); /* 12 */
2501 sth = st*cos_phi; /* 1 */
2502 nrij2 = iprod(r_ij, r_ij); /* 5 */
2503 nrkl2 = iprod(r_kl, r_kl); /* 5 */
2505 c = st*gmx_invsqrt(nrij2*nrkl2); /* 11 */
2506 cij = sth/nrij2; /* 10 */
2507 ckl = sth/nrkl2; /* 10 */
2509 for (m = 0; m < DIM; m++) /* 18+18 */
2511 f_i[m] = (c*r_kl[m]-cij*r_ij[m]);
2516 f_k[m] = (c*r_ij[m]-ckl*r_kl[m]);
2524 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
2527 rvec_inc(fshift[t1], f_i);
2528 rvec_dec(fshift[CENTRAL], f_i);
2533 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, al), dt);
2536 rvec_inc(fshift[t2], f_k);
2537 rvec_dec(fshift[CENTRAL], f_k);
2542 return vtot; /* 184 / 157 (bZAxis) total */
2545 real angres(int nbonds,
2546 const t_iatom forceatoms[], const t_iparams forceparams[],
2547 const rvec x[], rvec f[], rvec fshift[],
2548 const t_pbc *pbc, const t_graph *g,
2549 real lambda, real *dvdlambda,
2550 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2551 int gmx_unused *global_atom_index)
2553 return low_angres(nbonds, forceatoms, forceparams, x, f, fshift, pbc, g,
2554 lambda, dvdlambda, FALSE);
2557 real angresz(int nbonds,
2558 const t_iatom forceatoms[], const t_iparams forceparams[],
2559 const rvec x[], rvec f[], rvec fshift[],
2560 const t_pbc *pbc, const t_graph *g,
2561 real lambda, real *dvdlambda,
2562 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2563 int gmx_unused *global_atom_index)
2565 return low_angres(nbonds, forceatoms, forceparams, x, f, fshift, pbc, g,
2566 lambda, dvdlambda, TRUE);
2569 real dihres(int nbonds,
2570 const t_iatom forceatoms[], const t_iparams forceparams[],
2571 const rvec x[], rvec f[], rvec fshift[],
2572 const t_pbc *pbc, const t_graph *g,
2573 real lambda, real *dvdlambda,
2574 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2575 int gmx_unused *global_atom_index)
2578 int ai, aj, ak, al, i, k, type, t1, t2, t3;
2579 real phi0A, phi0B, dphiA, dphiB, kfacA, kfacB, phi0, dphi, kfac;
2580 real phi, ddphi, ddp, ddp2, dp, sign, d2r, fc, L1;
2581 rvec r_ij, r_kj, r_kl, m, n;
2588 for (i = 0; (i < nbonds); )
2590 type = forceatoms[i++];
2591 ai = forceatoms[i++];
2592 aj = forceatoms[i++];
2593 ak = forceatoms[i++];
2594 al = forceatoms[i++];
2596 phi0A = forceparams[type].dihres.phiA*d2r;
2597 dphiA = forceparams[type].dihres.dphiA*d2r;
2598 kfacA = forceparams[type].dihres.kfacA;
2600 phi0B = forceparams[type].dihres.phiB*d2r;
2601 dphiB = forceparams[type].dihres.dphiB*d2r;
2602 kfacB = forceparams[type].dihres.kfacB;
2604 phi0 = L1*phi0A + lambda*phi0B;
2605 dphi = L1*dphiA + lambda*dphiB;
2606 kfac = L1*kfacA + lambda*kfacB;
2608 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2609 &sign, &t1, &t2, &t3);
2614 fprintf(debug, "dihres[%d]: %d %d %d %d : phi=%f, dphi=%f, kfac=%f\n",
2615 k++, ai, aj, ak, al, phi0, dphi, kfac);
2617 /* phi can jump if phi0 is close to Pi/-Pi, which will cause huge
2618 * force changes if we just apply a normal harmonic.
2619 * Instead, we first calculate phi-phi0 and take it modulo (-Pi,Pi).
2620 * This means we will never have the periodicity problem, unless
2621 * the dihedral is Pi away from phiO, which is very unlikely due to
2625 make_dp_periodic(&dp);
2631 else if (dp < -dphi)
2643 vtot += 0.5*kfac*ddp2;
2646 *dvdlambda += 0.5*(kfacB - kfacA)*ddp2;
2647 /* lambda dependence from changing restraint distances */
2650 *dvdlambda -= kfac*ddp*((dphiB - dphiA)+(phi0B - phi0A));
2654 *dvdlambda += kfac*ddp*((dphiB - dphiA)-(phi0B - phi0A));
2656 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
2657 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2664 real unimplemented(int gmx_unused nbonds,
2665 const t_iatom gmx_unused forceatoms[], const t_iparams gmx_unused forceparams[],
2666 const rvec gmx_unused x[], rvec gmx_unused f[], rvec gmx_unused fshift[],
2667 const t_pbc gmx_unused *pbc, const t_graph gmx_unused *g,
2668 real gmx_unused lambda, real gmx_unused *dvdlambda,
2669 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2670 int gmx_unused *global_atom_index)
2672 gmx_impl("*** you are using a not implemented function");
2674 return 0.0; /* To make the compiler happy */
2677 real restrangles(int nbonds,
2678 const t_iatom forceatoms[], const t_iparams forceparams[],
2679 const rvec x[], rvec f[], rvec fshift[],
2680 const t_pbc *pbc, const t_graph *g,
2681 real gmx_unused lambda, real gmx_unused *dvdlambda,
2682 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2683 int gmx_unused *global_atom_index)
2685 int i, d, ai, aj, ak, type, m;
2689 ivec jt, dt_ij, dt_kj;
2691 real prefactor, ratio_ante, ratio_post;
2692 rvec delta_ante, delta_post, vec_temp;
2695 for (i = 0; (i < nbonds); )
2697 type = forceatoms[i++];
2698 ai = forceatoms[i++];
2699 aj = forceatoms[i++];
2700 ak = forceatoms[i++];
2702 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], vec_temp);
2703 pbc_rvec_sub(pbc, x[aj], x[ai], delta_ante);
2704 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], delta_post);
2707 /* This function computes factors needed for restricted angle potential.
2708 * The restricted angle potential is used in coarse-grained simulations to avoid singularities
2709 * when three particles align and the dihedral angle and dihedral potential
2710 * cannot be calculated. This potential is calculated using the formula:
2711 real restrangles(int nbonds,
2712 const t_iatom forceatoms[],const t_iparams forceparams[],
2713 const rvec x[],rvec f[],rvec fshift[],
2714 const t_pbc *pbc,const t_graph *g,
2715 real gmx_unused lambda,real gmx_unused *dvdlambda,
2716 const t_mdatoms gmx_unused *md,t_fcdata gmx_unused *fcd,
2717 int gmx_unused *global_atom_index)
2719 int i, d, ai, aj, ak, type, m;
2723 ivec jt,dt_ij,dt_kj;
2725 real prefactor, ratio_ante, ratio_post;
2726 rvec delta_ante, delta_post, vec_temp;
2729 for(i=0; (i<nbonds); )
2731 type = forceatoms[i++];
2732 ai = forceatoms[i++];
2733 aj = forceatoms[i++];
2734 ak = forceatoms[i++];
2736 * \f[V_{\rm ReB}(\theta_i) = \frac{1}{2} k_{\theta} \frac{(\cos\theta_i - \cos\theta_0)^2}
2737 * {\sin^2\theta_i}\f] ({eq:ReB} and ref \cite{MonicaGoga2013} from the manual).
2738 * For more explanations see comments file "restcbt.h". */
2740 compute_factors_restangles(type, forceparams, delta_ante, delta_post,
2741 &prefactor, &ratio_ante, &ratio_post, &v);
2743 /* Forces are computed per component */
2744 for (d = 0; d < DIM; d++)
2746 f_i[d] = prefactor * (ratio_ante * delta_ante[d] - delta_post[d]);
2747 f_j[d] = prefactor * ((ratio_post + 1.0) * delta_post[d] - (ratio_ante + 1.0) * delta_ante[d]);
2748 f_k[d] = prefactor * (delta_ante[d] - ratio_post * delta_post[d]);
2751 /* Computation of potential energy */
2757 for (m = 0; (m < DIM); m++)
2766 copy_ivec(SHIFT_IVEC(g, aj), jt);
2767 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
2768 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
2769 t1 = IVEC2IS(dt_ij);
2770 t2 = IVEC2IS(dt_kj);
2773 rvec_inc(fshift[t1], f_i);
2774 rvec_inc(fshift[CENTRAL], f_j);
2775 rvec_inc(fshift[t2], f_k);
2781 real restrdihs(int nbonds,
2782 const t_iatom forceatoms[], const t_iparams forceparams[],
2783 const rvec x[], rvec f[], rvec fshift[],
2784 const t_pbc *pbc, const t_graph *g,
2785 real gmx_unused lambda, real gmx_unused *dvlambda,
2786 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2787 int gmx_unused *global_atom_index)
2789 int i, d, type, ai, aj, ak, al;
2790 rvec f_i, f_j, f_k, f_l;
2792 ivec jt, dt_ij, dt_kj, dt_lj;
2795 rvec delta_ante, delta_crnt, delta_post, vec_temp;
2796 real factor_phi_ai_ante, factor_phi_ai_crnt, factor_phi_ai_post;
2797 real factor_phi_aj_ante, factor_phi_aj_crnt, factor_phi_aj_post;
2798 real factor_phi_ak_ante, factor_phi_ak_crnt, factor_phi_ak_post;
2799 real factor_phi_al_ante, factor_phi_al_crnt, factor_phi_al_post;
2804 for (i = 0; (i < nbonds); )
2806 type = forceatoms[i++];
2807 ai = forceatoms[i++];
2808 aj = forceatoms[i++];
2809 ak = forceatoms[i++];
2810 al = forceatoms[i++];
2812 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], vec_temp);
2813 pbc_rvec_sub(pbc, x[aj], x[ai], delta_ante);
2814 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], delta_crnt);
2815 t3 = pbc_rvec_sub(pbc, x[ak], x[al], vec_temp);
2816 pbc_rvec_sub(pbc, x[al], x[ak], delta_post);
2818 /* This function computes factors needed for restricted angle potential.
2819 * The restricted angle potential is used in coarse-grained simulations to avoid singularities
2820 * when three particles align and the dihedral angle and dihedral potential cannot be calculated.
2821 * This potential is calculated using the formula:
2822 * \f[V_{\rm ReB}(\theta_i) = \frac{1}{2} k_{\theta}
2823 * \frac{(\cos\theta_i - \cos\theta_0)^2}{\sin^2\theta_i}\f]
2824 * ({eq:ReB} and ref \cite{MonicaGoga2013} from the manual).
2825 * For more explanations see comments file "restcbt.h" */
2827 compute_factors_restrdihs(type, forceparams,
2828 delta_ante, delta_crnt, delta_post,
2829 &factor_phi_ai_ante, &factor_phi_ai_crnt, &factor_phi_ai_post,
2830 &factor_phi_aj_ante, &factor_phi_aj_crnt, &factor_phi_aj_post,
2831 &factor_phi_ak_ante, &factor_phi_ak_crnt, &factor_phi_ak_post,
2832 &factor_phi_al_ante, &factor_phi_al_crnt, &factor_phi_al_post,
2833 &prefactor_phi, &v);
2836 /* Computation of forces per component */
2837 for (d = 0; d < DIM; d++)
2839 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]);
2840 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]);
2841 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]);
2842 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]);
2844 /* Computation of the energy */
2850 /* Updating the forces */
2852 rvec_inc(f[ai], f_i);
2853 rvec_inc(f[aj], f_j);
2854 rvec_inc(f[ak], f_k);
2855 rvec_inc(f[al], f_l);
2858 /* Updating the fshift forces for the pressure coupling */
2861 copy_ivec(SHIFT_IVEC(g, aj), jt);
2862 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
2863 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
2864 ivec_sub(SHIFT_IVEC(g, al), jt, dt_lj);
2865 t1 = IVEC2IS(dt_ij);
2866 t2 = IVEC2IS(dt_kj);
2867 t3 = IVEC2IS(dt_lj);
2871 t3 = pbc_rvec_sub(pbc, x[al], x[aj], dx_jl);
2878 rvec_inc(fshift[t1], f_i);
2879 rvec_inc(fshift[CENTRAL], f_j);
2880 rvec_inc(fshift[t2], f_k);
2881 rvec_inc(fshift[t3], f_l);
2889 real cbtdihs(int nbonds,
2890 const t_iatom forceatoms[], const t_iparams forceparams[],
2891 const rvec x[], rvec f[], rvec fshift[],
2892 const t_pbc *pbc, const t_graph *g,
2893 real gmx_unused lambda, real gmx_unused *dvdlambda,
2894 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2895 int gmx_unused *global_atom_index)
2897 int type, ai, aj, ak, al, i, d;
2901 rvec f_i, f_j, f_k, f_l;
2902 ivec jt, dt_ij, dt_kj, dt_lj;
2904 rvec delta_ante, delta_crnt, delta_post;
2905 rvec f_phi_ai, f_phi_aj, f_phi_ak, f_phi_al;
2906 rvec f_theta_ante_ai, f_theta_ante_aj, f_theta_ante_ak;
2907 rvec f_theta_post_aj, f_theta_post_ak, f_theta_post_al;
2913 for (i = 0; (i < nbonds); )
2915 type = forceatoms[i++];
2916 ai = forceatoms[i++];
2917 aj = forceatoms[i++];
2918 ak = forceatoms[i++];
2919 al = forceatoms[i++];
2922 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], vec_temp);
2923 pbc_rvec_sub(pbc, x[aj], x[ai], delta_ante);
2924 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], vec_temp);
2925 pbc_rvec_sub(pbc, x[ak], x[aj], delta_crnt);
2926 t3 = pbc_rvec_sub(pbc, x[ak], x[al], vec_temp);
2927 pbc_rvec_sub(pbc, x[al], x[ak], delta_post);
2929 /* \brief Compute factors for CBT potential
2930 * The combined bending-torsion potential goes to zero in a very smooth manner, eliminating the numerical
2931 * instabilities, when three coarse-grained particles align and the dihedral angle and standard
2932 * dihedral potentials cannot be calculated. The CBT potential is calculated using the formula:
2933 * \f[V_{\rm CBT}(\theta_{i-1}, \theta_i, \phi_i) = k_{\phi} \sin^3\theta_{i-1} \sin^3\theta_{i}
2934 * \sum_{n=0}^4 { a_n \cos^n\phi_i}\f] ({eq:CBT} and ref \cite{MonicaGoga2013} from the manual).
2935 * It contains in its expression not only the dihedral angle \f$\phi\f$
2936 * but also \f[\theta_{i-1}\f] (theta_ante bellow) and \f[\theta_{i}\f] (theta_post bellow)
2937 * --- the adjacent bending angles.
2938 * For more explanations see comments file "restcbt.h". */
2940 compute_factors_cbtdihs(type, forceparams, delta_ante, delta_crnt, delta_post,
2941 f_phi_ai, f_phi_aj, f_phi_ak, f_phi_al,
2942 f_theta_ante_ai, f_theta_ante_aj, f_theta_ante_ak,
2943 f_theta_post_aj, f_theta_post_ak, f_theta_post_al,
2947 /* Acumulate the resuts per beads */
2948 for (d = 0; d < DIM; d++)
2950 f_i[d] = f_phi_ai[d] + f_theta_ante_ai[d];
2951 f_j[d] = f_phi_aj[d] + f_theta_ante_aj[d] + f_theta_post_aj[d];
2952 f_k[d] = f_phi_ak[d] + f_theta_ante_ak[d] + f_theta_post_ak[d];
2953 f_l[d] = f_phi_al[d] + f_theta_post_al[d];
2956 /* Compute the potential energy */
2961 /* Updating the forces */
2962 rvec_inc(f[ai], f_i);
2963 rvec_inc(f[aj], f_j);
2964 rvec_inc(f[ak], f_k);
2965 rvec_inc(f[al], f_l);
2968 /* Updating the fshift forces for the pressure coupling */
2971 copy_ivec(SHIFT_IVEC(g, aj), jt);
2972 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
2973 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
2974 ivec_sub(SHIFT_IVEC(g, al), jt, dt_lj);
2975 t1 = IVEC2IS(dt_ij);
2976 t2 = IVEC2IS(dt_kj);
2977 t3 = IVEC2IS(dt_lj);
2981 t3 = pbc_rvec_sub(pbc, x[al], x[aj], dx_jl);
2988 rvec_inc(fshift[t1], f_i);
2989 rvec_inc(fshift[CENTRAL], f_j);
2990 rvec_inc(fshift[t2], f_k);
2991 rvec_inc(fshift[t3], f_l);
2997 real rbdihs(int nbonds,
2998 const t_iatom forceatoms[], const t_iparams forceparams[],
2999 const rvec x[], rvec f[], rvec fshift[],
3000 const t_pbc *pbc, const t_graph *g,
3001 real lambda, real *dvdlambda,
3002 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3003 int gmx_unused *global_atom_index)
3005 const real c0 = 0.0, c1 = 1.0, c2 = 2.0, c3 = 3.0, c4 = 4.0, c5 = 5.0;
3006 int type, ai, aj, ak, al, i, j;
3008 rvec r_ij, r_kj, r_kl, m, n;
3009 real parmA[NR_RBDIHS];
3010 real parmB[NR_RBDIHS];
3011 real parm[NR_RBDIHS];
3012 real cos_phi, phi, rbp, rbpBA;
3013 real v, sign, ddphi, sin_phi;
3015 real L1 = 1.0-lambda;
3019 for (i = 0; (i < nbonds); )
3021 type = forceatoms[i++];
3022 ai = forceatoms[i++];
3023 aj = forceatoms[i++];
3024 ak = forceatoms[i++];
3025 al = forceatoms[i++];
3027 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
3028 &sign, &t1, &t2, &t3); /* 84 */
3030 /* Change to polymer convention */
3037 phi -= M_PI; /* 1 */
3041 /* Beware of accuracy loss, cannot use 1-sqrt(cos^2) ! */
3044 for (j = 0; (j < NR_RBDIHS); j++)
3046 parmA[j] = forceparams[type].rbdihs.rbcA[j];
3047 parmB[j] = forceparams[type].rbdihs.rbcB[j];
3048 parm[j] = L1*parmA[j]+lambda*parmB[j];
3050 /* Calculate cosine powers */
3051 /* Calculate the energy */
3052 /* Calculate the derivative */
3055 dvdl_term += (parmB[0]-parmA[0]);
3060 rbpBA = parmB[1]-parmA[1];
3061 ddphi += rbp*cosfac;
3064 dvdl_term += cosfac*rbpBA;
3066 rbpBA = parmB[2]-parmA[2];
3067 ddphi += c2*rbp*cosfac;
3070 dvdl_term += cosfac*rbpBA;
3072 rbpBA = parmB[3]-parmA[3];
3073 ddphi += c3*rbp*cosfac;
3076 dvdl_term += cosfac*rbpBA;
3078 rbpBA = parmB[4]-parmA[4];
3079 ddphi += c4*rbp*cosfac;
3082 dvdl_term += cosfac*rbpBA;
3084 rbpBA = parmB[5]-parmA[5];
3085 ddphi += c5*rbp*cosfac;
3088 dvdl_term += cosfac*rbpBA;
3090 ddphi = -ddphi*sin_phi; /* 11 */
3092 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
3093 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
3096 *dvdlambda += dvdl_term;
3101 int cmap_setup_grid_index(int ip, int grid_spacing, int *ipm1, int *ipp1, int *ipp2)
3107 ip = ip + grid_spacing - 1;
3109 else if (ip > grid_spacing)
3111 ip = ip - grid_spacing - 1;
3120 im1 = grid_spacing - 1;
3122 else if (ip == grid_spacing-2)
3126 else if (ip == grid_spacing-1)
3140 real cmap_dihs(int nbonds,
3141 const t_iatom forceatoms[], const t_iparams forceparams[],
3142 const gmx_cmap_t *cmap_grid,
3143 const rvec x[], rvec f[], rvec fshift[],
3144 const t_pbc *pbc, const t_graph *g,
3145 real gmx_unused lambda, real gmx_unused *dvdlambda,
3146 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3147 int gmx_unused *global_atom_index)
3149 int i, j, k, n, idx;
3150 int ai, aj, ak, al, am;
3151 int a1i, a1j, a1k, a1l, a2i, a2j, a2k, a2l;
3153 int t11, t21, t31, t12, t22, t32;
3154 int iphi1, ip1m1, ip1p1, ip1p2;
3155 int iphi2, ip2m1, ip2p1, ip2p2;
3157 int pos1, pos2, pos3, pos4, tmp;
3159 real ty[4], ty1[4], ty2[4], ty12[4], tc[16], tx[16];
3160 real phi1, psi1, cos_phi1, sin_phi1, sign1, xphi1;
3161 real phi2, psi2, cos_phi2, sin_phi2, sign2, xphi2;
3162 real dx, xx, tt, tu, e, df1, df2, ddf1, ddf2, ddf12, vtot;
3163 real ra21, rb21, rg21, rg1, rgr1, ra2r1, rb2r1, rabr1;
3164 real ra22, rb22, rg22, rg2, rgr2, ra2r2, rb2r2, rabr2;
3165 real fg1, hg1, fga1, hgb1, gaa1, gbb1;
3166 real fg2, hg2, fga2, hgb2, gaa2, gbb2;
3169 rvec r1_ij, r1_kj, r1_kl, m1, n1;
3170 rvec r2_ij, r2_kj, r2_kl, m2, n2;
3171 rvec f1_i, f1_j, f1_k, f1_l;
3172 rvec f2_i, f2_j, f2_k, f2_l;
3173 rvec a1, b1, a2, b2;
3174 rvec f1, g1, h1, f2, g2, h2;
3175 rvec dtf1, dtg1, dth1, dtf2, dtg2, dth2;
3176 ivec jt1, dt1_ij, dt1_kj, dt1_lj;
3177 ivec jt2, dt2_ij, dt2_kj, dt2_lj;
3181 int loop_index[4][4] = {
3188 /* Total CMAP energy */
3191 for (n = 0; n < nbonds; )
3193 /* Five atoms are involved in the two torsions */
3194 type = forceatoms[n++];
3195 ai = forceatoms[n++];
3196 aj = forceatoms[n++];
3197 ak = forceatoms[n++];
3198 al = forceatoms[n++];
3199 am = forceatoms[n++];
3201 /* Which CMAP type is this */
3202 cmapA = forceparams[type].cmap.cmapA;
3203 cmapd = cmap_grid->cmapdata[cmapA].cmap;
3211 phi1 = dih_angle(x[a1i], x[a1j], x[a1k], x[a1l], pbc, r1_ij, r1_kj, r1_kl, m1, n1,
3212 &sign1, &t11, &t21, &t31); /* 84 */
3214 cos_phi1 = cos(phi1);
3216 a1[0] = r1_ij[1]*r1_kj[2]-r1_ij[2]*r1_kj[1];
3217 a1[1] = r1_ij[2]*r1_kj[0]-r1_ij[0]*r1_kj[2];
3218 a1[2] = r1_ij[0]*r1_kj[1]-r1_ij[1]*r1_kj[0]; /* 9 */
3220 b1[0] = r1_kl[1]*r1_kj[2]-r1_kl[2]*r1_kj[1];
3221 b1[1] = r1_kl[2]*r1_kj[0]-r1_kl[0]*r1_kj[2];
3222 b1[2] = r1_kl[0]*r1_kj[1]-r1_kl[1]*r1_kj[0]; /* 9 */
3224 tmp = pbc_rvec_sub(pbc, x[a1l], x[a1k], h1);
3226 ra21 = iprod(a1, a1); /* 5 */
3227 rb21 = iprod(b1, b1); /* 5 */
3228 rg21 = iprod(r1_kj, r1_kj); /* 5 */
3234 rabr1 = sqrt(ra2r1*rb2r1);
3236 sin_phi1 = rg1 * rabr1 * iprod(a1, h1) * (-1);
3238 if (cos_phi1 < -0.5 || cos_phi1 > 0.5)
3240 phi1 = asin(sin_phi1);
3250 phi1 = -M_PI - phi1;
3256 phi1 = acos(cos_phi1);
3264 xphi1 = phi1 + M_PI; /* 1 */
3266 /* Second torsion */
3272 phi2 = dih_angle(x[a2i], x[a2j], x[a2k], x[a2l], pbc, r2_ij, r2_kj, r2_kl, m2, n2,
3273 &sign2, &t12, &t22, &t32); /* 84 */
3275 cos_phi2 = cos(phi2);
3277 a2[0] = r2_ij[1]*r2_kj[2]-r2_ij[2]*r2_kj[1];
3278 a2[1] = r2_ij[2]*r2_kj[0]-r2_ij[0]*r2_kj[2];
3279 a2[2] = r2_ij[0]*r2_kj[1]-r2_ij[1]*r2_kj[0]; /* 9 */
3281 b2[0] = r2_kl[1]*r2_kj[2]-r2_kl[2]*r2_kj[1];
3282 b2[1] = r2_kl[2]*r2_kj[0]-r2_kl[0]*r2_kj[2];
3283 b2[2] = r2_kl[0]*r2_kj[1]-r2_kl[1]*r2_kj[0]; /* 9 */
3285 tmp = pbc_rvec_sub(pbc, x[a2l], x[a2k], h2);
3287 ra22 = iprod(a2, a2); /* 5 */
3288 rb22 = iprod(b2, b2); /* 5 */
3289 rg22 = iprod(r2_kj, r2_kj); /* 5 */
3295 rabr2 = sqrt(ra2r2*rb2r2);
3297 sin_phi2 = rg2 * rabr2 * iprod(a2, h2) * (-1);
3299 if (cos_phi2 < -0.5 || cos_phi2 > 0.5)
3301 phi2 = asin(sin_phi2);
3311 phi2 = -M_PI - phi2;
3317 phi2 = acos(cos_phi2);
3325 xphi2 = phi2 + M_PI; /* 1 */
3327 /* Range mangling */
3330 xphi1 = xphi1 + 2*M_PI;
3332 else if (xphi1 >= 2*M_PI)
3334 xphi1 = xphi1 - 2*M_PI;
3339 xphi2 = xphi2 + 2*M_PI;
3341 else if (xphi2 >= 2*M_PI)
3343 xphi2 = xphi2 - 2*M_PI;
3346 /* Number of grid points */
3347 dx = 2*M_PI / cmap_grid->grid_spacing;
3349 /* Where on the grid are we */
3350 iphi1 = (int)(xphi1/dx);
3351 iphi2 = (int)(xphi2/dx);
3353 iphi1 = cmap_setup_grid_index(iphi1, cmap_grid->grid_spacing, &ip1m1, &ip1p1, &ip1p2);
3354 iphi2 = cmap_setup_grid_index(iphi2, cmap_grid->grid_spacing, &ip2m1, &ip2p1, &ip2p2);
3356 pos1 = iphi1*cmap_grid->grid_spacing+iphi2;
3357 pos2 = ip1p1*cmap_grid->grid_spacing+iphi2;
3358 pos3 = ip1p1*cmap_grid->grid_spacing+ip2p1;
3359 pos4 = iphi1*cmap_grid->grid_spacing+ip2p1;
3361 ty[0] = cmapd[pos1*4];
3362 ty[1] = cmapd[pos2*4];
3363 ty[2] = cmapd[pos3*4];
3364 ty[3] = cmapd[pos4*4];
3366 ty1[0] = cmapd[pos1*4+1];
3367 ty1[1] = cmapd[pos2*4+1];
3368 ty1[2] = cmapd[pos3*4+1];
3369 ty1[3] = cmapd[pos4*4+1];
3371 ty2[0] = cmapd[pos1*4+2];
3372 ty2[1] = cmapd[pos2*4+2];
3373 ty2[2] = cmapd[pos3*4+2];
3374 ty2[3] = cmapd[pos4*4+2];
3376 ty12[0] = cmapd[pos1*4+3];
3377 ty12[1] = cmapd[pos2*4+3];
3378 ty12[2] = cmapd[pos3*4+3];
3379 ty12[3] = cmapd[pos4*4+3];
3381 /* Switch to degrees */
3382 dx = 360.0 / cmap_grid->grid_spacing;
3383 xphi1 = xphi1 * RAD2DEG;
3384 xphi2 = xphi2 * RAD2DEG;
3386 for (i = 0; i < 4; i++) /* 16 */
3389 tx[i+4] = ty1[i]*dx;
3390 tx[i+8] = ty2[i]*dx;
3391 tx[i+12] = ty12[i]*dx*dx;
3395 for (i = 0; i < 4; i++) /* 1056 */
3397 for (j = 0; j < 4; j++)
3400 for (k = 0; k < 16; k++)
3402 xx = xx + cmap_coeff_matrix[k*16+idx]*tx[k];
3410 tt = (xphi1-iphi1*dx)/dx;
3411 tu = (xphi2-iphi2*dx)/dx;
3420 for (i = 3; i >= 0; i--)
3422 l1 = loop_index[i][3];
3423 l2 = loop_index[i][2];
3424 l3 = loop_index[i][1];
3426 e = tt * e + ((tc[i*4+3]*tu+tc[i*4+2])*tu + tc[i*4+1])*tu+tc[i*4];
3427 df1 = tu * df1 + (3.0*tc[l1]*tt+2.0*tc[l2])*tt+tc[l3];
3428 df2 = tt * df2 + (3.0*tc[i*4+3]*tu+2.0*tc[i*4+2])*tu+tc[i*4+1];
3429 ddf1 = tu * ddf1 + 2.0*3.0*tc[l1]*tt+2.0*tc[l2];
3430 ddf2 = tt * ddf2 + 2.0*3.0*tc[4*i+3]*tu+2.0*tc[4*i+2];
3433 ddf12 = tc[5] + 2.0*tc[9]*tt + 3.0*tc[13]*tt*tt + 2.0*tu*(tc[6]+2.0*tc[10]*tt+3.0*tc[14]*tt*tt) +
3434 3.0*tu*tu*(tc[7]+2.0*tc[11]*tt+3.0*tc[15]*tt*tt);
3439 ddf1 = ddf1 * fac * fac;
3440 ddf2 = ddf2 * fac * fac;
3441 ddf12 = ddf12 * fac * fac;
3446 /* Do forces - first torsion */
3447 fg1 = iprod(r1_ij, r1_kj);
3448 hg1 = iprod(r1_kl, r1_kj);
3449 fga1 = fg1*ra2r1*rgr1;
3450 hgb1 = hg1*rb2r1*rgr1;
3454 for (i = 0; i < DIM; i++)
3456 dtf1[i] = gaa1 * a1[i];
3457 dtg1[i] = fga1 * a1[i] - hgb1 * b1[i];
3458 dth1[i] = gbb1 * b1[i];
3460 f1[i] = df1 * dtf1[i];
3461 g1[i] = df1 * dtg1[i];
3462 h1[i] = df1 * dth1[i];
3465 f1_j[i] = -f1[i] - g1[i];
3466 f1_k[i] = h1[i] + g1[i];
3469 f[a1i][i] = f[a1i][i] + f1_i[i];
3470 f[a1j][i] = f[a1j][i] + f1_j[i]; /* - f1[i] - g1[i] */
3471 f[a1k][i] = f[a1k][i] + f1_k[i]; /* h1[i] + g1[i] */
3472 f[a1l][i] = f[a1l][i] + f1_l[i]; /* h1[i] */
3475 /* Do forces - second torsion */
3476 fg2 = iprod(r2_ij, r2_kj);
3477 hg2 = iprod(r2_kl, r2_kj);
3478 fga2 = fg2*ra2r2*rgr2;
3479 hgb2 = hg2*rb2r2*rgr2;
3483 for (i = 0; i < DIM; i++)
3485 dtf2[i] = gaa2 * a2[i];
3486 dtg2[i] = fga2 * a2[i] - hgb2 * b2[i];
3487 dth2[i] = gbb2 * b2[i];
3489 f2[i] = df2 * dtf2[i];
3490 g2[i] = df2 * dtg2[i];
3491 h2[i] = df2 * dth2[i];
3494 f2_j[i] = -f2[i] - g2[i];
3495 f2_k[i] = h2[i] + g2[i];
3498 f[a2i][i] = f[a2i][i] + f2_i[i]; /* f2[i] */
3499 f[a2j][i] = f[a2j][i] + f2_j[i]; /* - f2[i] - g2[i] */
3500 f[a2k][i] = f[a2k][i] + f2_k[i]; /* h2[i] + g2[i] */
3501 f[a2l][i] = f[a2l][i] + f2_l[i]; /* - h2[i] */
3507 copy_ivec(SHIFT_IVEC(g, a1j), jt1);
3508 ivec_sub(SHIFT_IVEC(g, a1i), jt1, dt1_ij);
3509 ivec_sub(SHIFT_IVEC(g, a1k), jt1, dt1_kj);
3510 ivec_sub(SHIFT_IVEC(g, a1l), jt1, dt1_lj);
3511 t11 = IVEC2IS(dt1_ij);
3512 t21 = IVEC2IS(dt1_kj);
3513 t31 = IVEC2IS(dt1_lj);
3515 copy_ivec(SHIFT_IVEC(g, a2j), jt2);
3516 ivec_sub(SHIFT_IVEC(g, a2i), jt2, dt2_ij);
3517 ivec_sub(SHIFT_IVEC(g, a2k), jt2, dt2_kj);
3518 ivec_sub(SHIFT_IVEC(g, a2l), jt2, dt2_lj);
3519 t12 = IVEC2IS(dt2_ij);
3520 t22 = IVEC2IS(dt2_kj);
3521 t32 = IVEC2IS(dt2_lj);
3525 t31 = pbc_rvec_sub(pbc, x[a1l], x[a1j], h1);
3526 t32 = pbc_rvec_sub(pbc, x[a2l], x[a2j], h2);
3534 rvec_inc(fshift[t11], f1_i);
3535 rvec_inc(fshift[CENTRAL], f1_j);
3536 rvec_inc(fshift[t21], f1_k);
3537 rvec_inc(fshift[t31], f1_l);
3539 rvec_inc(fshift[t21], f2_i);
3540 rvec_inc(fshift[CENTRAL], f2_j);
3541 rvec_inc(fshift[t22], f2_k);
3542 rvec_inc(fshift[t32], f2_l);
3549 /***********************************************************
3551 * G R O M O S 9 6 F U N C T I O N S
3553 ***********************************************************/
3554 real g96harmonic(real kA, real kB, real xA, real xB, real x, real lambda,
3557 const real half = 0.5;
3558 real L1, kk, x0, dx, dx2;
3559 real v, f, dvdlambda;
3562 kk = L1*kA+lambda*kB;
3563 x0 = L1*xA+lambda*xB;
3570 dvdlambda = half*(kB-kA)*dx2 + (xA-xB)*kk*dx;
3577 /* That was 21 flops */
3580 real g96bonds(int nbonds,
3581 const t_iatom forceatoms[], const t_iparams forceparams[],
3582 const rvec x[], rvec f[], rvec fshift[],
3583 const t_pbc *pbc, const t_graph *g,
3584 real lambda, real *dvdlambda,
3585 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3586 int gmx_unused *global_atom_index)
3588 int i, m, ki, ai, aj, type;
3589 real dr2, fbond, vbond, fij, vtot;
3594 for (i = 0; (i < nbonds); )
3596 type = forceatoms[i++];
3597 ai = forceatoms[i++];
3598 aj = forceatoms[i++];
3600 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
3601 dr2 = iprod(dx, dx); /* 5 */
3603 *dvdlambda += g96harmonic(forceparams[type].harmonic.krA,
3604 forceparams[type].harmonic.krB,
3605 forceparams[type].harmonic.rA,
3606 forceparams[type].harmonic.rB,
3607 dr2, lambda, &vbond, &fbond);
3609 vtot += 0.5*vbond; /* 1*/
3613 fprintf(debug, "G96-BONDS: dr = %10g vbond = %10g fbond = %10g\n",
3614 sqrt(dr2), vbond, fbond);
3620 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
3623 for (m = 0; (m < DIM); m++) /* 15 */
3628 fshift[ki][m] += fij;
3629 fshift[CENTRAL][m] -= fij;
3635 real g96bond_angle(const rvec xi, const rvec xj, const rvec xk, const t_pbc *pbc,
3636 rvec r_ij, rvec r_kj,
3638 /* Return value is the angle between the bonds i-j and j-k */
3642 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
3643 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
3645 costh = cos_angle(r_ij, r_kj); /* 25 */
3650 real g96angles(int nbonds,
3651 const t_iatom forceatoms[], const t_iparams forceparams[],
3652 const rvec x[], rvec f[], rvec fshift[],
3653 const t_pbc *pbc, const t_graph *g,
3654 real lambda, real *dvdlambda,
3655 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3656 int gmx_unused *global_atom_index)
3658 int i, ai, aj, ak, type, m, t1, t2;
3660 real cos_theta, dVdt, va, vtot;
3661 real rij_1, rij_2, rkj_1, rkj_2, rijrkj_1;
3663 ivec jt, dt_ij, dt_kj;
3666 for (i = 0; (i < nbonds); )
3668 type = forceatoms[i++];
3669 ai = forceatoms[i++];
3670 aj = forceatoms[i++];
3671 ak = forceatoms[i++];
3673 cos_theta = g96bond_angle(x[ai], x[aj], x[ak], pbc, r_ij, r_kj, &t1, &t2);
3675 *dvdlambda += g96harmonic(forceparams[type].harmonic.krA,
3676 forceparams[type].harmonic.krB,
3677 forceparams[type].harmonic.rA,
3678 forceparams[type].harmonic.rB,
3679 cos_theta, lambda, &va, &dVdt);
3682 rij_1 = gmx_invsqrt(iprod(r_ij, r_ij));
3683 rkj_1 = gmx_invsqrt(iprod(r_kj, r_kj));
3684 rij_2 = rij_1*rij_1;
3685 rkj_2 = rkj_1*rkj_1;
3686 rijrkj_1 = rij_1*rkj_1; /* 23 */
3691 fprintf(debug, "G96ANGLES: costheta = %10g vth = %10g dV/dct = %10g\n",
3692 cos_theta, va, dVdt);
3695 for (m = 0; (m < DIM); m++) /* 42 */
3697 f_i[m] = dVdt*(r_kj[m]*rijrkj_1 - r_ij[m]*rij_2*cos_theta);
3698 f_k[m] = dVdt*(r_ij[m]*rijrkj_1 - r_kj[m]*rkj_2*cos_theta);
3699 f_j[m] = -f_i[m]-f_k[m];
3707 copy_ivec(SHIFT_IVEC(g, aj), jt);
3709 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3710 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3711 t1 = IVEC2IS(dt_ij);
3712 t2 = IVEC2IS(dt_kj);
3714 rvec_inc(fshift[t1], f_i);
3715 rvec_inc(fshift[CENTRAL], f_j);
3716 rvec_inc(fshift[t2], f_k); /* 9 */
3722 real cross_bond_bond(int nbonds,
3723 const t_iatom forceatoms[], const t_iparams forceparams[],
3724 const rvec x[], rvec f[], rvec fshift[],
3725 const t_pbc *pbc, const t_graph *g,
3726 real gmx_unused lambda, real gmx_unused *dvdlambda,
3727 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3728 int gmx_unused *global_atom_index)
3730 /* Potential from Lawrence and Skimmer, Chem. Phys. Lett. 372 (2003)
3733 int i, ai, aj, ak, type, m, t1, t2;
3735 real vtot, vrr, s1, s2, r1, r2, r1e, r2e, krr;
3737 ivec jt, dt_ij, dt_kj;
3740 for (i = 0; (i < nbonds); )
3742 type = forceatoms[i++];
3743 ai = forceatoms[i++];
3744 aj = forceatoms[i++];
3745 ak = forceatoms[i++];
3746 r1e = forceparams[type].cross_bb.r1e;
3747 r2e = forceparams[type].cross_bb.r2e;
3748 krr = forceparams[type].cross_bb.krr;
3750 /* Compute distance vectors ... */
3751 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
3752 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
3754 /* ... and their lengths */
3758 /* Deviations from ideality */
3762 /* Energy (can be negative!) */
3767 svmul(-krr*s2/r1, r_ij, f_i);
3768 svmul(-krr*s1/r2, r_kj, f_k);
3770 for (m = 0; (m < DIM); m++) /* 12 */
3772 f_j[m] = -f_i[m] - f_k[m];
3781 copy_ivec(SHIFT_IVEC(g, aj), jt);
3783 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3784 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3785 t1 = IVEC2IS(dt_ij);
3786 t2 = IVEC2IS(dt_kj);
3788 rvec_inc(fshift[t1], f_i);
3789 rvec_inc(fshift[CENTRAL], f_j);
3790 rvec_inc(fshift[t2], f_k); /* 9 */
3796 real cross_bond_angle(int nbonds,
3797 const t_iatom forceatoms[], const t_iparams forceparams[],
3798 const rvec x[], rvec f[], rvec fshift[],
3799 const t_pbc *pbc, const t_graph *g,
3800 real gmx_unused lambda, real gmx_unused *dvdlambda,
3801 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3802 int gmx_unused *global_atom_index)
3804 /* Potential from Lawrence and Skimmer, Chem. Phys. Lett. 372 (2003)
3807 int i, ai, aj, ak, type, m, t1, t2, t3;
3808 rvec r_ij, r_kj, r_ik;
3809 real vtot, vrt, s1, s2, s3, r1, r2, r3, r1e, r2e, r3e, krt, k1, k2, k3;
3811 ivec jt, dt_ij, dt_kj;
3814 for (i = 0; (i < nbonds); )
3816 type = forceatoms[i++];
3817 ai = forceatoms[i++];
3818 aj = forceatoms[i++];
3819 ak = forceatoms[i++];
3820 r1e = forceparams[type].cross_ba.r1e;
3821 r2e = forceparams[type].cross_ba.r2e;
3822 r3e = forceparams[type].cross_ba.r3e;
3823 krt = forceparams[type].cross_ba.krt;
3825 /* Compute distance vectors ... */
3826 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
3827 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
3828 t3 = pbc_rvec_sub(pbc, x[ai], x[ak], r_ik);
3830 /* ... and their lengths */
3835 /* Deviations from ideality */
3840 /* Energy (can be negative!) */
3841 vrt = krt*s3*(s1+s2);
3847 k3 = -krt*(s1+s2)/r3;
3848 for (m = 0; (m < DIM); m++)
3850 f_i[m] = k1*r_ij[m] + k3*r_ik[m];
3851 f_k[m] = k2*r_kj[m] - k3*r_ik[m];
3852 f_j[m] = -f_i[m] - f_k[m];
3855 for (m = 0; (m < DIM); m++) /* 12 */
3865 copy_ivec(SHIFT_IVEC(g, aj), jt);
3867 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3868 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3869 t1 = IVEC2IS(dt_ij);
3870 t2 = IVEC2IS(dt_kj);
3872 rvec_inc(fshift[t1], f_i);
3873 rvec_inc(fshift[CENTRAL], f_j);
3874 rvec_inc(fshift[t2], f_k); /* 9 */
3880 static real bonded_tab(const char *type, int table_nr,
3881 const bondedtable_t *table, real kA, real kB, real r,
3882 real lambda, real *V, real *F)
3884 real k, tabscale, *VFtab, rt, eps, eps2, Yt, Ft, Geps, Heps2, Fp, VV, FF;
3886 real v, f, dvdlambda;
3888 k = (1.0 - lambda)*kA + lambda*kB;
3890 tabscale = table->scale;
3891 VFtab = table->data;
3897 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",
3898 type, table_nr, r, n0, n0+1, table->n);
3905 Geps = VFtab[nnn+2]*eps;
3906 Heps2 = VFtab[nnn+3]*eps2;
3907 Fp = Ft + Geps + Heps2;
3909 FF = Fp + Geps + 2.0*Heps2;
3911 *F = -k*FF*tabscale;
3913 dvdlambda = (kB - kA)*VV;
3917 /* That was 22 flops */
3920 real tab_bonds(int nbonds,
3921 const t_iatom forceatoms[], const t_iparams forceparams[],
3922 const rvec x[], rvec f[], rvec fshift[],
3923 const t_pbc *pbc, const t_graph *g,
3924 real lambda, real *dvdlambda,
3925 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
3926 int gmx_unused *global_atom_index)
3928 int i, m, ki, ai, aj, type, table;
3929 real dr, dr2, fbond, vbond, fij, vtot;
3934 for (i = 0; (i < nbonds); )
3936 type = forceatoms[i++];
3937 ai = forceatoms[i++];
3938 aj = forceatoms[i++];
3940 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
3941 dr2 = iprod(dx, dx); /* 5 */
3942 dr = dr2*gmx_invsqrt(dr2); /* 10 */
3944 table = forceparams[type].tab.table;
3946 *dvdlambda += bonded_tab("bond", table,
3947 &fcd->bondtab[table],
3948 forceparams[type].tab.kA,
3949 forceparams[type].tab.kB,
3950 dr, lambda, &vbond, &fbond); /* 22 */
3958 vtot += vbond; /* 1*/
3959 fbond *= gmx_invsqrt(dr2); /* 6 */
3963 fprintf(debug, "TABBONDS: dr = %10g vbond = %10g fbond = %10g\n",
3969 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
3972 for (m = 0; (m < DIM); m++) /* 15 */
3977 fshift[ki][m] += fij;
3978 fshift[CENTRAL][m] -= fij;
3984 real tab_angles(int nbonds,
3985 const t_iatom forceatoms[], const t_iparams forceparams[],
3986 const rvec x[], rvec f[], rvec fshift[],
3987 const t_pbc *pbc, const t_graph *g,
3988 real lambda, real *dvdlambda,
3989 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
3990 int gmx_unused *global_atom_index)
3992 int i, ai, aj, ak, t1, t2, type, table;
3994 real cos_theta, cos_theta2, theta, dVdt, va, vtot;
3995 ivec jt, dt_ij, dt_kj;
3998 for (i = 0; (i < nbonds); )
4000 type = forceatoms[i++];
4001 ai = forceatoms[i++];
4002 aj = forceatoms[i++];
4003 ak = forceatoms[i++];
4005 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
4006 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
4008 table = forceparams[type].tab.table;
4010 *dvdlambda += bonded_tab("angle", table,
4011 &fcd->angletab[table],
4012 forceparams[type].tab.kA,
4013 forceparams[type].tab.kB,
4014 theta, lambda, &va, &dVdt); /* 22 */
4017 cos_theta2 = sqr(cos_theta); /* 1 */
4026 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
4027 sth = st*cos_theta; /* 1 */
4031 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
4032 theta*RAD2DEG, va, dVdt);
4035 nrkj2 = iprod(r_kj, r_kj); /* 5 */
4036 nrij2 = iprod(r_ij, r_ij);
4038 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
4039 cii = sth/nrij2; /* 10 */
4040 ckk = sth/nrkj2; /* 10 */
4042 for (m = 0; (m < DIM); m++) /* 39 */
4044 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
4045 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
4046 f_j[m] = -f_i[m]-f_k[m];
4053 copy_ivec(SHIFT_IVEC(g, aj), jt);
4055 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
4056 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
4057 t1 = IVEC2IS(dt_ij);
4058 t2 = IVEC2IS(dt_kj);
4060 rvec_inc(fshift[t1], f_i);
4061 rvec_inc(fshift[CENTRAL], f_j);
4062 rvec_inc(fshift[t2], f_k);
4068 real tab_dihs(int nbonds,
4069 const t_iatom forceatoms[], const t_iparams forceparams[],
4070 const rvec x[], rvec f[], rvec fshift[],
4071 const t_pbc *pbc, const t_graph *g,
4072 real lambda, real *dvdlambda,
4073 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
4074 int gmx_unused *global_atom_index)
4076 int i, type, ai, aj, ak, al, table;
4078 rvec r_ij, r_kj, r_kl, m, n;
4079 real phi, sign, ddphi, vpd, vtot;
4082 for (i = 0; (i < nbonds); )
4084 type = forceatoms[i++];
4085 ai = forceatoms[i++];
4086 aj = forceatoms[i++];
4087 ak = forceatoms[i++];
4088 al = forceatoms[i++];
4090 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
4091 &sign, &t1, &t2, &t3); /* 84 */
4093 table = forceparams[type].tab.table;
4095 /* Hopefully phi+M_PI never results in values < 0 */
4096 *dvdlambda += bonded_tab("dihedral", table,
4097 &fcd->dihtab[table],
4098 forceparams[type].tab.kA,
4099 forceparams[type].tab.kB,
4100 phi+M_PI, lambda, &vpd, &ddphi);
4103 do_dih_fup(ai, aj, ak, al, -ddphi, r_ij, r_kj, r_kl, m, n,
4104 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
4107 fprintf(debug, "pdih: (%d,%d,%d,%d) phi=%g\n",
4108 ai, aj, ak, al, phi);
4115 /* Return if this is a potential calculated in bondfree.c,
4116 * i.e. an interaction that actually calculates a potential and
4117 * works on multiple atoms (not e.g. a connection or a position restraint).
4119 static gmx_inline gmx_bool ftype_is_bonded_potential(int ftype)
4122 (interaction_function[ftype].flags & IF_BOND) &&
4123 !(ftype == F_CONNBONDS || ftype == F_POSRES || ftype == F_FBPOSRES) &&
4124 (ftype < F_GB12 || ftype > F_GB14);
4127 static void divide_bondeds_over_threads(t_idef *idef, int nthreads)
4134 idef->nthreads = nthreads;
4136 if (F_NRE*(nthreads+1) > idef->il_thread_division_nalloc)
4138 idef->il_thread_division_nalloc = F_NRE*(nthreads+1);
4139 snew(idef->il_thread_division, idef->il_thread_division_nalloc);
4142 for (ftype = 0; ftype < F_NRE; ftype++)
4144 if (ftype_is_bonded_potential(ftype))
4146 nat1 = interaction_function[ftype].nratoms + 1;
4148 for (t = 0; t <= nthreads; t++)
4150 /* Divide the interactions equally over the threads.
4151 * When the different types of bonded interactions
4152 * are distributed roughly equally over the threads,
4153 * this should lead to well localized output into
4154 * the force buffer on each thread.
4155 * If this is not the case, a more advanced scheme
4156 * (not implemented yet) will do better.
4158 il_nr_thread = (((idef->il[ftype].nr/nat1)*t)/nthreads)*nat1;
4160 /* Ensure that distance restraint pairs with the same label
4161 * end up on the same thread.
4162 * This is slighlty tricky code, since the next for iteration
4163 * may have an initial il_nr_thread lower than the final value
4164 * in the previous iteration, but this will anyhow be increased
4165 * to the approriate value again by this while loop.
4167 while (ftype == F_DISRES &&
4169 il_nr_thread < idef->il[ftype].nr &&
4170 idef->iparams[idef->il[ftype].iatoms[il_nr_thread]].disres.label ==
4171 idef->iparams[idef->il[ftype].iatoms[il_nr_thread-nat1]].disres.label)
4173 il_nr_thread += nat1;
4176 idef->il_thread_division[ftype*(nthreads+1)+t] = il_nr_thread;
4183 calc_bonded_reduction_mask(const t_idef *idef,
4188 int ftype, nb, nat1, nb0, nb1, i, a;
4192 for (ftype = 0; ftype < F_NRE; ftype++)
4194 if (ftype_is_bonded_potential(ftype))
4196 nb = idef->il[ftype].nr;
4199 nat1 = interaction_function[ftype].nratoms + 1;
4201 /* Divide this interaction equally over the threads.
4202 * This is not stored: should match division in calc_bonds.
4204 nb0 = idef->il_thread_division[ftype*(nt+1)+t];
4205 nb1 = idef->il_thread_division[ftype*(nt+1)+t+1];
4207 for (i = nb0; i < nb1; i += nat1)
4209 for (a = 1; a < nat1; a++)
4211 mask |= (1U << (idef->il[ftype].iatoms[i+a]>>shift));
4221 void setup_bonded_threading(t_forcerec *fr, t_idef *idef)
4223 #define MAX_BLOCK_BITS 32
4227 assert(fr->nthreads >= 1);
4229 /* Divide the bonded interaction over the threads */
4230 divide_bondeds_over_threads(idef, fr->nthreads);
4232 if (fr->nthreads == 1)
4239 /* We divide the force array in a maximum of 32 blocks.
4240 * Minimum force block reduction size is 2^6=64.
4243 while (fr->natoms_force > (int)(MAX_BLOCK_BITS*(1U<<fr->red_ashift)))
4249 fprintf(debug, "bonded force buffer block atom shift %d bits\n",
4253 /* Determine to which blocks each thread's bonded force calculation
4254 * contributes. Store this is a mask for each thread.
4256 #pragma omp parallel for num_threads(fr->nthreads) schedule(static)
4257 for (t = 1; t < fr->nthreads; t++)
4259 fr->f_t[t].red_mask =
4260 calc_bonded_reduction_mask(idef, fr->red_ashift, t, fr->nthreads);
4263 /* Determine the maximum number of blocks we need to reduce over */
4266 for (t = 0; t < fr->nthreads; t++)
4269 for (b = 0; b < MAX_BLOCK_BITS; b++)
4271 if (fr->f_t[t].red_mask & (1U<<b))
4273 fr->red_nblock = max(fr->red_nblock, b+1);
4279 fprintf(debug, "thread %d flags %x count %d\n",
4280 t, fr->f_t[t].red_mask, c);
4286 fprintf(debug, "Number of blocks to reduce: %d of size %d\n",
4287 fr->red_nblock, 1<<fr->red_ashift);
4288 fprintf(debug, "Reduction density %.2f density/#thread %.2f\n",
4289 ctot*(1<<fr->red_ashift)/(double)fr->natoms_force,
4290 ctot*(1<<fr->red_ashift)/(double)(fr->natoms_force*fr->nthreads));
4294 static void zero_thread_forces(f_thread_t *f_t, int n,
4295 int nblock, int blocksize)
4297 int b, a0, a1, a, i, j;
4299 if (n > f_t->f_nalloc)
4301 f_t->f_nalloc = over_alloc_large(n);
4302 srenew(f_t->f, f_t->f_nalloc);
4305 if (f_t->red_mask != 0)
4307 for (b = 0; b < nblock; b++)
4309 if (f_t->red_mask && (1U<<b))
4312 a1 = min((b+1)*blocksize, n);
4313 for (a = a0; a < a1; a++)
4315 clear_rvec(f_t->f[a]);
4320 for (i = 0; i < SHIFTS; i++)
4322 clear_rvec(f_t->fshift[i]);
4324 for (i = 0; i < F_NRE; i++)
4328 for (i = 0; i < egNR; i++)
4330 for (j = 0; j < f_t->grpp.nener; j++)
4332 f_t->grpp.ener[i][j] = 0;
4335 for (i = 0; i < efptNR; i++)
4341 static void reduce_thread_force_buffer(int n, rvec *f,
4342 int nthreads, f_thread_t *f_t,
4343 int nblock, int block_size)
4345 /* The max thread number is arbitrary,
4346 * we used a fixed number to avoid memory management.
4347 * Using more than 16 threads is probably never useful performance wise.
4349 #define MAX_BONDED_THREADS 256
4352 if (nthreads > MAX_BONDED_THREADS)
4354 gmx_fatal(FARGS, "Can not reduce bonded forces on more than %d threads",
4355 MAX_BONDED_THREADS);
4358 /* This reduction can run on any number of threads,
4359 * independently of nthreads.
4361 #pragma omp parallel for num_threads(nthreads) schedule(static)
4362 for (b = 0; b < nblock; b++)
4364 rvec *fp[MAX_BONDED_THREADS];
4368 /* Determine which threads contribute to this block */
4370 for (ft = 1; ft < nthreads; ft++)
4372 if (f_t[ft].red_mask & (1U<<b))
4374 fp[nfb++] = f_t[ft].f;
4379 /* Reduce force buffers for threads that contribute */
4381 a1 = (b+1)*block_size;
4383 for (a = a0; a < a1; a++)
4385 for (fb = 0; fb < nfb; fb++)
4387 rvec_inc(f[a], fp[fb][a]);
4394 static void reduce_thread_forces(int n, rvec *f, rvec *fshift,
4395 real *ener, gmx_grppairener_t *grpp, real *dvdl,
4396 int nthreads, f_thread_t *f_t,
4397 int nblock, int block_size,
4398 gmx_bool bCalcEnerVir,
4403 /* Reduce the bonded force buffer */
4404 reduce_thread_force_buffer(n, f, nthreads, f_t, nblock, block_size);
4407 /* When necessary, reduce energy and virial using one thread only */
4412 for (i = 0; i < SHIFTS; i++)
4414 for (t = 1; t < nthreads; t++)
4416 rvec_inc(fshift[i], f_t[t].fshift[i]);
4419 for (i = 0; i < F_NRE; i++)
4421 for (t = 1; t < nthreads; t++)
4423 ener[i] += f_t[t].ener[i];
4426 for (i = 0; i < egNR; i++)
4428 for (j = 0; j < f_t[1].grpp.nener; j++)
4430 for (t = 1; t < nthreads; t++)
4433 grpp->ener[i][j] += f_t[t].grpp.ener[i][j];
4439 for (i = 0; i < efptNR; i++)
4442 for (t = 1; t < nthreads; t++)
4444 dvdl[i] += f_t[t].dvdl[i];
4451 static real calc_one_bond(int thread,
4452 int ftype, const t_idef *idef,
4453 rvec x[], rvec f[], rvec fshift[],
4455 const t_pbc *pbc, const t_graph *g,
4456 gmx_grppairener_t *grpp,
4458 real *lambda, real *dvdl,
4459 const t_mdatoms *md, t_fcdata *fcd,
4460 gmx_bool bCalcEnerVir,
4461 int *global_atom_index)
4463 int nat1, nbonds, efptFTYPE;
4468 if (IS_RESTRAINT_TYPE(ftype))
4470 efptFTYPE = efptRESTRAINT;
4474 efptFTYPE = efptBONDED;
4477 nat1 = interaction_function[ftype].nratoms + 1;
4478 nbonds = idef->il[ftype].nr/nat1;
4479 iatoms = idef->il[ftype].iatoms;
4481 nb0 = idef->il_thread_division[ftype*(idef->nthreads+1)+thread];
4482 nbn = idef->il_thread_division[ftype*(idef->nthreads+1)+thread+1] - nb0;
4484 if (!IS_LISTED_LJ_C(ftype))
4486 if (ftype == F_CMAP)
4488 v = cmap_dihs(nbn, iatoms+nb0,
4489 idef->iparams, &idef->cmap_grid,
4490 (const rvec*)x, f, fshift,
4491 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
4492 md, fcd, global_atom_index);
4494 #ifdef GMX_SIMD_HAVE_REAL
4495 else if (ftype == F_ANGLES &&
4496 !bCalcEnerVir && fr->efep == efepNO)
4498 /* No energies, shift forces, dvdl */
4499 angles_noener_simd(nbn, idef->il[ftype].iatoms+nb0,
4502 pbc, g, lambda[efptFTYPE], md, fcd,
4507 else if (ftype == F_PDIHS &&
4508 !bCalcEnerVir && fr->efep == efepNO)
4510 /* No energies, shift forces, dvdl */
4511 #ifdef GMX_SIMD_HAVE_REAL
4516 (nbn, idef->il[ftype].iatoms+nb0,
4519 pbc, g, lambda[efptFTYPE], md, fcd,
4525 v = interaction_function[ftype].ifunc(nbn, iatoms+nb0,
4527 (const rvec*)x, f, fshift,
4528 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
4529 md, fcd, global_atom_index);
4534 v = do_nonbonded_listed(ftype, nbn, iatoms+nb0, idef->iparams, (const rvec*)x, f, fshift,
4535 pbc, g, lambda, dvdl, md, fr, grpp, global_atom_index);
4540 inc_nrnb(nrnb, interaction_function[ftype].nrnb_ind, nbonds);
4546 void calc_bonds(const gmx_multisim_t *ms,
4548 rvec x[], history_t *hist,
4549 rvec f[], t_forcerec *fr,
4550 const t_pbc *pbc, const t_graph *g,
4551 gmx_enerdata_t *enerd, t_nrnb *nrnb,
4553 const t_mdatoms *md,
4554 t_fcdata *fcd, int *global_atom_index,
4555 t_atomtypes gmx_unused *atype, gmx_genborn_t gmx_unused *born,
4558 gmx_bool bCalcEnerVir;
4560 real v, dvdl[efptNR], dvdl_dum[efptNR]; /* The dummy array is to have a place to store the dhdl at other values
4561 of lambda, which will be thrown away in the end*/
4562 const t_pbc *pbc_null;
4566 assert(fr->nthreads == idef->nthreads);
4568 bCalcEnerVir = (force_flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY));
4570 for (i = 0; i < efptNR; i++)
4586 p_graph(debug, "Bondage is fun", g);
4590 /* Do pre force calculation stuff which might require communication */
4591 if (idef->il[F_ORIRES].nr)
4593 enerd->term[F_ORIRESDEV] =
4594 calc_orires_dev(ms, idef->il[F_ORIRES].nr,
4595 idef->il[F_ORIRES].iatoms,
4596 idef->iparams, md, (const rvec*)x,
4597 pbc_null, fcd, hist);
4599 if (idef->il[F_DISRES].nr)
4601 calc_disres_R_6(idef->il[F_DISRES].nr,
4602 idef->il[F_DISRES].iatoms,
4603 idef->iparams, (const rvec*)x, pbc_null,
4606 if (fcd->disres.nsystems > 1)
4608 gmx_sum_sim(2*fcd->disres.nres, fcd->disres.Rt_6, ms);
4613 #pragma omp parallel for num_threads(fr->nthreads) schedule(static)
4614 for (thread = 0; thread < fr->nthreads; thread++)
4621 gmx_grppairener_t *grpp;
4626 fshift = fr->fshift;
4628 grpp = &enerd->grpp;
4633 zero_thread_forces(&fr->f_t[thread], fr->natoms_force,
4634 fr->red_nblock, 1<<fr->red_ashift);
4636 ft = fr->f_t[thread].f;
4637 fshift = fr->f_t[thread].fshift;
4638 epot = fr->f_t[thread].ener;
4639 grpp = &fr->f_t[thread].grpp;
4640 dvdlt = fr->f_t[thread].dvdl;
4642 /* Loop over all bonded force types to calculate the bonded forces */
4643 for (ftype = 0; (ftype < F_NRE); ftype++)
4645 if (idef->il[ftype].nr > 0 && ftype_is_bonded_potential(ftype))
4647 v = calc_one_bond(thread, ftype, idef, x,
4648 ft, fshift, fr, pbc_null, g, grpp,
4649 nrnb, lambda, dvdlt,
4650 md, fcd, bCalcEnerVir,
4656 if (fr->nthreads > 1)
4658 reduce_thread_forces(fr->natoms_force, f, fr->fshift,
4659 enerd->term, &enerd->grpp, dvdl,
4660 fr->nthreads, fr->f_t,
4661 fr->red_nblock, 1<<fr->red_ashift,
4663 force_flags & GMX_FORCE_DHDL);
4665 if (force_flags & GMX_FORCE_DHDL)
4667 for (i = 0; i < efptNR; i++)
4669 enerd->dvdl_nonlin[i] += dvdl[i];
4673 /* Copy the sum of violations for the distance restraints from fcd */
4676 enerd->term[F_DISRESVIOL] = fcd->disres.sumviol;
4681 void calc_bonds_lambda(const t_idef *idef,
4684 const t_pbc *pbc, const t_graph *g,
4685 gmx_grppairener_t *grpp, real *epot, t_nrnb *nrnb,
4687 const t_mdatoms *md,
4689 int *global_atom_index)
4691 int i, ftype, nr_nonperturbed, nr;
4693 real dvdl_dum[efptNR];
4695 const t_pbc *pbc_null;
4707 /* Copy the whole idef, so we can modify the contents locally */
4709 idef_fe.nthreads = 1;
4710 snew(idef_fe.il_thread_division, F_NRE*(idef_fe.nthreads+1));
4712 /* We already have the forces, so we use temp buffers here */
4713 snew(f, fr->natoms_force);
4714 snew(fshift, SHIFTS);
4716 /* Loop over all bonded force types to calculate the bonded energies */
4717 for (ftype = 0; (ftype < F_NRE); ftype++)
4719 if (ftype_is_bonded_potential(ftype))
4721 /* Set the work range of thread 0 to the perturbed bondeds only */
4722 nr_nonperturbed = idef->il[ftype].nr_nonperturbed;
4723 nr = idef->il[ftype].nr;
4724 idef_fe.il_thread_division[ftype*2+0] = nr_nonperturbed;
4725 idef_fe.il_thread_division[ftype*2+1] = nr;
4727 /* This is only to get the flop count correct */
4728 idef_fe.il[ftype].nr = nr - nr_nonperturbed;
4730 if (nr - nr_nonperturbed > 0)
4732 v = calc_one_bond(0, ftype, &idef_fe,
4733 x, f, fshift, fr, pbc_null, g,
4734 grpp, nrnb, lambda, dvdl_dum,
4745 sfree(idef_fe.il_thread_division);