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.
45 #include "gromacs/math/utilities.h"
53 #include "gmx_fatal.h"
59 #include "nonbonded.h"
61 /* Include the SIMD macro file and then check for support */
62 #include "gromacs/simd/macros.h"
63 #if defined GMX_HAVE_SIMD_MACROS && defined GMX_SIMD_HAVE_TRIGONOMETRIC
65 #include "gromacs/simd/vector_operations.h"
68 /* Find a better place for this? */
69 const int cmap_coeff_matrix[] = {
70 1, 0, -3, 2, 0, 0, 0, 0, -3, 0, 9, -6, 2, 0, -6, 4,
71 0, 0, 0, 0, 0, 0, 0, 0, 3, 0, -9, 6, -2, 0, 6, -4,
72 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 9, -6, 0, 0, -6, 4,
73 0, 0, 3, -2, 0, 0, 0, 0, 0, 0, -9, 6, 0, 0, 6, -4,
74 0, 0, 0, 0, 1, 0, -3, 2, -2, 0, 6, -4, 1, 0, -3, 2,
75 0, 0, 0, 0, 0, 0, 0, 0, -1, 0, 3, -2, 1, 0, -3, 2,
76 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -3, 2, 0, 0, 3, -2,
77 0, 0, 0, 0, 0, 0, 3, -2, 0, 0, -6, 4, 0, 0, 3, -2,
78 0, 1, -2, 1, 0, 0, 0, 0, 0, -3, 6, -3, 0, 2, -4, 2,
79 0, 0, 0, 0, 0, 0, 0, 0, 0, 3, -6, 3, 0, -2, 4, -2,
80 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -3, 3, 0, 0, 2, -2,
81 0, 0, -1, 1, 0, 0, 0, 0, 0, 0, 3, -3, 0, 0, -2, 2,
82 0, 0, 0, 0, 0, 1, -2, 1, 0, -2, 4, -2, 0, 1, -2, 1,
83 0, 0, 0, 0, 0, 0, 0, 0, 0, -1, 2, -1, 0, 1, -2, 1,
84 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, -1, 0, 0, -1, 1,
85 0, 0, 0, 0, 0, 0, -1, 1, 0, 0, 2, -2, 0, 0, -1, 1
90 int glatnr(int *global_atom_index, int i)
94 if (global_atom_index == NULL)
100 atnr = global_atom_index[i] + 1;
106 static int pbc_rvec_sub(const t_pbc *pbc, const rvec xi, const rvec xj, rvec dx)
110 return pbc_dx_aiuc(pbc, xi, xj, dx);
114 rvec_sub(xi, xj, dx);
121 /* SIMD PBC data structure, containing 1/boxdiag and the box vectors */
134 /* Set the SIMD pbc data from a normal t_pbc struct */
135 static void set_pbc_simd(const t_pbc *pbc, pbc_simd_t *pbc_simd)
140 /* Setting inv_bdiag to 0 effectively turns off PBC */
141 clear_rvec(inv_bdiag);
144 for (d = 0; d < pbc->ndim_ePBC; d++)
146 inv_bdiag[d] = 1.0/pbc->box[d][d];
150 pbc_simd->inv_bzz = gmx_set1_pr(inv_bdiag[ZZ]);
151 pbc_simd->inv_byy = gmx_set1_pr(inv_bdiag[YY]);
152 pbc_simd->inv_bxx = gmx_set1_pr(inv_bdiag[XX]);
156 pbc_simd->bzx = gmx_set1_pr(pbc->box[ZZ][XX]);
157 pbc_simd->bzy = gmx_set1_pr(pbc->box[ZZ][YY]);
158 pbc_simd->bzz = gmx_set1_pr(pbc->box[ZZ][ZZ]);
159 pbc_simd->byx = gmx_set1_pr(pbc->box[YY][XX]);
160 pbc_simd->byy = gmx_set1_pr(pbc->box[YY][YY]);
161 pbc_simd->bxx = gmx_set1_pr(pbc->box[XX][XX]);
165 pbc_simd->bzx = gmx_setzero_pr();
166 pbc_simd->bzy = gmx_setzero_pr();
167 pbc_simd->bzz = gmx_setzero_pr();
168 pbc_simd->byx = gmx_setzero_pr();
169 pbc_simd->byy = gmx_setzero_pr();
170 pbc_simd->bxx = gmx_setzero_pr();
174 /* Correct distance vector *dx,*dy,*dz for PBC using SIMD */
175 static gmx_inline void
176 pbc_dx_simd(gmx_mm_pr *dx, gmx_mm_pr *dy, gmx_mm_pr *dz,
177 const pbc_simd_t *pbc)
181 sh = gmx_round_pr(gmx_mul_pr(*dz, pbc->inv_bzz));
182 *dx = gmx_nmsub_pr(sh, pbc->bzx, *dx);
183 *dy = gmx_nmsub_pr(sh, pbc->bzy, *dy);
184 *dz = gmx_nmsub_pr(sh, pbc->bzz, *dz);
186 sh = gmx_round_pr(gmx_mul_pr(*dy, pbc->inv_byy));
187 *dx = gmx_nmsub_pr(sh, pbc->byx, *dx);
188 *dy = gmx_nmsub_pr(sh, pbc->byy, *dy);
190 sh = gmx_round_pr(gmx_mul_pr(*dx, pbc->inv_bxx));
191 *dx = gmx_nmsub_pr(sh, pbc->bxx, *dx);
194 #endif /* SIMD_BONDEDS */
197 * Morse potential bond by Frank Everdij
199 * Three parameters needed:
201 * b0 = equilibrium distance in nm
202 * be = beta in nm^-1 (actually, it's nu_e*Sqrt(2*pi*pi*mu/D_e))
203 * cb = well depth in kJ/mol
205 * Note: the potential is referenced to be +cb at infinite separation
206 * and zero at the equilibrium distance!
209 real morse_bonds(int nbonds,
210 const t_iatom forceatoms[], const t_iparams forceparams[],
211 const rvec x[], rvec f[], rvec fshift[],
212 const t_pbc *pbc, const t_graph *g,
213 real lambda, real *dvdlambda,
214 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
215 int gmx_unused *global_atom_index)
217 const real one = 1.0;
218 const real two = 2.0;
219 real dr, dr2, temp, omtemp, cbomtemp, fbond, vbond, fij, vtot;
220 real b0, be, cb, b0A, beA, cbA, b0B, beB, cbB, L1;
222 int i, m, ki, type, ai, aj;
226 for (i = 0; (i < nbonds); )
228 type = forceatoms[i++];
229 ai = forceatoms[i++];
230 aj = forceatoms[i++];
232 b0A = forceparams[type].morse.b0A;
233 beA = forceparams[type].morse.betaA;
234 cbA = forceparams[type].morse.cbA;
236 b0B = forceparams[type].morse.b0B;
237 beB = forceparams[type].morse.betaB;
238 cbB = forceparams[type].morse.cbB;
240 L1 = one-lambda; /* 1 */
241 b0 = L1*b0A + lambda*b0B; /* 3 */
242 be = L1*beA + lambda*beB; /* 3 */
243 cb = L1*cbA + lambda*cbB; /* 3 */
245 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
246 dr2 = iprod(dx, dx); /* 5 */
247 dr = dr2*gmx_invsqrt(dr2); /* 10 */
248 temp = exp(-be*(dr-b0)); /* 12 */
252 /* bonds are constrainted. This may _not_ include bond constraints if they are lambda dependent */
253 *dvdlambda += cbB-cbA;
257 omtemp = one-temp; /* 1 */
258 cbomtemp = cb*omtemp; /* 1 */
259 vbond = cbomtemp*omtemp; /* 1 */
260 fbond = -two*be*temp*cbomtemp*gmx_invsqrt(dr2); /* 9 */
261 vtot += vbond; /* 1 */
263 *dvdlambda += (cbB - cbA) * omtemp * omtemp - (2-2*omtemp)*omtemp * cb * ((b0B-b0A)*be - (beB-beA)*(dr-b0)); /* 15 */
267 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
271 for (m = 0; (m < DIM); m++) /* 15 */
276 fshift[ki][m] += fij;
277 fshift[CENTRAL][m] -= fij;
283 real cubic_bonds(int nbonds,
284 const t_iatom forceatoms[], const t_iparams forceparams[],
285 const rvec x[], rvec f[], rvec fshift[],
286 const t_pbc *pbc, const t_graph *g,
287 real gmx_unused lambda, real gmx_unused *dvdlambda,
288 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
289 int gmx_unused *global_atom_index)
291 const real three = 3.0;
292 const real two = 2.0;
294 real dr, dr2, dist, kdist, kdist2, fbond, vbond, fij, vtot;
296 int i, m, ki, type, ai, aj;
300 for (i = 0; (i < nbonds); )
302 type = forceatoms[i++];
303 ai = forceatoms[i++];
304 aj = forceatoms[i++];
306 b0 = forceparams[type].cubic.b0;
307 kb = forceparams[type].cubic.kb;
308 kcub = forceparams[type].cubic.kcub;
310 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
311 dr2 = iprod(dx, dx); /* 5 */
318 dr = dr2*gmx_invsqrt(dr2); /* 10 */
323 vbond = kdist2 + kcub*kdist2*dist;
324 fbond = -(two*kdist + three*kdist2*kcub)/dr;
326 vtot += vbond; /* 21 */
330 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
333 for (m = 0; (m < DIM); m++) /* 15 */
338 fshift[ki][m] += fij;
339 fshift[CENTRAL][m] -= fij;
345 real FENE_bonds(int nbonds,
346 const t_iatom forceatoms[], const t_iparams forceparams[],
347 const rvec x[], rvec f[], rvec fshift[],
348 const t_pbc *pbc, const t_graph *g,
349 real gmx_unused lambda, real gmx_unused *dvdlambda,
350 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
351 int *global_atom_index)
353 const real half = 0.5;
354 const real one = 1.0;
356 real dr, dr2, bm2, omdr2obm2, fbond, vbond, fij, vtot;
358 int i, m, ki, type, ai, aj;
362 for (i = 0; (i < nbonds); )
364 type = forceatoms[i++];
365 ai = forceatoms[i++];
366 aj = forceatoms[i++];
368 bm = forceparams[type].fene.bm;
369 kb = forceparams[type].fene.kb;
371 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
372 dr2 = iprod(dx, dx); /* 5 */
384 "r^2 (%f) >= bm^2 (%f) in FENE bond between atoms %d and %d",
386 glatnr(global_atom_index, ai),
387 glatnr(global_atom_index, aj));
390 omdr2obm2 = one - dr2/bm2;
392 vbond = -half*kb*bm2*log(omdr2obm2);
393 fbond = -kb/omdr2obm2;
395 vtot += vbond; /* 35 */
399 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
402 for (m = 0; (m < DIM); m++) /* 15 */
407 fshift[ki][m] += fij;
408 fshift[CENTRAL][m] -= fij;
414 real harmonic(real kA, real kB, real xA, real xB, real x, real lambda,
417 const real half = 0.5;
418 real L1, kk, x0, dx, dx2;
419 real v, f, dvdlambda;
422 kk = L1*kA+lambda*kB;
423 x0 = L1*xA+lambda*xB;
430 dvdlambda = half*(kB-kA)*dx2 + (xA-xB)*kk*dx;
437 /* That was 19 flops */
441 real bonds(int nbonds,
442 const t_iatom forceatoms[], const t_iparams forceparams[],
443 const rvec x[], rvec f[], rvec fshift[],
444 const t_pbc *pbc, const t_graph *g,
445 real lambda, real *dvdlambda,
446 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
447 int gmx_unused *global_atom_index)
449 int i, m, ki, ai, aj, type;
450 real dr, dr2, fbond, vbond, fij, vtot;
455 for (i = 0; (i < nbonds); )
457 type = forceatoms[i++];
458 ai = forceatoms[i++];
459 aj = forceatoms[i++];
461 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
462 dr2 = iprod(dx, dx); /* 5 */
463 dr = dr2*gmx_invsqrt(dr2); /* 10 */
465 *dvdlambda += harmonic(forceparams[type].harmonic.krA,
466 forceparams[type].harmonic.krB,
467 forceparams[type].harmonic.rA,
468 forceparams[type].harmonic.rB,
469 dr, lambda, &vbond, &fbond); /* 19 */
477 vtot += vbond; /* 1*/
478 fbond *= gmx_invsqrt(dr2); /* 6 */
482 fprintf(debug, "BONDS: dr = %10g vbond = %10g fbond = %10g\n",
488 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
491 for (m = 0; (m < DIM); m++) /* 15 */
496 fshift[ki][m] += fij;
497 fshift[CENTRAL][m] -= fij;
503 real restraint_bonds(int nbonds,
504 const t_iatom forceatoms[], const t_iparams forceparams[],
505 const rvec x[], rvec f[], rvec fshift[],
506 const t_pbc *pbc, const t_graph *g,
507 real lambda, real *dvdlambda,
508 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
509 int gmx_unused *global_atom_index)
511 int i, m, ki, ai, aj, type;
512 real dr, dr2, fbond, vbond, fij, vtot;
514 real low, dlow, up1, dup1, up2, dup2, k, dk;
522 for (i = 0; (i < nbonds); )
524 type = forceatoms[i++];
525 ai = forceatoms[i++];
526 aj = forceatoms[i++];
528 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
529 dr2 = iprod(dx, dx); /* 5 */
530 dr = dr2*gmx_invsqrt(dr2); /* 10 */
532 low = L1*forceparams[type].restraint.lowA + lambda*forceparams[type].restraint.lowB;
533 dlow = -forceparams[type].restraint.lowA + forceparams[type].restraint.lowB;
534 up1 = L1*forceparams[type].restraint.up1A + lambda*forceparams[type].restraint.up1B;
535 dup1 = -forceparams[type].restraint.up1A + forceparams[type].restraint.up1B;
536 up2 = L1*forceparams[type].restraint.up2A + lambda*forceparams[type].restraint.up2B;
537 dup2 = -forceparams[type].restraint.up2A + forceparams[type].restraint.up2B;
538 k = L1*forceparams[type].restraint.kA + lambda*forceparams[type].restraint.kB;
539 dk = -forceparams[type].restraint.kA + forceparams[type].restraint.kB;
548 *dvdlambda += 0.5*dk*drh2 - k*dlow*drh;
561 *dvdlambda += 0.5*dk*drh2 - k*dup1*drh;
566 vbond = k*(up2 - up1)*(0.5*(up2 - up1) + drh);
567 fbond = -k*(up2 - up1);
568 *dvdlambda += dk*(up2 - up1)*(0.5*(up2 - up1) + drh)
569 + k*(dup2 - dup1)*(up2 - up1 + drh)
570 - k*(up2 - up1)*dup2;
578 vtot += vbond; /* 1*/
579 fbond *= gmx_invsqrt(dr2); /* 6 */
583 fprintf(debug, "BONDS: dr = %10g vbond = %10g fbond = %10g\n",
589 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
592 for (m = 0; (m < DIM); m++) /* 15 */
597 fshift[ki][m] += fij;
598 fshift[CENTRAL][m] -= fij;
605 real polarize(int nbonds,
606 const t_iatom forceatoms[], const t_iparams forceparams[],
607 const rvec x[], rvec f[], rvec fshift[],
608 const t_pbc *pbc, const t_graph *g,
609 real lambda, real *dvdlambda,
610 const t_mdatoms *md, t_fcdata gmx_unused *fcd,
611 int gmx_unused *global_atom_index)
613 int i, m, ki, ai, aj, type;
614 real dr, dr2, fbond, vbond, fij, vtot, ksh;
619 for (i = 0; (i < nbonds); )
621 type = forceatoms[i++];
622 ai = forceatoms[i++];
623 aj = forceatoms[i++];
624 ksh = sqr(md->chargeA[aj])*ONE_4PI_EPS0/forceparams[type].polarize.alpha;
627 fprintf(debug, "POL: local ai = %d aj = %d ksh = %.3f\n", ai, aj, ksh);
630 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
631 dr2 = iprod(dx, dx); /* 5 */
632 dr = dr2*gmx_invsqrt(dr2); /* 10 */
634 *dvdlambda += harmonic(ksh, ksh, 0, 0, dr, lambda, &vbond, &fbond); /* 19 */
641 vtot += vbond; /* 1*/
642 fbond *= gmx_invsqrt(dr2); /* 6 */
646 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
649 for (m = 0; (m < DIM); m++) /* 15 */
654 fshift[ki][m] += fij;
655 fshift[CENTRAL][m] -= fij;
661 real anharm_polarize(int nbonds,
662 const t_iatom forceatoms[], const t_iparams forceparams[],
663 const rvec x[], rvec f[], rvec fshift[],
664 const t_pbc *pbc, const t_graph *g,
665 real lambda, real *dvdlambda,
666 const t_mdatoms *md, t_fcdata gmx_unused *fcd,
667 int gmx_unused *global_atom_index)
669 int i, m, ki, ai, aj, type;
670 real dr, dr2, fbond, vbond, fij, vtot, ksh, khyp, drcut, ddr, ddr3;
675 for (i = 0; (i < nbonds); )
677 type = forceatoms[i++];
678 ai = forceatoms[i++];
679 aj = forceatoms[i++];
680 ksh = sqr(md->chargeA[aj])*ONE_4PI_EPS0/forceparams[type].anharm_polarize.alpha; /* 7*/
681 khyp = forceparams[type].anharm_polarize.khyp;
682 drcut = forceparams[type].anharm_polarize.drcut;
685 fprintf(debug, "POL: local ai = %d aj = %d ksh = %.3f\n", ai, aj, ksh);
688 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
689 dr2 = iprod(dx, dx); /* 5 */
690 dr = dr2*gmx_invsqrt(dr2); /* 10 */
692 *dvdlambda += harmonic(ksh, ksh, 0, 0, dr, lambda, &vbond, &fbond); /* 19 */
703 vbond += khyp*ddr*ddr3;
704 fbond -= 4*khyp*ddr3;
706 fbond *= gmx_invsqrt(dr2); /* 6 */
707 vtot += vbond; /* 1*/
711 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
714 for (m = 0; (m < DIM); m++) /* 15 */
719 fshift[ki][m] += fij;
720 fshift[CENTRAL][m] -= fij;
726 real water_pol(int nbonds,
727 const t_iatom forceatoms[], const t_iparams forceparams[],
728 const rvec x[], rvec f[], rvec gmx_unused fshift[],
729 const t_pbc gmx_unused *pbc, const t_graph gmx_unused *g,
730 real gmx_unused lambda, real gmx_unused *dvdlambda,
731 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
732 int gmx_unused *global_atom_index)
734 /* This routine implements anisotropic polarizibility for water, through
735 * a shell connected to a dummy with spring constant that differ in the
736 * three spatial dimensions in the molecular frame.
738 int i, m, aO, aH1, aH2, aD, aS, type, type0;
739 rvec dOH1, dOH2, dHH, dOD, dDS, nW, kk, dx, kdx, proj;
743 real vtot, fij, r_HH, r_OD, r_nW, tx, ty, tz, qS;
748 type0 = forceatoms[0];
750 qS = md->chargeA[aS];
751 kk[XX] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_x;
752 kk[YY] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_y;
753 kk[ZZ] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_z;
754 r_HH = 1.0/forceparams[type0].wpol.rHH;
755 r_OD = 1.0/forceparams[type0].wpol.rOD;
758 fprintf(debug, "WPOL: qS = %10.5f aS = %5d\n", qS, aS);
759 fprintf(debug, "WPOL: kk = %10.3f %10.3f %10.3f\n",
760 kk[XX], kk[YY], kk[ZZ]);
761 fprintf(debug, "WPOL: rOH = %10.3f rHH = %10.3f rOD = %10.3f\n",
762 forceparams[type0].wpol.rOH,
763 forceparams[type0].wpol.rHH,
764 forceparams[type0].wpol.rOD);
766 for (i = 0; (i < nbonds); i += 6)
768 type = forceatoms[i];
771 gmx_fatal(FARGS, "Sorry, type = %d, type0 = %d, file = %s, line = %d",
772 type, type0, __FILE__, __LINE__);
774 aO = forceatoms[i+1];
775 aH1 = forceatoms[i+2];
776 aH2 = forceatoms[i+3];
777 aD = forceatoms[i+4];
778 aS = forceatoms[i+5];
780 /* Compute vectors describing the water frame */
781 rvec_sub(x[aH1], x[aO], dOH1);
782 rvec_sub(x[aH2], x[aO], dOH2);
783 rvec_sub(x[aH2], x[aH1], dHH);
784 rvec_sub(x[aD], x[aO], dOD);
785 rvec_sub(x[aS], x[aD], dDS);
786 cprod(dOH1, dOH2, nW);
788 /* Compute inverse length of normal vector
789 * (this one could be precomputed, but I'm too lazy now)
791 r_nW = gmx_invsqrt(iprod(nW, nW));
792 /* This is for precision, but does not make a big difference,
795 r_OD = gmx_invsqrt(iprod(dOD, dOD));
797 /* Normalize the vectors in the water frame */
799 svmul(r_HH, dHH, dHH);
800 svmul(r_OD, dOD, dOD);
802 /* Compute displacement of shell along components of the vector */
803 dx[ZZ] = iprod(dDS, dOD);
804 /* Compute projection on the XY plane: dDS - dx[ZZ]*dOD */
805 for (m = 0; (m < DIM); m++)
807 proj[m] = dDS[m]-dx[ZZ]*dOD[m];
810 /*dx[XX] = iprod(dDS,nW);
811 dx[YY] = iprod(dDS,dHH);*/
812 dx[XX] = iprod(proj, nW);
813 for (m = 0; (m < DIM); m++)
815 proj[m] -= dx[XX]*nW[m];
817 dx[YY] = iprod(proj, dHH);
822 fprintf(debug, "WPOL: dx2=%10g dy2=%10g dz2=%10g sum=%10g dDS^2=%10g\n",
823 sqr(dx[XX]), sqr(dx[YY]), sqr(dx[ZZ]), iprod(dx, dx), iprod(dDS, dDS));
824 fprintf(debug, "WPOL: dHH=(%10g,%10g,%10g)\n", dHH[XX], dHH[YY], dHH[ZZ]);
825 fprintf(debug, "WPOL: dOD=(%10g,%10g,%10g), 1/r_OD = %10g\n",
826 dOD[XX], dOD[YY], dOD[ZZ], 1/r_OD);
827 fprintf(debug, "WPOL: nW =(%10g,%10g,%10g), 1/r_nW = %10g\n",
828 nW[XX], nW[YY], nW[ZZ], 1/r_nW);
829 fprintf(debug, "WPOL: dx =%10g, dy =%10g, dz =%10g\n",
830 dx[XX], dx[YY], dx[ZZ]);
831 fprintf(debug, "WPOL: dDSx=%10g, dDSy=%10g, dDSz=%10g\n",
832 dDS[XX], dDS[YY], dDS[ZZ]);
835 /* Now compute the forces and energy */
836 kdx[XX] = kk[XX]*dx[XX];
837 kdx[YY] = kk[YY]*dx[YY];
838 kdx[ZZ] = kk[ZZ]*dx[ZZ];
839 vtot += iprod(dx, kdx);
840 for (m = 0; (m < DIM); m++)
842 /* This is a tensor operation but written out for speed */
856 fprintf(debug, "WPOL: vwpol=%g\n", 0.5*iprod(dx, kdx));
857 fprintf(debug, "WPOL: df = (%10g, %10g, %10g)\n", df[XX], df[YY], df[ZZ]);
865 static real do_1_thole(const rvec xi, const rvec xj, rvec fi, rvec fj,
866 const t_pbc *pbc, real qq,
867 rvec fshift[], real afac)
870 real r12sq, r12_1, r12n, r12bar, v0, v1, fscal, ebar, fff;
873 t = pbc_rvec_sub(pbc, xi, xj, r12); /* 3 */
875 r12sq = iprod(r12, r12); /* 5 */
876 r12_1 = gmx_invsqrt(r12sq); /* 5 */
877 r12bar = afac/r12_1; /* 5 */
878 v0 = qq*ONE_4PI_EPS0*r12_1; /* 2 */
879 ebar = exp(-r12bar); /* 5 */
880 v1 = (1-(1+0.5*r12bar)*ebar); /* 4 */
881 fscal = ((v0*r12_1)*v1 - v0*0.5*afac*ebar*(r12bar+1))*r12_1; /* 9 */
884 fprintf(debug, "THOLE: v0 = %.3f v1 = %.3f r12= % .3f r12bar = %.3f fscal = %.3f ebar = %.3f\n", v0, v1, 1/r12_1, r12bar, fscal, ebar);
887 for (m = 0; (m < DIM); m++)
893 fshift[CENTRAL][m] -= fff;
896 return v0*v1; /* 1 */
900 real thole_pol(int nbonds,
901 const t_iatom forceatoms[], const t_iparams forceparams[],
902 const rvec x[], rvec f[], rvec fshift[],
903 const t_pbc *pbc, const t_graph gmx_unused *g,
904 real gmx_unused lambda, real gmx_unused *dvdlambda,
905 const t_mdatoms *md, t_fcdata gmx_unused *fcd,
906 int gmx_unused *global_atom_index)
908 /* Interaction between two pairs of particles with opposite charge */
909 int i, type, a1, da1, a2, da2;
910 real q1, q2, qq, a, al1, al2, afac;
913 for (i = 0; (i < nbonds); )
915 type = forceatoms[i++];
916 a1 = forceatoms[i++];
917 da1 = forceatoms[i++];
918 a2 = forceatoms[i++];
919 da2 = forceatoms[i++];
920 q1 = md->chargeA[da1];
921 q2 = md->chargeA[da2];
922 a = forceparams[type].thole.a;
923 al1 = forceparams[type].thole.alpha1;
924 al2 = forceparams[type].thole.alpha2;
926 afac = a*pow(al1*al2, -1.0/6.0);
927 V += do_1_thole(x[a1], x[a2], f[a1], f[a2], pbc, qq, fshift, afac);
928 V += do_1_thole(x[da1], x[a2], f[da1], f[a2], pbc, -qq, fshift, afac);
929 V += do_1_thole(x[a1], x[da2], f[a1], f[da2], pbc, -qq, fshift, afac);
930 V += do_1_thole(x[da1], x[da2], f[da1], f[da2], pbc, qq, fshift, afac);
936 real bond_angle(const rvec xi, const rvec xj, const rvec xk, const t_pbc *pbc,
937 rvec r_ij, rvec r_kj, real *costh,
939 /* Return value is the angle between the bonds i-j and j-k */
944 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
945 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
947 *costh = cos_angle(r_ij, r_kj); /* 25 */
948 th = acos(*costh); /* 10 */
953 real angles(int nbonds,
954 const t_iatom forceatoms[], const t_iparams forceparams[],
955 const rvec x[], rvec f[], rvec fshift[],
956 const t_pbc *pbc, const t_graph *g,
957 real lambda, real *dvdlambda,
958 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
959 int gmx_unused *global_atom_index)
961 int i, ai, aj, ak, t1, t2, type;
963 real cos_theta, cos_theta2, theta, dVdt, va, vtot;
964 ivec jt, dt_ij, dt_kj;
967 for (i = 0; i < nbonds; )
969 type = forceatoms[i++];
970 ai = forceatoms[i++];
971 aj = forceatoms[i++];
972 ak = forceatoms[i++];
974 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
975 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
977 *dvdlambda += harmonic(forceparams[type].harmonic.krA,
978 forceparams[type].harmonic.krB,
979 forceparams[type].harmonic.rA*DEG2RAD,
980 forceparams[type].harmonic.rB*DEG2RAD,
981 theta, lambda, &va, &dVdt); /* 21 */
984 cos_theta2 = sqr(cos_theta);
994 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
995 sth = st*cos_theta; /* 1 */
999 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1000 theta*RAD2DEG, va, dVdt);
1003 nrij2 = iprod(r_ij, r_ij); /* 5 */
1004 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1006 nrij_1 = gmx_invsqrt(nrij2); /* 10 */
1007 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1009 cik = st*nrij_1*nrkj_1; /* 2 */
1010 cii = sth*nrij_1*nrij_1; /* 2 */
1011 ckk = sth*nrkj_1*nrkj_1; /* 2 */
1013 for (m = 0; m < DIM; m++)
1015 f_i[m] = -(cik*r_kj[m] - cii*r_ij[m]);
1016 f_k[m] = -(cik*r_ij[m] - ckk*r_kj[m]);
1017 f_j[m] = -f_i[m] - f_k[m];
1024 copy_ivec(SHIFT_IVEC(g, aj), jt);
1026 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1027 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1028 t1 = IVEC2IS(dt_ij);
1029 t2 = IVEC2IS(dt_kj);
1031 rvec_inc(fshift[t1], f_i);
1032 rvec_inc(fshift[CENTRAL], f_j);
1033 rvec_inc(fshift[t2], f_k);
1042 /* As angles, but using SIMD to calculate many dihedrals at once.
1043 * This routines does not calculate energies and shift forces.
1045 static gmx_inline void
1046 angles_noener_simd(int nbonds,
1047 const t_iatom forceatoms[], const t_iparams forceparams[],
1048 const rvec x[], rvec f[],
1049 const t_pbc *pbc, const t_graph gmx_unused *g,
1050 real gmx_unused lambda,
1051 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1052 int gmx_unused *global_atom_index)
1054 #define UNROLL GMX_SIMD_WIDTH_HERE
1057 int type, ai[UNROLL], aj[UNROLL], ak[UNROLL];
1058 real coeff_array[2*UNROLL+UNROLL], *coeff;
1059 real dr_array[2*DIM*UNROLL+UNROLL], *dr;
1060 real f_buf_array[6*UNROLL+UNROLL], *f_buf;
1061 gmx_mm_pr k_S, theta0_S;
1062 gmx_mm_pr rijx_S, rijy_S, rijz_S;
1063 gmx_mm_pr rkjx_S, rkjy_S, rkjz_S;
1065 gmx_mm_pr min_one_plus_eps_S;
1066 gmx_mm_pr rij_rkj_S;
1067 gmx_mm_pr nrij2_S, nrij_1_S;
1068 gmx_mm_pr nrkj2_S, nrkj_1_S;
1069 gmx_mm_pr cos_S, invsin_S;
1071 gmx_mm_pr st_S, sth_S;
1072 gmx_mm_pr cik_S, cii_S, ckk_S;
1073 gmx_mm_pr f_ix_S, f_iy_S, f_iz_S;
1074 gmx_mm_pr f_kx_S, f_ky_S, f_kz_S;
1075 pbc_simd_t pbc_simd;
1077 /* Ensure register memory alignment */
1078 coeff = gmx_simd_align_real(coeff_array);
1079 dr = gmx_simd_align_real(dr_array);
1080 f_buf = gmx_simd_align_real(f_buf_array);
1082 set_pbc_simd(pbc, &pbc_simd);
1084 one_S = gmx_set1_pr(1.0);
1086 /* The smallest number > -1 */
1087 min_one_plus_eps_S = gmx_set1_pr(-1.0 + 2*GMX_REAL_EPS);
1089 /* nbonds is the number of angles times nfa1, here we step UNROLL angles */
1090 for (i = 0; (i < nbonds); i += UNROLL*nfa1)
1092 /* Collect atoms for UNROLL angles.
1093 * iu indexes into forceatoms, we should not let iu go beyond nbonds.
1096 for (s = 0; s < UNROLL; s++)
1098 type = forceatoms[iu];
1099 ai[s] = forceatoms[iu+1];
1100 aj[s] = forceatoms[iu+2];
1101 ak[s] = forceatoms[iu+3];
1103 coeff[s] = forceparams[type].harmonic.krA;
1104 coeff[UNROLL+s] = forceparams[type].harmonic.rA*DEG2RAD;
1106 /* If you can't use pbc_dx_simd below for PBC, e.g. because
1107 * you can't round in SIMD, use pbc_rvec_sub here.
1109 /* Store the non PBC corrected distances packed and aligned */
1110 for (m = 0; m < DIM; m++)
1112 dr[s + m *UNROLL] = x[ai[s]][m] - x[aj[s]][m];
1113 dr[s + (DIM+m)*UNROLL] = x[ak[s]][m] - x[aj[s]][m];
1116 /* At the end fill the arrays with identical entries */
1117 if (iu + nfa1 < nbonds)
1123 k_S = gmx_load_pr(coeff);
1124 theta0_S = gmx_load_pr(coeff+UNROLL);
1126 rijx_S = gmx_load_pr(dr + 0*UNROLL);
1127 rijy_S = gmx_load_pr(dr + 1*UNROLL);
1128 rijz_S = gmx_load_pr(dr + 2*UNROLL);
1129 rkjx_S = gmx_load_pr(dr + 3*UNROLL);
1130 rkjy_S = gmx_load_pr(dr + 4*UNROLL);
1131 rkjz_S = gmx_load_pr(dr + 5*UNROLL);
1133 pbc_dx_simd(&rijx_S, &rijy_S, &rijz_S, &pbc_simd);
1134 pbc_dx_simd(&rkjx_S, &rkjy_S, &rkjz_S, &pbc_simd);
1136 rij_rkj_S = gmx_iprod_pr(rijx_S, rijy_S, rijz_S,
1137 rkjx_S, rkjy_S, rkjz_S);
1139 nrij2_S = gmx_norm2_pr(rijx_S, rijy_S, rijz_S);
1140 nrkj2_S = gmx_norm2_pr(rkjx_S, rkjy_S, rkjz_S);
1142 nrij_1_S = gmx_invsqrt_pr(nrij2_S);
1143 nrkj_1_S = gmx_invsqrt_pr(nrkj2_S);
1145 cos_S = gmx_mul_pr(rij_rkj_S, gmx_mul_pr(nrij_1_S, nrkj_1_S));
1147 /* To allow for 180 degrees, we take the max of cos and -1 + 1bit,
1148 * so we can safely get the 1/sin from 1/sqrt(1 - cos^2).
1149 * This also ensures that rounding errors would cause the argument
1150 * of gmx_acos_pr to be < -1.
1151 * Note that we do not take precautions for cos(0)=1, so the outer
1152 * atoms in an angle should not be on top of each other.
1154 cos_S = gmx_max_pr(cos_S, min_one_plus_eps_S);
1156 theta_S = gmx_acos_pr(cos_S);
1158 invsin_S = gmx_invsqrt_pr(gmx_sub_pr(one_S, gmx_mul_pr(cos_S, cos_S)));
1160 st_S = gmx_mul_pr(gmx_mul_pr(k_S, gmx_sub_pr(theta0_S, theta_S)),
1162 sth_S = gmx_mul_pr(st_S, cos_S);
1164 cik_S = gmx_mul_pr(st_S, gmx_mul_pr(nrij_1_S, nrkj_1_S));
1165 cii_S = gmx_mul_pr(sth_S, gmx_mul_pr(nrij_1_S, nrij_1_S));
1166 ckk_S = gmx_mul_pr(sth_S, gmx_mul_pr(nrkj_1_S, nrkj_1_S));
1168 f_ix_S = gmx_mul_pr(cii_S, rijx_S);
1169 f_ix_S = gmx_nmsub_pr(cik_S, rkjx_S, f_ix_S);
1170 f_iy_S = gmx_mul_pr(cii_S, rijy_S);
1171 f_iy_S = gmx_nmsub_pr(cik_S, rkjy_S, f_iy_S);
1172 f_iz_S = gmx_mul_pr(cii_S, rijz_S);
1173 f_iz_S = gmx_nmsub_pr(cik_S, rkjz_S, f_iz_S);
1174 f_kx_S = gmx_mul_pr(ckk_S, rkjx_S);
1175 f_kx_S = gmx_nmsub_pr(cik_S, rijx_S, f_kx_S);
1176 f_ky_S = gmx_mul_pr(ckk_S, rkjy_S);
1177 f_ky_S = gmx_nmsub_pr(cik_S, rijy_S, f_ky_S);
1178 f_kz_S = gmx_mul_pr(ckk_S, rkjz_S);
1179 f_kz_S = gmx_nmsub_pr(cik_S, rijz_S, f_kz_S);
1181 gmx_store_pr(f_buf + 0*UNROLL, f_ix_S);
1182 gmx_store_pr(f_buf + 1*UNROLL, f_iy_S);
1183 gmx_store_pr(f_buf + 2*UNROLL, f_iz_S);
1184 gmx_store_pr(f_buf + 3*UNROLL, f_kx_S);
1185 gmx_store_pr(f_buf + 4*UNROLL, f_ky_S);
1186 gmx_store_pr(f_buf + 5*UNROLL, f_kz_S);
1192 for (m = 0; m < DIM; m++)
1194 f[ai[s]][m] += f_buf[s + m*UNROLL];
1195 f[aj[s]][m] -= f_buf[s + m*UNROLL] + f_buf[s + (DIM+m)*UNROLL];
1196 f[ak[s]][m] += f_buf[s + (DIM+m)*UNROLL];
1201 while (s < UNROLL && iu < nbonds);
1206 #endif /* SIMD_BONDEDS */
1208 real linear_angles(int nbonds,
1209 const t_iatom forceatoms[], const t_iparams forceparams[],
1210 const rvec x[], rvec f[], rvec fshift[],
1211 const t_pbc *pbc, const t_graph *g,
1212 real lambda, real *dvdlambda,
1213 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1214 int gmx_unused *global_atom_index)
1216 int i, m, ai, aj, ak, t1, t2, type;
1218 real L1, kA, kB, aA, aB, dr, dr2, va, vtot, a, b, klin;
1219 ivec jt, dt_ij, dt_kj;
1220 rvec r_ij, r_kj, r_ik, dx;
1224 for (i = 0; (i < nbonds); )
1226 type = forceatoms[i++];
1227 ai = forceatoms[i++];
1228 aj = forceatoms[i++];
1229 ak = forceatoms[i++];
1231 kA = forceparams[type].linangle.klinA;
1232 kB = forceparams[type].linangle.klinB;
1233 klin = L1*kA + lambda*kB;
1235 aA = forceparams[type].linangle.aA;
1236 aB = forceparams[type].linangle.aB;
1237 a = L1*aA+lambda*aB;
1240 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
1241 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
1242 rvec_sub(r_ij, r_kj, r_ik);
1245 for (m = 0; (m < DIM); m++)
1247 dr = -a * r_ij[m] - b * r_kj[m];
1252 f_j[m] = -(f_i[m]+f_k[m]);
1258 *dvdlambda += 0.5*(kB-kA)*dr2 + klin*(aB-aA)*iprod(dx, r_ik);
1264 copy_ivec(SHIFT_IVEC(g, aj), jt);
1266 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1267 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1268 t1 = IVEC2IS(dt_ij);
1269 t2 = IVEC2IS(dt_kj);
1271 rvec_inc(fshift[t1], f_i);
1272 rvec_inc(fshift[CENTRAL], f_j);
1273 rvec_inc(fshift[t2], f_k);
1278 real urey_bradley(int nbonds,
1279 const t_iatom forceatoms[], const t_iparams forceparams[],
1280 const rvec x[], rvec f[], rvec fshift[],
1281 const t_pbc *pbc, const t_graph *g,
1282 real lambda, real *dvdlambda,
1283 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1284 int gmx_unused *global_atom_index)
1286 int i, m, ai, aj, ak, t1, t2, type, ki;
1287 rvec r_ij, r_kj, r_ik;
1288 real cos_theta, cos_theta2, theta;
1289 real dVdt, va, vtot, dr, dr2, vbond, fbond, fik;
1290 real kthA, th0A, kUBA, r13A, kthB, th0B, kUBB, r13B;
1291 ivec jt, dt_ij, dt_kj, dt_ik;
1294 for (i = 0; (i < nbonds); )
1296 type = forceatoms[i++];
1297 ai = forceatoms[i++];
1298 aj = forceatoms[i++];
1299 ak = forceatoms[i++];
1300 th0A = forceparams[type].u_b.thetaA*DEG2RAD;
1301 kthA = forceparams[type].u_b.kthetaA;
1302 r13A = forceparams[type].u_b.r13A;
1303 kUBA = forceparams[type].u_b.kUBA;
1304 th0B = forceparams[type].u_b.thetaB*DEG2RAD;
1305 kthB = forceparams[type].u_b.kthetaB;
1306 r13B = forceparams[type].u_b.r13B;
1307 kUBB = forceparams[type].u_b.kUBB;
1309 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
1310 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
1312 *dvdlambda += harmonic(kthA, kthB, th0A, th0B, theta, lambda, &va, &dVdt); /* 21 */
1315 ki = pbc_rvec_sub(pbc, x[ai], x[ak], r_ik); /* 3 */
1316 dr2 = iprod(r_ik, r_ik); /* 5 */
1317 dr = dr2*gmx_invsqrt(dr2); /* 10 */
1319 *dvdlambda += harmonic(kUBA, kUBB, r13A, r13B, dr, lambda, &vbond, &fbond); /* 19 */
1321 cos_theta2 = sqr(cos_theta); /* 1 */
1329 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1330 sth = st*cos_theta; /* 1 */
1334 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1335 theta*RAD2DEG, va, dVdt);
1338 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1339 nrij2 = iprod(r_ij, r_ij);
1341 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
1342 cii = sth/nrij2; /* 10 */
1343 ckk = sth/nrkj2; /* 10 */
1345 for (m = 0; (m < DIM); m++) /* 39 */
1347 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
1348 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
1349 f_j[m] = -f_i[m]-f_k[m];
1356 copy_ivec(SHIFT_IVEC(g, aj), jt);
1358 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1359 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1360 t1 = IVEC2IS(dt_ij);
1361 t2 = IVEC2IS(dt_kj);
1363 rvec_inc(fshift[t1], f_i);
1364 rvec_inc(fshift[CENTRAL], f_j);
1365 rvec_inc(fshift[t2], f_k);
1367 /* Time for the bond calculations */
1373 vtot += vbond; /* 1*/
1374 fbond *= gmx_invsqrt(dr2); /* 6 */
1378 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ak), dt_ik);
1379 ki = IVEC2IS(dt_ik);
1381 for (m = 0; (m < DIM); m++) /* 15 */
1383 fik = fbond*r_ik[m];
1386 fshift[ki][m] += fik;
1387 fshift[CENTRAL][m] -= fik;
1393 real quartic_angles(int nbonds,
1394 const t_iatom forceatoms[], const t_iparams forceparams[],
1395 const rvec x[], rvec f[], rvec fshift[],
1396 const t_pbc *pbc, const t_graph *g,
1397 real gmx_unused lambda, real gmx_unused *dvdlambda,
1398 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1399 int gmx_unused *global_atom_index)
1401 int i, j, ai, aj, ak, t1, t2, type;
1403 real cos_theta, cos_theta2, theta, dt, dVdt, va, dtp, c, vtot;
1404 ivec jt, dt_ij, dt_kj;
1407 for (i = 0; (i < nbonds); )
1409 type = forceatoms[i++];
1410 ai = forceatoms[i++];
1411 aj = forceatoms[i++];
1412 ak = forceatoms[i++];
1414 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
1415 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
1417 dt = theta - forceparams[type].qangle.theta*DEG2RAD; /* 2 */
1420 va = forceparams[type].qangle.c[0];
1422 for (j = 1; j <= 4; j++)
1424 c = forceparams[type].qangle.c[j];
1433 cos_theta2 = sqr(cos_theta); /* 1 */
1442 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1443 sth = st*cos_theta; /* 1 */
1447 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1448 theta*RAD2DEG, va, dVdt);
1451 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1452 nrij2 = iprod(r_ij, r_ij);
1454 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
1455 cii = sth/nrij2; /* 10 */
1456 ckk = sth/nrkj2; /* 10 */
1458 for (m = 0; (m < DIM); m++) /* 39 */
1460 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
1461 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
1462 f_j[m] = -f_i[m]-f_k[m];
1469 copy_ivec(SHIFT_IVEC(g, aj), jt);
1471 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1472 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1473 t1 = IVEC2IS(dt_ij);
1474 t2 = IVEC2IS(dt_kj);
1476 rvec_inc(fshift[t1], f_i);
1477 rvec_inc(fshift[CENTRAL], f_j);
1478 rvec_inc(fshift[t2], f_k);
1484 real dih_angle(const rvec xi, const rvec xj, const rvec xk, const rvec xl,
1486 rvec r_ij, rvec r_kj, rvec r_kl, rvec m, rvec n,
1487 real *sign, int *t1, int *t2, int *t3)
1491 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
1492 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
1493 *t3 = pbc_rvec_sub(pbc, xk, xl, r_kl); /* 3 */
1495 cprod(r_ij, r_kj, m); /* 9 */
1496 cprod(r_kj, r_kl, n); /* 9 */
1497 phi = gmx_angle(m, n); /* 49 (assuming 25 for atan2) */
1498 ipr = iprod(r_ij, n); /* 5 */
1499 (*sign) = (ipr < 0.0) ? -1.0 : 1.0;
1500 phi = (*sign)*phi; /* 1 */
1508 /* As dih_angle above, but calculates 4 dihedral angles at once using SIMD,
1509 * also calculates the pre-factor required for the dihedral force update.
1510 * Note that bv and buf should be register aligned.
1512 static gmx_inline void
1513 dih_angle_simd(const rvec *x,
1514 const int *ai, const int *aj, const int *ak, const int *al,
1515 const pbc_simd_t *pbc,
1518 gmx_mm_pr *mx_S, gmx_mm_pr *my_S, gmx_mm_pr *mz_S,
1519 gmx_mm_pr *nx_S, gmx_mm_pr *ny_S, gmx_mm_pr *nz_S,
1520 gmx_mm_pr *nrkj_m2_S,
1521 gmx_mm_pr *nrkj_n2_S,
1525 #define UNROLL GMX_SIMD_WIDTH_HERE
1527 gmx_mm_pr rijx_S, rijy_S, rijz_S;
1528 gmx_mm_pr rkjx_S, rkjy_S, rkjz_S;
1529 gmx_mm_pr rklx_S, rkly_S, rklz_S;
1530 gmx_mm_pr cx_S, cy_S, cz_S;
1534 gmx_mm_pr iprm_S, iprn_S;
1535 gmx_mm_pr nrkj2_S, nrkj_1_S, nrkj_2_S, nrkj_S;
1538 gmx_mm_pr nrkj2_min_S;
1539 gmx_mm_pr real_eps_S;
1541 /* Used to avoid division by zero.
1542 * We take into acount that we multiply the result by real_eps_S.
1544 nrkj2_min_S = gmx_set1_pr(GMX_REAL_MIN/(2*GMX_REAL_EPS));
1546 /* The value of the last significant bit (GMX_REAL_EPS is half of that) */
1547 real_eps_S = gmx_set1_pr(2*GMX_REAL_EPS);
1549 for (s = 0; s < UNROLL; s++)
1551 /* If you can't use pbc_dx_simd below for PBC, e.g. because
1552 * you can't round in SIMD, use pbc_rvec_sub here.
1554 for (m = 0; m < DIM; m++)
1556 dr[s + (0*DIM + m)*UNROLL] = x[ai[s]][m] - x[aj[s]][m];
1557 dr[s + (1*DIM + m)*UNROLL] = x[ak[s]][m] - x[aj[s]][m];
1558 dr[s + (2*DIM + m)*UNROLL] = x[ak[s]][m] - x[al[s]][m];
1562 rijx_S = gmx_load_pr(dr + 0*UNROLL);
1563 rijy_S = gmx_load_pr(dr + 1*UNROLL);
1564 rijz_S = gmx_load_pr(dr + 2*UNROLL);
1565 rkjx_S = gmx_load_pr(dr + 3*UNROLL);
1566 rkjy_S = gmx_load_pr(dr + 4*UNROLL);
1567 rkjz_S = gmx_load_pr(dr + 5*UNROLL);
1568 rklx_S = gmx_load_pr(dr + 6*UNROLL);
1569 rkly_S = gmx_load_pr(dr + 7*UNROLL);
1570 rklz_S = gmx_load_pr(dr + 8*UNROLL);
1572 pbc_dx_simd(&rijx_S, &rijy_S, &rijz_S, pbc);
1573 pbc_dx_simd(&rkjx_S, &rkjy_S, &rkjz_S, pbc);
1574 pbc_dx_simd(&rklx_S, &rkly_S, &rklz_S, pbc);
1576 gmx_cprod_pr(rijx_S, rijy_S, rijz_S,
1577 rkjx_S, rkjy_S, rkjz_S,
1580 gmx_cprod_pr(rkjx_S, rkjy_S, rkjz_S,
1581 rklx_S, rkly_S, rklz_S,
1584 gmx_cprod_pr(*mx_S, *my_S, *mz_S,
1585 *nx_S, *ny_S, *nz_S,
1586 &cx_S, &cy_S, &cz_S);
1588 cn_S = gmx_sqrt_pr(gmx_norm2_pr(cx_S, cy_S, cz_S));
1590 s_S = gmx_iprod_pr(*mx_S, *my_S, *mz_S, *nx_S, *ny_S, *nz_S);
1592 /* Determine the dihedral angle, the sign might need correction */
1593 *phi_S = gmx_atan2_pr(cn_S, s_S);
1595 ipr_S = gmx_iprod_pr(rijx_S, rijy_S, rijz_S,
1596 *nx_S, *ny_S, *nz_S);
1598 iprm_S = gmx_norm2_pr(*mx_S, *my_S, *mz_S);
1599 iprn_S = gmx_norm2_pr(*nx_S, *ny_S, *nz_S);
1601 nrkj2_S = gmx_norm2_pr(rkjx_S, rkjy_S, rkjz_S);
1603 /* Avoid division by zero. When zero, the result is multiplied by 0
1604 * anyhow, so the 3 max below do not affect the final result.
1606 nrkj2_S = gmx_max_pr(nrkj2_S, nrkj2_min_S);
1607 nrkj_1_S = gmx_invsqrt_pr(nrkj2_S);
1608 nrkj_2_S = gmx_mul_pr(nrkj_1_S, nrkj_1_S);
1609 nrkj_S = gmx_mul_pr(nrkj2_S, nrkj_1_S);
1611 toler_S = gmx_mul_pr(nrkj2_S, real_eps_S);
1613 /* Here the plain-C code uses a conditional, but we can't do that in SIMD.
1614 * So we take a max with the tolerance instead. Since we multiply with
1615 * m or n later, the max does not affect the results.
1617 iprm_S = gmx_max_pr(iprm_S, toler_S);
1618 iprn_S = gmx_max_pr(iprn_S, toler_S);
1619 *nrkj_m2_S = gmx_mul_pr(nrkj_S, gmx_inv_pr(iprm_S));
1620 *nrkj_n2_S = gmx_mul_pr(nrkj_S, gmx_inv_pr(iprn_S));
1622 /* Set sign of phi_S with the sign of ipr_S; phi_S is currently positive */
1623 *phi_S = gmx_cpsgn_nonneg_pr(ipr_S, *phi_S);
1625 p_S = gmx_iprod_pr(rijx_S, rijy_S, rijz_S,
1626 rkjx_S, rkjy_S, rkjz_S);
1627 p_S = gmx_mul_pr(p_S, nrkj_2_S);
1629 q_S = gmx_iprod_pr(rklx_S, rkly_S, rklz_S,
1630 rkjx_S, rkjy_S, rkjz_S);
1631 q_S = gmx_mul_pr(q_S, nrkj_2_S);
1633 gmx_store_pr(p, p_S);
1634 gmx_store_pr(q, q_S);
1638 #endif /* SIMD_BONDEDS */
1641 void do_dih_fup(int i, int j, int k, int l, real ddphi,
1642 rvec r_ij, rvec r_kj, rvec r_kl,
1643 rvec m, rvec n, rvec f[], rvec fshift[],
1644 const t_pbc *pbc, const t_graph *g,
1645 const rvec x[], int t1, int t2, int t3)
1648 rvec f_i, f_j, f_k, f_l;
1649 rvec uvec, vvec, svec, dx_jl;
1650 real iprm, iprn, nrkj, nrkj2, nrkj_1, nrkj_2;
1651 real a, b, p, q, toler;
1652 ivec jt, dt_ij, dt_kj, dt_lj;
1654 iprm = iprod(m, m); /* 5 */
1655 iprn = iprod(n, n); /* 5 */
1656 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1657 toler = nrkj2*GMX_REAL_EPS;
1658 if ((iprm > toler) && (iprn > toler))
1660 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1661 nrkj_2 = nrkj_1*nrkj_1; /* 1 */
1662 nrkj = nrkj2*nrkj_1; /* 1 */
1663 a = -ddphi*nrkj/iprm; /* 11 */
1664 svmul(a, m, f_i); /* 3 */
1665 b = ddphi*nrkj/iprn; /* 11 */
1666 svmul(b, n, f_l); /* 3 */
1667 p = iprod(r_ij, r_kj); /* 5 */
1668 p *= nrkj_2; /* 1 */
1669 q = iprod(r_kl, r_kj); /* 5 */
1670 q *= nrkj_2; /* 1 */
1671 svmul(p, f_i, uvec); /* 3 */
1672 svmul(q, f_l, vvec); /* 3 */
1673 rvec_sub(uvec, vvec, svec); /* 3 */
1674 rvec_sub(f_i, svec, f_j); /* 3 */
1675 rvec_add(f_l, svec, f_k); /* 3 */
1676 rvec_inc(f[i], f_i); /* 3 */
1677 rvec_dec(f[j], f_j); /* 3 */
1678 rvec_dec(f[k], f_k); /* 3 */
1679 rvec_inc(f[l], f_l); /* 3 */
1683 copy_ivec(SHIFT_IVEC(g, j), jt);
1684 ivec_sub(SHIFT_IVEC(g, i), jt, dt_ij);
1685 ivec_sub(SHIFT_IVEC(g, k), jt, dt_kj);
1686 ivec_sub(SHIFT_IVEC(g, l), jt, dt_lj);
1687 t1 = IVEC2IS(dt_ij);
1688 t2 = IVEC2IS(dt_kj);
1689 t3 = IVEC2IS(dt_lj);
1693 t3 = pbc_rvec_sub(pbc, x[l], x[j], dx_jl);
1700 rvec_inc(fshift[t1], f_i);
1701 rvec_dec(fshift[CENTRAL], f_j);
1702 rvec_dec(fshift[t2], f_k);
1703 rvec_inc(fshift[t3], f_l);
1708 /* As do_dih_fup above, but without shift forces */
1710 do_dih_fup_noshiftf(int i, int j, int k, int l, real ddphi,
1711 rvec r_ij, rvec r_kj, rvec r_kl,
1712 rvec m, rvec n, rvec f[])
1714 rvec f_i, f_j, f_k, f_l;
1715 rvec uvec, vvec, svec, dx_jl;
1716 real iprm, iprn, nrkj, nrkj2, nrkj_1, nrkj_2;
1717 real a, b, p, q, toler;
1718 ivec jt, dt_ij, dt_kj, dt_lj;
1720 iprm = iprod(m, m); /* 5 */
1721 iprn = iprod(n, n); /* 5 */
1722 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1723 toler = nrkj2*GMX_REAL_EPS;
1724 if ((iprm > toler) && (iprn > toler))
1726 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1727 nrkj_2 = nrkj_1*nrkj_1; /* 1 */
1728 nrkj = nrkj2*nrkj_1; /* 1 */
1729 a = -ddphi*nrkj/iprm; /* 11 */
1730 svmul(a, m, f_i); /* 3 */
1731 b = ddphi*nrkj/iprn; /* 11 */
1732 svmul(b, n, f_l); /* 3 */
1733 p = iprod(r_ij, r_kj); /* 5 */
1734 p *= nrkj_2; /* 1 */
1735 q = iprod(r_kl, r_kj); /* 5 */
1736 q *= nrkj_2; /* 1 */
1737 svmul(p, f_i, uvec); /* 3 */
1738 svmul(q, f_l, vvec); /* 3 */
1739 rvec_sub(uvec, vvec, svec); /* 3 */
1740 rvec_sub(f_i, svec, f_j); /* 3 */
1741 rvec_add(f_l, svec, f_k); /* 3 */
1742 rvec_inc(f[i], f_i); /* 3 */
1743 rvec_dec(f[j], f_j); /* 3 */
1744 rvec_dec(f[k], f_k); /* 3 */
1745 rvec_inc(f[l], f_l); /* 3 */
1749 /* As do_dih_fup_noshiftf above, but with pre-calculated pre-factors */
1750 static gmx_inline void
1751 do_dih_fup_noshiftf_precalc(int i, int j, int k, int l,
1753 real f_i_x, real f_i_y, real f_i_z,
1754 real mf_l_x, real mf_l_y, real mf_l_z,
1757 rvec f_i, f_j, f_k, f_l;
1758 rvec uvec, vvec, svec;
1766 svmul(p, f_i, uvec);
1767 svmul(q, f_l, vvec);
1768 rvec_sub(uvec, vvec, svec);
1769 rvec_sub(f_i, svec, f_j);
1770 rvec_add(f_l, svec, f_k);
1771 rvec_inc(f[i], f_i);
1772 rvec_dec(f[j], f_j);
1773 rvec_dec(f[k], f_k);
1774 rvec_inc(f[l], f_l);
1778 real dopdihs(real cpA, real cpB, real phiA, real phiB, int mult,
1779 real phi, real lambda, real *V, real *F)
1781 real v, dvdlambda, mdphi, v1, sdphi, ddphi;
1782 real L1 = 1.0 - lambda;
1783 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1784 real dph0 = (phiB - phiA)*DEG2RAD;
1785 real cp = L1*cpA + lambda*cpB;
1787 mdphi = mult*phi - ph0;
1789 ddphi = -cp*mult*sdphi;
1790 v1 = 1.0 + cos(mdphi);
1793 dvdlambda = (cpB - cpA)*v1 + cp*dph0*sdphi;
1800 /* That was 40 flops */
1804 dopdihs_noener(real cpA, real cpB, real phiA, real phiB, int mult,
1805 real phi, real lambda, real *F)
1807 real mdphi, sdphi, ddphi;
1808 real L1 = 1.0 - lambda;
1809 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1810 real cp = L1*cpA + lambda*cpB;
1812 mdphi = mult*phi - ph0;
1814 ddphi = -cp*mult*sdphi;
1818 /* That was 20 flops */
1822 dopdihs_mdphi(real cpA, real cpB, real phiA, real phiB, int mult,
1823 real phi, real lambda, real *cp, real *mdphi)
1825 real L1 = 1.0 - lambda;
1826 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1828 *cp = L1*cpA + lambda*cpB;
1830 *mdphi = mult*phi - ph0;
1833 static real dopdihs_min(real cpA, real cpB, real phiA, real phiB, int mult,
1834 real phi, real lambda, real *V, real *F)
1835 /* similar to dopdihs, except for a minus sign *
1836 * and a different treatment of mult/phi0 */
1838 real v, dvdlambda, mdphi, v1, sdphi, ddphi;
1839 real L1 = 1.0 - lambda;
1840 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1841 real dph0 = (phiB - phiA)*DEG2RAD;
1842 real cp = L1*cpA + lambda*cpB;
1844 mdphi = mult*(phi-ph0);
1846 ddphi = cp*mult*sdphi;
1847 v1 = 1.0-cos(mdphi);
1850 dvdlambda = (cpB-cpA)*v1 + cp*dph0*sdphi;
1857 /* That was 40 flops */
1860 real pdihs(int nbonds,
1861 const t_iatom forceatoms[], const t_iparams forceparams[],
1862 const rvec x[], rvec f[], rvec fshift[],
1863 const t_pbc *pbc, const t_graph *g,
1864 real lambda, real *dvdlambda,
1865 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1866 int gmx_unused *global_atom_index)
1868 int i, type, ai, aj, ak, al;
1870 rvec r_ij, r_kj, r_kl, m, n;
1871 real phi, sign, ddphi, vpd, vtot;
1875 for (i = 0; (i < nbonds); )
1877 type = forceatoms[i++];
1878 ai = forceatoms[i++];
1879 aj = forceatoms[i++];
1880 ak = forceatoms[i++];
1881 al = forceatoms[i++];
1883 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
1884 &sign, &t1, &t2, &t3); /* 84 */
1885 *dvdlambda += dopdihs(forceparams[type].pdihs.cpA,
1886 forceparams[type].pdihs.cpB,
1887 forceparams[type].pdihs.phiA,
1888 forceparams[type].pdihs.phiB,
1889 forceparams[type].pdihs.mult,
1890 phi, lambda, &vpd, &ddphi);
1893 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
1894 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
1897 fprintf(debug, "pdih: (%d,%d,%d,%d) phi=%g\n",
1898 ai, aj, ak, al, phi);
1905 void make_dp_periodic(real *dp) /* 1 flop? */
1907 /* dp cannot be outside (-pi,pi) */
1912 else if (*dp < -M_PI)
1919 /* As pdihs above, but without calculating energies and shift forces */
1921 pdihs_noener(int nbonds,
1922 const t_iatom forceatoms[], const t_iparams forceparams[],
1923 const rvec x[], rvec f[],
1924 const t_pbc gmx_unused *pbc, const t_graph gmx_unused *g,
1926 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1927 int gmx_unused *global_atom_index)
1929 int i, type, ai, aj, ak, al;
1931 rvec r_ij, r_kj, r_kl, m, n;
1932 real phi, sign, ddphi_tot, ddphi;
1934 for (i = 0; (i < nbonds); )
1936 ai = forceatoms[i+1];
1937 aj = forceatoms[i+2];
1938 ak = forceatoms[i+3];
1939 al = forceatoms[i+4];
1941 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
1942 &sign, &t1, &t2, &t3);
1946 /* Loop over dihedrals working on the same atoms,
1947 * so we avoid recalculating angles and force distributions.
1951 type = forceatoms[i];
1952 dopdihs_noener(forceparams[type].pdihs.cpA,
1953 forceparams[type].pdihs.cpB,
1954 forceparams[type].pdihs.phiA,
1955 forceparams[type].pdihs.phiB,
1956 forceparams[type].pdihs.mult,
1957 phi, lambda, &ddphi);
1962 while (i < nbonds &&
1963 forceatoms[i+1] == ai &&
1964 forceatoms[i+2] == aj &&
1965 forceatoms[i+3] == ak &&
1966 forceatoms[i+4] == al);
1968 do_dih_fup_noshiftf(ai, aj, ak, al, ddphi_tot, r_ij, r_kj, r_kl, m, n, f);
1975 /* As pdihs_noner above, but using SIMD to calculate many dihedrals at once */
1977 pdihs_noener_simd(int nbonds,
1978 const t_iatom forceatoms[], const t_iparams forceparams[],
1979 const rvec x[], rvec f[],
1980 const t_pbc *pbc, const t_graph gmx_unused *g,
1981 real gmx_unused lambda,
1982 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1983 int gmx_unused *global_atom_index)
1985 #define UNROLL GMX_SIMD_WIDTH_HERE
1988 int type, ai[UNROLL], aj[UNROLL], ak[UNROLL], al[UNROLL];
1989 int t1[UNROLL], t2[UNROLL], t3[UNROLL];
1991 real dr_array[3*DIM*UNROLL+UNROLL], *dr;
1992 real buf_array[7*UNROLL+UNROLL], *buf;
1993 real *cp, *phi0, *mult, *phi, *p, *q, *sf_i, *msf_l;
1994 gmx_mm_pr phi0_S, phi_S;
1995 gmx_mm_pr mx_S, my_S, mz_S;
1996 gmx_mm_pr nx_S, ny_S, nz_S;
1997 gmx_mm_pr nrkj_m2_S, nrkj_n2_S;
1998 gmx_mm_pr cp_S, mdphi_S, mult_S;
1999 gmx_mm_pr sin_S, cos_S;
2001 gmx_mm_pr sf_i_S, msf_l_S;
2002 pbc_simd_t pbc_simd;
2004 /* Ensure SIMD register alignment */
2005 dr = gmx_simd_align_real(dr_array);
2006 buf = gmx_simd_align_real(buf_array);
2008 /* Extract aligned pointer for parameters and variables */
2009 cp = buf + 0*UNROLL;
2010 phi0 = buf + 1*UNROLL;
2011 mult = buf + 2*UNROLL;
2014 sf_i = buf + 5*UNROLL;
2015 msf_l = buf + 6*UNROLL;
2017 set_pbc_simd(pbc, &pbc_simd);
2019 /* nbonds is the number of dihedrals times nfa1, here we step UNROLL dihs */
2020 for (i = 0; (i < nbonds); i += UNROLL*nfa1)
2022 /* Collect atoms quadruplets for UNROLL dihedrals.
2023 * iu indexes into forceatoms, we should not let iu go beyond nbonds.
2026 for (s = 0; s < UNROLL; s++)
2028 type = forceatoms[iu];
2029 ai[s] = forceatoms[iu+1];
2030 aj[s] = forceatoms[iu+2];
2031 ak[s] = forceatoms[iu+3];
2032 al[s] = forceatoms[iu+4];
2034 cp[s] = forceparams[type].pdihs.cpA;
2035 phi0[s] = forceparams[type].pdihs.phiA*DEG2RAD;
2036 mult[s] = forceparams[type].pdihs.mult;
2038 /* At the end fill the arrays with identical entries */
2039 if (iu + nfa1 < nbonds)
2045 /* Caclulate UNROLL dihedral angles at once */
2046 dih_angle_simd(x, ai, aj, ak, al, &pbc_simd,
2049 &mx_S, &my_S, &mz_S,
2050 &nx_S, &ny_S, &nz_S,
2055 cp_S = gmx_load_pr(cp);
2056 phi0_S = gmx_load_pr(phi0);
2057 mult_S = gmx_load_pr(mult);
2059 mdphi_S = gmx_sub_pr(gmx_mul_pr(mult_S, phi_S), phi0_S);
2061 /* Calculate UNROLL sines at once */
2062 gmx_sincos_pr(mdphi_S, &sin_S, &cos_S);
2063 mddphi_S = gmx_mul_pr(gmx_mul_pr(cp_S, mult_S), sin_S);
2064 sf_i_S = gmx_mul_pr(mddphi_S, nrkj_m2_S);
2065 msf_l_S = gmx_mul_pr(mddphi_S, nrkj_n2_S);
2067 /* After this m?_S will contain f[i] */
2068 mx_S = gmx_mul_pr(sf_i_S, mx_S);
2069 my_S = gmx_mul_pr(sf_i_S, my_S);
2070 mz_S = gmx_mul_pr(sf_i_S, mz_S);
2072 /* After this m?_S will contain -f[l] */
2073 nx_S = gmx_mul_pr(msf_l_S, nx_S);
2074 ny_S = gmx_mul_pr(msf_l_S, ny_S);
2075 nz_S = gmx_mul_pr(msf_l_S, nz_S);
2077 gmx_store_pr(dr + 0*UNROLL, mx_S);
2078 gmx_store_pr(dr + 1*UNROLL, my_S);
2079 gmx_store_pr(dr + 2*UNROLL, mz_S);
2080 gmx_store_pr(dr + 3*UNROLL, nx_S);
2081 gmx_store_pr(dr + 4*UNROLL, ny_S);
2082 gmx_store_pr(dr + 5*UNROLL, nz_S);
2088 do_dih_fup_noshiftf_precalc(ai[s], aj[s], ak[s], al[s],
2093 dr[(DIM+XX)*UNROLL+s],
2094 dr[(DIM+YY)*UNROLL+s],
2095 dr[(DIM+ZZ)*UNROLL+s],
2100 while (s < UNROLL && iu < nbonds);
2105 #endif /* SIMD_BONDEDS */
2108 real idihs(int nbonds,
2109 const t_iatom forceatoms[], const t_iparams forceparams[],
2110 const rvec x[], rvec f[], rvec fshift[],
2111 const t_pbc *pbc, const t_graph *g,
2112 real lambda, real *dvdlambda,
2113 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2114 int gmx_unused *global_atom_index)
2116 int i, type, ai, aj, ak, al;
2118 real phi, phi0, dphi0, ddphi, sign, vtot;
2119 rvec r_ij, r_kj, r_kl, m, n;
2120 real L1, kk, dp, dp2, kA, kB, pA, pB, dvdl_term;
2125 for (i = 0; (i < nbonds); )
2127 type = forceatoms[i++];
2128 ai = forceatoms[i++];
2129 aj = forceatoms[i++];
2130 ak = forceatoms[i++];
2131 al = forceatoms[i++];
2133 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2134 &sign, &t1, &t2, &t3); /* 84 */
2136 /* phi can jump if phi0 is close to Pi/-Pi, which will cause huge
2137 * force changes if we just apply a normal harmonic.
2138 * Instead, we first calculate phi-phi0 and take it modulo (-Pi,Pi).
2139 * This means we will never have the periodicity problem, unless
2140 * the dihedral is Pi away from phiO, which is very unlikely due to
2143 kA = forceparams[type].harmonic.krA;
2144 kB = forceparams[type].harmonic.krB;
2145 pA = forceparams[type].harmonic.rA;
2146 pB = forceparams[type].harmonic.rB;
2148 kk = L1*kA + lambda*kB;
2149 phi0 = (L1*pA + lambda*pB)*DEG2RAD;
2150 dphi0 = (pB - pA)*DEG2RAD;
2154 make_dp_periodic(&dp);
2161 dvdl_term += 0.5*(kB - kA)*dp2 - kk*dphi0*dp;
2163 do_dih_fup(ai, aj, ak, al, (real)(-ddphi), r_ij, r_kj, r_kl, m, n,
2164 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2169 fprintf(debug, "idih: (%d,%d,%d,%d) phi=%g\n",
2170 ai, aj, ak, al, phi);
2175 *dvdlambda += dvdl_term;
2180 /*! \brief returns dx, rdist, and dpdl for functions posres() and fbposres()
2182 static void posres_dx(const rvec x, const rvec pos0A, const rvec pos0B,
2183 const rvec comA_sc, const rvec comB_sc,
2185 t_pbc *pbc, int refcoord_scaling, int npbcdim,
2186 rvec dx, rvec rdist, rvec dpdl)
2189 real posA, posB, L1, ref = 0.;
2194 for (m = 0; m < DIM; m++)
2200 switch (refcoord_scaling)
2204 rdist[m] = L1*posA + lambda*posB;
2205 dpdl[m] = posB - posA;
2208 /* Box relative coordinates are stored for dimensions with pbc */
2209 posA *= pbc->box[m][m];
2210 posB *= pbc->box[m][m];
2211 for (d = m+1; d < npbcdim; d++)
2213 posA += pos0A[d]*pbc->box[d][m];
2214 posB += pos0B[d]*pbc->box[d][m];
2216 ref = L1*posA + lambda*posB;
2218 dpdl[m] = posB - posA;
2221 ref = L1*comA_sc[m] + lambda*comB_sc[m];
2222 rdist[m] = L1*posA + lambda*posB;
2223 dpdl[m] = comB_sc[m] - comA_sc[m] + posB - posA;
2226 gmx_fatal(FARGS, "No such scaling method implemented");
2231 ref = L1*posA + lambda*posB;
2233 dpdl[m] = posB - posA;
2236 /* We do pbc_dx with ref+rdist,
2237 * since with only ref we can be up to half a box vector wrong.
2239 pos[m] = ref + rdist[m];
2244 pbc_dx(pbc, x, pos, dx);
2248 rvec_sub(x, pos, dx);
2252 /*! \brief Adds forces of flat-bottomed positions restraints to f[]
2253 * and fixes vir_diag. Returns the flat-bottomed potential. */
2254 real fbposres(int nbonds,
2255 const t_iatom forceatoms[], const t_iparams forceparams[],
2256 const rvec x[], rvec f[], rvec vir_diag,
2258 int refcoord_scaling, int ePBC, rvec com)
2259 /* compute flat-bottomed positions restraints */
2261 int i, ai, m, d, type, npbcdim = 0, fbdim;
2262 const t_iparams *pr;
2264 real ref = 0, dr, dr2, rpot, rfb, rfb2, fact, invdr;
2265 rvec com_sc, rdist, pos, dx, dpdl, fm;
2268 npbcdim = ePBC2npbcdim(ePBC);
2270 if (refcoord_scaling == erscCOM)
2273 for (m = 0; m < npbcdim; m++)
2275 for (d = m; d < npbcdim; d++)
2277 com_sc[m] += com[d]*pbc->box[d][m];
2283 for (i = 0; (i < nbonds); )
2285 type = forceatoms[i++];
2286 ai = forceatoms[i++];
2287 pr = &forceparams[type];
2289 /* same calculation as for normal posres, but with identical A and B states, and lambda==0 */
2290 posres_dx(x[ai], forceparams[type].fbposres.pos0, forceparams[type].fbposres.pos0,
2291 com_sc, com_sc, 0.0,
2292 pbc, refcoord_scaling, npbcdim,
2298 kk = pr->fbposres.k;
2299 rfb = pr->fbposres.r;
2302 /* with rfb<0, push particle out of the sphere/cylinder/layer */
2310 switch (pr->fbposres.geom)
2312 case efbposresSPHERE:
2313 /* spherical flat-bottom posres */
2316 ( (dr2 > rfb2 && bInvert == FALSE ) || (dr2 < rfb2 && bInvert == TRUE ) )
2320 v = 0.5*kk*sqr(dr - rfb);
2321 fact = -kk*(dr-rfb)/dr; /* Force pointing to the center pos0 */
2322 svmul(fact, dx, fm);
2325 case efbposresCYLINDER:
2326 /* cylidrical flat-bottom posres in x-y plane. fm[ZZ] = 0. */
2327 dr2 = sqr(dx[XX])+sqr(dx[YY]);
2329 ( (dr2 > rfb2 && bInvert == FALSE ) || (dr2 < rfb2 && bInvert == TRUE ) )
2334 v = 0.5*kk*sqr(dr - rfb);
2335 fm[XX] = -kk*(dr-rfb)*dx[XX]*invdr; /* Force pointing to the center */
2336 fm[YY] = -kk*(dr-rfb)*dx[YY]*invdr;
2339 case efbposresX: /* fbdim=XX */
2340 case efbposresY: /* fbdim=YY */
2341 case efbposresZ: /* fbdim=ZZ */
2342 /* 1D flat-bottom potential */
2343 fbdim = pr->fbposres.geom - efbposresX;
2345 if ( ( dr > rfb && bInvert == FALSE ) || ( 0 < dr && dr < rfb && bInvert == TRUE ) )
2347 v = 0.5*kk*sqr(dr - rfb);
2348 fm[fbdim] = -kk*(dr - rfb);
2350 else if ( (dr < (-rfb) && bInvert == FALSE ) || ( (-rfb) < dr && dr < 0 && bInvert == TRUE ))
2352 v = 0.5*kk*sqr(dr + rfb);
2353 fm[fbdim] = -kk*(dr + rfb);
2360 for (m = 0; (m < DIM); m++)
2363 /* Here we correct for the pbc_dx which included rdist */
2364 vir_diag[m] -= 0.5*(dx[m] + rdist[m])*fm[m];
2372 real posres(int nbonds,
2373 const t_iatom forceatoms[], const t_iparams forceparams[],
2374 const rvec x[], rvec f[], rvec vir_diag,
2376 real lambda, real *dvdlambda,
2377 int refcoord_scaling, int ePBC, rvec comA, rvec comB)
2379 int i, ai, m, d, type, ki, npbcdim = 0;
2380 const t_iparams *pr;
2383 real posA, posB, ref = 0;
2384 rvec comA_sc, comB_sc, rdist, dpdl, pos, dx;
2385 gmx_bool bForceValid = TRUE;
2387 if ((f == NULL) || (vir_diag == NULL)) /* should both be null together! */
2389 bForceValid = FALSE;
2392 npbcdim = ePBC2npbcdim(ePBC);
2394 if (refcoord_scaling == erscCOM)
2396 clear_rvec(comA_sc);
2397 clear_rvec(comB_sc);
2398 for (m = 0; m < npbcdim; m++)
2400 for (d = m; d < npbcdim; d++)
2402 comA_sc[m] += comA[d]*pbc->box[d][m];
2403 comB_sc[m] += comB[d]*pbc->box[d][m];
2411 for (i = 0; (i < nbonds); )
2413 type = forceatoms[i++];
2414 ai = forceatoms[i++];
2415 pr = &forceparams[type];
2417 /* return dx, rdist, and dpdl */
2418 posres_dx(x[ai], forceparams[type].posres.pos0A, forceparams[type].posres.pos0B,
2419 comA_sc, comB_sc, lambda,
2420 pbc, refcoord_scaling, npbcdim,
2423 for (m = 0; (m < DIM); m++)
2425 kk = L1*pr->posres.fcA[m] + lambda*pr->posres.fcB[m];
2427 vtot += 0.5*kk*dx[m]*dx[m];
2429 0.5*(pr->posres.fcB[m] - pr->posres.fcA[m])*dx[m]*dx[m]
2432 /* Here we correct for the pbc_dx which included rdist */
2436 vir_diag[m] -= 0.5*(dx[m] + rdist[m])*fm;
2444 static real low_angres(int nbonds,
2445 const t_iatom forceatoms[], const t_iparams forceparams[],
2446 const rvec x[], rvec f[], rvec fshift[],
2447 const t_pbc *pbc, const t_graph *g,
2448 real lambda, real *dvdlambda,
2451 int i, m, type, ai, aj, ak, al;
2453 real phi, cos_phi, cos_phi2, vid, vtot, dVdphi;
2454 rvec r_ij, r_kl, f_i, f_k = {0, 0, 0};
2455 real st, sth, nrij2, nrkl2, c, cij, ckl;
2458 t2 = 0; /* avoid warning with gcc-3.3. It is never used uninitialized */
2461 ak = al = 0; /* to avoid warnings */
2462 for (i = 0; i < nbonds; )
2464 type = forceatoms[i++];
2465 ai = forceatoms[i++];
2466 aj = forceatoms[i++];
2467 t1 = pbc_rvec_sub(pbc, x[aj], x[ai], r_ij); /* 3 */
2470 ak = forceatoms[i++];
2471 al = forceatoms[i++];
2472 t2 = pbc_rvec_sub(pbc, x[al], x[ak], r_kl); /* 3 */
2481 cos_phi = cos_angle(r_ij, r_kl); /* 25 */
2482 phi = acos(cos_phi); /* 10 */
2484 *dvdlambda += dopdihs_min(forceparams[type].pdihs.cpA,
2485 forceparams[type].pdihs.cpB,
2486 forceparams[type].pdihs.phiA,
2487 forceparams[type].pdihs.phiB,
2488 forceparams[type].pdihs.mult,
2489 phi, lambda, &vid, &dVdphi); /* 40 */
2493 cos_phi2 = sqr(cos_phi); /* 1 */
2496 st = -dVdphi*gmx_invsqrt(1 - cos_phi2); /* 12 */
2497 sth = st*cos_phi; /* 1 */
2498 nrij2 = iprod(r_ij, r_ij); /* 5 */
2499 nrkl2 = iprod(r_kl, r_kl); /* 5 */
2501 c = st*gmx_invsqrt(nrij2*nrkl2); /* 11 */
2502 cij = sth/nrij2; /* 10 */
2503 ckl = sth/nrkl2; /* 10 */
2505 for (m = 0; m < DIM; m++) /* 18+18 */
2507 f_i[m] = (c*r_kl[m]-cij*r_ij[m]);
2512 f_k[m] = (c*r_ij[m]-ckl*r_kl[m]);
2520 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
2523 rvec_inc(fshift[t1], f_i);
2524 rvec_dec(fshift[CENTRAL], f_i);
2529 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, al), dt);
2532 rvec_inc(fshift[t2], f_k);
2533 rvec_dec(fshift[CENTRAL], f_k);
2538 return vtot; /* 184 / 157 (bZAxis) total */
2541 real angres(int nbonds,
2542 const t_iatom forceatoms[], const t_iparams forceparams[],
2543 const rvec x[], rvec f[], rvec fshift[],
2544 const t_pbc *pbc, const t_graph *g,
2545 real lambda, real *dvdlambda,
2546 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2547 int gmx_unused *global_atom_index)
2549 return low_angres(nbonds, forceatoms, forceparams, x, f, fshift, pbc, g,
2550 lambda, dvdlambda, FALSE);
2553 real angresz(int nbonds,
2554 const t_iatom forceatoms[], const t_iparams forceparams[],
2555 const rvec x[], rvec f[], rvec fshift[],
2556 const t_pbc *pbc, const t_graph *g,
2557 real lambda, real *dvdlambda,
2558 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2559 int gmx_unused *global_atom_index)
2561 return low_angres(nbonds, forceatoms, forceparams, x, f, fshift, pbc, g,
2562 lambda, dvdlambda, TRUE);
2565 real dihres(int nbonds,
2566 const t_iatom forceatoms[], const t_iparams forceparams[],
2567 const rvec x[], rvec f[], rvec fshift[],
2568 const t_pbc *pbc, const t_graph *g,
2569 real lambda, real *dvdlambda,
2570 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2571 int gmx_unused *global_atom_index)
2574 int ai, aj, ak, al, i, k, type, t1, t2, t3;
2575 real phi0A, phi0B, dphiA, dphiB, kfacA, kfacB, phi0, dphi, kfac;
2576 real phi, ddphi, ddp, ddp2, dp, sign, d2r, fc, L1;
2577 rvec r_ij, r_kj, r_kl, m, n;
2584 for (i = 0; (i < nbonds); )
2586 type = forceatoms[i++];
2587 ai = forceatoms[i++];
2588 aj = forceatoms[i++];
2589 ak = forceatoms[i++];
2590 al = forceatoms[i++];
2592 phi0A = forceparams[type].dihres.phiA*d2r;
2593 dphiA = forceparams[type].dihres.dphiA*d2r;
2594 kfacA = forceparams[type].dihres.kfacA;
2596 phi0B = forceparams[type].dihres.phiB*d2r;
2597 dphiB = forceparams[type].dihres.dphiB*d2r;
2598 kfacB = forceparams[type].dihres.kfacB;
2600 phi0 = L1*phi0A + lambda*phi0B;
2601 dphi = L1*dphiA + lambda*dphiB;
2602 kfac = L1*kfacA + lambda*kfacB;
2604 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2605 &sign, &t1, &t2, &t3);
2610 fprintf(debug, "dihres[%d]: %d %d %d %d : phi=%f, dphi=%f, kfac=%f\n",
2611 k++, ai, aj, ak, al, phi0, dphi, kfac);
2613 /* phi can jump if phi0 is close to Pi/-Pi, which will cause huge
2614 * force changes if we just apply a normal harmonic.
2615 * Instead, we first calculate phi-phi0 and take it modulo (-Pi,Pi).
2616 * This means we will never have the periodicity problem, unless
2617 * the dihedral is Pi away from phiO, which is very unlikely due to
2621 make_dp_periodic(&dp);
2627 else if (dp < -dphi)
2639 vtot += 0.5*kfac*ddp2;
2642 *dvdlambda += 0.5*(kfacB - kfacA)*ddp2;
2643 /* lambda dependence from changing restraint distances */
2646 *dvdlambda -= kfac*ddp*((dphiB - dphiA)+(phi0B - phi0A));
2650 *dvdlambda += kfac*ddp*((dphiB - dphiA)-(phi0B - phi0A));
2652 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
2653 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2660 real unimplemented(int gmx_unused nbonds,
2661 const t_iatom gmx_unused forceatoms[], const t_iparams gmx_unused forceparams[],
2662 const rvec gmx_unused x[], rvec gmx_unused f[], rvec gmx_unused fshift[],
2663 const t_pbc gmx_unused *pbc, const t_graph gmx_unused *g,
2664 real gmx_unused lambda, real gmx_unused *dvdlambda,
2665 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2666 int gmx_unused *global_atom_index)
2668 gmx_impl("*** you are using a not implemented function");
2670 return 0.0; /* To make the compiler happy */
2673 real rbdihs(int nbonds,
2674 const t_iatom forceatoms[], const t_iparams forceparams[],
2675 const rvec x[], rvec f[], rvec fshift[],
2676 const t_pbc *pbc, const t_graph *g,
2677 real lambda, real *dvdlambda,
2678 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2679 int gmx_unused *global_atom_index)
2681 const real c0 = 0.0, c1 = 1.0, c2 = 2.0, c3 = 3.0, c4 = 4.0, c5 = 5.0;
2682 int type, ai, aj, ak, al, i, j;
2684 rvec r_ij, r_kj, r_kl, m, n;
2685 real parmA[NR_RBDIHS];
2686 real parmB[NR_RBDIHS];
2687 real parm[NR_RBDIHS];
2688 real cos_phi, phi, rbp, rbpBA;
2689 real v, sign, ddphi, sin_phi;
2691 real L1 = 1.0-lambda;
2695 for (i = 0; (i < nbonds); )
2697 type = forceatoms[i++];
2698 ai = forceatoms[i++];
2699 aj = forceatoms[i++];
2700 ak = forceatoms[i++];
2701 al = forceatoms[i++];
2703 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2704 &sign, &t1, &t2, &t3); /* 84 */
2706 /* Change to polymer convention */
2713 phi -= M_PI; /* 1 */
2717 /* Beware of accuracy loss, cannot use 1-sqrt(cos^2) ! */
2720 for (j = 0; (j < NR_RBDIHS); j++)
2722 parmA[j] = forceparams[type].rbdihs.rbcA[j];
2723 parmB[j] = forceparams[type].rbdihs.rbcB[j];
2724 parm[j] = L1*parmA[j]+lambda*parmB[j];
2726 /* Calculate cosine powers */
2727 /* Calculate the energy */
2728 /* Calculate the derivative */
2731 dvdl_term += (parmB[0]-parmA[0]);
2736 rbpBA = parmB[1]-parmA[1];
2737 ddphi += rbp*cosfac;
2740 dvdl_term += cosfac*rbpBA;
2742 rbpBA = parmB[2]-parmA[2];
2743 ddphi += c2*rbp*cosfac;
2746 dvdl_term += cosfac*rbpBA;
2748 rbpBA = parmB[3]-parmA[3];
2749 ddphi += c3*rbp*cosfac;
2752 dvdl_term += cosfac*rbpBA;
2754 rbpBA = parmB[4]-parmA[4];
2755 ddphi += c4*rbp*cosfac;
2758 dvdl_term += cosfac*rbpBA;
2760 rbpBA = parmB[5]-parmA[5];
2761 ddphi += c5*rbp*cosfac;
2764 dvdl_term += cosfac*rbpBA;
2766 ddphi = -ddphi*sin_phi; /* 11 */
2768 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
2769 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2772 *dvdlambda += dvdl_term;
2777 int cmap_setup_grid_index(int ip, int grid_spacing, int *ipm1, int *ipp1, int *ipp2)
2783 ip = ip + grid_spacing - 1;
2785 else if (ip > grid_spacing)
2787 ip = ip - grid_spacing - 1;
2796 im1 = grid_spacing - 1;
2798 else if (ip == grid_spacing-2)
2802 else if (ip == grid_spacing-1)
2816 real cmap_dihs(int nbonds,
2817 const t_iatom forceatoms[], const t_iparams forceparams[],
2818 const gmx_cmap_t *cmap_grid,
2819 const rvec x[], rvec f[], rvec fshift[],
2820 const t_pbc *pbc, const t_graph *g,
2821 real gmx_unused lambda, real gmx_unused *dvdlambda,
2822 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2823 int gmx_unused *global_atom_index)
2825 int i, j, k, n, idx;
2826 int ai, aj, ak, al, am;
2827 int a1i, a1j, a1k, a1l, a2i, a2j, a2k, a2l;
2829 int t11, t21, t31, t12, t22, t32;
2830 int iphi1, ip1m1, ip1p1, ip1p2;
2831 int iphi2, ip2m1, ip2p1, ip2p2;
2833 int pos1, pos2, pos3, pos4, tmp;
2835 real ty[4], ty1[4], ty2[4], ty12[4], tc[16], tx[16];
2836 real phi1, psi1, cos_phi1, sin_phi1, sign1, xphi1;
2837 real phi2, psi2, cos_phi2, sin_phi2, sign2, xphi2;
2838 real dx, xx, tt, tu, e, df1, df2, ddf1, ddf2, ddf12, vtot;
2839 real ra21, rb21, rg21, rg1, rgr1, ra2r1, rb2r1, rabr1;
2840 real ra22, rb22, rg22, rg2, rgr2, ra2r2, rb2r2, rabr2;
2841 real fg1, hg1, fga1, hgb1, gaa1, gbb1;
2842 real fg2, hg2, fga2, hgb2, gaa2, gbb2;
2845 rvec r1_ij, r1_kj, r1_kl, m1, n1;
2846 rvec r2_ij, r2_kj, r2_kl, m2, n2;
2847 rvec f1_i, f1_j, f1_k, f1_l;
2848 rvec f2_i, f2_j, f2_k, f2_l;
2849 rvec a1, b1, a2, b2;
2850 rvec f1, g1, h1, f2, g2, h2;
2851 rvec dtf1, dtg1, dth1, dtf2, dtg2, dth2;
2852 ivec jt1, dt1_ij, dt1_kj, dt1_lj;
2853 ivec jt2, dt2_ij, dt2_kj, dt2_lj;
2857 int loop_index[4][4] = {
2864 /* Total CMAP energy */
2867 for (n = 0; n < nbonds; )
2869 /* Five atoms are involved in the two torsions */
2870 type = forceatoms[n++];
2871 ai = forceatoms[n++];
2872 aj = forceatoms[n++];
2873 ak = forceatoms[n++];
2874 al = forceatoms[n++];
2875 am = forceatoms[n++];
2877 /* Which CMAP type is this */
2878 cmapA = forceparams[type].cmap.cmapA;
2879 cmapd = cmap_grid->cmapdata[cmapA].cmap;
2887 phi1 = dih_angle(x[a1i], x[a1j], x[a1k], x[a1l], pbc, r1_ij, r1_kj, r1_kl, m1, n1,
2888 &sign1, &t11, &t21, &t31); /* 84 */
2890 cos_phi1 = cos(phi1);
2892 a1[0] = r1_ij[1]*r1_kj[2]-r1_ij[2]*r1_kj[1];
2893 a1[1] = r1_ij[2]*r1_kj[0]-r1_ij[0]*r1_kj[2];
2894 a1[2] = r1_ij[0]*r1_kj[1]-r1_ij[1]*r1_kj[0]; /* 9 */
2896 b1[0] = r1_kl[1]*r1_kj[2]-r1_kl[2]*r1_kj[1];
2897 b1[1] = r1_kl[2]*r1_kj[0]-r1_kl[0]*r1_kj[2];
2898 b1[2] = r1_kl[0]*r1_kj[1]-r1_kl[1]*r1_kj[0]; /* 9 */
2900 tmp = pbc_rvec_sub(pbc, x[a1l], x[a1k], h1);
2902 ra21 = iprod(a1, a1); /* 5 */
2903 rb21 = iprod(b1, b1); /* 5 */
2904 rg21 = iprod(r1_kj, r1_kj); /* 5 */
2910 rabr1 = sqrt(ra2r1*rb2r1);
2912 sin_phi1 = rg1 * rabr1 * iprod(a1, h1) * (-1);
2914 if (cos_phi1 < -0.5 || cos_phi1 > 0.5)
2916 phi1 = asin(sin_phi1);
2926 phi1 = -M_PI - phi1;
2932 phi1 = acos(cos_phi1);
2940 xphi1 = phi1 + M_PI; /* 1 */
2942 /* Second torsion */
2948 phi2 = dih_angle(x[a2i], x[a2j], x[a2k], x[a2l], pbc, r2_ij, r2_kj, r2_kl, m2, n2,
2949 &sign2, &t12, &t22, &t32); /* 84 */
2951 cos_phi2 = cos(phi2);
2953 a2[0] = r2_ij[1]*r2_kj[2]-r2_ij[2]*r2_kj[1];
2954 a2[1] = r2_ij[2]*r2_kj[0]-r2_ij[0]*r2_kj[2];
2955 a2[2] = r2_ij[0]*r2_kj[1]-r2_ij[1]*r2_kj[0]; /* 9 */
2957 b2[0] = r2_kl[1]*r2_kj[2]-r2_kl[2]*r2_kj[1];
2958 b2[1] = r2_kl[2]*r2_kj[0]-r2_kl[0]*r2_kj[2];
2959 b2[2] = r2_kl[0]*r2_kj[1]-r2_kl[1]*r2_kj[0]; /* 9 */
2961 tmp = pbc_rvec_sub(pbc, x[a2l], x[a2k], h2);
2963 ra22 = iprod(a2, a2); /* 5 */
2964 rb22 = iprod(b2, b2); /* 5 */
2965 rg22 = iprod(r2_kj, r2_kj); /* 5 */
2971 rabr2 = sqrt(ra2r2*rb2r2);
2973 sin_phi2 = rg2 * rabr2 * iprod(a2, h2) * (-1);
2975 if (cos_phi2 < -0.5 || cos_phi2 > 0.5)
2977 phi2 = asin(sin_phi2);
2987 phi2 = -M_PI - phi2;
2993 phi2 = acos(cos_phi2);
3001 xphi2 = phi2 + M_PI; /* 1 */
3003 /* Range mangling */
3006 xphi1 = xphi1 + 2*M_PI;
3008 else if (xphi1 >= 2*M_PI)
3010 xphi1 = xphi1 - 2*M_PI;
3015 xphi2 = xphi2 + 2*M_PI;
3017 else if (xphi2 >= 2*M_PI)
3019 xphi2 = xphi2 - 2*M_PI;
3022 /* Number of grid points */
3023 dx = 2*M_PI / cmap_grid->grid_spacing;
3025 /* Where on the grid are we */
3026 iphi1 = (int)(xphi1/dx);
3027 iphi2 = (int)(xphi2/dx);
3029 iphi1 = cmap_setup_grid_index(iphi1, cmap_grid->grid_spacing, &ip1m1, &ip1p1, &ip1p2);
3030 iphi2 = cmap_setup_grid_index(iphi2, cmap_grid->grid_spacing, &ip2m1, &ip2p1, &ip2p2);
3032 pos1 = iphi1*cmap_grid->grid_spacing+iphi2;
3033 pos2 = ip1p1*cmap_grid->grid_spacing+iphi2;
3034 pos3 = ip1p1*cmap_grid->grid_spacing+ip2p1;
3035 pos4 = iphi1*cmap_grid->grid_spacing+ip2p1;
3037 ty[0] = cmapd[pos1*4];
3038 ty[1] = cmapd[pos2*4];
3039 ty[2] = cmapd[pos3*4];
3040 ty[3] = cmapd[pos4*4];
3042 ty1[0] = cmapd[pos1*4+1];
3043 ty1[1] = cmapd[pos2*4+1];
3044 ty1[2] = cmapd[pos3*4+1];
3045 ty1[3] = cmapd[pos4*4+1];
3047 ty2[0] = cmapd[pos1*4+2];
3048 ty2[1] = cmapd[pos2*4+2];
3049 ty2[2] = cmapd[pos3*4+2];
3050 ty2[3] = cmapd[pos4*4+2];
3052 ty12[0] = cmapd[pos1*4+3];
3053 ty12[1] = cmapd[pos2*4+3];
3054 ty12[2] = cmapd[pos3*4+3];
3055 ty12[3] = cmapd[pos4*4+3];
3057 /* Switch to degrees */
3058 dx = 360.0 / cmap_grid->grid_spacing;
3059 xphi1 = xphi1 * RAD2DEG;
3060 xphi2 = xphi2 * RAD2DEG;
3062 for (i = 0; i < 4; i++) /* 16 */
3065 tx[i+4] = ty1[i]*dx;
3066 tx[i+8] = ty2[i]*dx;
3067 tx[i+12] = ty12[i]*dx*dx;
3071 for (i = 0; i < 4; i++) /* 1056 */
3073 for (j = 0; j < 4; j++)
3076 for (k = 0; k < 16; k++)
3078 xx = xx + cmap_coeff_matrix[k*16+idx]*tx[k];
3086 tt = (xphi1-iphi1*dx)/dx;
3087 tu = (xphi2-iphi2*dx)/dx;
3096 for (i = 3; i >= 0; i--)
3098 l1 = loop_index[i][3];
3099 l2 = loop_index[i][2];
3100 l3 = loop_index[i][1];
3102 e = tt * e + ((tc[i*4+3]*tu+tc[i*4+2])*tu + tc[i*4+1])*tu+tc[i*4];
3103 df1 = tu * df1 + (3.0*tc[l1]*tt+2.0*tc[l2])*tt+tc[l3];
3104 df2 = tt * df2 + (3.0*tc[i*4+3]*tu+2.0*tc[i*4+2])*tu+tc[i*4+1];
3105 ddf1 = tu * ddf1 + 2.0*3.0*tc[l1]*tt+2.0*tc[l2];
3106 ddf2 = tt * ddf2 + 2.0*3.0*tc[4*i+3]*tu+2.0*tc[4*i+2];
3109 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) +
3110 3.0*tu*tu*(tc[7]+2.0*tc[11]*tt+3.0*tc[15]*tt*tt);
3115 ddf1 = ddf1 * fac * fac;
3116 ddf2 = ddf2 * fac * fac;
3117 ddf12 = ddf12 * fac * fac;
3122 /* Do forces - first torsion */
3123 fg1 = iprod(r1_ij, r1_kj);
3124 hg1 = iprod(r1_kl, r1_kj);
3125 fga1 = fg1*ra2r1*rgr1;
3126 hgb1 = hg1*rb2r1*rgr1;
3130 for (i = 0; i < DIM; i++)
3132 dtf1[i] = gaa1 * a1[i];
3133 dtg1[i] = fga1 * a1[i] - hgb1 * b1[i];
3134 dth1[i] = gbb1 * b1[i];
3136 f1[i] = df1 * dtf1[i];
3137 g1[i] = df1 * dtg1[i];
3138 h1[i] = df1 * dth1[i];
3141 f1_j[i] = -f1[i] - g1[i];
3142 f1_k[i] = h1[i] + g1[i];
3145 f[a1i][i] = f[a1i][i] + f1_i[i];
3146 f[a1j][i] = f[a1j][i] + f1_j[i]; /* - f1[i] - g1[i] */
3147 f[a1k][i] = f[a1k][i] + f1_k[i]; /* h1[i] + g1[i] */
3148 f[a1l][i] = f[a1l][i] + f1_l[i]; /* h1[i] */
3151 /* Do forces - second torsion */
3152 fg2 = iprod(r2_ij, r2_kj);
3153 hg2 = iprod(r2_kl, r2_kj);
3154 fga2 = fg2*ra2r2*rgr2;
3155 hgb2 = hg2*rb2r2*rgr2;
3159 for (i = 0; i < DIM; i++)
3161 dtf2[i] = gaa2 * a2[i];
3162 dtg2[i] = fga2 * a2[i] - hgb2 * b2[i];
3163 dth2[i] = gbb2 * b2[i];
3165 f2[i] = df2 * dtf2[i];
3166 g2[i] = df2 * dtg2[i];
3167 h2[i] = df2 * dth2[i];
3170 f2_j[i] = -f2[i] - g2[i];
3171 f2_k[i] = h2[i] + g2[i];
3174 f[a2i][i] = f[a2i][i] + f2_i[i]; /* f2[i] */
3175 f[a2j][i] = f[a2j][i] + f2_j[i]; /* - f2[i] - g2[i] */
3176 f[a2k][i] = f[a2k][i] + f2_k[i]; /* h2[i] + g2[i] */
3177 f[a2l][i] = f[a2l][i] + f2_l[i]; /* - h2[i] */
3183 copy_ivec(SHIFT_IVEC(g, a1j), jt1);
3184 ivec_sub(SHIFT_IVEC(g, a1i), jt1, dt1_ij);
3185 ivec_sub(SHIFT_IVEC(g, a1k), jt1, dt1_kj);
3186 ivec_sub(SHIFT_IVEC(g, a1l), jt1, dt1_lj);
3187 t11 = IVEC2IS(dt1_ij);
3188 t21 = IVEC2IS(dt1_kj);
3189 t31 = IVEC2IS(dt1_lj);
3191 copy_ivec(SHIFT_IVEC(g, a2j), jt2);
3192 ivec_sub(SHIFT_IVEC(g, a2i), jt2, dt2_ij);
3193 ivec_sub(SHIFT_IVEC(g, a2k), jt2, dt2_kj);
3194 ivec_sub(SHIFT_IVEC(g, a2l), jt2, dt2_lj);
3195 t12 = IVEC2IS(dt2_ij);
3196 t22 = IVEC2IS(dt2_kj);
3197 t32 = IVEC2IS(dt2_lj);
3201 t31 = pbc_rvec_sub(pbc, x[a1l], x[a1j], h1);
3202 t32 = pbc_rvec_sub(pbc, x[a2l], x[a2j], h2);
3210 rvec_inc(fshift[t11], f1_i);
3211 rvec_inc(fshift[CENTRAL], f1_j);
3212 rvec_inc(fshift[t21], f1_k);
3213 rvec_inc(fshift[t31], f1_l);
3215 rvec_inc(fshift[t21], f2_i);
3216 rvec_inc(fshift[CENTRAL], f2_j);
3217 rvec_inc(fshift[t22], f2_k);
3218 rvec_inc(fshift[t32], f2_l);
3225 /***********************************************************
3227 * G R O M O S 9 6 F U N C T I O N S
3229 ***********************************************************/
3230 real g96harmonic(real kA, real kB, real xA, real xB, real x, real lambda,
3233 const real half = 0.5;
3234 real L1, kk, x0, dx, dx2;
3235 real v, f, dvdlambda;
3238 kk = L1*kA+lambda*kB;
3239 x0 = L1*xA+lambda*xB;
3246 dvdlambda = half*(kB-kA)*dx2 + (xA-xB)*kk*dx;
3253 /* That was 21 flops */
3256 real g96bonds(int nbonds,
3257 const t_iatom forceatoms[], const t_iparams forceparams[],
3258 const rvec x[], rvec f[], rvec fshift[],
3259 const t_pbc *pbc, const t_graph *g,
3260 real lambda, real *dvdlambda,
3261 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3262 int gmx_unused *global_atom_index)
3264 int i, m, ki, ai, aj, type;
3265 real dr2, fbond, vbond, fij, vtot;
3270 for (i = 0; (i < nbonds); )
3272 type = forceatoms[i++];
3273 ai = forceatoms[i++];
3274 aj = forceatoms[i++];
3276 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
3277 dr2 = iprod(dx, dx); /* 5 */
3279 *dvdlambda += g96harmonic(forceparams[type].harmonic.krA,
3280 forceparams[type].harmonic.krB,
3281 forceparams[type].harmonic.rA,
3282 forceparams[type].harmonic.rB,
3283 dr2, lambda, &vbond, &fbond);
3285 vtot += 0.5*vbond; /* 1*/
3289 fprintf(debug, "G96-BONDS: dr = %10g vbond = %10g fbond = %10g\n",
3290 sqrt(dr2), vbond, fbond);
3296 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
3299 for (m = 0; (m < DIM); m++) /* 15 */
3304 fshift[ki][m] += fij;
3305 fshift[CENTRAL][m] -= fij;
3311 real g96bond_angle(const rvec xi, const rvec xj, const rvec xk, const t_pbc *pbc,
3312 rvec r_ij, rvec r_kj,
3314 /* Return value is the angle between the bonds i-j and j-k */
3318 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
3319 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
3321 costh = cos_angle(r_ij, r_kj); /* 25 */
3326 real g96angles(int nbonds,
3327 const t_iatom forceatoms[], const t_iparams forceparams[],
3328 const rvec x[], rvec f[], rvec fshift[],
3329 const t_pbc *pbc, const t_graph *g,
3330 real lambda, real *dvdlambda,
3331 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3332 int gmx_unused *global_atom_index)
3334 int i, ai, aj, ak, type, m, t1, t2;
3336 real cos_theta, dVdt, va, vtot;
3337 real rij_1, rij_2, rkj_1, rkj_2, rijrkj_1;
3339 ivec jt, dt_ij, dt_kj;
3342 for (i = 0; (i < nbonds); )
3344 type = forceatoms[i++];
3345 ai = forceatoms[i++];
3346 aj = forceatoms[i++];
3347 ak = forceatoms[i++];
3349 cos_theta = g96bond_angle(x[ai], x[aj], x[ak], pbc, r_ij, r_kj, &t1, &t2);
3351 *dvdlambda += g96harmonic(forceparams[type].harmonic.krA,
3352 forceparams[type].harmonic.krB,
3353 forceparams[type].harmonic.rA,
3354 forceparams[type].harmonic.rB,
3355 cos_theta, lambda, &va, &dVdt);
3358 rij_1 = gmx_invsqrt(iprod(r_ij, r_ij));
3359 rkj_1 = gmx_invsqrt(iprod(r_kj, r_kj));
3360 rij_2 = rij_1*rij_1;
3361 rkj_2 = rkj_1*rkj_1;
3362 rijrkj_1 = rij_1*rkj_1; /* 23 */
3367 fprintf(debug, "G96ANGLES: costheta = %10g vth = %10g dV/dct = %10g\n",
3368 cos_theta, va, dVdt);
3371 for (m = 0; (m < DIM); m++) /* 42 */
3373 f_i[m] = dVdt*(r_kj[m]*rijrkj_1 - r_ij[m]*rij_2*cos_theta);
3374 f_k[m] = dVdt*(r_ij[m]*rijrkj_1 - r_kj[m]*rkj_2*cos_theta);
3375 f_j[m] = -f_i[m]-f_k[m];
3383 copy_ivec(SHIFT_IVEC(g, aj), jt);
3385 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3386 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3387 t1 = IVEC2IS(dt_ij);
3388 t2 = IVEC2IS(dt_kj);
3390 rvec_inc(fshift[t1], f_i);
3391 rvec_inc(fshift[CENTRAL], f_j);
3392 rvec_inc(fshift[t2], f_k); /* 9 */
3398 real cross_bond_bond(int nbonds,
3399 const t_iatom forceatoms[], const t_iparams forceparams[],
3400 const rvec x[], rvec f[], rvec fshift[],
3401 const t_pbc *pbc, const t_graph *g,
3402 real gmx_unused lambda, real gmx_unused *dvdlambda,
3403 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3404 int gmx_unused *global_atom_index)
3406 /* Potential from Lawrence and Skimmer, Chem. Phys. Lett. 372 (2003)
3409 int i, ai, aj, ak, type, m, t1, t2;
3411 real vtot, vrr, s1, s2, r1, r2, r1e, r2e, krr;
3413 ivec jt, dt_ij, dt_kj;
3416 for (i = 0; (i < nbonds); )
3418 type = forceatoms[i++];
3419 ai = forceatoms[i++];
3420 aj = forceatoms[i++];
3421 ak = forceatoms[i++];
3422 r1e = forceparams[type].cross_bb.r1e;
3423 r2e = forceparams[type].cross_bb.r2e;
3424 krr = forceparams[type].cross_bb.krr;
3426 /* Compute distance vectors ... */
3427 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
3428 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
3430 /* ... and their lengths */
3434 /* Deviations from ideality */
3438 /* Energy (can be negative!) */
3443 svmul(-krr*s2/r1, r_ij, f_i);
3444 svmul(-krr*s1/r2, r_kj, f_k);
3446 for (m = 0; (m < DIM); m++) /* 12 */
3448 f_j[m] = -f_i[m] - f_k[m];
3457 copy_ivec(SHIFT_IVEC(g, aj), jt);
3459 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3460 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3461 t1 = IVEC2IS(dt_ij);
3462 t2 = IVEC2IS(dt_kj);
3464 rvec_inc(fshift[t1], f_i);
3465 rvec_inc(fshift[CENTRAL], f_j);
3466 rvec_inc(fshift[t2], f_k); /* 9 */
3472 real cross_bond_angle(int nbonds,
3473 const t_iatom forceatoms[], const t_iparams forceparams[],
3474 const rvec x[], rvec f[], rvec fshift[],
3475 const t_pbc *pbc, const t_graph *g,
3476 real gmx_unused lambda, real gmx_unused *dvdlambda,
3477 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3478 int gmx_unused *global_atom_index)
3480 /* Potential from Lawrence and Skimmer, Chem. Phys. Lett. 372 (2003)
3483 int i, ai, aj, ak, type, m, t1, t2, t3;
3484 rvec r_ij, r_kj, r_ik;
3485 real vtot, vrt, s1, s2, s3, r1, r2, r3, r1e, r2e, r3e, krt, k1, k2, k3;
3487 ivec jt, dt_ij, dt_kj;
3490 for (i = 0; (i < nbonds); )
3492 type = forceatoms[i++];
3493 ai = forceatoms[i++];
3494 aj = forceatoms[i++];
3495 ak = forceatoms[i++];
3496 r1e = forceparams[type].cross_ba.r1e;
3497 r2e = forceparams[type].cross_ba.r2e;
3498 r3e = forceparams[type].cross_ba.r3e;
3499 krt = forceparams[type].cross_ba.krt;
3501 /* Compute distance vectors ... */
3502 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
3503 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
3504 t3 = pbc_rvec_sub(pbc, x[ai], x[ak], r_ik);
3506 /* ... and their lengths */
3511 /* Deviations from ideality */
3516 /* Energy (can be negative!) */
3517 vrt = krt*s3*(s1+s2);
3523 k3 = -krt*(s1+s2)/r3;
3524 for (m = 0; (m < DIM); m++)
3526 f_i[m] = k1*r_ij[m] + k3*r_ik[m];
3527 f_k[m] = k2*r_kj[m] - k3*r_ik[m];
3528 f_j[m] = -f_i[m] - f_k[m];
3531 for (m = 0; (m < DIM); m++) /* 12 */
3541 copy_ivec(SHIFT_IVEC(g, aj), jt);
3543 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3544 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3545 t1 = IVEC2IS(dt_ij);
3546 t2 = IVEC2IS(dt_kj);
3548 rvec_inc(fshift[t1], f_i);
3549 rvec_inc(fshift[CENTRAL], f_j);
3550 rvec_inc(fshift[t2], f_k); /* 9 */
3556 static real bonded_tab(const char *type, int table_nr,
3557 const bondedtable_t *table, real kA, real kB, real r,
3558 real lambda, real *V, real *F)
3560 real k, tabscale, *VFtab, rt, eps, eps2, Yt, Ft, Geps, Heps2, Fp, VV, FF;
3562 real v, f, dvdlambda;
3564 k = (1.0 - lambda)*kA + lambda*kB;
3566 tabscale = table->scale;
3567 VFtab = table->data;
3573 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",
3574 type, table_nr, r, n0, n0+1, table->n);
3581 Geps = VFtab[nnn+2]*eps;
3582 Heps2 = VFtab[nnn+3]*eps2;
3583 Fp = Ft + Geps + Heps2;
3585 FF = Fp + Geps + 2.0*Heps2;
3587 *F = -k*FF*tabscale;
3589 dvdlambda = (kB - kA)*VV;
3593 /* That was 22 flops */
3596 real tab_bonds(int nbonds,
3597 const t_iatom forceatoms[], const t_iparams forceparams[],
3598 const rvec x[], rvec f[], rvec fshift[],
3599 const t_pbc *pbc, const t_graph *g,
3600 real lambda, real *dvdlambda,
3601 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
3602 int gmx_unused *global_atom_index)
3604 int i, m, ki, ai, aj, type, table;
3605 real dr, dr2, fbond, vbond, fij, vtot;
3610 for (i = 0; (i < nbonds); )
3612 type = forceatoms[i++];
3613 ai = forceatoms[i++];
3614 aj = forceatoms[i++];
3616 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
3617 dr2 = iprod(dx, dx); /* 5 */
3618 dr = dr2*gmx_invsqrt(dr2); /* 10 */
3620 table = forceparams[type].tab.table;
3622 *dvdlambda += bonded_tab("bond", table,
3623 &fcd->bondtab[table],
3624 forceparams[type].tab.kA,
3625 forceparams[type].tab.kB,
3626 dr, lambda, &vbond, &fbond); /* 22 */
3634 vtot += vbond; /* 1*/
3635 fbond *= gmx_invsqrt(dr2); /* 6 */
3639 fprintf(debug, "TABBONDS: dr = %10g vbond = %10g fbond = %10g\n",
3645 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
3648 for (m = 0; (m < DIM); m++) /* 15 */
3653 fshift[ki][m] += fij;
3654 fshift[CENTRAL][m] -= fij;
3660 real tab_angles(int nbonds,
3661 const t_iatom forceatoms[], const t_iparams forceparams[],
3662 const rvec x[], rvec f[], rvec fshift[],
3663 const t_pbc *pbc, const t_graph *g,
3664 real lambda, real *dvdlambda,
3665 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
3666 int gmx_unused *global_atom_index)
3668 int i, ai, aj, ak, t1, t2, type, table;
3670 real cos_theta, cos_theta2, theta, dVdt, va, vtot;
3671 ivec jt, dt_ij, dt_kj;
3674 for (i = 0; (i < nbonds); )
3676 type = forceatoms[i++];
3677 ai = forceatoms[i++];
3678 aj = forceatoms[i++];
3679 ak = forceatoms[i++];
3681 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
3682 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
3684 table = forceparams[type].tab.table;
3686 *dvdlambda += bonded_tab("angle", table,
3687 &fcd->angletab[table],
3688 forceparams[type].tab.kA,
3689 forceparams[type].tab.kB,
3690 theta, lambda, &va, &dVdt); /* 22 */
3693 cos_theta2 = sqr(cos_theta); /* 1 */
3702 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
3703 sth = st*cos_theta; /* 1 */
3707 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
3708 theta*RAD2DEG, va, dVdt);
3711 nrkj2 = iprod(r_kj, r_kj); /* 5 */
3712 nrij2 = iprod(r_ij, r_ij);
3714 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
3715 cii = sth/nrij2; /* 10 */
3716 ckk = sth/nrkj2; /* 10 */
3718 for (m = 0; (m < DIM); m++) /* 39 */
3720 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
3721 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
3722 f_j[m] = -f_i[m]-f_k[m];
3729 copy_ivec(SHIFT_IVEC(g, aj), jt);
3731 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3732 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3733 t1 = IVEC2IS(dt_ij);
3734 t2 = IVEC2IS(dt_kj);
3736 rvec_inc(fshift[t1], f_i);
3737 rvec_inc(fshift[CENTRAL], f_j);
3738 rvec_inc(fshift[t2], f_k);
3744 real tab_dihs(int nbonds,
3745 const t_iatom forceatoms[], const t_iparams forceparams[],
3746 const rvec x[], rvec f[], rvec fshift[],
3747 const t_pbc *pbc, const t_graph *g,
3748 real lambda, real *dvdlambda,
3749 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
3750 int gmx_unused *global_atom_index)
3752 int i, type, ai, aj, ak, al, table;
3754 rvec r_ij, r_kj, r_kl, m, n;
3755 real phi, sign, ddphi, vpd, vtot;
3758 for (i = 0; (i < nbonds); )
3760 type = forceatoms[i++];
3761 ai = forceatoms[i++];
3762 aj = forceatoms[i++];
3763 ak = forceatoms[i++];
3764 al = forceatoms[i++];
3766 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
3767 &sign, &t1, &t2, &t3); /* 84 */
3769 table = forceparams[type].tab.table;
3771 /* Hopefully phi+M_PI never results in values < 0 */
3772 *dvdlambda += bonded_tab("dihedral", table,
3773 &fcd->dihtab[table],
3774 forceparams[type].tab.kA,
3775 forceparams[type].tab.kB,
3776 phi+M_PI, lambda, &vpd, &ddphi);
3779 do_dih_fup(ai, aj, ak, al, -ddphi, r_ij, r_kj, r_kl, m, n,
3780 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
3783 fprintf(debug, "pdih: (%d,%d,%d,%d) phi=%g\n",
3784 ai, aj, ak, al, phi);
3791 /* Return if this is a potential calculated in bondfree.c,
3792 * i.e. an interaction that actually calculates a potential and
3793 * works on multiple atoms (not e.g. a connection or a position restraint).
3795 static gmx_inline gmx_bool ftype_is_bonded_potential(int ftype)
3798 (interaction_function[ftype].flags & IF_BOND) &&
3799 !(ftype == F_CONNBONDS || ftype == F_POSRES || ftype == F_FBPOSRES) &&
3800 (ftype < F_GB12 || ftype > F_GB14);
3803 static void divide_bondeds_over_threads(t_idef *idef, int nthreads)
3810 idef->nthreads = nthreads;
3812 if (F_NRE*(nthreads+1) > idef->il_thread_division_nalloc)
3814 idef->il_thread_division_nalloc = F_NRE*(nthreads+1);
3815 snew(idef->il_thread_division, idef->il_thread_division_nalloc);
3818 for (ftype = 0; ftype < F_NRE; ftype++)
3820 if (ftype_is_bonded_potential(ftype))
3822 nat1 = interaction_function[ftype].nratoms + 1;
3824 for (t = 0; t <= nthreads; t++)
3826 /* Divide the interactions equally over the threads.
3827 * When the different types of bonded interactions
3828 * are distributed roughly equally over the threads,
3829 * this should lead to well localized output into
3830 * the force buffer on each thread.
3831 * If this is not the case, a more advanced scheme
3832 * (not implemented yet) will do better.
3834 il_nr_thread = (((idef->il[ftype].nr/nat1)*t)/nthreads)*nat1;
3836 /* Ensure that distance restraint pairs with the same label
3837 * end up on the same thread.
3838 * This is slighlty tricky code, since the next for iteration
3839 * may have an initial il_nr_thread lower than the final value
3840 * in the previous iteration, but this will anyhow be increased
3841 * to the approriate value again by this while loop.
3843 while (ftype == F_DISRES &&
3845 il_nr_thread < idef->il[ftype].nr &&
3846 idef->iparams[idef->il[ftype].iatoms[il_nr_thread]].disres.label ==
3847 idef->iparams[idef->il[ftype].iatoms[il_nr_thread-nat1]].disres.label)
3849 il_nr_thread += nat1;
3852 idef->il_thread_division[ftype*(nthreads+1)+t] = il_nr_thread;
3859 calc_bonded_reduction_mask(const t_idef *idef,
3864 int ftype, nb, nat1, nb0, nb1, i, a;
3868 for (ftype = 0; ftype < F_NRE; ftype++)
3870 if (ftype_is_bonded_potential(ftype))
3872 nb = idef->il[ftype].nr;
3875 nat1 = interaction_function[ftype].nratoms + 1;
3877 /* Divide this interaction equally over the threads.
3878 * This is not stored: should match division in calc_bonds.
3880 nb0 = idef->il_thread_division[ftype*(nt+1)+t];
3881 nb1 = idef->il_thread_division[ftype*(nt+1)+t+1];
3883 for (i = nb0; i < nb1; i += nat1)
3885 for (a = 1; a < nat1; a++)
3887 mask |= (1U << (idef->il[ftype].iatoms[i+a]>>shift));
3897 void setup_bonded_threading(t_forcerec *fr, t_idef *idef)
3899 #define MAX_BLOCK_BITS 32
3903 assert(fr->nthreads >= 1);
3905 /* Divide the bonded interaction over the threads */
3906 divide_bondeds_over_threads(idef, fr->nthreads);
3908 if (fr->nthreads == 1)
3915 /* We divide the force array in a maximum of 32 blocks.
3916 * Minimum force block reduction size is 2^6=64.
3919 while (fr->natoms_force > (int)(MAX_BLOCK_BITS*(1U<<fr->red_ashift)))
3925 fprintf(debug, "bonded force buffer block atom shift %d bits\n",
3929 /* Determine to which blocks each thread's bonded force calculation
3930 * contributes. Store this is a mask for each thread.
3932 #pragma omp parallel for num_threads(fr->nthreads) schedule(static)
3933 for (t = 1; t < fr->nthreads; t++)
3935 fr->f_t[t].red_mask =
3936 calc_bonded_reduction_mask(idef, fr->red_ashift, t, fr->nthreads);
3939 /* Determine the maximum number of blocks we need to reduce over */
3942 for (t = 0; t < fr->nthreads; t++)
3945 for (b = 0; b < MAX_BLOCK_BITS; b++)
3947 if (fr->f_t[t].red_mask & (1U<<b))
3949 fr->red_nblock = max(fr->red_nblock, b+1);
3955 fprintf(debug, "thread %d flags %x count %d\n",
3956 t, fr->f_t[t].red_mask, c);
3962 fprintf(debug, "Number of blocks to reduce: %d of size %d\n",
3963 fr->red_nblock, 1<<fr->red_ashift);
3964 fprintf(debug, "Reduction density %.2f density/#thread %.2f\n",
3965 ctot*(1<<fr->red_ashift)/(double)fr->natoms_force,
3966 ctot*(1<<fr->red_ashift)/(double)(fr->natoms_force*fr->nthreads));
3970 static void zero_thread_forces(f_thread_t *f_t, int n,
3971 int nblock, int blocksize)
3973 int b, a0, a1, a, i, j;
3975 if (n > f_t->f_nalloc)
3977 f_t->f_nalloc = over_alloc_large(n);
3978 srenew(f_t->f, f_t->f_nalloc);
3981 if (f_t->red_mask != 0)
3983 for (b = 0; b < nblock; b++)
3985 if (f_t->red_mask && (1U<<b))
3988 a1 = min((b+1)*blocksize, n);
3989 for (a = a0; a < a1; a++)
3991 clear_rvec(f_t->f[a]);
3996 for (i = 0; i < SHIFTS; i++)
3998 clear_rvec(f_t->fshift[i]);
4000 for (i = 0; i < F_NRE; i++)
4004 for (i = 0; i < egNR; i++)
4006 for (j = 0; j < f_t->grpp.nener; j++)
4008 f_t->grpp.ener[i][j] = 0;
4011 for (i = 0; i < efptNR; i++)
4017 static void reduce_thread_force_buffer(int n, rvec *f,
4018 int nthreads, f_thread_t *f_t,
4019 int nblock, int block_size)
4021 /* The max thread number is arbitrary,
4022 * we used a fixed number to avoid memory management.
4023 * Using more than 16 threads is probably never useful performance wise.
4025 #define MAX_BONDED_THREADS 256
4028 if (nthreads > MAX_BONDED_THREADS)
4030 gmx_fatal(FARGS, "Can not reduce bonded forces on more than %d threads",
4031 MAX_BONDED_THREADS);
4034 /* This reduction can run on any number of threads,
4035 * independently of nthreads.
4037 #pragma omp parallel for num_threads(nthreads) schedule(static)
4038 for (b = 0; b < nblock; b++)
4040 rvec *fp[MAX_BONDED_THREADS];
4044 /* Determine which threads contribute to this block */
4046 for (ft = 1; ft < nthreads; ft++)
4048 if (f_t[ft].red_mask & (1U<<b))
4050 fp[nfb++] = f_t[ft].f;
4055 /* Reduce force buffers for threads that contribute */
4057 a1 = (b+1)*block_size;
4059 for (a = a0; a < a1; a++)
4061 for (fb = 0; fb < nfb; fb++)
4063 rvec_inc(f[a], fp[fb][a]);
4070 static void reduce_thread_forces(int n, rvec *f, rvec *fshift,
4071 real *ener, gmx_grppairener_t *grpp, real *dvdl,
4072 int nthreads, f_thread_t *f_t,
4073 int nblock, int block_size,
4074 gmx_bool bCalcEnerVir,
4079 /* Reduce the bonded force buffer */
4080 reduce_thread_force_buffer(n, f, nthreads, f_t, nblock, block_size);
4083 /* When necessary, reduce energy and virial using one thread only */
4088 for (i = 0; i < SHIFTS; i++)
4090 for (t = 1; t < nthreads; t++)
4092 rvec_inc(fshift[i], f_t[t].fshift[i]);
4095 for (i = 0; i < F_NRE; i++)
4097 for (t = 1; t < nthreads; t++)
4099 ener[i] += f_t[t].ener[i];
4102 for (i = 0; i < egNR; i++)
4104 for (j = 0; j < f_t[1].grpp.nener; j++)
4106 for (t = 1; t < nthreads; t++)
4109 grpp->ener[i][j] += f_t[t].grpp.ener[i][j];
4115 for (i = 0; i < efptNR; i++)
4118 for (t = 1; t < nthreads; t++)
4120 dvdl[i] += f_t[t].dvdl[i];
4127 static real calc_one_bond(FILE *fplog, int thread,
4128 int ftype, const t_idef *idef,
4129 rvec x[], rvec f[], rvec fshift[],
4131 const t_pbc *pbc, const t_graph *g,
4132 gmx_grppairener_t *grpp,
4134 real *lambda, real *dvdl,
4135 const t_mdatoms *md, t_fcdata *fcd,
4136 gmx_bool bCalcEnerVir,
4137 int *global_atom_index, gmx_bool bPrintSepPot)
4139 int nat1, nbonds, efptFTYPE;
4144 if (IS_RESTRAINT_TYPE(ftype))
4146 efptFTYPE = efptRESTRAINT;
4150 efptFTYPE = efptBONDED;
4153 nat1 = interaction_function[ftype].nratoms + 1;
4154 nbonds = idef->il[ftype].nr/nat1;
4155 iatoms = idef->il[ftype].iatoms;
4157 nb0 = idef->il_thread_division[ftype*(idef->nthreads+1)+thread];
4158 nbn = idef->il_thread_division[ftype*(idef->nthreads+1)+thread+1] - nb0;
4160 if (!IS_LISTED_LJ_C(ftype))
4162 if (ftype == F_CMAP)
4164 v = cmap_dihs(nbn, iatoms+nb0,
4165 idef->iparams, &idef->cmap_grid,
4166 (const rvec*)x, f, fshift,
4167 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
4168 md, fcd, global_atom_index);
4171 else if (ftype == F_ANGLES &&
4172 !bCalcEnerVir && fr->efep == efepNO)
4174 /* No energies, shift forces, dvdl */
4175 angles_noener_simd(nbn, idef->il[ftype].iatoms+nb0,
4178 pbc, g, lambda[efptFTYPE], md, fcd,
4183 else if (ftype == F_PDIHS &&
4184 !bCalcEnerVir && fr->efep == efepNO)
4186 /* No energies, shift forces, dvdl */
4187 #ifndef SIMD_BONDEDS
4192 (nbn, idef->il[ftype].iatoms+nb0,
4195 pbc, g, lambda[efptFTYPE], md, fcd,
4201 v = interaction_function[ftype].ifunc(nbn, iatoms+nb0,
4203 (const rvec*)x, f, fshift,
4204 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
4205 md, fcd, global_atom_index);
4209 fprintf(fplog, " %-23s #%4d V %12.5e dVdl %12.5e\n",
4210 interaction_function[ftype].longname,
4211 nbonds, v, lambda[efptFTYPE]);
4216 v = do_nonbonded_listed(ftype, nbn, iatoms+nb0, idef->iparams, (const rvec*)x, f, fshift,
4217 pbc, g, lambda, dvdl, md, fr, grpp, global_atom_index);
4221 fprintf(fplog, " %-5s + %-15s #%4d dVdl %12.5e\n",
4222 interaction_function[ftype].longname,
4223 interaction_function[F_LJ14].longname, nbonds, dvdl[efptVDW]);
4224 fprintf(fplog, " %-5s + %-15s #%4d dVdl %12.5e\n",
4225 interaction_function[ftype].longname,
4226 interaction_function[F_COUL14].longname, nbonds, dvdl[efptCOUL]);
4232 inc_nrnb(nrnb, interaction_function[ftype].nrnb_ind, nbonds);
4238 void calc_bonds(FILE *fplog, const gmx_multisim_t *ms,
4240 rvec x[], history_t *hist,
4241 rvec f[], t_forcerec *fr,
4242 const t_pbc *pbc, const t_graph *g,
4243 gmx_enerdata_t *enerd, t_nrnb *nrnb,
4245 const t_mdatoms *md,
4246 t_fcdata *fcd, int *global_atom_index,
4247 t_atomtypes gmx_unused *atype, gmx_genborn_t gmx_unused *born,
4249 gmx_bool bPrintSepPot, gmx_int64_t step)
4251 gmx_bool bCalcEnerVir;
4253 real v, dvdl[efptNR], dvdl_dum[efptNR]; /* The dummy array is to have a place to store the dhdl at other values
4254 of lambda, which will be thrown away in the end*/
4255 const t_pbc *pbc_null;
4259 assert(fr->nthreads == idef->nthreads);
4261 bCalcEnerVir = (force_flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY));
4263 for (i = 0; i < efptNR; i++)
4277 fprintf(fplog, "Step %s: bonded V and dVdl for this node\n",
4278 gmx_step_str(step, buf));
4284 p_graph(debug, "Bondage is fun", g);
4288 /* Do pre force calculation stuff which might require communication */
4289 if (idef->il[F_ORIRES].nr)
4291 enerd->term[F_ORIRESDEV] =
4292 calc_orires_dev(ms, idef->il[F_ORIRES].nr,
4293 idef->il[F_ORIRES].iatoms,
4294 idef->iparams, md, (const rvec*)x,
4295 pbc_null, fcd, hist);
4297 if (idef->il[F_DISRES].nr)
4299 calc_disres_R_6(idef->il[F_DISRES].nr,
4300 idef->il[F_DISRES].iatoms,
4301 idef->iparams, (const rvec*)x, pbc_null,
4304 if (fcd->disres.nsystems > 1)
4306 gmx_sum_sim(2*fcd->disres.nres, fcd->disres.Rt_6, ms);
4311 #pragma omp parallel for num_threads(fr->nthreads) schedule(static)
4312 for (thread = 0; thread < fr->nthreads; thread++)
4319 gmx_grppairener_t *grpp;
4324 fshift = fr->fshift;
4326 grpp = &enerd->grpp;
4331 zero_thread_forces(&fr->f_t[thread], fr->natoms_force,
4332 fr->red_nblock, 1<<fr->red_ashift);
4334 ft = fr->f_t[thread].f;
4335 fshift = fr->f_t[thread].fshift;
4336 epot = fr->f_t[thread].ener;
4337 grpp = &fr->f_t[thread].grpp;
4338 dvdlt = fr->f_t[thread].dvdl;
4340 /* Loop over all bonded force types to calculate the bonded forces */
4341 for (ftype = 0; (ftype < F_NRE); ftype++)
4343 if (idef->il[ftype].nr > 0 && ftype_is_bonded_potential(ftype))
4345 v = calc_one_bond(fplog, thread, ftype, idef, x,
4346 ft, fshift, fr, pbc_null, g, grpp,
4347 nrnb, lambda, dvdlt,
4348 md, fcd, bCalcEnerVir,
4349 global_atom_index, bPrintSepPot);
4354 if (fr->nthreads > 1)
4356 reduce_thread_forces(fr->natoms_force, f, fr->fshift,
4357 enerd->term, &enerd->grpp, dvdl,
4358 fr->nthreads, fr->f_t,
4359 fr->red_nblock, 1<<fr->red_ashift,
4361 force_flags & GMX_FORCE_DHDL);
4363 if (force_flags & GMX_FORCE_DHDL)
4365 for (i = 0; i < efptNR; i++)
4367 enerd->dvdl_nonlin[i] += dvdl[i];
4371 /* Copy the sum of violations for the distance restraints from fcd */
4374 enerd->term[F_DISRESVIOL] = fcd->disres.sumviol;
4379 void calc_bonds_lambda(FILE *fplog,
4383 const t_pbc *pbc, const t_graph *g,
4384 gmx_grppairener_t *grpp, real *epot, t_nrnb *nrnb,
4386 const t_mdatoms *md,
4388 int *global_atom_index)
4390 int i, ftype, nr_nonperturbed, nr;
4392 real dvdl_dum[efptNR];
4394 const t_pbc *pbc_null;
4406 /* Copy the whole idef, so we can modify the contents locally */
4408 idef_fe.nthreads = 1;
4409 snew(idef_fe.il_thread_division, F_NRE*(idef_fe.nthreads+1));
4411 /* We already have the forces, so we use temp buffers here */
4412 snew(f, fr->natoms_force);
4413 snew(fshift, SHIFTS);
4415 /* Loop over all bonded force types to calculate the bonded energies */
4416 for (ftype = 0; (ftype < F_NRE); ftype++)
4418 if (ftype_is_bonded_potential(ftype))
4420 /* Set the work range of thread 0 to the perturbed bondeds only */
4421 nr_nonperturbed = idef->il[ftype].nr_nonperturbed;
4422 nr = idef->il[ftype].nr;
4423 idef_fe.il_thread_division[ftype*2+0] = nr_nonperturbed;
4424 idef_fe.il_thread_division[ftype*2+1] = nr;
4426 /* This is only to get the flop count correct */
4427 idef_fe.il[ftype].nr = nr - nr_nonperturbed;
4429 if (nr - nr_nonperturbed > 0)
4431 v = calc_one_bond(fplog, 0, ftype, &idef_fe,
4432 x, f, fshift, fr, pbc_null, g,
4433 grpp, nrnb, lambda, dvdl_dum,
4435 global_atom_index, FALSE);
4444 sfree(idef_fe.il_thread_division);