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 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
54 #include "gmx_fatal.h"
60 #include "nonbonded.h"
62 /* Include the SIMD macro file and then check for support */
63 #include "gmx_simd_macros.h"
64 #if defined GMX_HAVE_SIMD_MACROS && defined GMX_SIMD_HAVE_TRIGONOMETRIC
66 #include "gmx_simd_vec.h"
69 /* Find a better place for this? */
70 const int cmap_coeff_matrix[] = {
71 1, 0, -3, 2, 0, 0, 0, 0, -3, 0, 9, -6, 2, 0, -6, 4,
72 0, 0, 0, 0, 0, 0, 0, 0, 3, 0, -9, 6, -2, 0, 6, -4,
73 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 9, -6, 0, 0, -6, 4,
74 0, 0, 3, -2, 0, 0, 0, 0, 0, 0, -9, 6, 0, 0, 6, -4,
75 0, 0, 0, 0, 1, 0, -3, 2, -2, 0, 6, -4, 1, 0, -3, 2,
76 0, 0, 0, 0, 0, 0, 0, 0, -1, 0, 3, -2, 1, 0, -3, 2,
77 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -3, 2, 0, 0, 3, -2,
78 0, 0, 0, 0, 0, 0, 3, -2, 0, 0, -6, 4, 0, 0, 3, -2,
79 0, 1, -2, 1, 0, 0, 0, 0, 0, -3, 6, -3, 0, 2, -4, 2,
80 0, 0, 0, 0, 0, 0, 0, 0, 0, 3, -6, 3, 0, -2, 4, -2,
81 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -3, 3, 0, 0, 2, -2,
82 0, 0, -1, 1, 0, 0, 0, 0, 0, 0, 3, -3, 0, 0, -2, 2,
83 0, 0, 0, 0, 0, 1, -2, 1, 0, -2, 4, -2, 0, 1, -2, 1,
84 0, 0, 0, 0, 0, 0, 0, 0, 0, -1, 2, -1, 0, 1, -2, 1,
85 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, -1, 0, 0, -1, 1,
86 0, 0, 0, 0, 0, 0, -1, 1, 0, 0, 2, -2, 0, 0, -1, 1
91 int glatnr(int *global_atom_index, int i)
95 if (global_atom_index == NULL)
101 atnr = global_atom_index[i] + 1;
107 static int pbc_rvec_sub(const t_pbc *pbc, const rvec xi, const rvec xj, rvec dx)
111 return pbc_dx_aiuc(pbc, xi, xj, dx);
115 rvec_sub(xi, xj, dx);
122 /* SIMD PBC data structure, containing 1/boxdiag and the box vectors */
135 /* Set the SIMD pbc data from a normal t_pbc struct */
136 static void set_pbc_simd(const t_pbc *pbc, pbc_simd_t *pbc_simd)
141 /* Setting inv_bdiag to 0 effectively turns off PBC */
142 clear_rvec(inv_bdiag);
145 for (d = 0; d < pbc->ndim_ePBC; d++)
147 inv_bdiag[d] = 1.0/pbc->box[d][d];
151 pbc_simd->inv_bzz = gmx_set1_pr(inv_bdiag[ZZ]);
152 pbc_simd->inv_byy = gmx_set1_pr(inv_bdiag[YY]);
153 pbc_simd->inv_bxx = gmx_set1_pr(inv_bdiag[XX]);
157 pbc_simd->bzx = gmx_set1_pr(pbc->box[ZZ][XX]);
158 pbc_simd->bzy = gmx_set1_pr(pbc->box[ZZ][YY]);
159 pbc_simd->bzz = gmx_set1_pr(pbc->box[ZZ][ZZ]);
160 pbc_simd->byx = gmx_set1_pr(pbc->box[YY][XX]);
161 pbc_simd->byy = gmx_set1_pr(pbc->box[YY][YY]);
162 pbc_simd->bxx = gmx_set1_pr(pbc->box[XX][XX]);
166 pbc_simd->bzx = gmx_setzero_pr();
167 pbc_simd->bzy = gmx_setzero_pr();
168 pbc_simd->bzz = gmx_setzero_pr();
169 pbc_simd->byx = gmx_setzero_pr();
170 pbc_simd->byy = gmx_setzero_pr();
171 pbc_simd->bxx = gmx_setzero_pr();
175 /* Correct distance vector *dx,*dy,*dz for PBC using SIMD */
176 static gmx_inline void
177 pbc_dx_simd(gmx_mm_pr *dx, gmx_mm_pr *dy, gmx_mm_pr *dz,
178 const pbc_simd_t *pbc)
182 sh = gmx_round_pr(gmx_mul_pr(*dz, pbc->inv_bzz));
183 *dx = gmx_nmsub_pr(sh, pbc->bzx, *dx);
184 *dy = gmx_nmsub_pr(sh, pbc->bzy, *dy);
185 *dz = gmx_nmsub_pr(sh, pbc->bzz, *dz);
187 sh = gmx_round_pr(gmx_mul_pr(*dy, pbc->inv_byy));
188 *dx = gmx_nmsub_pr(sh, pbc->byx, *dx);
189 *dy = gmx_nmsub_pr(sh, pbc->byy, *dy);
191 sh = gmx_round_pr(gmx_mul_pr(*dx, pbc->inv_bxx));
192 *dx = gmx_nmsub_pr(sh, pbc->bxx, *dx);
195 #endif /* SIMD_BONDEDS */
198 * Morse potential bond by Frank Everdij
200 * Three parameters needed:
202 * b0 = equilibrium distance in nm
203 * be = beta in nm^-1 (actually, it's nu_e*Sqrt(2*pi*pi*mu/D_e))
204 * cb = well depth in kJ/mol
206 * Note: the potential is referenced to be +cb at infinite separation
207 * and zero at the equilibrium distance!
210 real morse_bonds(int nbonds,
211 const t_iatom forceatoms[], const t_iparams forceparams[],
212 const rvec x[], rvec f[], rvec fshift[],
213 const t_pbc *pbc, const t_graph *g,
214 real lambda, real *dvdlambda,
215 const t_mdatoms *md, t_fcdata *fcd,
216 int *global_atom_index)
218 const real one = 1.0;
219 const real two = 2.0;
220 real dr, dr2, temp, omtemp, cbomtemp, fbond, vbond, fij, vtot;
221 real b0, be, cb, b0A, beA, cbA, b0B, beB, cbB, L1;
223 int i, m, ki, type, ai, aj;
227 for (i = 0; (i < nbonds); )
229 type = forceatoms[i++];
230 ai = forceatoms[i++];
231 aj = forceatoms[i++];
233 b0A = forceparams[type].morse.b0A;
234 beA = forceparams[type].morse.betaA;
235 cbA = forceparams[type].morse.cbA;
237 b0B = forceparams[type].morse.b0B;
238 beB = forceparams[type].morse.betaB;
239 cbB = forceparams[type].morse.cbB;
241 L1 = one-lambda; /* 1 */
242 b0 = L1*b0A + lambda*b0B; /* 3 */
243 be = L1*beA + lambda*beB; /* 3 */
244 cb = L1*cbA + lambda*cbB; /* 3 */
246 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
247 dr2 = iprod(dx, dx); /* 5 */
248 dr = dr2*gmx_invsqrt(dr2); /* 10 */
249 temp = exp(-be*(dr-b0)); /* 12 */
253 /* bonds are constrainted. This may _not_ include bond constraints if they are lambda dependent */
254 *dvdlambda += cbB-cbA;
258 omtemp = one-temp; /* 1 */
259 cbomtemp = cb*omtemp; /* 1 */
260 vbond = cbomtemp*omtemp; /* 1 */
261 fbond = -two*be*temp*cbomtemp*gmx_invsqrt(dr2); /* 9 */
262 vtot += vbond; /* 1 */
264 *dvdlambda += (cbB - cbA) * omtemp * omtemp - (2-2*omtemp)*omtemp * cb * ((b0B-b0A)*be - (beB-beA)*(dr-b0)); /* 15 */
268 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
272 for (m = 0; (m < DIM); m++) /* 15 */
277 fshift[ki][m] += fij;
278 fshift[CENTRAL][m] -= fij;
284 real cubic_bonds(int nbonds,
285 const t_iatom forceatoms[], const t_iparams forceparams[],
286 const rvec x[], rvec f[], rvec fshift[],
287 const t_pbc *pbc, const t_graph *g,
288 real lambda, real *dvdlambda,
289 const t_mdatoms *md, t_fcdata *fcd,
290 int *global_atom_index)
292 const real three = 3.0;
293 const real two = 2.0;
295 real dr, dr2, dist, kdist, kdist2, fbond, vbond, fij, vtot;
297 int i, m, ki, type, ai, aj;
301 for (i = 0; (i < nbonds); )
303 type = forceatoms[i++];
304 ai = forceatoms[i++];
305 aj = forceatoms[i++];
307 b0 = forceparams[type].cubic.b0;
308 kb = forceparams[type].cubic.kb;
309 kcub = forceparams[type].cubic.kcub;
311 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
312 dr2 = iprod(dx, dx); /* 5 */
319 dr = dr2*gmx_invsqrt(dr2); /* 10 */
324 vbond = kdist2 + kcub*kdist2*dist;
325 fbond = -(two*kdist + three*kdist2*kcub)/dr;
327 vtot += vbond; /* 21 */
331 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
334 for (m = 0; (m < DIM); m++) /* 15 */
339 fshift[ki][m] += fij;
340 fshift[CENTRAL][m] -= fij;
346 real FENE_bonds(int nbonds,
347 const t_iatom forceatoms[], const t_iparams forceparams[],
348 const rvec x[], rvec f[], rvec fshift[],
349 const t_pbc *pbc, const t_graph *g,
350 real lambda, real *dvdlambda,
351 const t_mdatoms *md, t_fcdata *fcd,
352 int *global_atom_index)
354 const real half = 0.5;
355 const real one = 1.0;
357 real dr, dr2, bm2, omdr2obm2, fbond, vbond, fij, vtot;
359 int i, m, ki, type, ai, aj;
363 for (i = 0; (i < nbonds); )
365 type = forceatoms[i++];
366 ai = forceatoms[i++];
367 aj = forceatoms[i++];
369 bm = forceparams[type].fene.bm;
370 kb = forceparams[type].fene.kb;
372 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
373 dr2 = iprod(dx, dx); /* 5 */
385 "r^2 (%f) >= bm^2 (%f) in FENE bond between atoms %d and %d",
387 glatnr(global_atom_index, ai),
388 glatnr(global_atom_index, aj));
391 omdr2obm2 = one - dr2/bm2;
393 vbond = -half*kb*bm2*log(omdr2obm2);
394 fbond = -kb/omdr2obm2;
396 vtot += vbond; /* 35 */
400 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
403 for (m = 0; (m < DIM); m++) /* 15 */
408 fshift[ki][m] += fij;
409 fshift[CENTRAL][m] -= fij;
415 real harmonic(real kA, real kB, real xA, real xB, real x, real lambda,
418 const real half = 0.5;
419 real L1, kk, x0, dx, dx2;
420 real v, f, dvdlambda;
423 kk = L1*kA+lambda*kB;
424 x0 = L1*xA+lambda*xB;
431 dvdlambda = half*(kB-kA)*dx2 + (xA-xB)*kk*dx;
438 /* That was 19 flops */
442 real bonds(int nbonds,
443 const t_iatom forceatoms[], const t_iparams forceparams[],
444 const rvec x[], rvec f[], rvec fshift[],
445 const t_pbc *pbc, const t_graph *g,
446 real lambda, real *dvdlambda,
447 const t_mdatoms *md, t_fcdata *fcd,
448 int *global_atom_index)
450 int i, m, ki, ai, aj, type;
451 real dr, dr2, fbond, vbond, fij, vtot;
456 for (i = 0; (i < nbonds); )
458 type = forceatoms[i++];
459 ai = forceatoms[i++];
460 aj = forceatoms[i++];
462 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
463 dr2 = iprod(dx, dx); /* 5 */
464 dr = dr2*gmx_invsqrt(dr2); /* 10 */
466 *dvdlambda += harmonic(forceparams[type].harmonic.krA,
467 forceparams[type].harmonic.krB,
468 forceparams[type].harmonic.rA,
469 forceparams[type].harmonic.rB,
470 dr, lambda, &vbond, &fbond); /* 19 */
478 vtot += vbond; /* 1*/
479 fbond *= gmx_invsqrt(dr2); /* 6 */
483 fprintf(debug, "BONDS: dr = %10g vbond = %10g fbond = %10g\n",
489 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
492 for (m = 0; (m < DIM); m++) /* 15 */
497 fshift[ki][m] += fij;
498 fshift[CENTRAL][m] -= fij;
504 real restraint_bonds(int nbonds,
505 const t_iatom forceatoms[], const t_iparams forceparams[],
506 const rvec x[], rvec f[], rvec fshift[],
507 const t_pbc *pbc, const t_graph *g,
508 real lambda, real *dvdlambda,
509 const t_mdatoms *md, t_fcdata *fcd,
510 int *global_atom_index)
512 int i, m, ki, ai, aj, type;
513 real dr, dr2, fbond, vbond, fij, vtot;
515 real low, dlow, up1, dup1, up2, dup2, k, dk;
523 for (i = 0; (i < nbonds); )
525 type = forceatoms[i++];
526 ai = forceatoms[i++];
527 aj = forceatoms[i++];
529 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
530 dr2 = iprod(dx, dx); /* 5 */
531 dr = dr2*gmx_invsqrt(dr2); /* 10 */
533 low = L1*forceparams[type].restraint.lowA + lambda*forceparams[type].restraint.lowB;
534 dlow = -forceparams[type].restraint.lowA + forceparams[type].restraint.lowB;
535 up1 = L1*forceparams[type].restraint.up1A + lambda*forceparams[type].restraint.up1B;
536 dup1 = -forceparams[type].restraint.up1A + forceparams[type].restraint.up1B;
537 up2 = L1*forceparams[type].restraint.up2A + lambda*forceparams[type].restraint.up2B;
538 dup2 = -forceparams[type].restraint.up2A + forceparams[type].restraint.up2B;
539 k = L1*forceparams[type].restraint.kA + lambda*forceparams[type].restraint.kB;
540 dk = -forceparams[type].restraint.kA + forceparams[type].restraint.kB;
549 *dvdlambda += 0.5*dk*drh2 - k*dlow*drh;
562 *dvdlambda += 0.5*dk*drh2 - k*dup1*drh;
567 vbond = k*(up2 - up1)*(0.5*(up2 - up1) + drh);
568 fbond = -k*(up2 - up1);
569 *dvdlambda += dk*(up2 - up1)*(0.5*(up2 - up1) + drh)
570 + k*(dup2 - dup1)*(up2 - up1 + drh)
571 - k*(up2 - up1)*dup2;
579 vtot += vbond; /* 1*/
580 fbond *= gmx_invsqrt(dr2); /* 6 */
584 fprintf(debug, "BONDS: dr = %10g vbond = %10g fbond = %10g\n",
590 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
593 for (m = 0; (m < DIM); m++) /* 15 */
598 fshift[ki][m] += fij;
599 fshift[CENTRAL][m] -= fij;
606 real polarize(int nbonds,
607 const t_iatom forceatoms[], const t_iparams forceparams[],
608 const rvec x[], rvec f[], rvec fshift[],
609 const t_pbc *pbc, const t_graph *g,
610 real lambda, real *dvdlambda,
611 const t_mdatoms *md, t_fcdata *fcd,
612 int *global_atom_index)
614 int i, m, ki, ai, aj, type;
615 real dr, dr2, fbond, vbond, fij, vtot, ksh;
620 for (i = 0; (i < nbonds); )
622 type = forceatoms[i++];
623 ai = forceatoms[i++];
624 aj = forceatoms[i++];
625 ksh = sqr(md->chargeA[aj])*ONE_4PI_EPS0/forceparams[type].polarize.alpha;
628 fprintf(debug, "POL: local ai = %d aj = %d ksh = %.3f\n", ai, aj, ksh);
631 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
632 dr2 = iprod(dx, dx); /* 5 */
633 dr = dr2*gmx_invsqrt(dr2); /* 10 */
635 *dvdlambda += harmonic(ksh, ksh, 0, 0, dr, lambda, &vbond, &fbond); /* 19 */
642 vtot += vbond; /* 1*/
643 fbond *= gmx_invsqrt(dr2); /* 6 */
647 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
650 for (m = 0; (m < DIM); m++) /* 15 */
655 fshift[ki][m] += fij;
656 fshift[CENTRAL][m] -= fij;
662 real anharm_polarize(int nbonds,
663 const t_iatom forceatoms[], const t_iparams forceparams[],
664 const rvec x[], rvec f[], rvec fshift[],
665 const t_pbc *pbc, const t_graph *g,
666 real lambda, real *dvdlambda,
667 const t_mdatoms *md, t_fcdata *fcd,
668 int *global_atom_index)
670 int i, m, ki, ai, aj, type;
671 real dr, dr2, fbond, vbond, fij, vtot, ksh, khyp, drcut, ddr, ddr3;
676 for (i = 0; (i < nbonds); )
678 type = forceatoms[i++];
679 ai = forceatoms[i++];
680 aj = forceatoms[i++];
681 ksh = sqr(md->chargeA[aj])*ONE_4PI_EPS0/forceparams[type].anharm_polarize.alpha; /* 7*/
682 khyp = forceparams[type].anharm_polarize.khyp;
683 drcut = forceparams[type].anharm_polarize.drcut;
686 fprintf(debug, "POL: local ai = %d aj = %d ksh = %.3f\n", ai, aj, ksh);
689 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
690 dr2 = iprod(dx, dx); /* 5 */
691 dr = dr2*gmx_invsqrt(dr2); /* 10 */
693 *dvdlambda += harmonic(ksh, ksh, 0, 0, dr, lambda, &vbond, &fbond); /* 19 */
704 vbond += khyp*ddr*ddr3;
705 fbond -= 4*khyp*ddr3;
707 fbond *= gmx_invsqrt(dr2); /* 6 */
708 vtot += vbond; /* 1*/
712 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
715 for (m = 0; (m < DIM); m++) /* 15 */
720 fshift[ki][m] += fij;
721 fshift[CENTRAL][m] -= fij;
727 real water_pol(int nbonds,
728 const t_iatom forceatoms[], const t_iparams forceparams[],
729 const rvec x[], rvec f[], rvec fshift[],
730 const t_pbc *pbc, const t_graph *g,
731 real lambda, real *dvdlambda,
732 const t_mdatoms *md, t_fcdata *fcd,
733 int *global_atom_index)
735 /* This routine implements anisotropic polarizibility for water, through
736 * a shell connected to a dummy with spring constant that differ in the
737 * three spatial dimensions in the molecular frame.
739 int i, m, aO, aH1, aH2, aD, aS, type, type0;
740 rvec dOH1, dOH2, dHH, dOD, dDS, nW, kk, dx, kdx, proj;
744 real vtot, fij, r_HH, r_OD, r_nW, tx, ty, tz, qS;
749 type0 = forceatoms[0];
751 qS = md->chargeA[aS];
752 kk[XX] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_x;
753 kk[YY] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_y;
754 kk[ZZ] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_z;
755 r_HH = 1.0/forceparams[type0].wpol.rHH;
756 r_OD = 1.0/forceparams[type0].wpol.rOD;
759 fprintf(debug, "WPOL: qS = %10.5f aS = %5d\n", qS, aS);
760 fprintf(debug, "WPOL: kk = %10.3f %10.3f %10.3f\n",
761 kk[XX], kk[YY], kk[ZZ]);
762 fprintf(debug, "WPOL: rOH = %10.3f rHH = %10.3f rOD = %10.3f\n",
763 forceparams[type0].wpol.rOH,
764 forceparams[type0].wpol.rHH,
765 forceparams[type0].wpol.rOD);
767 for (i = 0; (i < nbonds); i += 6)
769 type = forceatoms[i];
772 gmx_fatal(FARGS, "Sorry, type = %d, type0 = %d, file = %s, line = %d",
773 type, type0, __FILE__, __LINE__);
775 aO = forceatoms[i+1];
776 aH1 = forceatoms[i+2];
777 aH2 = forceatoms[i+3];
778 aD = forceatoms[i+4];
779 aS = forceatoms[i+5];
781 /* Compute vectors describing the water frame */
782 rvec_sub(x[aH1], x[aO], dOH1);
783 rvec_sub(x[aH2], x[aO], dOH2);
784 rvec_sub(x[aH2], x[aH1], dHH);
785 rvec_sub(x[aD], x[aO], dOD);
786 rvec_sub(x[aS], x[aD], dDS);
787 cprod(dOH1, dOH2, nW);
789 /* Compute inverse length of normal vector
790 * (this one could be precomputed, but I'm too lazy now)
792 r_nW = gmx_invsqrt(iprod(nW, nW));
793 /* This is for precision, but does not make a big difference,
796 r_OD = gmx_invsqrt(iprod(dOD, dOD));
798 /* Normalize the vectors in the water frame */
800 svmul(r_HH, dHH, dHH);
801 svmul(r_OD, dOD, dOD);
803 /* Compute displacement of shell along components of the vector */
804 dx[ZZ] = iprod(dDS, dOD);
805 /* Compute projection on the XY plane: dDS - dx[ZZ]*dOD */
806 for (m = 0; (m < DIM); m++)
808 proj[m] = dDS[m]-dx[ZZ]*dOD[m];
811 /*dx[XX] = iprod(dDS,nW);
812 dx[YY] = iprod(dDS,dHH);*/
813 dx[XX] = iprod(proj, nW);
814 for (m = 0; (m < DIM); m++)
816 proj[m] -= dx[XX]*nW[m];
818 dx[YY] = iprod(proj, dHH);
823 fprintf(debug, "WPOL: dx2=%10g dy2=%10g dz2=%10g sum=%10g dDS^2=%10g\n",
824 sqr(dx[XX]), sqr(dx[YY]), sqr(dx[ZZ]), iprod(dx, dx), iprod(dDS, dDS));
825 fprintf(debug, "WPOL: dHH=(%10g,%10g,%10g)\n", dHH[XX], dHH[YY], dHH[ZZ]);
826 fprintf(debug, "WPOL: dOD=(%10g,%10g,%10g), 1/r_OD = %10g\n",
827 dOD[XX], dOD[YY], dOD[ZZ], 1/r_OD);
828 fprintf(debug, "WPOL: nW =(%10g,%10g,%10g), 1/r_nW = %10g\n",
829 nW[XX], nW[YY], nW[ZZ], 1/r_nW);
830 fprintf(debug, "WPOL: dx =%10g, dy =%10g, dz =%10g\n",
831 dx[XX], dx[YY], dx[ZZ]);
832 fprintf(debug, "WPOL: dDSx=%10g, dDSy=%10g, dDSz=%10g\n",
833 dDS[XX], dDS[YY], dDS[ZZ]);
836 /* Now compute the forces and energy */
837 kdx[XX] = kk[XX]*dx[XX];
838 kdx[YY] = kk[YY]*dx[YY];
839 kdx[ZZ] = kk[ZZ]*dx[ZZ];
840 vtot += iprod(dx, kdx);
841 for (m = 0; (m < DIM); m++)
843 /* This is a tensor operation but written out for speed */
857 fprintf(debug, "WPOL: vwpol=%g\n", 0.5*iprod(dx, kdx));
858 fprintf(debug, "WPOL: df = (%10g, %10g, %10g)\n", df[XX], df[YY], df[ZZ]);
866 static real do_1_thole(const rvec xi, const rvec xj, rvec fi, rvec fj,
867 const t_pbc *pbc, real qq,
868 rvec fshift[], real afac)
871 real r12sq, r12_1, r12n, r12bar, v0, v1, fscal, ebar, fff;
874 t = pbc_rvec_sub(pbc, xi, xj, r12); /* 3 */
876 r12sq = iprod(r12, r12); /* 5 */
877 r12_1 = gmx_invsqrt(r12sq); /* 5 */
878 r12bar = afac/r12_1; /* 5 */
879 v0 = qq*ONE_4PI_EPS0*r12_1; /* 2 */
880 ebar = exp(-r12bar); /* 5 */
881 v1 = (1-(1+0.5*r12bar)*ebar); /* 4 */
882 fscal = ((v0*r12_1)*v1 - v0*0.5*afac*ebar*(r12bar+1))*r12_1; /* 9 */
885 fprintf(debug, "THOLE: v0 = %.3f v1 = %.3f r12= % .3f r12bar = %.3f fscal = %.3f ebar = %.3f\n", v0, v1, 1/r12_1, r12bar, fscal, ebar);
888 for (m = 0; (m < DIM); m++)
894 fshift[CENTRAL][m] -= fff;
897 return v0*v1; /* 1 */
901 real thole_pol(int nbonds,
902 const t_iatom forceatoms[], const t_iparams forceparams[],
903 const rvec x[], rvec f[], rvec fshift[],
904 const t_pbc *pbc, const t_graph *g,
905 real lambda, real *dvdlambda,
906 const t_mdatoms *md, t_fcdata *fcd,
907 int *global_atom_index)
909 /* Interaction between two pairs of particles with opposite charge */
910 int i, type, a1, da1, a2, da2;
911 real q1, q2, qq, a, al1, al2, afac;
914 for (i = 0; (i < nbonds); )
916 type = forceatoms[i++];
917 a1 = forceatoms[i++];
918 da1 = forceatoms[i++];
919 a2 = forceatoms[i++];
920 da2 = forceatoms[i++];
921 q1 = md->chargeA[da1];
922 q2 = md->chargeA[da2];
923 a = forceparams[type].thole.a;
924 al1 = forceparams[type].thole.alpha1;
925 al2 = forceparams[type].thole.alpha2;
927 afac = a*pow(al1*al2, -1.0/6.0);
928 V += do_1_thole(x[a1], x[a2], f[a1], f[a2], pbc, qq, fshift, afac);
929 V += do_1_thole(x[da1], x[a2], f[da1], f[a2], pbc, -qq, fshift, afac);
930 V += do_1_thole(x[a1], x[da2], f[a1], f[da2], pbc, -qq, fshift, afac);
931 V += do_1_thole(x[da1], x[da2], f[da1], f[da2], pbc, qq, fshift, afac);
937 real bond_angle(const rvec xi, const rvec xj, const rvec xk, const t_pbc *pbc,
938 rvec r_ij, rvec r_kj, real *costh,
940 /* Return value is the angle between the bonds i-j and j-k */
945 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
946 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
948 *costh = cos_angle(r_ij, r_kj); /* 25 */
949 th = acos(*costh); /* 10 */
954 real angles(int nbonds,
955 const t_iatom forceatoms[], const t_iparams forceparams[],
956 const rvec x[], rvec f[], rvec fshift[],
957 const t_pbc *pbc, const t_graph *g,
958 real lambda, real *dvdlambda,
959 const t_mdatoms *md, t_fcdata *fcd,
960 int *global_atom_index)
962 int i, ai, aj, ak, t1, t2, type;
964 real cos_theta, cos_theta2, theta, dVdt, va, vtot;
965 ivec jt, dt_ij, dt_kj;
968 for (i = 0; i < nbonds; )
970 type = forceatoms[i++];
971 ai = forceatoms[i++];
972 aj = forceatoms[i++];
973 ak = forceatoms[i++];
975 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
976 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
978 *dvdlambda += harmonic(forceparams[type].harmonic.krA,
979 forceparams[type].harmonic.krB,
980 forceparams[type].harmonic.rA*DEG2RAD,
981 forceparams[type].harmonic.rB*DEG2RAD,
982 theta, lambda, &va, &dVdt); /* 21 */
985 cos_theta2 = sqr(cos_theta);
995 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
996 sth = st*cos_theta; /* 1 */
1000 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1001 theta*RAD2DEG, va, dVdt);
1004 nrij2 = iprod(r_ij, r_ij); /* 5 */
1005 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1007 nrij_1 = gmx_invsqrt(nrij2); /* 10 */
1008 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1010 cik = st*nrij_1*nrkj_1; /* 2 */
1011 cii = sth*nrij_1*nrij_1; /* 2 */
1012 ckk = sth*nrkj_1*nrkj_1; /* 2 */
1014 for (m = 0; m < DIM; m++)
1016 f_i[m] = -(cik*r_kj[m] - cii*r_ij[m]);
1017 f_k[m] = -(cik*r_ij[m] - ckk*r_kj[m]);
1018 f_j[m] = -f_i[m] - f_k[m];
1025 copy_ivec(SHIFT_IVEC(g, aj), jt);
1027 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1028 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1029 t1 = IVEC2IS(dt_ij);
1030 t2 = IVEC2IS(dt_kj);
1032 rvec_inc(fshift[t1], f_i);
1033 rvec_inc(fshift[CENTRAL], f_j);
1034 rvec_inc(fshift[t2], f_k);
1043 /* As angles, but using SIMD to calculate many dihedrals at once.
1044 * This routines does not calculate energies and shift forces.
1046 static gmx_inline void
1047 angles_noener_simd(int nbonds,
1048 const t_iatom forceatoms[], const t_iparams forceparams[],
1049 const rvec x[], rvec f[],
1050 const t_pbc *pbc, const t_graph *g,
1052 const t_mdatoms *md, t_fcdata *fcd,
1053 int *global_atom_index)
1055 #define UNROLL GMX_SIMD_WIDTH_HERE
1058 int type, ai[UNROLL], aj[UNROLL], ak[UNROLL];
1059 real coeff_array[2*UNROLL+UNROLL], *coeff;
1060 real dr_array[2*DIM*UNROLL+UNROLL], *dr;
1061 real f_buf_array[6*UNROLL+UNROLL], *f_buf;
1062 gmx_mm_pr k_S, theta0_S;
1063 gmx_mm_pr rijx_S, rijy_S, rijz_S;
1064 gmx_mm_pr rkjx_S, rkjy_S, rkjz_S;
1066 gmx_mm_pr min_one_plus_eps_S;
1067 gmx_mm_pr rij_rkj_S;
1068 gmx_mm_pr nrij2_S, nrij_1_S;
1069 gmx_mm_pr nrkj2_S, nrkj_1_S;
1070 gmx_mm_pr cos_S, invsin_S;
1072 gmx_mm_pr st_S, sth_S;
1073 gmx_mm_pr cik_S, cii_S, ckk_S;
1074 gmx_mm_pr f_ix_S, f_iy_S, f_iz_S;
1075 gmx_mm_pr f_kx_S, f_ky_S, f_kz_S;
1076 pbc_simd_t pbc_simd;
1078 /* Ensure register memory alignment */
1079 coeff = gmx_simd_align_real(coeff_array);
1080 dr = gmx_simd_align_real(dr_array);
1081 f_buf = gmx_simd_align_real(f_buf_array);
1083 set_pbc_simd(pbc,&pbc_simd);
1085 one_S = gmx_set1_pr(1.0);
1087 /* The smallest number > -1 */
1088 min_one_plus_eps_S = gmx_set1_pr(-1.0 + 2*GMX_REAL_EPS);
1090 /* nbonds is the number of angles times nfa1, here we step UNROLL angles */
1091 for (i = 0; (i < nbonds); i += UNROLL*nfa1)
1093 /* Collect atoms for UNROLL angles.
1094 * iu indexes into forceatoms, we should not let iu go beyond nbonds.
1097 for (s = 0; s < UNROLL; s++)
1099 type = forceatoms[iu];
1100 ai[s] = forceatoms[iu+1];
1101 aj[s] = forceatoms[iu+2];
1102 ak[s] = forceatoms[iu+3];
1104 coeff[s] = forceparams[type].harmonic.krA;
1105 coeff[UNROLL+s] = forceparams[type].harmonic.rA*DEG2RAD;
1107 /* If you can't use pbc_dx_simd below for PBC, e.g. because
1108 * you can't round in SIMD, use pbc_rvec_sub here.
1110 /* Store the non PBC corrected distances packed and aligned */
1111 for (m = 0; m < DIM; m++)
1113 dr[s + m *UNROLL] = x[ai[s]][m] - x[aj[s]][m];
1114 dr[s + (DIM+m)*UNROLL] = x[ak[s]][m] - x[aj[s]][m];
1117 /* At the end fill the arrays with identical entries */
1118 if (iu + nfa1 < nbonds)
1124 k_S = gmx_load_pr(coeff);
1125 theta0_S = gmx_load_pr(coeff+UNROLL);
1127 rijx_S = gmx_load_pr(dr + 0*UNROLL);
1128 rijy_S = gmx_load_pr(dr + 1*UNROLL);
1129 rijz_S = gmx_load_pr(dr + 2*UNROLL);
1130 rkjx_S = gmx_load_pr(dr + 3*UNROLL);
1131 rkjy_S = gmx_load_pr(dr + 4*UNROLL);
1132 rkjz_S = gmx_load_pr(dr + 5*UNROLL);
1134 pbc_dx_simd(&rijx_S, &rijy_S, &rijz_S, &pbc_simd);
1135 pbc_dx_simd(&rkjx_S, &rkjy_S, &rkjz_S, &pbc_simd);
1137 rij_rkj_S = gmx_iprod_pr(rijx_S, rijy_S, rijz_S,
1138 rkjx_S, rkjy_S, rkjz_S);
1140 nrij2_S = gmx_norm2_pr(rijx_S, rijy_S, rijz_S);
1141 nrkj2_S = gmx_norm2_pr(rkjx_S, rkjy_S, rkjz_S);
1143 nrij_1_S = gmx_invsqrt_pr(nrij2_S);
1144 nrkj_1_S = gmx_invsqrt_pr(nrkj2_S);
1146 cos_S = gmx_mul_pr(rij_rkj_S, gmx_mul_pr(nrij_1_S, nrkj_1_S));
1148 /* To allow for 180 degrees, we take the max of cos and -1 + 1bit,
1149 * so we can safely get the 1/sin from 1/sqrt(1 - cos^2).
1150 * This also ensures that rounding errors would cause the argument
1151 * of gmx_acos_pr to be < -1.
1152 * Note that we do not take precautions for cos(0)=1, so the outer
1153 * atoms in an angle should not be on top of each other.
1155 cos_S = gmx_max_pr(cos_S, min_one_plus_eps_S);
1157 theta_S = gmx_acos_pr(cos_S);
1159 invsin_S = gmx_invsqrt_pr(gmx_sub_pr(one_S, gmx_mul_pr(cos_S, cos_S)));
1161 st_S = gmx_mul_pr(gmx_mul_pr(k_S, gmx_sub_pr(theta0_S, theta_S)),
1163 sth_S = gmx_mul_pr(st_S, cos_S);
1165 cik_S = gmx_mul_pr(st_S, gmx_mul_pr(nrij_1_S, nrkj_1_S));
1166 cii_S = gmx_mul_pr(sth_S, gmx_mul_pr(nrij_1_S, nrij_1_S));
1167 ckk_S = gmx_mul_pr(sth_S, gmx_mul_pr(nrkj_1_S, nrkj_1_S));
1169 f_ix_S = gmx_mul_pr(cii_S, rijx_S);
1170 f_ix_S = gmx_nmsub_pr(cik_S, rkjx_S, f_ix_S);
1171 f_iy_S = gmx_mul_pr(cii_S, rijy_S);
1172 f_iy_S = gmx_nmsub_pr(cik_S, rkjy_S, f_iy_S);
1173 f_iz_S = gmx_mul_pr(cii_S, rijz_S);
1174 f_iz_S = gmx_nmsub_pr(cik_S, rkjz_S, f_iz_S);
1175 f_kx_S = gmx_mul_pr(ckk_S, rkjx_S);
1176 f_kx_S = gmx_nmsub_pr(cik_S, rijx_S, f_kx_S);
1177 f_ky_S = gmx_mul_pr(ckk_S, rkjy_S);
1178 f_ky_S = gmx_nmsub_pr(cik_S, rijy_S, f_ky_S);
1179 f_kz_S = gmx_mul_pr(ckk_S, rkjz_S);
1180 f_kz_S = gmx_nmsub_pr(cik_S, rijz_S, f_kz_S);
1182 gmx_store_pr(f_buf + 0*UNROLL, f_ix_S);
1183 gmx_store_pr(f_buf + 1*UNROLL, f_iy_S);
1184 gmx_store_pr(f_buf + 2*UNROLL, f_iz_S);
1185 gmx_store_pr(f_buf + 3*UNROLL, f_kx_S);
1186 gmx_store_pr(f_buf + 4*UNROLL, f_ky_S);
1187 gmx_store_pr(f_buf + 5*UNROLL, f_kz_S);
1193 for (m = 0; m < DIM; m++)
1195 f[ai[s]][m] += f_buf[s + m*UNROLL];
1196 f[aj[s]][m] -= f_buf[s + m*UNROLL] + f_buf[s + (DIM+m)*UNROLL];
1197 f[ak[s]][m] += f_buf[s + (DIM+m)*UNROLL];
1202 while (s < UNROLL && iu < nbonds);
1207 #endif /* SIMD_BONDEDS */
1209 real linear_angles(int nbonds,
1210 const t_iatom forceatoms[], const t_iparams forceparams[],
1211 const rvec x[], rvec f[], rvec fshift[],
1212 const t_pbc *pbc, const t_graph *g,
1213 real lambda, real *dvdlambda,
1214 const t_mdatoms *md, t_fcdata *fcd,
1215 int *global_atom_index)
1217 int i, m, ai, aj, ak, t1, t2, type;
1219 real L1, kA, kB, aA, aB, dr, dr2, va, vtot, a, b, klin;
1220 ivec jt, dt_ij, dt_kj;
1221 rvec r_ij, r_kj, r_ik, dx;
1225 for (i = 0; (i < nbonds); )
1227 type = forceatoms[i++];
1228 ai = forceatoms[i++];
1229 aj = forceatoms[i++];
1230 ak = forceatoms[i++];
1232 kA = forceparams[type].linangle.klinA;
1233 kB = forceparams[type].linangle.klinB;
1234 klin = L1*kA + lambda*kB;
1236 aA = forceparams[type].linangle.aA;
1237 aB = forceparams[type].linangle.aB;
1238 a = L1*aA+lambda*aB;
1241 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
1242 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
1243 rvec_sub(r_ij, r_kj, r_ik);
1246 for (m = 0; (m < DIM); m++)
1248 dr = -a * r_ij[m] - b * r_kj[m];
1253 f_j[m] = -(f_i[m]+f_k[m]);
1259 *dvdlambda += 0.5*(kB-kA)*dr2 + klin*(aB-aA)*iprod(dx, r_ik);
1265 copy_ivec(SHIFT_IVEC(g, aj), jt);
1267 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1268 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1269 t1 = IVEC2IS(dt_ij);
1270 t2 = IVEC2IS(dt_kj);
1272 rvec_inc(fshift[t1], f_i);
1273 rvec_inc(fshift[CENTRAL], f_j);
1274 rvec_inc(fshift[t2], f_k);
1279 real urey_bradley(int nbonds,
1280 const t_iatom forceatoms[], const t_iparams forceparams[],
1281 const rvec x[], rvec f[], rvec fshift[],
1282 const t_pbc *pbc, const t_graph *g,
1283 real lambda, real *dvdlambda,
1284 const t_mdatoms *md, t_fcdata *fcd,
1285 int *global_atom_index)
1287 int i, m, ai, aj, ak, t1, t2, type, ki;
1288 rvec r_ij, r_kj, r_ik;
1289 real cos_theta, cos_theta2, theta;
1290 real dVdt, va, vtot, dr, dr2, vbond, fbond, fik;
1291 real kthA, th0A, kUBA, r13A, kthB, th0B, kUBB, r13B;
1292 ivec jt, dt_ij, dt_kj, dt_ik;
1295 for (i = 0; (i < nbonds); )
1297 type = forceatoms[i++];
1298 ai = forceatoms[i++];
1299 aj = forceatoms[i++];
1300 ak = forceatoms[i++];
1301 th0A = forceparams[type].u_b.thetaA*DEG2RAD;
1302 kthA = forceparams[type].u_b.kthetaA;
1303 r13A = forceparams[type].u_b.r13A;
1304 kUBA = forceparams[type].u_b.kUBA;
1305 th0B = forceparams[type].u_b.thetaB*DEG2RAD;
1306 kthB = forceparams[type].u_b.kthetaB;
1307 r13B = forceparams[type].u_b.r13B;
1308 kUBB = forceparams[type].u_b.kUBB;
1310 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
1311 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
1313 *dvdlambda += harmonic(kthA, kthB, th0A, th0B, theta, lambda, &va, &dVdt); /* 21 */
1316 ki = pbc_rvec_sub(pbc, x[ai], x[ak], r_ik); /* 3 */
1317 dr2 = iprod(r_ik, r_ik); /* 5 */
1318 dr = dr2*gmx_invsqrt(dr2); /* 10 */
1320 *dvdlambda += harmonic(kUBA, kUBB, r13A, r13B, dr, lambda, &vbond, &fbond); /* 19 */
1322 cos_theta2 = sqr(cos_theta); /* 1 */
1330 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1331 sth = st*cos_theta; /* 1 */
1335 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1336 theta*RAD2DEG, va, dVdt);
1339 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1340 nrij2 = iprod(r_ij, r_ij);
1342 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
1343 cii = sth/nrij2; /* 10 */
1344 ckk = sth/nrkj2; /* 10 */
1346 for (m = 0; (m < DIM); m++) /* 39 */
1348 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
1349 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
1350 f_j[m] = -f_i[m]-f_k[m];
1357 copy_ivec(SHIFT_IVEC(g, aj), jt);
1359 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1360 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1361 t1 = IVEC2IS(dt_ij);
1362 t2 = IVEC2IS(dt_kj);
1364 rvec_inc(fshift[t1], f_i);
1365 rvec_inc(fshift[CENTRAL], f_j);
1366 rvec_inc(fshift[t2], f_k);
1368 /* Time for the bond calculations */
1374 vtot += vbond; /* 1*/
1375 fbond *= gmx_invsqrt(dr2); /* 6 */
1379 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ak), dt_ik);
1380 ki = IVEC2IS(dt_ik);
1382 for (m = 0; (m < DIM); m++) /* 15 */
1384 fik = fbond*r_ik[m];
1387 fshift[ki][m] += fik;
1388 fshift[CENTRAL][m] -= fik;
1394 real quartic_angles(int nbonds,
1395 const t_iatom forceatoms[], const t_iparams forceparams[],
1396 const rvec x[], rvec f[], rvec fshift[],
1397 const t_pbc *pbc, const t_graph *g,
1398 real lambda, real *dvdlambda,
1399 const t_mdatoms *md, t_fcdata *fcd,
1400 int *global_atom_index)
1402 int i, j, ai, aj, ak, t1, t2, type;
1404 real cos_theta, cos_theta2, theta, dt, dVdt, va, dtp, c, vtot;
1405 ivec jt, dt_ij, dt_kj;
1408 for (i = 0; (i < nbonds); )
1410 type = forceatoms[i++];
1411 ai = forceatoms[i++];
1412 aj = forceatoms[i++];
1413 ak = forceatoms[i++];
1415 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
1416 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
1418 dt = theta - forceparams[type].qangle.theta*DEG2RAD; /* 2 */
1421 va = forceparams[type].qangle.c[0];
1423 for (j = 1; j <= 4; j++)
1425 c = forceparams[type].qangle.c[j];
1434 cos_theta2 = sqr(cos_theta); /* 1 */
1443 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1444 sth = st*cos_theta; /* 1 */
1448 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1449 theta*RAD2DEG, va, dVdt);
1452 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1453 nrij2 = iprod(r_ij, r_ij);
1455 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
1456 cii = sth/nrij2; /* 10 */
1457 ckk = sth/nrkj2; /* 10 */
1459 for (m = 0; (m < DIM); m++) /* 39 */
1461 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
1462 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
1463 f_j[m] = -f_i[m]-f_k[m];
1470 copy_ivec(SHIFT_IVEC(g, aj), jt);
1472 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1473 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1474 t1 = IVEC2IS(dt_ij);
1475 t2 = IVEC2IS(dt_kj);
1477 rvec_inc(fshift[t1], f_i);
1478 rvec_inc(fshift[CENTRAL], f_j);
1479 rvec_inc(fshift[t2], f_k);
1485 real dih_angle(const rvec xi, const rvec xj, const rvec xk, const rvec xl,
1487 rvec r_ij, rvec r_kj, rvec r_kl, rvec m, rvec n,
1488 real *sign, int *t1, int *t2, int *t3)
1492 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
1493 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
1494 *t3 = pbc_rvec_sub(pbc, xk, xl, r_kl); /* 3 */
1496 cprod(r_ij, r_kj, m); /* 9 */
1497 cprod(r_kj, r_kl, n); /* 9 */
1498 phi = gmx_angle(m, n); /* 49 (assuming 25 for atan2) */
1499 ipr = iprod(r_ij, n); /* 5 */
1500 (*sign) = (ipr < 0.0) ? -1.0 : 1.0;
1501 phi = (*sign)*phi; /* 1 */
1509 /* As dih_angle above, but calculates 4 dihedral angles at once using SIMD,
1510 * also calculates the pre-factor required for the dihedral force update.
1511 * Note that bv and buf should be register aligned.
1513 static gmx_inline void
1514 dih_angle_simd(const rvec *x,
1515 const int *ai, const int *aj, const int *ak, const int *al,
1516 const pbc_simd_t *pbc,
1519 gmx_mm_pr *mx_S, gmx_mm_pr *my_S, gmx_mm_pr *mz_S,
1520 gmx_mm_pr *nx_S, gmx_mm_pr *ny_S, gmx_mm_pr *nz_S,
1521 gmx_mm_pr *nrkj_m2_S,
1522 gmx_mm_pr *nrkj_n2_S,
1526 #define UNROLL GMX_SIMD_WIDTH_HERE
1528 gmx_mm_pr rijx_S, rijy_S, rijz_S;
1529 gmx_mm_pr rkjx_S, rkjy_S, rkjz_S;
1530 gmx_mm_pr rklx_S, rkly_S, rklz_S;
1531 gmx_mm_pr cx_S, cy_S, cz_S;
1535 gmx_mm_pr iprm_S, iprn_S;
1536 gmx_mm_pr nrkj2_S, nrkj_1_S, nrkj_2_S, nrkj_S;
1539 gmx_mm_pr nrkj2_min_S;
1540 gmx_mm_pr real_eps_S;
1542 /* Used to avoid division by zero.
1543 * We take into acount that we multiply the result by real_eps_S.
1545 nrkj2_min_S = gmx_set1_pr(GMX_REAL_MIN/(2*GMX_REAL_EPS));
1547 /* The value of the last significant bit (GMX_REAL_EPS is half of that) */
1548 real_eps_S = gmx_set1_pr(2*GMX_REAL_EPS);
1550 for (s = 0; s < UNROLL; s++)
1552 /* If you can't use pbc_dx_simd below for PBC, e.g. because
1553 * you can't round in SIMD, use pbc_rvec_sub here.
1555 for (m = 0; m < DIM; m++)
1557 dr[s + (0*DIM + m)*UNROLL] = x[ai[s]][m] - x[aj[s]][m];
1558 dr[s + (1*DIM + m)*UNROLL] = x[ak[s]][m] - x[aj[s]][m];
1559 dr[s + (2*DIM + m)*UNROLL] = x[ak[s]][m] - x[al[s]][m];
1563 rijx_S = gmx_load_pr(dr + 0*UNROLL);
1564 rijy_S = gmx_load_pr(dr + 1*UNROLL);
1565 rijz_S = gmx_load_pr(dr + 2*UNROLL);
1566 rkjx_S = gmx_load_pr(dr + 3*UNROLL);
1567 rkjy_S = gmx_load_pr(dr + 4*UNROLL);
1568 rkjz_S = gmx_load_pr(dr + 5*UNROLL);
1569 rklx_S = gmx_load_pr(dr + 6*UNROLL);
1570 rkly_S = gmx_load_pr(dr + 7*UNROLL);
1571 rklz_S = gmx_load_pr(dr + 8*UNROLL);
1573 pbc_dx_simd(&rijx_S, &rijy_S, &rijz_S, pbc);
1574 pbc_dx_simd(&rkjx_S, &rkjy_S, &rkjz_S, pbc);
1575 pbc_dx_simd(&rklx_S, &rkly_S, &rklz_S, pbc);
1577 gmx_cprod_pr(rijx_S, rijy_S, rijz_S,
1578 rkjx_S, rkjy_S, rkjz_S,
1581 gmx_cprod_pr(rkjx_S, rkjy_S, rkjz_S,
1582 rklx_S, rkly_S, rklz_S,
1585 gmx_cprod_pr(*mx_S, *my_S, *mz_S,
1586 *nx_S, *ny_S, *nz_S,
1587 &cx_S, &cy_S, &cz_S);
1589 cn_S = gmx_sqrt_pr(gmx_norm2_pr(cx_S, cy_S, cz_S));
1591 s_S = gmx_iprod_pr(*mx_S, *my_S, *mz_S, *nx_S, *ny_S, *nz_S);
1593 /* Determine the dihedral angle, the sign might need correction */
1594 *phi_S = gmx_atan2_pr(cn_S, s_S);
1596 ipr_S = gmx_iprod_pr(rijx_S, rijy_S, rijz_S,
1597 *nx_S, *ny_S, *nz_S);
1599 iprm_S = gmx_norm2_pr(*mx_S, *my_S, *mz_S);
1600 iprn_S = gmx_norm2_pr(*nx_S, *ny_S, *nz_S);
1602 nrkj2_S = gmx_norm2_pr(rkjx_S, rkjy_S, rkjz_S);
1604 /* Avoid division by zero. When zero, the result is multiplied by 0
1605 * anyhow, so the 3 max below do not affect the final result.
1607 nrkj2_S = gmx_max_pr(nrkj2_S, nrkj2_min_S);
1608 nrkj_1_S = gmx_invsqrt_pr(nrkj2_S);
1609 nrkj_2_S = gmx_mul_pr(nrkj_1_S, nrkj_1_S);
1610 nrkj_S = gmx_mul_pr(nrkj2_S, nrkj_1_S);
1612 toler_S = gmx_mul_pr(nrkj2_S, real_eps_S);
1614 /* Here the plain-C code uses a conditional, but we can't do that in SIMD.
1615 * So we take a max with the tolerance instead. Since we multiply with
1616 * m or n later, the max does not affect the results.
1618 iprm_S = gmx_max_pr(iprm_S, toler_S);
1619 iprn_S = gmx_max_pr(iprn_S, toler_S);
1620 *nrkj_m2_S = gmx_mul_pr(nrkj_S, gmx_inv_pr(iprm_S));
1621 *nrkj_n2_S = gmx_mul_pr(nrkj_S, gmx_inv_pr(iprn_S));
1623 /* Set sign of phi_S with the sign of ipr_S; phi_S is currently positive */
1624 *phi_S = gmx_cpsgn_nonneg_pr(ipr_S, *phi_S);
1626 p_S = gmx_iprod_pr(rijx_S, rijy_S, rijz_S,
1627 rkjx_S, rkjy_S, rkjz_S);
1628 p_S = gmx_mul_pr(p_S, nrkj_2_S);
1630 q_S = gmx_iprod_pr(rklx_S, rkly_S, rklz_S,
1631 rkjx_S, rkjy_S, rkjz_S);
1632 q_S = gmx_mul_pr(q_S, nrkj_2_S);
1634 gmx_store_pr(p, p_S);
1635 gmx_store_pr(q, q_S);
1639 #endif /* SIMD_BONDEDS */
1642 void do_dih_fup(int i, int j, int k, int l, real ddphi,
1643 rvec r_ij, rvec r_kj, rvec r_kl,
1644 rvec m, rvec n, rvec f[], rvec fshift[],
1645 const t_pbc *pbc, const t_graph *g,
1646 const rvec x[], int t1, int t2, int t3)
1649 rvec f_i, f_j, f_k, f_l;
1650 rvec uvec, vvec, svec, dx_jl;
1651 real iprm, iprn, nrkj, nrkj2, nrkj_1, nrkj_2;
1652 real a, b, p, q, toler;
1653 ivec jt, dt_ij, dt_kj, dt_lj;
1655 iprm = iprod(m, m); /* 5 */
1656 iprn = iprod(n, n); /* 5 */
1657 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1658 toler = nrkj2*GMX_REAL_EPS;
1659 if ((iprm > toler) && (iprn > toler))
1661 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1662 nrkj_2 = nrkj_1*nrkj_1; /* 1 */
1663 nrkj = nrkj2*nrkj_1; /* 1 */
1664 a = -ddphi*nrkj/iprm; /* 11 */
1665 svmul(a, m, f_i); /* 3 */
1666 b = ddphi*nrkj/iprn; /* 11 */
1667 svmul(b, n, f_l); /* 3 */
1668 p = iprod(r_ij, r_kj); /* 5 */
1669 p *= nrkj_2; /* 1 */
1670 q = iprod(r_kl, r_kj); /* 5 */
1671 q *= nrkj_2; /* 1 */
1672 svmul(p, f_i, uvec); /* 3 */
1673 svmul(q, f_l, vvec); /* 3 */
1674 rvec_sub(uvec, vvec, svec); /* 3 */
1675 rvec_sub(f_i, svec, f_j); /* 3 */
1676 rvec_add(f_l, svec, f_k); /* 3 */
1677 rvec_inc(f[i], f_i); /* 3 */
1678 rvec_dec(f[j], f_j); /* 3 */
1679 rvec_dec(f[k], f_k); /* 3 */
1680 rvec_inc(f[l], f_l); /* 3 */
1684 copy_ivec(SHIFT_IVEC(g, j), jt);
1685 ivec_sub(SHIFT_IVEC(g, i), jt, dt_ij);
1686 ivec_sub(SHIFT_IVEC(g, k), jt, dt_kj);
1687 ivec_sub(SHIFT_IVEC(g, l), jt, dt_lj);
1688 t1 = IVEC2IS(dt_ij);
1689 t2 = IVEC2IS(dt_kj);
1690 t3 = IVEC2IS(dt_lj);
1694 t3 = pbc_rvec_sub(pbc, x[l], x[j], dx_jl);
1701 rvec_inc(fshift[t1], f_i);
1702 rvec_dec(fshift[CENTRAL], f_j);
1703 rvec_dec(fshift[t2], f_k);
1704 rvec_inc(fshift[t3], f_l);
1709 /* As do_dih_fup above, but without shift forces */
1711 do_dih_fup_noshiftf(int i, int j, int k, int l, real ddphi,
1712 rvec r_ij, rvec r_kj, rvec r_kl,
1713 rvec m, rvec n, rvec f[])
1715 rvec f_i, f_j, f_k, f_l;
1716 rvec uvec, vvec, svec, dx_jl;
1717 real iprm, iprn, nrkj, nrkj2, nrkj_1, nrkj_2;
1718 real a, b, p, q, toler;
1719 ivec jt, dt_ij, dt_kj, dt_lj;
1721 iprm = iprod(m, m); /* 5 */
1722 iprn = iprod(n, n); /* 5 */
1723 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1724 toler = nrkj2*GMX_REAL_EPS;
1725 if ((iprm > toler) && (iprn > toler))
1727 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1728 nrkj_2 = nrkj_1*nrkj_1; /* 1 */
1729 nrkj = nrkj2*nrkj_1; /* 1 */
1730 a = -ddphi*nrkj/iprm; /* 11 */
1731 svmul(a, m, f_i); /* 3 */
1732 b = ddphi*nrkj/iprn; /* 11 */
1733 svmul(b, n, f_l); /* 3 */
1734 p = iprod(r_ij, r_kj); /* 5 */
1735 p *= nrkj_2; /* 1 */
1736 q = iprod(r_kl, r_kj); /* 5 */
1737 q *= nrkj_2; /* 1 */
1738 svmul(p, f_i, uvec); /* 3 */
1739 svmul(q, f_l, vvec); /* 3 */
1740 rvec_sub(uvec, vvec, svec); /* 3 */
1741 rvec_sub(f_i, svec, f_j); /* 3 */
1742 rvec_add(f_l, svec, f_k); /* 3 */
1743 rvec_inc(f[i], f_i); /* 3 */
1744 rvec_dec(f[j], f_j); /* 3 */
1745 rvec_dec(f[k], f_k); /* 3 */
1746 rvec_inc(f[l], f_l); /* 3 */
1750 /* As do_dih_fup_noshiftf above, but with pre-calculated pre-factors */
1751 static gmx_inline void
1752 do_dih_fup_noshiftf_precalc(int i, int j, int k, int l,
1754 real f_i_x, real f_i_y, real f_i_z,
1755 real mf_l_x, real mf_l_y, real mf_l_z,
1758 rvec f_i, f_j, f_k, f_l;
1759 rvec uvec, vvec, svec;
1767 svmul(p, f_i, uvec);
1768 svmul(q, f_l, vvec);
1769 rvec_sub(uvec, vvec, svec);
1770 rvec_sub(f_i, svec, f_j);
1771 rvec_add(f_l, svec, f_k);
1772 rvec_inc(f[i], f_i);
1773 rvec_dec(f[j], f_j);
1774 rvec_dec(f[k], f_k);
1775 rvec_inc(f[l], f_l);
1779 real dopdihs(real cpA, real cpB, real phiA, real phiB, int mult,
1780 real phi, real lambda, real *V, real *F)
1782 real v, dvdlambda, mdphi, v1, sdphi, ddphi;
1783 real L1 = 1.0 - lambda;
1784 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1785 real dph0 = (phiB - phiA)*DEG2RAD;
1786 real cp = L1*cpA + lambda*cpB;
1788 mdphi = mult*phi - ph0;
1790 ddphi = -cp*mult*sdphi;
1791 v1 = 1.0 + cos(mdphi);
1794 dvdlambda = (cpB - cpA)*v1 + cp*dph0*sdphi;
1801 /* That was 40 flops */
1805 dopdihs_noener(real cpA, real cpB, real phiA, real phiB, int mult,
1806 real phi, real lambda, real *F)
1808 real mdphi, sdphi, ddphi;
1809 real L1 = 1.0 - lambda;
1810 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1811 real cp = L1*cpA + lambda*cpB;
1813 mdphi = mult*phi - ph0;
1815 ddphi = -cp*mult*sdphi;
1819 /* That was 20 flops */
1823 dopdihs_mdphi(real cpA, real cpB, real phiA, real phiB, int mult,
1824 real phi, real lambda, real *cp, real *mdphi)
1826 real L1 = 1.0 - lambda;
1827 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1829 *cp = L1*cpA + lambda*cpB;
1831 *mdphi = mult*phi - ph0;
1834 static real dopdihs_min(real cpA, real cpB, real phiA, real phiB, int mult,
1835 real phi, real lambda, real *V, real *F)
1836 /* similar to dopdihs, except for a minus sign *
1837 * and a different treatment of mult/phi0 */
1839 real v, dvdlambda, mdphi, v1, sdphi, ddphi;
1840 real L1 = 1.0 - lambda;
1841 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1842 real dph0 = (phiB - phiA)*DEG2RAD;
1843 real cp = L1*cpA + lambda*cpB;
1845 mdphi = mult*(phi-ph0);
1847 ddphi = cp*mult*sdphi;
1848 v1 = 1.0-cos(mdphi);
1851 dvdlambda = (cpB-cpA)*v1 + cp*dph0*sdphi;
1858 /* That was 40 flops */
1861 real pdihs(int nbonds,
1862 const t_iatom forceatoms[], const t_iparams forceparams[],
1863 const rvec x[], rvec f[], rvec fshift[],
1864 const t_pbc *pbc, const t_graph *g,
1865 real lambda, real *dvdlambda,
1866 const t_mdatoms *md, t_fcdata *fcd,
1867 int *global_atom_index)
1869 int i, type, ai, aj, ak, al;
1871 rvec r_ij, r_kj, r_kl, m, n;
1872 real phi, sign, ddphi, vpd, vtot;
1876 for (i = 0; (i < nbonds); )
1878 type = forceatoms[i++];
1879 ai = forceatoms[i++];
1880 aj = forceatoms[i++];
1881 ak = forceatoms[i++];
1882 al = forceatoms[i++];
1884 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
1885 &sign, &t1, &t2, &t3); /* 84 */
1886 *dvdlambda += dopdihs(forceparams[type].pdihs.cpA,
1887 forceparams[type].pdihs.cpB,
1888 forceparams[type].pdihs.phiA,
1889 forceparams[type].pdihs.phiB,
1890 forceparams[type].pdihs.mult,
1891 phi, lambda, &vpd, &ddphi);
1894 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
1895 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
1898 fprintf(debug, "pdih: (%d,%d,%d,%d) phi=%g\n",
1899 ai, aj, ak, al, phi);
1906 void make_dp_periodic(real *dp) /* 1 flop? */
1908 /* dp cannot be outside (-pi,pi) */
1913 else if (*dp < -M_PI)
1920 /* As pdihs above, but without calculating energies and shift forces */
1922 pdihs_noener(int nbonds,
1923 const t_iatom forceatoms[], const t_iparams forceparams[],
1924 const rvec x[], rvec f[],
1925 const t_pbc *pbc, const t_graph *g,
1927 const t_mdatoms *md, t_fcdata *fcd,
1928 int *global_atom_index)
1930 int i, type, ai, aj, ak, al;
1932 rvec r_ij, r_kj, r_kl, m, n;
1933 real phi, sign, ddphi_tot, ddphi;
1935 for (i = 0; (i < nbonds); )
1937 ai = forceatoms[i+1];
1938 aj = forceatoms[i+2];
1939 ak = forceatoms[i+3];
1940 al = forceatoms[i+4];
1942 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
1943 &sign, &t1, &t2, &t3);
1947 /* Loop over dihedrals working on the same atoms,
1948 * so we avoid recalculating angles and force distributions.
1952 type = forceatoms[i];
1953 dopdihs_noener(forceparams[type].pdihs.cpA,
1954 forceparams[type].pdihs.cpB,
1955 forceparams[type].pdihs.phiA,
1956 forceparams[type].pdihs.phiB,
1957 forceparams[type].pdihs.mult,
1958 phi, lambda, &ddphi);
1963 while (i < nbonds &&
1964 forceatoms[i+1] == ai &&
1965 forceatoms[i+2] == aj &&
1966 forceatoms[i+3] == ak &&
1967 forceatoms[i+4] == al);
1969 do_dih_fup_noshiftf(ai, aj, ak, al, ddphi_tot, r_ij, r_kj, r_kl, m, n, f);
1976 /* As pdihs_noner above, but using SIMD to calculate many dihedrals at once */
1978 pdihs_noener_simd(int nbonds,
1979 const t_iatom forceatoms[], const t_iparams forceparams[],
1980 const rvec x[], rvec f[],
1981 const t_pbc *pbc, const t_graph *g,
1983 const t_mdatoms *md, t_fcdata *fcd,
1984 int *global_atom_index)
1986 #define UNROLL GMX_SIMD_WIDTH_HERE
1989 int type, ai[UNROLL], aj[UNROLL], ak[UNROLL], al[UNROLL];
1990 int t1[UNROLL], t2[UNROLL], t3[UNROLL];
1992 real dr_array[3*DIM*UNROLL+UNROLL], *dr;
1993 real buf_array[7*UNROLL+UNROLL], *buf;
1994 real *cp, *phi0, *mult, *phi, *p, *q, *sf_i, *msf_l;
1995 gmx_mm_pr phi0_S, phi_S;
1996 gmx_mm_pr mx_S, my_S, mz_S;
1997 gmx_mm_pr nx_S, ny_S, nz_S;
1998 gmx_mm_pr nrkj_m2_S, nrkj_n2_S;
1999 gmx_mm_pr cp_S, mdphi_S, mult_S;
2000 gmx_mm_pr sin_S, cos_S;
2002 gmx_mm_pr sf_i_S, msf_l_S;
2003 pbc_simd_t pbc_simd;
2005 /* Ensure SIMD register alignment */
2006 dr = gmx_simd_align_real(dr_array);
2007 buf = gmx_simd_align_real(buf_array);
2009 /* Extract aligned pointer for parameters and variables */
2010 cp = buf + 0*UNROLL;
2011 phi0 = buf + 1*UNROLL;
2012 mult = buf + 2*UNROLL;
2015 sf_i = buf + 5*UNROLL;
2016 msf_l = buf + 6*UNROLL;
2018 set_pbc_simd(pbc, &pbc_simd);
2020 /* nbonds is the number of dihedrals times nfa1, here we step UNROLL dihs */
2021 for (i = 0; (i < nbonds); i += UNROLL*nfa1)
2023 /* Collect atoms quadruplets for UNROLL dihedrals.
2024 * iu indexes into forceatoms, we should not let iu go beyond nbonds.
2027 for (s = 0; s < UNROLL; s++)
2029 type = forceatoms[iu];
2030 ai[s] = forceatoms[iu+1];
2031 aj[s] = forceatoms[iu+2];
2032 ak[s] = forceatoms[iu+3];
2033 al[s] = forceatoms[iu+4];
2035 cp[s] = forceparams[type].pdihs.cpA;
2036 phi0[s] = forceparams[type].pdihs.phiA*DEG2RAD;
2037 mult[s] = forceparams[type].pdihs.mult;
2039 /* At the end fill the arrays with identical entries */
2040 if (iu + nfa1 < nbonds)
2046 /* Caclulate UNROLL dihedral angles at once */
2047 dih_angle_simd(x, ai, aj, ak, al, &pbc_simd,
2050 &mx_S, &my_S, &mz_S,
2051 &nx_S, &ny_S, &nz_S,
2056 cp_S = gmx_load_pr(cp);
2057 phi0_S = gmx_load_pr(phi0);
2058 mult_S = gmx_load_pr(mult);
2060 mdphi_S = gmx_sub_pr(gmx_mul_pr(mult_S, phi_S), phi0_S);
2062 /* Calculate UNROLL sines at once */
2063 gmx_sincos_pr(mdphi_S, &sin_S, &cos_S);
2064 mddphi_S = gmx_mul_pr(gmx_mul_pr(cp_S, mult_S), sin_S);
2065 sf_i_S = gmx_mul_pr(mddphi_S, nrkj_m2_S);
2066 msf_l_S = gmx_mul_pr(mddphi_S, nrkj_n2_S);
2068 /* After this m?_S will contain f[i] */
2069 mx_S = gmx_mul_pr(sf_i_S, mx_S);
2070 my_S = gmx_mul_pr(sf_i_S, my_S);
2071 mz_S = gmx_mul_pr(sf_i_S, mz_S);
2073 /* After this m?_S will contain -f[l] */
2074 nx_S = gmx_mul_pr(msf_l_S, nx_S);
2075 ny_S = gmx_mul_pr(msf_l_S, ny_S);
2076 nz_S = gmx_mul_pr(msf_l_S, nz_S);
2078 gmx_store_pr(dr + 0*UNROLL, mx_S);
2079 gmx_store_pr(dr + 1*UNROLL, my_S);
2080 gmx_store_pr(dr + 2*UNROLL, mz_S);
2081 gmx_store_pr(dr + 3*UNROLL, nx_S);
2082 gmx_store_pr(dr + 4*UNROLL, ny_S);
2083 gmx_store_pr(dr + 5*UNROLL, nz_S);
2089 do_dih_fup_noshiftf_precalc(ai[s], aj[s], ak[s], al[s],
2094 dr[(DIM+XX)*UNROLL+s],
2095 dr[(DIM+YY)*UNROLL+s],
2096 dr[(DIM+ZZ)*UNROLL+s],
2101 while (s < UNROLL && iu < nbonds);
2106 #endif /* SIMD_BONDEDS */
2109 real idihs(int nbonds,
2110 const t_iatom forceatoms[], const t_iparams forceparams[],
2111 const rvec x[], rvec f[], rvec fshift[],
2112 const t_pbc *pbc, const t_graph *g,
2113 real lambda, real *dvdlambda,
2114 const t_mdatoms *md, t_fcdata *fcd,
2115 int *global_atom_index)
2117 int i, type, ai, aj, ak, al;
2119 real phi, phi0, dphi0, ddphi, sign, vtot;
2120 rvec r_ij, r_kj, r_kl, m, n;
2121 real L1, kk, dp, dp2, kA, kB, pA, pB, dvdl_term;
2126 for (i = 0; (i < nbonds); )
2128 type = forceatoms[i++];
2129 ai = forceatoms[i++];
2130 aj = forceatoms[i++];
2131 ak = forceatoms[i++];
2132 al = forceatoms[i++];
2134 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2135 &sign, &t1, &t2, &t3); /* 84 */
2137 /* phi can jump if phi0 is close to Pi/-Pi, which will cause huge
2138 * force changes if we just apply a normal harmonic.
2139 * Instead, we first calculate phi-phi0 and take it modulo (-Pi,Pi).
2140 * This means we will never have the periodicity problem, unless
2141 * the dihedral is Pi away from phiO, which is very unlikely due to
2144 kA = forceparams[type].harmonic.krA;
2145 kB = forceparams[type].harmonic.krB;
2146 pA = forceparams[type].harmonic.rA;
2147 pB = forceparams[type].harmonic.rB;
2149 kk = L1*kA + lambda*kB;
2150 phi0 = (L1*pA + lambda*pB)*DEG2RAD;
2151 dphi0 = (pB - pA)*DEG2RAD;
2155 make_dp_periodic(&dp);
2162 dvdl_term += 0.5*(kB - kA)*dp2 - kk*dphi0*dp;
2164 do_dih_fup(ai, aj, ak, al, (real)(-ddphi), r_ij, r_kj, r_kl, m, n,
2165 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2170 fprintf(debug, "idih: (%d,%d,%d,%d) phi=%g\n",
2171 ai, aj, ak, al, phi);
2176 *dvdlambda += dvdl_term;
2181 real posres(int nbonds,
2182 const t_iatom forceatoms[], const t_iparams forceparams[],
2183 const rvec x[], rvec f[], rvec vir_diag,
2185 real lambda, real *dvdlambda,
2186 int refcoord_scaling, int ePBC, rvec comA, rvec comB)
2188 int i, ai, m, d, type, ki, npbcdim = 0;
2189 const t_iparams *pr;
2192 real posA, posB, ref = 0;
2193 rvec comA_sc, comB_sc, rdist, dpdl, pos, dx;
2194 gmx_bool bForceValid = TRUE;
2196 if ((f == NULL) || (vir_diag == NULL)) /* should both be null together! */
2198 bForceValid = FALSE;
2201 npbcdim = ePBC2npbcdim(ePBC);
2203 if (refcoord_scaling == erscCOM)
2205 clear_rvec(comA_sc);
2206 clear_rvec(comB_sc);
2207 for (m = 0; m < npbcdim; m++)
2209 for (d = m; d < npbcdim; d++)
2211 comA_sc[m] += comA[d]*pbc->box[d][m];
2212 comB_sc[m] += comB[d]*pbc->box[d][m];
2220 for (i = 0; (i < nbonds); )
2222 type = forceatoms[i++];
2223 ai = forceatoms[i++];
2224 pr = &forceparams[type];
2226 for (m = 0; m < DIM; m++)
2228 posA = forceparams[type].posres.pos0A[m];
2229 posB = forceparams[type].posres.pos0B[m];
2232 switch (refcoord_scaling)
2236 rdist[m] = L1*posA + lambda*posB;
2237 dpdl[m] = posB - posA;
2240 /* Box relative coordinates are stored for dimensions with pbc */
2241 posA *= pbc->box[m][m];
2242 posB *= pbc->box[m][m];
2243 for (d = m+1; d < npbcdim; d++)
2245 posA += forceparams[type].posres.pos0A[d]*pbc->box[d][m];
2246 posB += forceparams[type].posres.pos0B[d]*pbc->box[d][m];
2248 ref = L1*posA + lambda*posB;
2250 dpdl[m] = posB - posA;
2253 ref = L1*comA_sc[m] + lambda*comB_sc[m];
2254 rdist[m] = L1*posA + lambda*posB;
2255 dpdl[m] = comB_sc[m] - comA_sc[m] + posB - posA;
2258 gmx_fatal(FARGS, "No such scaling method implemented");
2263 ref = L1*posA + lambda*posB;
2265 dpdl[m] = posB - posA;
2268 /* We do pbc_dx with ref+rdist,
2269 * since with only ref we can be up to half a box vector wrong.
2271 pos[m] = ref + rdist[m];
2276 pbc_dx(pbc, x[ai], pos, dx);
2280 rvec_sub(x[ai], pos, dx);
2283 for (m = 0; (m < DIM); m++)
2285 kk = L1*pr->posres.fcA[m] + lambda*pr->posres.fcB[m];
2287 vtot += 0.5*kk*dx[m]*dx[m];
2289 0.5*(pr->posres.fcB[m] - pr->posres.fcA[m])*dx[m]*dx[m]
2292 /* Here we correct for the pbc_dx which included rdist */
2296 vir_diag[m] -= 0.5*(dx[m] + rdist[m])*fm;
2304 static real low_angres(int nbonds,
2305 const t_iatom forceatoms[], const t_iparams forceparams[],
2306 const rvec x[], rvec f[], rvec fshift[],
2307 const t_pbc *pbc, const t_graph *g,
2308 real lambda, real *dvdlambda,
2311 int i, m, type, ai, aj, ak, al;
2313 real phi, cos_phi, cos_phi2, vid, vtot, dVdphi;
2314 rvec r_ij, r_kl, f_i, f_k = {0, 0, 0};
2315 real st, sth, nrij2, nrkl2, c, cij, ckl;
2318 t2 = 0; /* avoid warning with gcc-3.3. It is never used uninitialized */
2321 ak = al = 0; /* to avoid warnings */
2322 for (i = 0; i < nbonds; )
2324 type = forceatoms[i++];
2325 ai = forceatoms[i++];
2326 aj = forceatoms[i++];
2327 t1 = pbc_rvec_sub(pbc, x[aj], x[ai], r_ij); /* 3 */
2330 ak = forceatoms[i++];
2331 al = forceatoms[i++];
2332 t2 = pbc_rvec_sub(pbc, x[al], x[ak], r_kl); /* 3 */
2341 cos_phi = cos_angle(r_ij, r_kl); /* 25 */
2342 phi = acos(cos_phi); /* 10 */
2344 *dvdlambda += dopdihs_min(forceparams[type].pdihs.cpA,
2345 forceparams[type].pdihs.cpB,
2346 forceparams[type].pdihs.phiA,
2347 forceparams[type].pdihs.phiB,
2348 forceparams[type].pdihs.mult,
2349 phi, lambda, &vid, &dVdphi); /* 40 */
2353 cos_phi2 = sqr(cos_phi); /* 1 */
2356 st = -dVdphi*gmx_invsqrt(1 - cos_phi2); /* 12 */
2357 sth = st*cos_phi; /* 1 */
2358 nrij2 = iprod(r_ij, r_ij); /* 5 */
2359 nrkl2 = iprod(r_kl, r_kl); /* 5 */
2361 c = st*gmx_invsqrt(nrij2*nrkl2); /* 11 */
2362 cij = sth/nrij2; /* 10 */
2363 ckl = sth/nrkl2; /* 10 */
2365 for (m = 0; m < DIM; m++) /* 18+18 */
2367 f_i[m] = (c*r_kl[m]-cij*r_ij[m]);
2372 f_k[m] = (c*r_ij[m]-ckl*r_kl[m]);
2380 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
2383 rvec_inc(fshift[t1], f_i);
2384 rvec_dec(fshift[CENTRAL], f_i);
2389 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, al), dt);
2392 rvec_inc(fshift[t2], f_k);
2393 rvec_dec(fshift[CENTRAL], f_k);
2398 return vtot; /* 184 / 157 (bZAxis) total */
2401 real angres(int nbonds,
2402 const t_iatom forceatoms[], const t_iparams forceparams[],
2403 const rvec x[], rvec f[], rvec fshift[],
2404 const t_pbc *pbc, const t_graph *g,
2405 real lambda, real *dvdlambda,
2406 const t_mdatoms *md, t_fcdata *fcd,
2407 int *global_atom_index)
2409 return low_angres(nbonds, forceatoms, forceparams, x, f, fshift, pbc, g,
2410 lambda, dvdlambda, FALSE);
2413 real angresz(int nbonds,
2414 const t_iatom forceatoms[], const t_iparams forceparams[],
2415 const rvec x[], rvec f[], rvec fshift[],
2416 const t_pbc *pbc, const t_graph *g,
2417 real lambda, real *dvdlambda,
2418 const t_mdatoms *md, t_fcdata *fcd,
2419 int *global_atom_index)
2421 return low_angres(nbonds, forceatoms, forceparams, x, f, fshift, pbc, g,
2422 lambda, dvdlambda, TRUE);
2425 real dihres(int nbonds,
2426 const t_iatom forceatoms[], const t_iparams forceparams[],
2427 const rvec x[], rvec f[], rvec fshift[],
2428 const t_pbc *pbc, const t_graph *g,
2429 real lambda, real *dvdlambda,
2430 const t_mdatoms *md, t_fcdata *fcd,
2431 int *global_atom_index)
2434 int ai, aj, ak, al, i, k, type, t1, t2, t3;
2435 real phi0A, phi0B, dphiA, dphiB, kfacA, kfacB, phi0, dphi, kfac;
2436 real phi, ddphi, ddp, ddp2, dp, sign, d2r, fc, L1;
2437 rvec r_ij, r_kj, r_kl, m, n;
2444 for (i = 0; (i < nbonds); )
2446 type = forceatoms[i++];
2447 ai = forceatoms[i++];
2448 aj = forceatoms[i++];
2449 ak = forceatoms[i++];
2450 al = forceatoms[i++];
2452 phi0A = forceparams[type].dihres.phiA*d2r;
2453 dphiA = forceparams[type].dihres.dphiA*d2r;
2454 kfacA = forceparams[type].dihres.kfacA;
2456 phi0B = forceparams[type].dihres.phiB*d2r;
2457 dphiB = forceparams[type].dihres.dphiB*d2r;
2458 kfacB = forceparams[type].dihres.kfacB;
2460 phi0 = L1*phi0A + lambda*phi0B;
2461 dphi = L1*dphiA + lambda*dphiB;
2462 kfac = L1*kfacA + lambda*kfacB;
2464 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2465 &sign, &t1, &t2, &t3);
2470 fprintf(debug, "dihres[%d]: %d %d %d %d : phi=%f, dphi=%f, kfac=%f\n",
2471 k++, ai, aj, ak, al, phi0, dphi, kfac);
2473 /* phi can jump if phi0 is close to Pi/-Pi, which will cause huge
2474 * force changes if we just apply a normal harmonic.
2475 * Instead, we first calculate phi-phi0 and take it modulo (-Pi,Pi).
2476 * This means we will never have the periodicity problem, unless
2477 * the dihedral is Pi away from phiO, which is very unlikely due to
2481 make_dp_periodic(&dp);
2487 else if (dp < -dphi)
2499 vtot += 0.5*kfac*ddp2;
2502 *dvdlambda += 0.5*(kfacB - kfacA)*ddp2;
2503 /* lambda dependence from changing restraint distances */
2506 *dvdlambda -= kfac*ddp*((dphiB - dphiA)+(phi0B - phi0A));
2510 *dvdlambda += kfac*ddp*((dphiB - dphiA)-(phi0B - phi0A));
2512 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
2513 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2520 real unimplemented(int nbonds,
2521 const t_iatom forceatoms[], const t_iparams forceparams[],
2522 const rvec x[], rvec f[], rvec fshift[],
2523 const t_pbc *pbc, const t_graph *g,
2524 real lambda, real *dvdlambda,
2525 const t_mdatoms *md, t_fcdata *fcd,
2526 int *global_atom_index)
2528 gmx_impl("*** you are using a not implemented function");
2530 return 0.0; /* To make the compiler happy */
2533 real rbdihs(int nbonds,
2534 const t_iatom forceatoms[], const t_iparams forceparams[],
2535 const rvec x[], rvec f[], rvec fshift[],
2536 const t_pbc *pbc, const t_graph *g,
2537 real lambda, real *dvdlambda,
2538 const t_mdatoms *md, t_fcdata *fcd,
2539 int *global_atom_index)
2541 const real c0 = 0.0, c1 = 1.0, c2 = 2.0, c3 = 3.0, c4 = 4.0, c5 = 5.0;
2542 int type, ai, aj, ak, al, i, j;
2544 rvec r_ij, r_kj, r_kl, m, n;
2545 real parmA[NR_RBDIHS];
2546 real parmB[NR_RBDIHS];
2547 real parm[NR_RBDIHS];
2548 real cos_phi, phi, rbp, rbpBA;
2549 real v, sign, ddphi, sin_phi;
2551 real L1 = 1.0-lambda;
2555 for (i = 0; (i < nbonds); )
2557 type = forceatoms[i++];
2558 ai = forceatoms[i++];
2559 aj = forceatoms[i++];
2560 ak = forceatoms[i++];
2561 al = forceatoms[i++];
2563 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2564 &sign, &t1, &t2, &t3); /* 84 */
2566 /* Change to polymer convention */
2573 phi -= M_PI; /* 1 */
2577 /* Beware of accuracy loss, cannot use 1-sqrt(cos^2) ! */
2580 for (j = 0; (j < NR_RBDIHS); j++)
2582 parmA[j] = forceparams[type].rbdihs.rbcA[j];
2583 parmB[j] = forceparams[type].rbdihs.rbcB[j];
2584 parm[j] = L1*parmA[j]+lambda*parmB[j];
2586 /* Calculate cosine powers */
2587 /* Calculate the energy */
2588 /* Calculate the derivative */
2591 dvdl_term += (parmB[0]-parmA[0]);
2596 rbpBA = parmB[1]-parmA[1];
2597 ddphi += rbp*cosfac;
2600 dvdl_term += cosfac*rbpBA;
2602 rbpBA = parmB[2]-parmA[2];
2603 ddphi += c2*rbp*cosfac;
2606 dvdl_term += cosfac*rbpBA;
2608 rbpBA = parmB[3]-parmA[3];
2609 ddphi += c3*rbp*cosfac;
2612 dvdl_term += cosfac*rbpBA;
2614 rbpBA = parmB[4]-parmA[4];
2615 ddphi += c4*rbp*cosfac;
2618 dvdl_term += cosfac*rbpBA;
2620 rbpBA = parmB[5]-parmA[5];
2621 ddphi += c5*rbp*cosfac;
2624 dvdl_term += cosfac*rbpBA;
2626 ddphi = -ddphi*sin_phi; /* 11 */
2628 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
2629 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2632 *dvdlambda += dvdl_term;
2637 int cmap_setup_grid_index(int ip, int grid_spacing, int *ipm1, int *ipp1, int *ipp2)
2643 ip = ip + grid_spacing - 1;
2645 else if (ip > grid_spacing)
2647 ip = ip - grid_spacing - 1;
2656 im1 = grid_spacing - 1;
2658 else if (ip == grid_spacing-2)
2662 else if (ip == grid_spacing-1)
2676 real cmap_dihs(int nbonds,
2677 const t_iatom forceatoms[], const t_iparams forceparams[],
2678 const gmx_cmap_t *cmap_grid,
2679 const rvec x[], rvec f[], rvec fshift[],
2680 const t_pbc *pbc, const t_graph *g,
2681 real lambda, real *dvdlambda,
2682 const t_mdatoms *md, t_fcdata *fcd,
2683 int *global_atom_index)
2685 int i, j, k, n, idx;
2686 int ai, aj, ak, al, am;
2687 int a1i, a1j, a1k, a1l, a2i, a2j, a2k, a2l;
2689 int t11, t21, t31, t12, t22, t32;
2690 int iphi1, ip1m1, ip1p1, ip1p2;
2691 int iphi2, ip2m1, ip2p1, ip2p2;
2693 int pos1, pos2, pos3, pos4, tmp;
2695 real ty[4], ty1[4], ty2[4], ty12[4], tc[16], tx[16];
2696 real phi1, psi1, cos_phi1, sin_phi1, sign1, xphi1;
2697 real phi2, psi2, cos_phi2, sin_phi2, sign2, xphi2;
2698 real dx, xx, tt, tu, e, df1, df2, ddf1, ddf2, ddf12, vtot;
2699 real ra21, rb21, rg21, rg1, rgr1, ra2r1, rb2r1, rabr1;
2700 real ra22, rb22, rg22, rg2, rgr2, ra2r2, rb2r2, rabr2;
2701 real fg1, hg1, fga1, hgb1, gaa1, gbb1;
2702 real fg2, hg2, fga2, hgb2, gaa2, gbb2;
2705 rvec r1_ij, r1_kj, r1_kl, m1, n1;
2706 rvec r2_ij, r2_kj, r2_kl, m2, n2;
2707 rvec f1_i, f1_j, f1_k, f1_l;
2708 rvec f2_i, f2_j, f2_k, f2_l;
2709 rvec a1, b1, a2, b2;
2710 rvec f1, g1, h1, f2, g2, h2;
2711 rvec dtf1, dtg1, dth1, dtf2, dtg2, dth2;
2712 ivec jt1, dt1_ij, dt1_kj, dt1_lj;
2713 ivec jt2, dt2_ij, dt2_kj, dt2_lj;
2717 int loop_index[4][4] = {
2724 /* Total CMAP energy */
2727 for (n = 0; n < nbonds; )
2729 /* Five atoms are involved in the two torsions */
2730 type = forceatoms[n++];
2731 ai = forceatoms[n++];
2732 aj = forceatoms[n++];
2733 ak = forceatoms[n++];
2734 al = forceatoms[n++];
2735 am = forceatoms[n++];
2737 /* Which CMAP type is this */
2738 cmapA = forceparams[type].cmap.cmapA;
2739 cmapd = cmap_grid->cmapdata[cmapA].cmap;
2747 phi1 = dih_angle(x[a1i], x[a1j], x[a1k], x[a1l], pbc, r1_ij, r1_kj, r1_kl, m1, n1,
2748 &sign1, &t11, &t21, &t31); /* 84 */
2750 cos_phi1 = cos(phi1);
2752 a1[0] = r1_ij[1]*r1_kj[2]-r1_ij[2]*r1_kj[1];
2753 a1[1] = r1_ij[2]*r1_kj[0]-r1_ij[0]*r1_kj[2];
2754 a1[2] = r1_ij[0]*r1_kj[1]-r1_ij[1]*r1_kj[0]; /* 9 */
2756 b1[0] = r1_kl[1]*r1_kj[2]-r1_kl[2]*r1_kj[1];
2757 b1[1] = r1_kl[2]*r1_kj[0]-r1_kl[0]*r1_kj[2];
2758 b1[2] = r1_kl[0]*r1_kj[1]-r1_kl[1]*r1_kj[0]; /* 9 */
2760 tmp = pbc_rvec_sub(pbc, x[a1l], x[a1k], h1);
2762 ra21 = iprod(a1, a1); /* 5 */
2763 rb21 = iprod(b1, b1); /* 5 */
2764 rg21 = iprod(r1_kj, r1_kj); /* 5 */
2770 rabr1 = sqrt(ra2r1*rb2r1);
2772 sin_phi1 = rg1 * rabr1 * iprod(a1, h1) * (-1);
2774 if (cos_phi1 < -0.5 || cos_phi1 > 0.5)
2776 phi1 = asin(sin_phi1);
2786 phi1 = -M_PI - phi1;
2792 phi1 = acos(cos_phi1);
2800 xphi1 = phi1 + M_PI; /* 1 */
2802 /* Second torsion */
2808 phi2 = dih_angle(x[a2i], x[a2j], x[a2k], x[a2l], pbc, r2_ij, r2_kj, r2_kl, m2, n2,
2809 &sign2, &t12, &t22, &t32); /* 84 */
2811 cos_phi2 = cos(phi2);
2813 a2[0] = r2_ij[1]*r2_kj[2]-r2_ij[2]*r2_kj[1];
2814 a2[1] = r2_ij[2]*r2_kj[0]-r2_ij[0]*r2_kj[2];
2815 a2[2] = r2_ij[0]*r2_kj[1]-r2_ij[1]*r2_kj[0]; /* 9 */
2817 b2[0] = r2_kl[1]*r2_kj[2]-r2_kl[2]*r2_kj[1];
2818 b2[1] = r2_kl[2]*r2_kj[0]-r2_kl[0]*r2_kj[2];
2819 b2[2] = r2_kl[0]*r2_kj[1]-r2_kl[1]*r2_kj[0]; /* 9 */
2821 tmp = pbc_rvec_sub(pbc, x[a2l], x[a2k], h2);
2823 ra22 = iprod(a2, a2); /* 5 */
2824 rb22 = iprod(b2, b2); /* 5 */
2825 rg22 = iprod(r2_kj, r2_kj); /* 5 */
2831 rabr2 = sqrt(ra2r2*rb2r2);
2833 sin_phi2 = rg2 * rabr2 * iprod(a2, h2) * (-1);
2835 if (cos_phi2 < -0.5 || cos_phi2 > 0.5)
2837 phi2 = asin(sin_phi2);
2847 phi2 = -M_PI - phi2;
2853 phi2 = acos(cos_phi2);
2861 xphi2 = phi2 + M_PI; /* 1 */
2863 /* Range mangling */
2866 xphi1 = xphi1 + 2*M_PI;
2868 else if (xphi1 >= 2*M_PI)
2870 xphi1 = xphi1 - 2*M_PI;
2875 xphi2 = xphi2 + 2*M_PI;
2877 else if (xphi2 >= 2*M_PI)
2879 xphi2 = xphi2 - 2*M_PI;
2882 /* Number of grid points */
2883 dx = 2*M_PI / cmap_grid->grid_spacing;
2885 /* Where on the grid are we */
2886 iphi1 = (int)(xphi1/dx);
2887 iphi2 = (int)(xphi2/dx);
2889 iphi1 = cmap_setup_grid_index(iphi1, cmap_grid->grid_spacing, &ip1m1, &ip1p1, &ip1p2);
2890 iphi2 = cmap_setup_grid_index(iphi2, cmap_grid->grid_spacing, &ip2m1, &ip2p1, &ip2p2);
2892 pos1 = iphi1*cmap_grid->grid_spacing+iphi2;
2893 pos2 = ip1p1*cmap_grid->grid_spacing+iphi2;
2894 pos3 = ip1p1*cmap_grid->grid_spacing+ip2p1;
2895 pos4 = iphi1*cmap_grid->grid_spacing+ip2p1;
2897 ty[0] = cmapd[pos1*4];
2898 ty[1] = cmapd[pos2*4];
2899 ty[2] = cmapd[pos3*4];
2900 ty[3] = cmapd[pos4*4];
2902 ty1[0] = cmapd[pos1*4+1];
2903 ty1[1] = cmapd[pos2*4+1];
2904 ty1[2] = cmapd[pos3*4+1];
2905 ty1[3] = cmapd[pos4*4+1];
2907 ty2[0] = cmapd[pos1*4+2];
2908 ty2[1] = cmapd[pos2*4+2];
2909 ty2[2] = cmapd[pos3*4+2];
2910 ty2[3] = cmapd[pos4*4+2];
2912 ty12[0] = cmapd[pos1*4+3];
2913 ty12[1] = cmapd[pos2*4+3];
2914 ty12[2] = cmapd[pos3*4+3];
2915 ty12[3] = cmapd[pos4*4+3];
2917 /* Switch to degrees */
2918 dx = 360.0 / cmap_grid->grid_spacing;
2919 xphi1 = xphi1 * RAD2DEG;
2920 xphi2 = xphi2 * RAD2DEG;
2922 for (i = 0; i < 4; i++) /* 16 */
2925 tx[i+4] = ty1[i]*dx;
2926 tx[i+8] = ty2[i]*dx;
2927 tx[i+12] = ty12[i]*dx*dx;
2931 for (i = 0; i < 4; i++) /* 1056 */
2933 for (j = 0; j < 4; j++)
2936 for (k = 0; k < 16; k++)
2938 xx = xx + cmap_coeff_matrix[k*16+idx]*tx[k];
2946 tt = (xphi1-iphi1*dx)/dx;
2947 tu = (xphi2-iphi2*dx)/dx;
2956 for (i = 3; i >= 0; i--)
2958 l1 = loop_index[i][3];
2959 l2 = loop_index[i][2];
2960 l3 = loop_index[i][1];
2962 e = tt * e + ((tc[i*4+3]*tu+tc[i*4+2])*tu + tc[i*4+1])*tu+tc[i*4];
2963 df1 = tu * df1 + (3.0*tc[l1]*tt+2.0*tc[l2])*tt+tc[l3];
2964 df2 = tt * df2 + (3.0*tc[i*4+3]*tu+2.0*tc[i*4+2])*tu+tc[i*4+1];
2965 ddf1 = tu * ddf1 + 2.0*3.0*tc[l1]*tt+2.0*tc[l2];
2966 ddf2 = tt * ddf2 + 2.0*3.0*tc[4*i+3]*tu+2.0*tc[4*i+2];
2969 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) +
2970 3.0*tu*tu*(tc[7]+2.0*tc[11]*tt+3.0*tc[15]*tt*tt);
2975 ddf1 = ddf1 * fac * fac;
2976 ddf2 = ddf2 * fac * fac;
2977 ddf12 = ddf12 * fac * fac;
2982 /* Do forces - first torsion */
2983 fg1 = iprod(r1_ij, r1_kj);
2984 hg1 = iprod(r1_kl, r1_kj);
2985 fga1 = fg1*ra2r1*rgr1;
2986 hgb1 = hg1*rb2r1*rgr1;
2990 for (i = 0; i < DIM; i++)
2992 dtf1[i] = gaa1 * a1[i];
2993 dtg1[i] = fga1 * a1[i] - hgb1 * b1[i];
2994 dth1[i] = gbb1 * b1[i];
2996 f1[i] = df1 * dtf1[i];
2997 g1[i] = df1 * dtg1[i];
2998 h1[i] = df1 * dth1[i];
3001 f1_j[i] = -f1[i] - g1[i];
3002 f1_k[i] = h1[i] + g1[i];
3005 f[a1i][i] = f[a1i][i] + f1_i[i];
3006 f[a1j][i] = f[a1j][i] + f1_j[i]; /* - f1[i] - g1[i] */
3007 f[a1k][i] = f[a1k][i] + f1_k[i]; /* h1[i] + g1[i] */
3008 f[a1l][i] = f[a1l][i] + f1_l[i]; /* h1[i] */
3011 /* Do forces - second torsion */
3012 fg2 = iprod(r2_ij, r2_kj);
3013 hg2 = iprod(r2_kl, r2_kj);
3014 fga2 = fg2*ra2r2*rgr2;
3015 hgb2 = hg2*rb2r2*rgr2;
3019 for (i = 0; i < DIM; i++)
3021 dtf2[i] = gaa2 * a2[i];
3022 dtg2[i] = fga2 * a2[i] - hgb2 * b2[i];
3023 dth2[i] = gbb2 * b2[i];
3025 f2[i] = df2 * dtf2[i];
3026 g2[i] = df2 * dtg2[i];
3027 h2[i] = df2 * dth2[i];
3030 f2_j[i] = -f2[i] - g2[i];
3031 f2_k[i] = h2[i] + g2[i];
3034 f[a2i][i] = f[a2i][i] + f2_i[i]; /* f2[i] */
3035 f[a2j][i] = f[a2j][i] + f2_j[i]; /* - f2[i] - g2[i] */
3036 f[a2k][i] = f[a2k][i] + f2_k[i]; /* h2[i] + g2[i] */
3037 f[a2l][i] = f[a2l][i] + f2_l[i]; /* - h2[i] */
3043 copy_ivec(SHIFT_IVEC(g, a1j), jt1);
3044 ivec_sub(SHIFT_IVEC(g, a1i), jt1, dt1_ij);
3045 ivec_sub(SHIFT_IVEC(g, a1k), jt1, dt1_kj);
3046 ivec_sub(SHIFT_IVEC(g, a1l), jt1, dt1_lj);
3047 t11 = IVEC2IS(dt1_ij);
3048 t21 = IVEC2IS(dt1_kj);
3049 t31 = IVEC2IS(dt1_lj);
3051 copy_ivec(SHIFT_IVEC(g, a2j), jt2);
3052 ivec_sub(SHIFT_IVEC(g, a2i), jt2, dt2_ij);
3053 ivec_sub(SHIFT_IVEC(g, a2k), jt2, dt2_kj);
3054 ivec_sub(SHIFT_IVEC(g, a2l), jt2, dt2_lj);
3055 t12 = IVEC2IS(dt2_ij);
3056 t22 = IVEC2IS(dt2_kj);
3057 t32 = IVEC2IS(dt2_lj);
3061 t31 = pbc_rvec_sub(pbc, x[a1l], x[a1j], h1);
3062 t32 = pbc_rvec_sub(pbc, x[a2l], x[a2j], h2);
3070 rvec_inc(fshift[t11], f1_i);
3071 rvec_inc(fshift[CENTRAL], f1_j);
3072 rvec_inc(fshift[t21], f1_k);
3073 rvec_inc(fshift[t31], f1_l);
3075 rvec_inc(fshift[t21], f2_i);
3076 rvec_inc(fshift[CENTRAL], f2_j);
3077 rvec_inc(fshift[t22], f2_k);
3078 rvec_inc(fshift[t32], f2_l);
3085 /***********************************************************
3087 * G R O M O S 9 6 F U N C T I O N S
3089 ***********************************************************/
3090 real g96harmonic(real kA, real kB, real xA, real xB, real x, real lambda,
3093 const real half = 0.5;
3094 real L1, kk, x0, dx, dx2;
3095 real v, f, dvdlambda;
3098 kk = L1*kA+lambda*kB;
3099 x0 = L1*xA+lambda*xB;
3106 dvdlambda = half*(kB-kA)*dx2 + (xA-xB)*kk*dx;
3113 /* That was 21 flops */
3116 real g96bonds(int nbonds,
3117 const t_iatom forceatoms[], const t_iparams forceparams[],
3118 const rvec x[], rvec f[], rvec fshift[],
3119 const t_pbc *pbc, const t_graph *g,
3120 real lambda, real *dvdlambda,
3121 const t_mdatoms *md, t_fcdata *fcd,
3122 int *global_atom_index)
3124 int i, m, ki, ai, aj, type;
3125 real dr2, fbond, vbond, fij, vtot;
3130 for (i = 0; (i < nbonds); )
3132 type = forceatoms[i++];
3133 ai = forceatoms[i++];
3134 aj = forceatoms[i++];
3136 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
3137 dr2 = iprod(dx, dx); /* 5 */
3139 *dvdlambda += g96harmonic(forceparams[type].harmonic.krA,
3140 forceparams[type].harmonic.krB,
3141 forceparams[type].harmonic.rA,
3142 forceparams[type].harmonic.rB,
3143 dr2, lambda, &vbond, &fbond);
3145 vtot += 0.5*vbond; /* 1*/
3149 fprintf(debug, "G96-BONDS: dr = %10g vbond = %10g fbond = %10g\n",
3150 sqrt(dr2), vbond, fbond);
3156 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
3159 for (m = 0; (m < DIM); m++) /* 15 */
3164 fshift[ki][m] += fij;
3165 fshift[CENTRAL][m] -= fij;
3171 real g96bond_angle(const rvec xi, const rvec xj, const rvec xk, const t_pbc *pbc,
3172 rvec r_ij, rvec r_kj,
3174 /* Return value is the angle between the bonds i-j and j-k */
3178 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
3179 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
3181 costh = cos_angle(r_ij, r_kj); /* 25 */
3186 real g96angles(int nbonds,
3187 const t_iatom forceatoms[], const t_iparams forceparams[],
3188 const rvec x[], rvec f[], rvec fshift[],
3189 const t_pbc *pbc, const t_graph *g,
3190 real lambda, real *dvdlambda,
3191 const t_mdatoms *md, t_fcdata *fcd,
3192 int *global_atom_index)
3194 int i, ai, aj, ak, type, m, t1, t2;
3196 real cos_theta, dVdt, va, vtot;
3197 real rij_1, rij_2, rkj_1, rkj_2, rijrkj_1;
3199 ivec jt, dt_ij, dt_kj;
3202 for (i = 0; (i < nbonds); )
3204 type = forceatoms[i++];
3205 ai = forceatoms[i++];
3206 aj = forceatoms[i++];
3207 ak = forceatoms[i++];
3209 cos_theta = g96bond_angle(x[ai], x[aj], x[ak], pbc, r_ij, r_kj, &t1, &t2);
3211 *dvdlambda += g96harmonic(forceparams[type].harmonic.krA,
3212 forceparams[type].harmonic.krB,
3213 forceparams[type].harmonic.rA,
3214 forceparams[type].harmonic.rB,
3215 cos_theta, lambda, &va, &dVdt);
3218 rij_1 = gmx_invsqrt(iprod(r_ij, r_ij));
3219 rkj_1 = gmx_invsqrt(iprod(r_kj, r_kj));
3220 rij_2 = rij_1*rij_1;
3221 rkj_2 = rkj_1*rkj_1;
3222 rijrkj_1 = rij_1*rkj_1; /* 23 */
3227 fprintf(debug, "G96ANGLES: costheta = %10g vth = %10g dV/dct = %10g\n",
3228 cos_theta, va, dVdt);
3231 for (m = 0; (m < DIM); m++) /* 42 */
3233 f_i[m] = dVdt*(r_kj[m]*rijrkj_1 - r_ij[m]*rij_2*cos_theta);
3234 f_k[m] = dVdt*(r_ij[m]*rijrkj_1 - r_kj[m]*rkj_2*cos_theta);
3235 f_j[m] = -f_i[m]-f_k[m];
3243 copy_ivec(SHIFT_IVEC(g, aj), jt);
3245 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3246 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3247 t1 = IVEC2IS(dt_ij);
3248 t2 = IVEC2IS(dt_kj);
3250 rvec_inc(fshift[t1], f_i);
3251 rvec_inc(fshift[CENTRAL], f_j);
3252 rvec_inc(fshift[t2], f_k); /* 9 */
3258 real cross_bond_bond(int nbonds,
3259 const t_iatom forceatoms[], const t_iparams forceparams[],
3260 const rvec x[], rvec f[], rvec fshift[],
3261 const t_pbc *pbc, const t_graph *g,
3262 real lambda, real *dvdlambda,
3263 const t_mdatoms *md, t_fcdata *fcd,
3264 int *global_atom_index)
3266 /* Potential from Lawrence and Skimmer, Chem. Phys. Lett. 372 (2003)
3269 int i, ai, aj, ak, type, m, t1, t2;
3271 real vtot, vrr, s1, s2, r1, r2, r1e, r2e, krr;
3273 ivec jt, dt_ij, dt_kj;
3276 for (i = 0; (i < nbonds); )
3278 type = forceatoms[i++];
3279 ai = forceatoms[i++];
3280 aj = forceatoms[i++];
3281 ak = forceatoms[i++];
3282 r1e = forceparams[type].cross_bb.r1e;
3283 r2e = forceparams[type].cross_bb.r2e;
3284 krr = forceparams[type].cross_bb.krr;
3286 /* Compute distance vectors ... */
3287 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
3288 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
3290 /* ... and their lengths */
3294 /* Deviations from ideality */
3298 /* Energy (can be negative!) */
3303 svmul(-krr*s2/r1, r_ij, f_i);
3304 svmul(-krr*s1/r2, r_kj, f_k);
3306 for (m = 0; (m < DIM); m++) /* 12 */
3308 f_j[m] = -f_i[m] - f_k[m];
3317 copy_ivec(SHIFT_IVEC(g, aj), jt);
3319 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3320 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3321 t1 = IVEC2IS(dt_ij);
3322 t2 = IVEC2IS(dt_kj);
3324 rvec_inc(fshift[t1], f_i);
3325 rvec_inc(fshift[CENTRAL], f_j);
3326 rvec_inc(fshift[t2], f_k); /* 9 */
3332 real cross_bond_angle(int nbonds,
3333 const t_iatom forceatoms[], const t_iparams forceparams[],
3334 const rvec x[], rvec f[], rvec fshift[],
3335 const t_pbc *pbc, const t_graph *g,
3336 real lambda, real *dvdlambda,
3337 const t_mdatoms *md, t_fcdata *fcd,
3338 int *global_atom_index)
3340 /* Potential from Lawrence and Skimmer, Chem. Phys. Lett. 372 (2003)
3343 int i, ai, aj, ak, type, m, t1, t2, t3;
3344 rvec r_ij, r_kj, r_ik;
3345 real vtot, vrt, s1, s2, s3, r1, r2, r3, r1e, r2e, r3e, krt, k1, k2, k3;
3347 ivec jt, dt_ij, dt_kj;
3350 for (i = 0; (i < nbonds); )
3352 type = forceatoms[i++];
3353 ai = forceatoms[i++];
3354 aj = forceatoms[i++];
3355 ak = forceatoms[i++];
3356 r1e = forceparams[type].cross_ba.r1e;
3357 r2e = forceparams[type].cross_ba.r2e;
3358 r3e = forceparams[type].cross_ba.r3e;
3359 krt = forceparams[type].cross_ba.krt;
3361 /* Compute distance vectors ... */
3362 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
3363 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
3364 t3 = pbc_rvec_sub(pbc, x[ai], x[ak], r_ik);
3366 /* ... and their lengths */
3371 /* Deviations from ideality */
3376 /* Energy (can be negative!) */
3377 vrt = krt*s3*(s1+s2);
3383 k3 = -krt*(s1+s2)/r3;
3384 for (m = 0; (m < DIM); m++)
3386 f_i[m] = k1*r_ij[m] + k3*r_ik[m];
3387 f_k[m] = k2*r_kj[m] - k3*r_ik[m];
3388 f_j[m] = -f_i[m] - f_k[m];
3391 for (m = 0; (m < DIM); m++) /* 12 */
3401 copy_ivec(SHIFT_IVEC(g, aj), jt);
3403 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3404 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3405 t1 = IVEC2IS(dt_ij);
3406 t2 = IVEC2IS(dt_kj);
3408 rvec_inc(fshift[t1], f_i);
3409 rvec_inc(fshift[CENTRAL], f_j);
3410 rvec_inc(fshift[t2], f_k); /* 9 */
3416 static real bonded_tab(const char *type, int table_nr,
3417 const bondedtable_t *table, real kA, real kB, real r,
3418 real lambda, real *V, real *F)
3420 real k, tabscale, *VFtab, rt, eps, eps2, Yt, Ft, Geps, Heps2, Fp, VV, FF;
3422 real v, f, dvdlambda;
3424 k = (1.0 - lambda)*kA + lambda*kB;
3426 tabscale = table->scale;
3427 VFtab = table->data;
3433 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",
3434 type, table_nr, r, n0, n0+1, table->n);
3441 Geps = VFtab[nnn+2]*eps;
3442 Heps2 = VFtab[nnn+3]*eps2;
3443 Fp = Ft + Geps + Heps2;
3445 FF = Fp + Geps + 2.0*Heps2;
3447 *F = -k*FF*tabscale;
3449 dvdlambda = (kB - kA)*VV;
3453 /* That was 22 flops */
3456 real tab_bonds(int nbonds,
3457 const t_iatom forceatoms[], const t_iparams forceparams[],
3458 const rvec x[], rvec f[], rvec fshift[],
3459 const t_pbc *pbc, const t_graph *g,
3460 real lambda, real *dvdlambda,
3461 const t_mdatoms *md, t_fcdata *fcd,
3462 int *global_atom_index)
3464 int i, m, ki, ai, aj, type, table;
3465 real dr, dr2, fbond, vbond, fij, vtot;
3470 for (i = 0; (i < nbonds); )
3472 type = forceatoms[i++];
3473 ai = forceatoms[i++];
3474 aj = forceatoms[i++];
3476 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
3477 dr2 = iprod(dx, dx); /* 5 */
3478 dr = dr2*gmx_invsqrt(dr2); /* 10 */
3480 table = forceparams[type].tab.table;
3482 *dvdlambda += bonded_tab("bond", table,
3483 &fcd->bondtab[table],
3484 forceparams[type].tab.kA,
3485 forceparams[type].tab.kB,
3486 dr, lambda, &vbond, &fbond); /* 22 */
3494 vtot += vbond; /* 1*/
3495 fbond *= gmx_invsqrt(dr2); /* 6 */
3499 fprintf(debug, "TABBONDS: dr = %10g vbond = %10g fbond = %10g\n",
3505 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
3508 for (m = 0; (m < DIM); m++) /* 15 */
3513 fshift[ki][m] += fij;
3514 fshift[CENTRAL][m] -= fij;
3520 real tab_angles(int nbonds,
3521 const t_iatom forceatoms[], const t_iparams forceparams[],
3522 const rvec x[], rvec f[], rvec fshift[],
3523 const t_pbc *pbc, const t_graph *g,
3524 real lambda, real *dvdlambda,
3525 const t_mdatoms *md, t_fcdata *fcd,
3526 int *global_atom_index)
3528 int i, ai, aj, ak, t1, t2, type, table;
3530 real cos_theta, cos_theta2, theta, dVdt, va, vtot;
3531 ivec jt, dt_ij, dt_kj;
3534 for (i = 0; (i < nbonds); )
3536 type = forceatoms[i++];
3537 ai = forceatoms[i++];
3538 aj = forceatoms[i++];
3539 ak = forceatoms[i++];
3541 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
3542 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
3544 table = forceparams[type].tab.table;
3546 *dvdlambda += bonded_tab("angle", table,
3547 &fcd->angletab[table],
3548 forceparams[type].tab.kA,
3549 forceparams[type].tab.kB,
3550 theta, lambda, &va, &dVdt); /* 22 */
3553 cos_theta2 = sqr(cos_theta); /* 1 */
3562 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
3563 sth = st*cos_theta; /* 1 */
3567 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
3568 theta*RAD2DEG, va, dVdt);
3571 nrkj2 = iprod(r_kj, r_kj); /* 5 */
3572 nrij2 = iprod(r_ij, r_ij);
3574 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
3575 cii = sth/nrij2; /* 10 */
3576 ckk = sth/nrkj2; /* 10 */
3578 for (m = 0; (m < DIM); m++) /* 39 */
3580 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
3581 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
3582 f_j[m] = -f_i[m]-f_k[m];
3589 copy_ivec(SHIFT_IVEC(g, aj), jt);
3591 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3592 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3593 t1 = IVEC2IS(dt_ij);
3594 t2 = IVEC2IS(dt_kj);
3596 rvec_inc(fshift[t1], f_i);
3597 rvec_inc(fshift[CENTRAL], f_j);
3598 rvec_inc(fshift[t2], f_k);
3604 real tab_dihs(int nbonds,
3605 const t_iatom forceatoms[], const t_iparams forceparams[],
3606 const rvec x[], rvec f[], rvec fshift[],
3607 const t_pbc *pbc, const t_graph *g,
3608 real lambda, real *dvdlambda,
3609 const t_mdatoms *md, t_fcdata *fcd,
3610 int *global_atom_index)
3612 int i, type, ai, aj, ak, al, table;
3614 rvec r_ij, r_kj, r_kl, m, n;
3615 real phi, sign, ddphi, vpd, vtot;
3618 for (i = 0; (i < nbonds); )
3620 type = forceatoms[i++];
3621 ai = forceatoms[i++];
3622 aj = forceatoms[i++];
3623 ak = forceatoms[i++];
3624 al = forceatoms[i++];
3626 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
3627 &sign, &t1, &t2, &t3); /* 84 */
3629 table = forceparams[type].tab.table;
3631 /* Hopefully phi+M_PI never results in values < 0 */
3632 *dvdlambda += bonded_tab("dihedral", table,
3633 &fcd->dihtab[table],
3634 forceparams[type].tab.kA,
3635 forceparams[type].tab.kB,
3636 phi+M_PI, lambda, &vpd, &ddphi);
3639 do_dih_fup(ai, aj, ak, al, -ddphi, r_ij, r_kj, r_kl, m, n,
3640 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
3643 fprintf(debug, "pdih: (%d,%d,%d,%d) phi=%g\n",
3644 ai, aj, ak, al, phi);
3651 /* Return if this is a potential calculated in bondfree.c,
3652 * i.e. an interaction that actually calculates a potential and
3653 * works on multiple atoms (not e.g. a connection or a position restraint).
3655 static gmx_inline gmx_bool ftype_is_bonded_potential(int ftype)
3658 (interaction_function[ftype].flags & IF_BOND) &&
3659 !(ftype == F_CONNBONDS || ftype == F_POSRES) &&
3660 (ftype < F_GB12 || ftype > F_GB14);
3663 static void divide_bondeds_over_threads(t_idef *idef, int nthreads)
3670 idef->nthreads = nthreads;
3672 if (F_NRE*(nthreads+1) > idef->il_thread_division_nalloc)
3674 idef->il_thread_division_nalloc = F_NRE*(nthreads+1);
3675 snew(idef->il_thread_division, idef->il_thread_division_nalloc);
3678 for (ftype = 0; ftype < F_NRE; ftype++)
3680 if (ftype_is_bonded_potential(ftype))
3682 nat1 = interaction_function[ftype].nratoms + 1;
3684 for (t = 0; t <= nthreads; t++)
3686 /* Divide the interactions equally over the threads.
3687 * When the different types of bonded interactions
3688 * are distributed roughly equally over the threads,
3689 * this should lead to well localized output into
3690 * the force buffer on each thread.
3691 * If this is not the case, a more advanced scheme
3692 * (not implemented yet) will do better.
3694 il_nr_thread = (((idef->il[ftype].nr/nat1)*t)/nthreads)*nat1;
3696 /* Ensure that distance restraint pairs with the same label
3697 * end up on the same thread.
3698 * This is slighlty tricky code, since the next for iteration
3699 * may have an initial il_nr_thread lower than the final value
3700 * in the previous iteration, but this will anyhow be increased
3701 * to the approriate value again by this while loop.
3703 while (ftype == F_DISRES &&
3705 il_nr_thread < idef->il[ftype].nr &&
3706 idef->iparams[idef->il[ftype].iatoms[il_nr_thread]].disres.label ==
3707 idef->iparams[idef->il[ftype].iatoms[il_nr_thread-nat1]].disres.label)
3709 il_nr_thread += nat1;
3712 idef->il_thread_division[ftype*(nthreads+1)+t] = il_nr_thread;
3719 calc_bonded_reduction_mask(const t_idef *idef,
3724 int ftype, nb, nat1, nb0, nb1, i, a;
3728 for (ftype = 0; ftype < F_NRE; ftype++)
3730 if (ftype_is_bonded_potential(ftype))
3732 nb = idef->il[ftype].nr;
3735 nat1 = interaction_function[ftype].nratoms + 1;
3737 /* Divide this interaction equally over the threads.
3738 * This is not stored: should match division in calc_bonds.
3740 nb0 = idef->il_thread_division[ftype*(nt+1)+t];
3741 nb1 = idef->il_thread_division[ftype*(nt+1)+t+1];
3743 for (i = nb0; i < nb1; i += nat1)
3745 for (a = 1; a < nat1; a++)
3747 mask |= (1U << (idef->il[ftype].iatoms[i+a]>>shift));
3757 void setup_bonded_threading(t_forcerec *fr, t_idef *idef)
3759 #define MAX_BLOCK_BITS 32
3764 assert(fr->nthreads >= 1);
3767 /* Divide the bonded interaction over the threads */
3768 divide_bondeds_over_threads(idef, fr->nthreads);
3770 if (fr->nthreads == 1)
3777 /* We divide the force array in a maximum of 32 blocks.
3778 * Minimum force block reduction size is 2^6=64.
3781 while (fr->natoms_force > (int)(MAX_BLOCK_BITS*(1U<<fr->red_ashift)))
3787 fprintf(debug, "bonded force buffer block atom shift %d bits\n",
3791 /* Determine to which blocks each thread's bonded force calculation
3792 * contributes. Store this is a mask for each thread.
3794 #pragma omp parallel for num_threads(fr->nthreads) schedule(static)
3795 for (t = 1; t < fr->nthreads; t++)
3797 fr->f_t[t].red_mask =
3798 calc_bonded_reduction_mask(idef, fr->red_ashift, t, fr->nthreads);
3801 /* Determine the maximum number of blocks we need to reduce over */
3804 for (t = 0; t < fr->nthreads; t++)
3807 for (b = 0; b < MAX_BLOCK_BITS; b++)
3809 if (fr->f_t[t].red_mask & (1U<<b))
3811 fr->red_nblock = max(fr->red_nblock, b+1);
3817 fprintf(debug, "thread %d flags %x count %d\n",
3818 t, fr->f_t[t].red_mask, c);
3824 fprintf(debug, "Number of blocks to reduce: %d of size %d\n",
3825 fr->red_nblock, 1<<fr->red_ashift);
3826 fprintf(debug, "Reduction density %.2f density/#thread %.2f\n",
3827 ctot*(1<<fr->red_ashift)/(double)fr->natoms_force,
3828 ctot*(1<<fr->red_ashift)/(double)(fr->natoms_force*fr->nthreads));
3832 static void zero_thread_forces(f_thread_t *f_t, int n,
3833 int nblock, int blocksize)
3835 int b, a0, a1, a, i, j;
3837 if (n > f_t->f_nalloc)
3839 f_t->f_nalloc = over_alloc_large(n);
3840 srenew(f_t->f, f_t->f_nalloc);
3843 if (f_t->red_mask != 0)
3845 for (b = 0; b < nblock; b++)
3847 if (f_t->red_mask && (1U<<b))
3850 a1 = min((b+1)*blocksize, n);
3851 for (a = a0; a < a1; a++)
3853 clear_rvec(f_t->f[a]);
3858 for (i = 0; i < SHIFTS; i++)
3860 clear_rvec(f_t->fshift[i]);
3862 for (i = 0; i < F_NRE; i++)
3866 for (i = 0; i < egNR; i++)
3868 for (j = 0; j < f_t->grpp.nener; j++)
3870 f_t->grpp.ener[i][j] = 0;
3873 for (i = 0; i < efptNR; i++)
3879 static void reduce_thread_force_buffer(int n, rvec *f,
3880 int nthreads, f_thread_t *f_t,
3881 int nblock, int block_size)
3883 /* The max thread number is arbitrary,
3884 * we used a fixed number to avoid memory management.
3885 * Using more than 16 threads is probably never useful performance wise.
3887 #define MAX_BONDED_THREADS 256
3890 if (nthreads > MAX_BONDED_THREADS)
3892 gmx_fatal(FARGS, "Can not reduce bonded forces on more than %d threads",
3893 MAX_BONDED_THREADS);
3896 /* This reduction can run on any number of threads,
3897 * independently of nthreads.
3899 #pragma omp parallel for num_threads(nthreads) schedule(static)
3900 for (b = 0; b < nblock; b++)
3902 rvec *fp[MAX_BONDED_THREADS];
3906 /* Determine which threads contribute to this block */
3908 for (ft = 1; ft < nthreads; ft++)
3910 if (f_t[ft].red_mask & (1U<<b))
3912 fp[nfb++] = f_t[ft].f;
3917 /* Reduce force buffers for threads that contribute */
3919 a1 = (b+1)*block_size;
3921 for (a = a0; a < a1; a++)
3923 for (fb = 0; fb < nfb; fb++)
3925 rvec_inc(f[a], fp[fb][a]);
3932 static void reduce_thread_forces(int n, rvec *f, rvec *fshift,
3933 real *ener, gmx_grppairener_t *grpp, real *dvdl,
3934 int nthreads, f_thread_t *f_t,
3935 int nblock, int block_size,
3936 gmx_bool bCalcEnerVir,
3941 /* Reduce the bonded force buffer */
3942 reduce_thread_force_buffer(n, f, nthreads, f_t, nblock, block_size);
3945 /* When necessary, reduce energy and virial using one thread only */
3950 for (i = 0; i < SHIFTS; i++)
3952 for (t = 1; t < nthreads; t++)
3954 rvec_inc(fshift[i], f_t[t].fshift[i]);
3957 for (i = 0; i < F_NRE; i++)
3959 for (t = 1; t < nthreads; t++)
3961 ener[i] += f_t[t].ener[i];
3964 for (i = 0; i < egNR; i++)
3966 for (j = 0; j < f_t[1].grpp.nener; j++)
3968 for (t = 1; t < nthreads; t++)
3971 grpp->ener[i][j] += f_t[t].grpp.ener[i][j];
3977 for (i = 0; i < efptNR; i++)
3980 for (t = 1; t < nthreads; t++)
3982 dvdl[i] += f_t[t].dvdl[i];
3989 static real calc_one_bond(FILE *fplog, int thread,
3990 int ftype, const t_idef *idef,
3991 rvec x[], rvec f[], rvec fshift[],
3993 const t_pbc *pbc, const t_graph *g,
3994 gmx_grppairener_t *grpp,
3996 real *lambda, real *dvdl,
3997 const t_mdatoms *md, t_fcdata *fcd,
3998 gmx_bool bCalcEnerVir,
3999 int *global_atom_index, gmx_bool bPrintSepPot)
4001 int nat1, nbonds, efptFTYPE;
4006 if (IS_RESTRAINT_TYPE(ftype))
4008 efptFTYPE = efptRESTRAINT;
4012 efptFTYPE = efptBONDED;
4015 nat1 = interaction_function[ftype].nratoms + 1;
4016 nbonds = idef->il[ftype].nr/nat1;
4017 iatoms = idef->il[ftype].iatoms;
4019 nb0 = idef->il_thread_division[ftype*(idef->nthreads+1)+thread];
4020 nbn = idef->il_thread_division[ftype*(idef->nthreads+1)+thread+1] - nb0;
4022 if (!IS_LISTED_LJ_C(ftype))
4024 if (ftype == F_CMAP)
4026 v = cmap_dihs(nbn, iatoms+nb0,
4027 idef->iparams, &idef->cmap_grid,
4028 (const rvec*)x, f, fshift,
4029 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
4030 md, fcd, global_atom_index);
4033 else if (ftype == F_ANGLES &&
4034 !bCalcEnerVir && fr->efep == efepNO)
4036 /* No energies, shift forces, dvdl */
4037 angles_noener_simd(nbn, idef->il[ftype].iatoms+nb0,
4040 pbc, g, lambda[efptFTYPE], md, fcd,
4045 else if (ftype == F_PDIHS &&
4046 !bCalcEnerVir && fr->efep == efepNO)
4048 /* No energies, shift forces, dvdl */
4049 #ifndef SIMD_BONDEDS
4054 (nbn, idef->il[ftype].iatoms+nb0,
4057 pbc, g, lambda[efptFTYPE], md, fcd,
4063 v = interaction_function[ftype].ifunc(nbn, iatoms+nb0,
4065 (const rvec*)x, f, fshift,
4066 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
4067 md, fcd, global_atom_index);
4071 fprintf(fplog, " %-23s #%4d V %12.5e dVdl %12.5e\n",
4072 interaction_function[ftype].longname,
4073 nbonds, v, lambda[efptFTYPE]);
4078 v = do_nonbonded_listed(ftype, nbn, iatoms+nb0, idef->iparams, (const rvec*)x, f, fshift,
4079 pbc, g, lambda, dvdl, md, fr, grpp, global_atom_index);
4083 fprintf(fplog, " %-5s + %-15s #%4d dVdl %12.5e\n",
4084 interaction_function[ftype].longname,
4085 interaction_function[F_LJ14].longname, nbonds, dvdl[efptVDW]);
4086 fprintf(fplog, " %-5s + %-15s #%4d dVdl %12.5e\n",
4087 interaction_function[ftype].longname,
4088 interaction_function[F_COUL14].longname, nbonds, dvdl[efptCOUL]);
4094 inc_nrnb(nrnb, interaction_function[ftype].nrnb_ind, nbonds);
4100 void calc_bonds(FILE *fplog, const gmx_multisim_t *ms,
4102 rvec x[], history_t *hist,
4103 rvec f[], t_forcerec *fr,
4104 const t_pbc *pbc, const t_graph *g,
4105 gmx_enerdata_t *enerd, t_nrnb *nrnb,
4107 const t_mdatoms *md,
4108 t_fcdata *fcd, int *global_atom_index,
4109 t_atomtypes *atype, gmx_genborn_t *born,
4111 gmx_bool bPrintSepPot, gmx_large_int_t step)
4113 gmx_bool bCalcEnerVir;
4115 real v, dvdl[efptNR], dvdl_dum[efptNR]; /* The dummy array is to have a place to store the dhdl at other values
4116 of lambda, which will be thrown away in the end*/
4117 const t_pbc *pbc_null;
4122 assert(fr->nthreads == idef->nthreads);
4125 bCalcEnerVir = (force_flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY));
4127 for (i = 0; i < efptNR; i++)
4141 fprintf(fplog, "Step %s: bonded V and dVdl for this node\n",
4142 gmx_step_str(step, buf));
4148 p_graph(debug, "Bondage is fun", g);
4152 /* Do pre force calculation stuff which might require communication */
4153 if (idef->il[F_ORIRES].nr)
4155 enerd->term[F_ORIRESDEV] =
4156 calc_orires_dev(ms, idef->il[F_ORIRES].nr,
4157 idef->il[F_ORIRES].iatoms,
4158 idef->iparams, md, (const rvec*)x,
4159 pbc_null, fcd, hist);
4161 if (idef->il[F_DISRES].nr)
4163 calc_disres_R_6(ms, idef->il[F_DISRES].nr,
4164 idef->il[F_DISRES].iatoms,
4165 idef->iparams, (const rvec*)x, pbc_null,
4169 #pragma omp parallel for num_threads(fr->nthreads) schedule(static)
4170 for (thread = 0; thread < fr->nthreads; thread++)
4177 gmx_grppairener_t *grpp;
4182 fshift = fr->fshift;
4184 grpp = &enerd->grpp;
4189 zero_thread_forces(&fr->f_t[thread], fr->natoms_force,
4190 fr->red_nblock, 1<<fr->red_ashift);
4192 ft = fr->f_t[thread].f;
4193 fshift = fr->f_t[thread].fshift;
4194 epot = fr->f_t[thread].ener;
4195 grpp = &fr->f_t[thread].grpp;
4196 dvdlt = fr->f_t[thread].dvdl;
4198 /* Loop over all bonded force types to calculate the bonded forces */
4199 for (ftype = 0; (ftype < F_NRE); ftype++)
4201 if (idef->il[ftype].nr > 0 && ftype_is_bonded_potential(ftype))
4203 v = calc_one_bond(fplog, thread, ftype, idef, x,
4204 ft, fshift, fr, pbc_null, g, grpp,
4205 nrnb, lambda, dvdlt,
4206 md, fcd, bCalcEnerVir,
4207 global_atom_index, bPrintSepPot);
4212 if (fr->nthreads > 1)
4214 reduce_thread_forces(fr->natoms_force, f, fr->fshift,
4215 enerd->term, &enerd->grpp, dvdl,
4216 fr->nthreads, fr->f_t,
4217 fr->red_nblock, 1<<fr->red_ashift,
4219 force_flags & GMX_FORCE_DHDL);
4221 if (force_flags & GMX_FORCE_DHDL)
4223 for (i = 0; i < efptNR; i++)
4225 enerd->dvdl_nonlin[i] += dvdl[i];
4229 /* Copy the sum of violations for the distance restraints from fcd */
4232 enerd->term[F_DISRESVIOL] = fcd->disres.sumviol;
4237 void calc_bonds_lambda(FILE *fplog,
4241 const t_pbc *pbc, const t_graph *g,
4242 gmx_grppairener_t *grpp, real *epot, t_nrnb *nrnb,
4244 const t_mdatoms *md,
4246 int *global_atom_index)
4248 int i, ftype, nr_nonperturbed, nr;
4250 real dvdl_dum[efptNR];
4252 const t_pbc *pbc_null;
4264 /* Copy the whole idef, so we can modify the contents locally */
4266 idef_fe.nthreads = 1;
4267 snew(idef_fe.il_thread_division, F_NRE*(idef_fe.nthreads+1));
4269 /* We already have the forces, so we use temp buffers here */
4270 snew(f, fr->natoms_force);
4271 snew(fshift, SHIFTS);
4273 /* Loop over all bonded force types to calculate the bonded energies */
4274 for (ftype = 0; (ftype < F_NRE); ftype++)
4276 if (ftype_is_bonded_potential(ftype))
4278 /* Set the work range of thread 0 to the perturbed bondeds only */
4279 nr_nonperturbed = idef->il[ftype].nr_nonperturbed;
4280 nr = idef->il[ftype].nr;
4281 idef_fe.il_thread_division[ftype*2+0] = nr_nonperturbed;
4282 idef_fe.il_thread_division[ftype*2+1] = nr;
4284 /* This is only to get the flop count correct */
4285 idef_fe.il[ftype].nr = nr - nr_nonperturbed;
4287 if (nr - nr_nonperturbed > 0)
4289 v = calc_one_bond(fplog, 0, ftype, &idef_fe,
4290 x, f, fshift, fr, pbc_null, g,
4291 grpp, nrnb, lambda, dvdl_dum,
4293 global_atom_index, FALSE);
4302 sfree(idef_fe.il_thread_division);