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, ki;
741 rvec dOH1, dOH2, dHH, dOD, dDS, nW, kk, dx, kdx, proj;
745 real vtot, fij, r_HH, r_OD, r_nW, tx, ty, tz, qS;
750 type0 = forceatoms[0];
752 qS = md->chargeA[aS];
753 kk[XX] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_x;
754 kk[YY] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_y;
755 kk[ZZ] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_z;
756 r_HH = 1.0/forceparams[type0].wpol.rHH;
757 r_OD = 1.0/forceparams[type0].wpol.rOD;
760 fprintf(debug, "WPOL: qS = %10.5f aS = %5d\n", qS, aS);
761 fprintf(debug, "WPOL: kk = %10.3f %10.3f %10.3f\n",
762 kk[XX], kk[YY], kk[ZZ]);
763 fprintf(debug, "WPOL: rOH = %10.3f rHH = %10.3f rOD = %10.3f\n",
764 forceparams[type0].wpol.rOH,
765 forceparams[type0].wpol.rHH,
766 forceparams[type0].wpol.rOD);
768 for (i = 0; (i < nbonds); i += 6)
770 type = forceatoms[i];
773 gmx_fatal(FARGS, "Sorry, type = %d, type0 = %d, file = %s, line = %d",
774 type, type0, __FILE__, __LINE__);
776 aO = forceatoms[i+1];
777 aH1 = forceatoms[i+2];
778 aH2 = forceatoms[i+3];
779 aD = forceatoms[i+4];
780 aS = forceatoms[i+5];
782 /* Compute vectors describing the water frame */
783 pbc_rvec_sub(pbc, x[aH1], x[aO], dOH1);
784 pbc_rvec_sub(pbc, x[aH2], x[aO], dOH2);
785 pbc_rvec_sub(pbc, x[aH2], x[aH1], dHH);
786 pbc_rvec_sub(pbc, x[aD], x[aO], dOD);
787 ki = pbc_rvec_sub(pbc, x[aS], x[aD], dDS);
788 cprod(dOH1, dOH2, nW);
790 /* Compute inverse length of normal vector
791 * (this one could be precomputed, but I'm too lazy now)
793 r_nW = gmx_invsqrt(iprod(nW, nW));
794 /* This is for precision, but does not make a big difference,
797 r_OD = gmx_invsqrt(iprod(dOD, dOD));
799 /* Normalize the vectors in the water frame */
801 svmul(r_HH, dHH, dHH);
802 svmul(r_OD, dOD, dOD);
804 /* Compute displacement of shell along components of the vector */
805 dx[ZZ] = iprod(dDS, dOD);
806 /* Compute projection on the XY plane: dDS - dx[ZZ]*dOD */
807 for (m = 0; (m < DIM); m++)
809 proj[m] = dDS[m]-dx[ZZ]*dOD[m];
812 /*dx[XX] = iprod(dDS,nW);
813 dx[YY] = iprod(dDS,dHH);*/
814 dx[XX] = iprod(proj, nW);
815 for (m = 0; (m < DIM); m++)
817 proj[m] -= dx[XX]*nW[m];
819 dx[YY] = iprod(proj, dHH);
824 fprintf(debug, "WPOL: dx2=%10g dy2=%10g dz2=%10g sum=%10g dDS^2=%10g\n",
825 sqr(dx[XX]), sqr(dx[YY]), sqr(dx[ZZ]), iprod(dx, dx), iprod(dDS, dDS));
826 fprintf(debug, "WPOL: dHH=(%10g,%10g,%10g)\n", dHH[XX], dHH[YY], dHH[ZZ]);
827 fprintf(debug, "WPOL: dOD=(%10g,%10g,%10g), 1/r_OD = %10g\n",
828 dOD[XX], dOD[YY], dOD[ZZ], 1/r_OD);
829 fprintf(debug, "WPOL: nW =(%10g,%10g,%10g), 1/r_nW = %10g\n",
830 nW[XX], nW[YY], nW[ZZ], 1/r_nW);
831 fprintf(debug, "WPOL: dx =%10g, dy =%10g, dz =%10g\n",
832 dx[XX], dx[YY], dx[ZZ]);
833 fprintf(debug, "WPOL: dDSx=%10g, dDSy=%10g, dDSz=%10g\n",
834 dDS[XX], dDS[YY], dDS[ZZ]);
837 /* Now compute the forces and energy */
838 kdx[XX] = kk[XX]*dx[XX];
839 kdx[YY] = kk[YY]*dx[YY];
840 kdx[ZZ] = kk[ZZ]*dx[ZZ];
841 vtot += iprod(dx, kdx);
845 ivec_sub(SHIFT_IVEC(g, aS), SHIFT_IVEC(g, aD), dt);
849 for (m = 0; (m < DIM); m++)
851 /* This is a tensor operation but written out for speed */
861 fshift[ki][m] += fij;
862 fshift[CENTRAL][m] -= fij;
867 fprintf(debug, "WPOL: vwpol=%g\n", 0.5*iprod(dx, kdx));
868 fprintf(debug, "WPOL: df = (%10g, %10g, %10g)\n", df[XX], df[YY], df[ZZ]);
876 static real do_1_thole(const rvec xi, const rvec xj, rvec fi, rvec fj,
877 const t_pbc *pbc, real qq,
878 rvec fshift[], real afac)
881 real r12sq, r12_1, r12n, r12bar, v0, v1, fscal, ebar, fff;
884 t = pbc_rvec_sub(pbc, xi, xj, r12); /* 3 */
886 r12sq = iprod(r12, r12); /* 5 */
887 r12_1 = gmx_invsqrt(r12sq); /* 5 */
888 r12bar = afac/r12_1; /* 5 */
889 v0 = qq*ONE_4PI_EPS0*r12_1; /* 2 */
890 ebar = exp(-r12bar); /* 5 */
891 v1 = (1-(1+0.5*r12bar)*ebar); /* 4 */
892 fscal = ((v0*r12_1)*v1 - v0*0.5*afac*ebar*(r12bar+1))*r12_1; /* 9 */
895 fprintf(debug, "THOLE: v0 = %.3f v1 = %.3f r12= % .3f r12bar = %.3f fscal = %.3f ebar = %.3f\n", v0, v1, 1/r12_1, r12bar, fscal, ebar);
898 for (m = 0; (m < DIM); m++)
904 fshift[CENTRAL][m] -= fff;
907 return v0*v1; /* 1 */
911 real thole_pol(int nbonds,
912 const t_iatom forceatoms[], const t_iparams forceparams[],
913 const rvec x[], rvec f[], rvec fshift[],
914 const t_pbc *pbc, const t_graph *g,
915 real lambda, real *dvdlambda,
916 const t_mdatoms *md, t_fcdata *fcd,
917 int *global_atom_index)
919 /* Interaction between two pairs of particles with opposite charge */
920 int i, type, a1, da1, a2, da2;
921 real q1, q2, qq, a, al1, al2, afac;
924 for (i = 0; (i < nbonds); )
926 type = forceatoms[i++];
927 a1 = forceatoms[i++];
928 da1 = forceatoms[i++];
929 a2 = forceatoms[i++];
930 da2 = forceatoms[i++];
931 q1 = md->chargeA[da1];
932 q2 = md->chargeA[da2];
933 a = forceparams[type].thole.a;
934 al1 = forceparams[type].thole.alpha1;
935 al2 = forceparams[type].thole.alpha2;
937 afac = a*pow(al1*al2, -1.0/6.0);
938 V += do_1_thole(x[a1], x[a2], f[a1], f[a2], pbc, qq, fshift, afac);
939 V += do_1_thole(x[da1], x[a2], f[da1], f[a2], pbc, -qq, fshift, afac);
940 V += do_1_thole(x[a1], x[da2], f[a1], f[da2], pbc, -qq, fshift, afac);
941 V += do_1_thole(x[da1], x[da2], f[da1], f[da2], pbc, qq, fshift, afac);
947 real bond_angle(const rvec xi, const rvec xj, const rvec xk, const t_pbc *pbc,
948 rvec r_ij, rvec r_kj, real *costh,
950 /* Return value is the angle between the bonds i-j and j-k */
955 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
956 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
958 *costh = cos_angle(r_ij, r_kj); /* 25 */
959 th = acos(*costh); /* 10 */
964 real angles(int nbonds,
965 const t_iatom forceatoms[], const t_iparams forceparams[],
966 const rvec x[], rvec f[], rvec fshift[],
967 const t_pbc *pbc, const t_graph *g,
968 real lambda, real *dvdlambda,
969 const t_mdatoms *md, t_fcdata *fcd,
970 int *global_atom_index)
972 int i, ai, aj, ak, t1, t2, type;
974 real cos_theta, cos_theta2, theta, dVdt, va, vtot;
975 ivec jt, dt_ij, dt_kj;
978 for (i = 0; i < nbonds; )
980 type = forceatoms[i++];
981 ai = forceatoms[i++];
982 aj = forceatoms[i++];
983 ak = forceatoms[i++];
985 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
986 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
988 *dvdlambda += harmonic(forceparams[type].harmonic.krA,
989 forceparams[type].harmonic.krB,
990 forceparams[type].harmonic.rA*DEG2RAD,
991 forceparams[type].harmonic.rB*DEG2RAD,
992 theta, lambda, &va, &dVdt); /* 21 */
995 cos_theta2 = sqr(cos_theta);
1002 real nrkj_1, nrij_1;
1005 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1006 sth = st*cos_theta; /* 1 */
1010 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1011 theta*RAD2DEG, va, dVdt);
1014 nrij2 = iprod(r_ij, r_ij); /* 5 */
1015 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1017 nrij_1 = gmx_invsqrt(nrij2); /* 10 */
1018 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1020 cik = st*nrij_1*nrkj_1; /* 2 */
1021 cii = sth*nrij_1*nrij_1; /* 2 */
1022 ckk = sth*nrkj_1*nrkj_1; /* 2 */
1024 for (m = 0; m < DIM; m++)
1026 f_i[m] = -(cik*r_kj[m] - cii*r_ij[m]);
1027 f_k[m] = -(cik*r_ij[m] - ckk*r_kj[m]);
1028 f_j[m] = -f_i[m] - f_k[m];
1035 copy_ivec(SHIFT_IVEC(g, aj), jt);
1037 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1038 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1039 t1 = IVEC2IS(dt_ij);
1040 t2 = IVEC2IS(dt_kj);
1042 rvec_inc(fshift[t1], f_i);
1043 rvec_inc(fshift[CENTRAL], f_j);
1044 rvec_inc(fshift[t2], f_k);
1053 /* As angles, but using SIMD to calculate many dihedrals at once.
1054 * This routines does not calculate energies and shift forces.
1056 static gmx_inline void
1057 angles_noener_simd(int nbonds,
1058 const t_iatom forceatoms[], const t_iparams forceparams[],
1059 const rvec x[], rvec f[],
1060 const t_pbc *pbc, const t_graph *g,
1062 const t_mdatoms *md, t_fcdata *fcd,
1063 int *global_atom_index)
1065 #define UNROLL GMX_SIMD_WIDTH_HERE
1068 int type, ai[UNROLL], aj[UNROLL], ak[UNROLL];
1069 real coeff_array[2*UNROLL+UNROLL], *coeff;
1070 real dr_array[2*DIM*UNROLL+UNROLL], *dr;
1071 real f_buf_array[6*UNROLL+UNROLL], *f_buf;
1072 gmx_mm_pr k_S, theta0_S;
1073 gmx_mm_pr rijx_S, rijy_S, rijz_S;
1074 gmx_mm_pr rkjx_S, rkjy_S, rkjz_S;
1076 gmx_mm_pr min_one_plus_eps_S;
1077 gmx_mm_pr rij_rkj_S;
1078 gmx_mm_pr nrij2_S, nrij_1_S;
1079 gmx_mm_pr nrkj2_S, nrkj_1_S;
1080 gmx_mm_pr cos_S, invsin_S;
1082 gmx_mm_pr st_S, sth_S;
1083 gmx_mm_pr cik_S, cii_S, ckk_S;
1084 gmx_mm_pr f_ix_S, f_iy_S, f_iz_S;
1085 gmx_mm_pr f_kx_S, f_ky_S, f_kz_S;
1086 pbc_simd_t pbc_simd;
1088 /* Ensure register memory alignment */
1089 coeff = gmx_simd_align_real(coeff_array);
1090 dr = gmx_simd_align_real(dr_array);
1091 f_buf = gmx_simd_align_real(f_buf_array);
1093 set_pbc_simd(pbc,&pbc_simd);
1095 one_S = gmx_set1_pr(1.0);
1097 /* The smallest number > -1 */
1098 min_one_plus_eps_S = gmx_set1_pr(-1.0 + 2*GMX_REAL_EPS);
1100 /* nbonds is the number of angles times nfa1, here we step UNROLL angles */
1101 for (i = 0; (i < nbonds); i += UNROLL*nfa1)
1103 /* Collect atoms for UNROLL angles.
1104 * iu indexes into forceatoms, we should not let iu go beyond nbonds.
1107 for (s = 0; s < UNROLL; s++)
1109 type = forceatoms[iu];
1110 ai[s] = forceatoms[iu+1];
1111 aj[s] = forceatoms[iu+2];
1112 ak[s] = forceatoms[iu+3];
1114 coeff[s] = forceparams[type].harmonic.krA;
1115 coeff[UNROLL+s] = forceparams[type].harmonic.rA*DEG2RAD;
1117 /* If you can't use pbc_dx_simd below for PBC, e.g. because
1118 * you can't round in SIMD, use pbc_rvec_sub here.
1120 /* Store the non PBC corrected distances packed and aligned */
1121 for (m = 0; m < DIM; m++)
1123 dr[s + m *UNROLL] = x[ai[s]][m] - x[aj[s]][m];
1124 dr[s + (DIM+m)*UNROLL] = x[ak[s]][m] - x[aj[s]][m];
1127 /* At the end fill the arrays with identical entries */
1128 if (iu + nfa1 < nbonds)
1134 k_S = gmx_load_pr(coeff);
1135 theta0_S = gmx_load_pr(coeff+UNROLL);
1137 rijx_S = gmx_load_pr(dr + 0*UNROLL);
1138 rijy_S = gmx_load_pr(dr + 1*UNROLL);
1139 rijz_S = gmx_load_pr(dr + 2*UNROLL);
1140 rkjx_S = gmx_load_pr(dr + 3*UNROLL);
1141 rkjy_S = gmx_load_pr(dr + 4*UNROLL);
1142 rkjz_S = gmx_load_pr(dr + 5*UNROLL);
1144 pbc_dx_simd(&rijx_S, &rijy_S, &rijz_S, &pbc_simd);
1145 pbc_dx_simd(&rkjx_S, &rkjy_S, &rkjz_S, &pbc_simd);
1147 rij_rkj_S = gmx_iprod_pr(rijx_S, rijy_S, rijz_S,
1148 rkjx_S, rkjy_S, rkjz_S);
1150 nrij2_S = gmx_norm2_pr(rijx_S, rijy_S, rijz_S);
1151 nrkj2_S = gmx_norm2_pr(rkjx_S, rkjy_S, rkjz_S);
1153 nrij_1_S = gmx_invsqrt_pr(nrij2_S);
1154 nrkj_1_S = gmx_invsqrt_pr(nrkj2_S);
1156 cos_S = gmx_mul_pr(rij_rkj_S, gmx_mul_pr(nrij_1_S, nrkj_1_S));
1158 /* To allow for 180 degrees, we take the max of cos and -1 + 1bit,
1159 * so we can safely get the 1/sin from 1/sqrt(1 - cos^2).
1160 * This also ensures that rounding errors would cause the argument
1161 * of gmx_acos_pr to be < -1.
1162 * Note that we do not take precautions for cos(0)=1, so the outer
1163 * atoms in an angle should not be on top of each other.
1165 cos_S = gmx_max_pr(cos_S, min_one_plus_eps_S);
1167 theta_S = gmx_acos_pr(cos_S);
1169 invsin_S = gmx_invsqrt_pr(gmx_sub_pr(one_S, gmx_mul_pr(cos_S, cos_S)));
1171 st_S = gmx_mul_pr(gmx_mul_pr(k_S, gmx_sub_pr(theta0_S, theta_S)),
1173 sth_S = gmx_mul_pr(st_S, cos_S);
1175 cik_S = gmx_mul_pr(st_S, gmx_mul_pr(nrij_1_S, nrkj_1_S));
1176 cii_S = gmx_mul_pr(sth_S, gmx_mul_pr(nrij_1_S, nrij_1_S));
1177 ckk_S = gmx_mul_pr(sth_S, gmx_mul_pr(nrkj_1_S, nrkj_1_S));
1179 f_ix_S = gmx_mul_pr(cii_S, rijx_S);
1180 f_ix_S = gmx_nmsub_pr(cik_S, rkjx_S, f_ix_S);
1181 f_iy_S = gmx_mul_pr(cii_S, rijy_S);
1182 f_iy_S = gmx_nmsub_pr(cik_S, rkjy_S, f_iy_S);
1183 f_iz_S = gmx_mul_pr(cii_S, rijz_S);
1184 f_iz_S = gmx_nmsub_pr(cik_S, rkjz_S, f_iz_S);
1185 f_kx_S = gmx_mul_pr(ckk_S, rkjx_S);
1186 f_kx_S = gmx_nmsub_pr(cik_S, rijx_S, f_kx_S);
1187 f_ky_S = gmx_mul_pr(ckk_S, rkjy_S);
1188 f_ky_S = gmx_nmsub_pr(cik_S, rijy_S, f_ky_S);
1189 f_kz_S = gmx_mul_pr(ckk_S, rkjz_S);
1190 f_kz_S = gmx_nmsub_pr(cik_S, rijz_S, f_kz_S);
1192 gmx_store_pr(f_buf + 0*UNROLL, f_ix_S);
1193 gmx_store_pr(f_buf + 1*UNROLL, f_iy_S);
1194 gmx_store_pr(f_buf + 2*UNROLL, f_iz_S);
1195 gmx_store_pr(f_buf + 3*UNROLL, f_kx_S);
1196 gmx_store_pr(f_buf + 4*UNROLL, f_ky_S);
1197 gmx_store_pr(f_buf + 5*UNROLL, f_kz_S);
1203 for (m = 0; m < DIM; m++)
1205 f[ai[s]][m] += f_buf[s + m*UNROLL];
1206 f[aj[s]][m] -= f_buf[s + m*UNROLL] + f_buf[s + (DIM+m)*UNROLL];
1207 f[ak[s]][m] += f_buf[s + (DIM+m)*UNROLL];
1212 while (s < UNROLL && iu < nbonds);
1217 #endif /* SIMD_BONDEDS */
1219 real linear_angles(int nbonds,
1220 const t_iatom forceatoms[], const t_iparams forceparams[],
1221 const rvec x[], rvec f[], rvec fshift[],
1222 const t_pbc *pbc, const t_graph *g,
1223 real lambda, real *dvdlambda,
1224 const t_mdatoms *md, t_fcdata *fcd,
1225 int *global_atom_index)
1227 int i, m, ai, aj, ak, t1, t2, type;
1229 real L1, kA, kB, aA, aB, dr, dr2, va, vtot, a, b, klin;
1230 ivec jt, dt_ij, dt_kj;
1231 rvec r_ij, r_kj, r_ik, dx;
1235 for (i = 0; (i < nbonds); )
1237 type = forceatoms[i++];
1238 ai = forceatoms[i++];
1239 aj = forceatoms[i++];
1240 ak = forceatoms[i++];
1242 kA = forceparams[type].linangle.klinA;
1243 kB = forceparams[type].linangle.klinB;
1244 klin = L1*kA + lambda*kB;
1246 aA = forceparams[type].linangle.aA;
1247 aB = forceparams[type].linangle.aB;
1248 a = L1*aA+lambda*aB;
1251 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
1252 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
1253 rvec_sub(r_ij, r_kj, r_ik);
1256 for (m = 0; (m < DIM); m++)
1258 dr = -a * r_ij[m] - b * r_kj[m];
1263 f_j[m] = -(f_i[m]+f_k[m]);
1269 *dvdlambda += 0.5*(kB-kA)*dr2 + klin*(aB-aA)*iprod(dx, r_ik);
1275 copy_ivec(SHIFT_IVEC(g, aj), jt);
1277 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1278 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1279 t1 = IVEC2IS(dt_ij);
1280 t2 = IVEC2IS(dt_kj);
1282 rvec_inc(fshift[t1], f_i);
1283 rvec_inc(fshift[CENTRAL], f_j);
1284 rvec_inc(fshift[t2], f_k);
1289 real urey_bradley(int nbonds,
1290 const t_iatom forceatoms[], const t_iparams forceparams[],
1291 const rvec x[], rvec f[], rvec fshift[],
1292 const t_pbc *pbc, const t_graph *g,
1293 real lambda, real *dvdlambda,
1294 const t_mdatoms *md, t_fcdata *fcd,
1295 int *global_atom_index)
1297 int i, m, ai, aj, ak, t1, t2, type, ki;
1298 rvec r_ij, r_kj, r_ik;
1299 real cos_theta, cos_theta2, theta;
1300 real dVdt, va, vtot, dr, dr2, vbond, fbond, fik;
1301 real kthA, th0A, kUBA, r13A, kthB, th0B, kUBB, r13B;
1302 ivec jt, dt_ij, dt_kj, dt_ik;
1305 for (i = 0; (i < nbonds); )
1307 type = forceatoms[i++];
1308 ai = forceatoms[i++];
1309 aj = forceatoms[i++];
1310 ak = forceatoms[i++];
1311 th0A = forceparams[type].u_b.thetaA*DEG2RAD;
1312 kthA = forceparams[type].u_b.kthetaA;
1313 r13A = forceparams[type].u_b.r13A;
1314 kUBA = forceparams[type].u_b.kUBA;
1315 th0B = forceparams[type].u_b.thetaB*DEG2RAD;
1316 kthB = forceparams[type].u_b.kthetaB;
1317 r13B = forceparams[type].u_b.r13B;
1318 kUBB = forceparams[type].u_b.kUBB;
1320 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
1321 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
1323 *dvdlambda += harmonic(kthA, kthB, th0A, th0B, theta, lambda, &va, &dVdt); /* 21 */
1326 ki = pbc_rvec_sub(pbc, x[ai], x[ak], r_ik); /* 3 */
1327 dr2 = iprod(r_ik, r_ik); /* 5 */
1328 dr = dr2*gmx_invsqrt(dr2); /* 10 */
1330 *dvdlambda += harmonic(kUBA, kUBB, r13A, r13B, dr, lambda, &vbond, &fbond); /* 19 */
1332 cos_theta2 = sqr(cos_theta); /* 1 */
1340 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1341 sth = st*cos_theta; /* 1 */
1345 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1346 theta*RAD2DEG, va, dVdt);
1349 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1350 nrij2 = iprod(r_ij, r_ij);
1352 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
1353 cii = sth/nrij2; /* 10 */
1354 ckk = sth/nrkj2; /* 10 */
1356 for (m = 0; (m < DIM); m++) /* 39 */
1358 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
1359 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
1360 f_j[m] = -f_i[m]-f_k[m];
1367 copy_ivec(SHIFT_IVEC(g, aj), jt);
1369 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1370 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1371 t1 = IVEC2IS(dt_ij);
1372 t2 = IVEC2IS(dt_kj);
1374 rvec_inc(fshift[t1], f_i);
1375 rvec_inc(fshift[CENTRAL], f_j);
1376 rvec_inc(fshift[t2], f_k);
1378 /* Time for the bond calculations */
1384 vtot += vbond; /* 1*/
1385 fbond *= gmx_invsqrt(dr2); /* 6 */
1389 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ak), dt_ik);
1390 ki = IVEC2IS(dt_ik);
1392 for (m = 0; (m < DIM); m++) /* 15 */
1394 fik = fbond*r_ik[m];
1397 fshift[ki][m] += fik;
1398 fshift[CENTRAL][m] -= fik;
1404 real quartic_angles(int nbonds,
1405 const t_iatom forceatoms[], const t_iparams forceparams[],
1406 const rvec x[], rvec f[], rvec fshift[],
1407 const t_pbc *pbc, const t_graph *g,
1408 real lambda, real *dvdlambda,
1409 const t_mdatoms *md, t_fcdata *fcd,
1410 int *global_atom_index)
1412 int i, j, ai, aj, ak, t1, t2, type;
1414 real cos_theta, cos_theta2, theta, dt, dVdt, va, dtp, c, vtot;
1415 ivec jt, dt_ij, dt_kj;
1418 for (i = 0; (i < nbonds); )
1420 type = forceatoms[i++];
1421 ai = forceatoms[i++];
1422 aj = forceatoms[i++];
1423 ak = forceatoms[i++];
1425 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
1426 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
1428 dt = theta - forceparams[type].qangle.theta*DEG2RAD; /* 2 */
1431 va = forceparams[type].qangle.c[0];
1433 for (j = 1; j <= 4; j++)
1435 c = forceparams[type].qangle.c[j];
1444 cos_theta2 = sqr(cos_theta); /* 1 */
1453 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1454 sth = st*cos_theta; /* 1 */
1458 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1459 theta*RAD2DEG, va, dVdt);
1462 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1463 nrij2 = iprod(r_ij, r_ij);
1465 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
1466 cii = sth/nrij2; /* 10 */
1467 ckk = sth/nrkj2; /* 10 */
1469 for (m = 0; (m < DIM); m++) /* 39 */
1471 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
1472 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
1473 f_j[m] = -f_i[m]-f_k[m];
1480 copy_ivec(SHIFT_IVEC(g, aj), jt);
1482 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1483 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1484 t1 = IVEC2IS(dt_ij);
1485 t2 = IVEC2IS(dt_kj);
1487 rvec_inc(fshift[t1], f_i);
1488 rvec_inc(fshift[CENTRAL], f_j);
1489 rvec_inc(fshift[t2], f_k);
1495 real dih_angle(const rvec xi, const rvec xj, const rvec xk, const rvec xl,
1497 rvec r_ij, rvec r_kj, rvec r_kl, rvec m, rvec n,
1498 real *sign, int *t1, int *t2, int *t3)
1502 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
1503 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
1504 *t3 = pbc_rvec_sub(pbc, xk, xl, r_kl); /* 3 */
1506 cprod(r_ij, r_kj, m); /* 9 */
1507 cprod(r_kj, r_kl, n); /* 9 */
1508 phi = gmx_angle(m, n); /* 49 (assuming 25 for atan2) */
1509 ipr = iprod(r_ij, n); /* 5 */
1510 (*sign) = (ipr < 0.0) ? -1.0 : 1.0;
1511 phi = (*sign)*phi; /* 1 */
1519 /* As dih_angle above, but calculates 4 dihedral angles at once using SIMD,
1520 * also calculates the pre-factor required for the dihedral force update.
1521 * Note that bv and buf should be register aligned.
1523 static gmx_inline void
1524 dih_angle_simd(const rvec *x,
1525 const int *ai, const int *aj, const int *ak, const int *al,
1526 const pbc_simd_t *pbc,
1529 gmx_mm_pr *mx_S, gmx_mm_pr *my_S, gmx_mm_pr *mz_S,
1530 gmx_mm_pr *nx_S, gmx_mm_pr *ny_S, gmx_mm_pr *nz_S,
1531 gmx_mm_pr *nrkj_m2_S,
1532 gmx_mm_pr *nrkj_n2_S,
1536 #define UNROLL GMX_SIMD_WIDTH_HERE
1538 gmx_mm_pr rijx_S, rijy_S, rijz_S;
1539 gmx_mm_pr rkjx_S, rkjy_S, rkjz_S;
1540 gmx_mm_pr rklx_S, rkly_S, rklz_S;
1541 gmx_mm_pr cx_S, cy_S, cz_S;
1545 gmx_mm_pr iprm_S, iprn_S;
1546 gmx_mm_pr nrkj2_S, nrkj_1_S, nrkj_2_S, nrkj_S;
1549 gmx_mm_pr nrkj2_min_S;
1550 gmx_mm_pr real_eps_S;
1552 /* Used to avoid division by zero.
1553 * We take into acount that we multiply the result by real_eps_S.
1555 nrkj2_min_S = gmx_set1_pr(GMX_REAL_MIN/(2*GMX_REAL_EPS));
1557 /* The value of the last significant bit (GMX_REAL_EPS is half of that) */
1558 real_eps_S = gmx_set1_pr(2*GMX_REAL_EPS);
1560 for (s = 0; s < UNROLL; s++)
1562 /* If you can't use pbc_dx_simd below for PBC, e.g. because
1563 * you can't round in SIMD, use pbc_rvec_sub here.
1565 for (m = 0; m < DIM; m++)
1567 dr[s + (0*DIM + m)*UNROLL] = x[ai[s]][m] - x[aj[s]][m];
1568 dr[s + (1*DIM + m)*UNROLL] = x[ak[s]][m] - x[aj[s]][m];
1569 dr[s + (2*DIM + m)*UNROLL] = x[ak[s]][m] - x[al[s]][m];
1573 rijx_S = gmx_load_pr(dr + 0*UNROLL);
1574 rijy_S = gmx_load_pr(dr + 1*UNROLL);
1575 rijz_S = gmx_load_pr(dr + 2*UNROLL);
1576 rkjx_S = gmx_load_pr(dr + 3*UNROLL);
1577 rkjy_S = gmx_load_pr(dr + 4*UNROLL);
1578 rkjz_S = gmx_load_pr(dr + 5*UNROLL);
1579 rklx_S = gmx_load_pr(dr + 6*UNROLL);
1580 rkly_S = gmx_load_pr(dr + 7*UNROLL);
1581 rklz_S = gmx_load_pr(dr + 8*UNROLL);
1583 pbc_dx_simd(&rijx_S, &rijy_S, &rijz_S, pbc);
1584 pbc_dx_simd(&rkjx_S, &rkjy_S, &rkjz_S, pbc);
1585 pbc_dx_simd(&rklx_S, &rkly_S, &rklz_S, pbc);
1587 gmx_cprod_pr(rijx_S, rijy_S, rijz_S,
1588 rkjx_S, rkjy_S, rkjz_S,
1591 gmx_cprod_pr(rkjx_S, rkjy_S, rkjz_S,
1592 rklx_S, rkly_S, rklz_S,
1595 gmx_cprod_pr(*mx_S, *my_S, *mz_S,
1596 *nx_S, *ny_S, *nz_S,
1597 &cx_S, &cy_S, &cz_S);
1599 cn_S = gmx_sqrt_pr(gmx_norm2_pr(cx_S, cy_S, cz_S));
1601 s_S = gmx_iprod_pr(*mx_S, *my_S, *mz_S, *nx_S, *ny_S, *nz_S);
1603 /* Determine the dihedral angle, the sign might need correction */
1604 *phi_S = gmx_atan2_pr(cn_S, s_S);
1606 ipr_S = gmx_iprod_pr(rijx_S, rijy_S, rijz_S,
1607 *nx_S, *ny_S, *nz_S);
1609 iprm_S = gmx_norm2_pr(*mx_S, *my_S, *mz_S);
1610 iprn_S = gmx_norm2_pr(*nx_S, *ny_S, *nz_S);
1612 nrkj2_S = gmx_norm2_pr(rkjx_S, rkjy_S, rkjz_S);
1614 /* Avoid division by zero. When zero, the result is multiplied by 0
1615 * anyhow, so the 3 max below do not affect the final result.
1617 nrkj2_S = gmx_max_pr(nrkj2_S, nrkj2_min_S);
1618 nrkj_1_S = gmx_invsqrt_pr(nrkj2_S);
1619 nrkj_2_S = gmx_mul_pr(nrkj_1_S, nrkj_1_S);
1620 nrkj_S = gmx_mul_pr(nrkj2_S, nrkj_1_S);
1622 toler_S = gmx_mul_pr(nrkj2_S, real_eps_S);
1624 /* Here the plain-C code uses a conditional, but we can't do that in SIMD.
1625 * So we take a max with the tolerance instead. Since we multiply with
1626 * m or n later, the max does not affect the results.
1628 iprm_S = gmx_max_pr(iprm_S, toler_S);
1629 iprn_S = gmx_max_pr(iprn_S, toler_S);
1630 *nrkj_m2_S = gmx_mul_pr(nrkj_S, gmx_inv_pr(iprm_S));
1631 *nrkj_n2_S = gmx_mul_pr(nrkj_S, gmx_inv_pr(iprn_S));
1633 /* Set sign of phi_S with the sign of ipr_S; phi_S is currently positive */
1634 *phi_S = gmx_cpsgn_nonneg_pr(ipr_S, *phi_S);
1636 p_S = gmx_iprod_pr(rijx_S, rijy_S, rijz_S,
1637 rkjx_S, rkjy_S, rkjz_S);
1638 p_S = gmx_mul_pr(p_S, nrkj_2_S);
1640 q_S = gmx_iprod_pr(rklx_S, rkly_S, rklz_S,
1641 rkjx_S, rkjy_S, rkjz_S);
1642 q_S = gmx_mul_pr(q_S, nrkj_2_S);
1644 gmx_store_pr(p, p_S);
1645 gmx_store_pr(q, q_S);
1649 #endif /* SIMD_BONDEDS */
1652 void do_dih_fup(int i, int j, int k, int l, real ddphi,
1653 rvec r_ij, rvec r_kj, rvec r_kl,
1654 rvec m, rvec n, rvec f[], rvec fshift[],
1655 const t_pbc *pbc, const t_graph *g,
1656 const rvec x[], int t1, int t2, int t3)
1659 rvec f_i, f_j, f_k, f_l;
1660 rvec uvec, vvec, svec, dx_jl;
1661 real iprm, iprn, nrkj, nrkj2, nrkj_1, nrkj_2;
1662 real a, b, p, q, toler;
1663 ivec jt, dt_ij, dt_kj, dt_lj;
1665 iprm = iprod(m, m); /* 5 */
1666 iprn = iprod(n, n); /* 5 */
1667 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1668 toler = nrkj2*GMX_REAL_EPS;
1669 if ((iprm > toler) && (iprn > toler))
1671 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1672 nrkj_2 = nrkj_1*nrkj_1; /* 1 */
1673 nrkj = nrkj2*nrkj_1; /* 1 */
1674 a = -ddphi*nrkj/iprm; /* 11 */
1675 svmul(a, m, f_i); /* 3 */
1676 b = ddphi*nrkj/iprn; /* 11 */
1677 svmul(b, n, f_l); /* 3 */
1678 p = iprod(r_ij, r_kj); /* 5 */
1679 p *= nrkj_2; /* 1 */
1680 q = iprod(r_kl, r_kj); /* 5 */
1681 q *= nrkj_2; /* 1 */
1682 svmul(p, f_i, uvec); /* 3 */
1683 svmul(q, f_l, vvec); /* 3 */
1684 rvec_sub(uvec, vvec, svec); /* 3 */
1685 rvec_sub(f_i, svec, f_j); /* 3 */
1686 rvec_add(f_l, svec, f_k); /* 3 */
1687 rvec_inc(f[i], f_i); /* 3 */
1688 rvec_dec(f[j], f_j); /* 3 */
1689 rvec_dec(f[k], f_k); /* 3 */
1690 rvec_inc(f[l], f_l); /* 3 */
1694 copy_ivec(SHIFT_IVEC(g, j), jt);
1695 ivec_sub(SHIFT_IVEC(g, i), jt, dt_ij);
1696 ivec_sub(SHIFT_IVEC(g, k), jt, dt_kj);
1697 ivec_sub(SHIFT_IVEC(g, l), jt, dt_lj);
1698 t1 = IVEC2IS(dt_ij);
1699 t2 = IVEC2IS(dt_kj);
1700 t3 = IVEC2IS(dt_lj);
1704 t3 = pbc_rvec_sub(pbc, x[l], x[j], dx_jl);
1711 rvec_inc(fshift[t1], f_i);
1712 rvec_dec(fshift[CENTRAL], f_j);
1713 rvec_dec(fshift[t2], f_k);
1714 rvec_inc(fshift[t3], f_l);
1719 /* As do_dih_fup above, but without shift forces */
1721 do_dih_fup_noshiftf(int i, int j, int k, int l, real ddphi,
1722 rvec r_ij, rvec r_kj, rvec r_kl,
1723 rvec m, rvec n, rvec f[])
1725 rvec f_i, f_j, f_k, f_l;
1726 rvec uvec, vvec, svec, dx_jl;
1727 real iprm, iprn, nrkj, nrkj2, nrkj_1, nrkj_2;
1728 real a, b, p, q, toler;
1729 ivec jt, dt_ij, dt_kj, dt_lj;
1731 iprm = iprod(m, m); /* 5 */
1732 iprn = iprod(n, n); /* 5 */
1733 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1734 toler = nrkj2*GMX_REAL_EPS;
1735 if ((iprm > toler) && (iprn > toler))
1737 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1738 nrkj_2 = nrkj_1*nrkj_1; /* 1 */
1739 nrkj = nrkj2*nrkj_1; /* 1 */
1740 a = -ddphi*nrkj/iprm; /* 11 */
1741 svmul(a, m, f_i); /* 3 */
1742 b = ddphi*nrkj/iprn; /* 11 */
1743 svmul(b, n, f_l); /* 3 */
1744 p = iprod(r_ij, r_kj); /* 5 */
1745 p *= nrkj_2; /* 1 */
1746 q = iprod(r_kl, r_kj); /* 5 */
1747 q *= nrkj_2; /* 1 */
1748 svmul(p, f_i, uvec); /* 3 */
1749 svmul(q, f_l, vvec); /* 3 */
1750 rvec_sub(uvec, vvec, svec); /* 3 */
1751 rvec_sub(f_i, svec, f_j); /* 3 */
1752 rvec_add(f_l, svec, f_k); /* 3 */
1753 rvec_inc(f[i], f_i); /* 3 */
1754 rvec_dec(f[j], f_j); /* 3 */
1755 rvec_dec(f[k], f_k); /* 3 */
1756 rvec_inc(f[l], f_l); /* 3 */
1760 /* As do_dih_fup_noshiftf above, but with pre-calculated pre-factors */
1761 static gmx_inline void
1762 do_dih_fup_noshiftf_precalc(int i, int j, int k, int l,
1764 real f_i_x, real f_i_y, real f_i_z,
1765 real mf_l_x, real mf_l_y, real mf_l_z,
1768 rvec f_i, f_j, f_k, f_l;
1769 rvec uvec, vvec, svec;
1777 svmul(p, f_i, uvec);
1778 svmul(q, f_l, vvec);
1779 rvec_sub(uvec, vvec, svec);
1780 rvec_sub(f_i, svec, f_j);
1781 rvec_add(f_l, svec, f_k);
1782 rvec_inc(f[i], f_i);
1783 rvec_dec(f[j], f_j);
1784 rvec_dec(f[k], f_k);
1785 rvec_inc(f[l], f_l);
1789 real dopdihs(real cpA, real cpB, real phiA, real phiB, int mult,
1790 real phi, real lambda, real *V, real *F)
1792 real v, dvdlambda, mdphi, v1, sdphi, ddphi;
1793 real L1 = 1.0 - lambda;
1794 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1795 real dph0 = (phiB - phiA)*DEG2RAD;
1796 real cp = L1*cpA + lambda*cpB;
1798 mdphi = mult*phi - ph0;
1800 ddphi = -cp*mult*sdphi;
1801 v1 = 1.0 + cos(mdphi);
1804 dvdlambda = (cpB - cpA)*v1 + cp*dph0*sdphi;
1811 /* That was 40 flops */
1815 dopdihs_noener(real cpA, real cpB, real phiA, real phiB, int mult,
1816 real phi, real lambda, real *F)
1818 real mdphi, sdphi, ddphi;
1819 real L1 = 1.0 - lambda;
1820 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1821 real cp = L1*cpA + lambda*cpB;
1823 mdphi = mult*phi - ph0;
1825 ddphi = -cp*mult*sdphi;
1829 /* That was 20 flops */
1833 dopdihs_mdphi(real cpA, real cpB, real phiA, real phiB, int mult,
1834 real phi, real lambda, real *cp, real *mdphi)
1836 real L1 = 1.0 - lambda;
1837 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1839 *cp = L1*cpA + lambda*cpB;
1841 *mdphi = mult*phi - ph0;
1844 static real dopdihs_min(real cpA, real cpB, real phiA, real phiB, int mult,
1845 real phi, real lambda, real *V, real *F)
1846 /* similar to dopdihs, except for a minus sign *
1847 * and a different treatment of mult/phi0 */
1849 real v, dvdlambda, mdphi, v1, sdphi, ddphi;
1850 real L1 = 1.0 - lambda;
1851 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1852 real dph0 = (phiB - phiA)*DEG2RAD;
1853 real cp = L1*cpA + lambda*cpB;
1855 mdphi = mult*(phi-ph0);
1857 ddphi = cp*mult*sdphi;
1858 v1 = 1.0-cos(mdphi);
1861 dvdlambda = (cpB-cpA)*v1 + cp*dph0*sdphi;
1868 /* That was 40 flops */
1871 real pdihs(int nbonds,
1872 const t_iatom forceatoms[], const t_iparams forceparams[],
1873 const rvec x[], rvec f[], rvec fshift[],
1874 const t_pbc *pbc, const t_graph *g,
1875 real lambda, real *dvdlambda,
1876 const t_mdatoms *md, t_fcdata *fcd,
1877 int *global_atom_index)
1879 int i, type, ai, aj, ak, al;
1881 rvec r_ij, r_kj, r_kl, m, n;
1882 real phi, sign, ddphi, vpd, vtot;
1886 for (i = 0; (i < nbonds); )
1888 type = forceatoms[i++];
1889 ai = forceatoms[i++];
1890 aj = forceatoms[i++];
1891 ak = forceatoms[i++];
1892 al = forceatoms[i++];
1894 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
1895 &sign, &t1, &t2, &t3); /* 84 */
1896 *dvdlambda += dopdihs(forceparams[type].pdihs.cpA,
1897 forceparams[type].pdihs.cpB,
1898 forceparams[type].pdihs.phiA,
1899 forceparams[type].pdihs.phiB,
1900 forceparams[type].pdihs.mult,
1901 phi, lambda, &vpd, &ddphi);
1904 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
1905 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
1908 fprintf(debug, "pdih: (%d,%d,%d,%d) phi=%g\n",
1909 ai, aj, ak, al, phi);
1916 void make_dp_periodic(real *dp) /* 1 flop? */
1918 /* dp cannot be outside (-pi,pi) */
1923 else if (*dp < -M_PI)
1930 /* As pdihs above, but without calculating energies and shift forces */
1932 pdihs_noener(int nbonds,
1933 const t_iatom forceatoms[], const t_iparams forceparams[],
1934 const rvec x[], rvec f[],
1935 const t_pbc *pbc, const t_graph *g,
1937 const t_mdatoms *md, t_fcdata *fcd,
1938 int *global_atom_index)
1940 int i, type, ai, aj, ak, al;
1942 rvec r_ij, r_kj, r_kl, m, n;
1943 real phi, sign, ddphi_tot, ddphi;
1945 for (i = 0; (i < nbonds); )
1947 ai = forceatoms[i+1];
1948 aj = forceatoms[i+2];
1949 ak = forceatoms[i+3];
1950 al = forceatoms[i+4];
1952 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
1953 &sign, &t1, &t2, &t3);
1957 /* Loop over dihedrals working on the same atoms,
1958 * so we avoid recalculating angles and force distributions.
1962 type = forceatoms[i];
1963 dopdihs_noener(forceparams[type].pdihs.cpA,
1964 forceparams[type].pdihs.cpB,
1965 forceparams[type].pdihs.phiA,
1966 forceparams[type].pdihs.phiB,
1967 forceparams[type].pdihs.mult,
1968 phi, lambda, &ddphi);
1973 while (i < nbonds &&
1974 forceatoms[i+1] == ai &&
1975 forceatoms[i+2] == aj &&
1976 forceatoms[i+3] == ak &&
1977 forceatoms[i+4] == al);
1979 do_dih_fup_noshiftf(ai, aj, ak, al, ddphi_tot, r_ij, r_kj, r_kl, m, n, f);
1986 /* As pdihs_noner above, but using SIMD to calculate many dihedrals at once */
1988 pdihs_noener_simd(int nbonds,
1989 const t_iatom forceatoms[], const t_iparams forceparams[],
1990 const rvec x[], rvec f[],
1991 const t_pbc *pbc, const t_graph *g,
1993 const t_mdatoms *md, t_fcdata *fcd,
1994 int *global_atom_index)
1996 #define UNROLL GMX_SIMD_WIDTH_HERE
1999 int type, ai[UNROLL], aj[UNROLL], ak[UNROLL], al[UNROLL];
2000 int t1[UNROLL], t2[UNROLL], t3[UNROLL];
2002 real dr_array[3*DIM*UNROLL+UNROLL], *dr;
2003 real buf_array[7*UNROLL+UNROLL], *buf;
2004 real *cp, *phi0, *mult, *phi, *p, *q, *sf_i, *msf_l;
2005 gmx_mm_pr phi0_S, phi_S;
2006 gmx_mm_pr mx_S, my_S, mz_S;
2007 gmx_mm_pr nx_S, ny_S, nz_S;
2008 gmx_mm_pr nrkj_m2_S, nrkj_n2_S;
2009 gmx_mm_pr cp_S, mdphi_S, mult_S;
2010 gmx_mm_pr sin_S, cos_S;
2012 gmx_mm_pr sf_i_S, msf_l_S;
2013 pbc_simd_t pbc_simd;
2015 /* Ensure SIMD register alignment */
2016 dr = gmx_simd_align_real(dr_array);
2017 buf = gmx_simd_align_real(buf_array);
2019 /* Extract aligned pointer for parameters and variables */
2020 cp = buf + 0*UNROLL;
2021 phi0 = buf + 1*UNROLL;
2022 mult = buf + 2*UNROLL;
2025 sf_i = buf + 5*UNROLL;
2026 msf_l = buf + 6*UNROLL;
2028 set_pbc_simd(pbc, &pbc_simd);
2030 /* nbonds is the number of dihedrals times nfa1, here we step UNROLL dihs */
2031 for (i = 0; (i < nbonds); i += UNROLL*nfa1)
2033 /* Collect atoms quadruplets for UNROLL dihedrals.
2034 * iu indexes into forceatoms, we should not let iu go beyond nbonds.
2037 for (s = 0; s < UNROLL; s++)
2039 type = forceatoms[iu];
2040 ai[s] = forceatoms[iu+1];
2041 aj[s] = forceatoms[iu+2];
2042 ak[s] = forceatoms[iu+3];
2043 al[s] = forceatoms[iu+4];
2045 cp[s] = forceparams[type].pdihs.cpA;
2046 phi0[s] = forceparams[type].pdihs.phiA*DEG2RAD;
2047 mult[s] = forceparams[type].pdihs.mult;
2049 /* At the end fill the arrays with identical entries */
2050 if (iu + nfa1 < nbonds)
2056 /* Caclulate UNROLL dihedral angles at once */
2057 dih_angle_simd(x, ai, aj, ak, al, &pbc_simd,
2060 &mx_S, &my_S, &mz_S,
2061 &nx_S, &ny_S, &nz_S,
2066 cp_S = gmx_load_pr(cp);
2067 phi0_S = gmx_load_pr(phi0);
2068 mult_S = gmx_load_pr(mult);
2070 mdphi_S = gmx_sub_pr(gmx_mul_pr(mult_S, phi_S), phi0_S);
2072 /* Calculate UNROLL sines at once */
2073 gmx_sincos_pr(mdphi_S, &sin_S, &cos_S);
2074 mddphi_S = gmx_mul_pr(gmx_mul_pr(cp_S, mult_S), sin_S);
2075 sf_i_S = gmx_mul_pr(mddphi_S, nrkj_m2_S);
2076 msf_l_S = gmx_mul_pr(mddphi_S, nrkj_n2_S);
2078 /* After this m?_S will contain f[i] */
2079 mx_S = gmx_mul_pr(sf_i_S, mx_S);
2080 my_S = gmx_mul_pr(sf_i_S, my_S);
2081 mz_S = gmx_mul_pr(sf_i_S, mz_S);
2083 /* After this m?_S will contain -f[l] */
2084 nx_S = gmx_mul_pr(msf_l_S, nx_S);
2085 ny_S = gmx_mul_pr(msf_l_S, ny_S);
2086 nz_S = gmx_mul_pr(msf_l_S, nz_S);
2088 gmx_store_pr(dr + 0*UNROLL, mx_S);
2089 gmx_store_pr(dr + 1*UNROLL, my_S);
2090 gmx_store_pr(dr + 2*UNROLL, mz_S);
2091 gmx_store_pr(dr + 3*UNROLL, nx_S);
2092 gmx_store_pr(dr + 4*UNROLL, ny_S);
2093 gmx_store_pr(dr + 5*UNROLL, nz_S);
2099 do_dih_fup_noshiftf_precalc(ai[s], aj[s], ak[s], al[s],
2104 dr[(DIM+XX)*UNROLL+s],
2105 dr[(DIM+YY)*UNROLL+s],
2106 dr[(DIM+ZZ)*UNROLL+s],
2111 while (s < UNROLL && iu < nbonds);
2116 #endif /* SIMD_BONDEDS */
2119 real idihs(int nbonds,
2120 const t_iatom forceatoms[], const t_iparams forceparams[],
2121 const rvec x[], rvec f[], rvec fshift[],
2122 const t_pbc *pbc, const t_graph *g,
2123 real lambda, real *dvdlambda,
2124 const t_mdatoms *md, t_fcdata *fcd,
2125 int *global_atom_index)
2127 int i, type, ai, aj, ak, al;
2129 real phi, phi0, dphi0, ddphi, sign, vtot;
2130 rvec r_ij, r_kj, r_kl, m, n;
2131 real L1, kk, dp, dp2, kA, kB, pA, pB, dvdl_term;
2136 for (i = 0; (i < nbonds); )
2138 type = forceatoms[i++];
2139 ai = forceatoms[i++];
2140 aj = forceatoms[i++];
2141 ak = forceatoms[i++];
2142 al = forceatoms[i++];
2144 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2145 &sign, &t1, &t2, &t3); /* 84 */
2147 /* phi can jump if phi0 is close to Pi/-Pi, which will cause huge
2148 * force changes if we just apply a normal harmonic.
2149 * Instead, we first calculate phi-phi0 and take it modulo (-Pi,Pi).
2150 * This means we will never have the periodicity problem, unless
2151 * the dihedral is Pi away from phiO, which is very unlikely due to
2154 kA = forceparams[type].harmonic.krA;
2155 kB = forceparams[type].harmonic.krB;
2156 pA = forceparams[type].harmonic.rA;
2157 pB = forceparams[type].harmonic.rB;
2159 kk = L1*kA + lambda*kB;
2160 phi0 = (L1*pA + lambda*pB)*DEG2RAD;
2161 dphi0 = (pB - pA)*DEG2RAD;
2165 make_dp_periodic(&dp);
2172 dvdl_term += 0.5*(kB - kA)*dp2 - kk*dphi0*dp;
2174 do_dih_fup(ai, aj, ak, al, (real)(-ddphi), r_ij, r_kj, r_kl, m, n,
2175 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2180 fprintf(debug, "idih: (%d,%d,%d,%d) phi=%g\n",
2181 ai, aj, ak, al, phi);
2186 *dvdlambda += dvdl_term;
2191 real posres(int nbonds,
2192 const t_iatom forceatoms[], const t_iparams forceparams[],
2193 const rvec x[], rvec f[], rvec vir_diag,
2195 real lambda, real *dvdlambda,
2196 int refcoord_scaling, int ePBC, rvec comA, rvec comB)
2198 int i, ai, m, d, type, ki, npbcdim = 0;
2199 const t_iparams *pr;
2202 real posA, posB, ref = 0;
2203 rvec comA_sc, comB_sc, rdist, dpdl, pos, dx;
2204 gmx_bool bForceValid = TRUE;
2206 if ((f == NULL) || (vir_diag == NULL)) /* should both be null together! */
2208 bForceValid = FALSE;
2211 npbcdim = ePBC2npbcdim(ePBC);
2213 if (refcoord_scaling == erscCOM)
2215 clear_rvec(comA_sc);
2216 clear_rvec(comB_sc);
2217 for (m = 0; m < npbcdim; m++)
2219 for (d = m; d < npbcdim; d++)
2221 comA_sc[m] += comA[d]*pbc->box[d][m];
2222 comB_sc[m] += comB[d]*pbc->box[d][m];
2230 for (i = 0; (i < nbonds); )
2232 type = forceatoms[i++];
2233 ai = forceatoms[i++];
2234 pr = &forceparams[type];
2236 for (m = 0; m < DIM; m++)
2238 posA = forceparams[type].posres.pos0A[m];
2239 posB = forceparams[type].posres.pos0B[m];
2242 switch (refcoord_scaling)
2246 rdist[m] = L1*posA + lambda*posB;
2247 dpdl[m] = posB - posA;
2250 /* Box relative coordinates are stored for dimensions with pbc */
2251 posA *= pbc->box[m][m];
2252 posB *= pbc->box[m][m];
2253 for (d = m+1; d < npbcdim; d++)
2255 posA += forceparams[type].posres.pos0A[d]*pbc->box[d][m];
2256 posB += forceparams[type].posres.pos0B[d]*pbc->box[d][m];
2258 ref = L1*posA + lambda*posB;
2260 dpdl[m] = posB - posA;
2263 ref = L1*comA_sc[m] + lambda*comB_sc[m];
2264 rdist[m] = L1*posA + lambda*posB;
2265 dpdl[m] = comB_sc[m] - comA_sc[m] + posB - posA;
2268 gmx_fatal(FARGS, "No such scaling method implemented");
2273 ref = L1*posA + lambda*posB;
2275 dpdl[m] = posB - posA;
2278 /* We do pbc_dx with ref+rdist,
2279 * since with only ref we can be up to half a box vector wrong.
2281 pos[m] = ref + rdist[m];
2286 pbc_dx(pbc, x[ai], pos, dx);
2290 rvec_sub(x[ai], pos, dx);
2293 for (m = 0; (m < DIM); m++)
2295 kk = L1*pr->posres.fcA[m] + lambda*pr->posres.fcB[m];
2297 vtot += 0.5*kk*dx[m]*dx[m];
2299 0.5*(pr->posres.fcB[m] - pr->posres.fcA[m])*dx[m]*dx[m]
2302 /* Here we correct for the pbc_dx which included rdist */
2306 vir_diag[m] -= 0.5*(dx[m] + rdist[m])*fm;
2314 static real low_angres(int nbonds,
2315 const t_iatom forceatoms[], const t_iparams forceparams[],
2316 const rvec x[], rvec f[], rvec fshift[],
2317 const t_pbc *pbc, const t_graph *g,
2318 real lambda, real *dvdlambda,
2321 int i, m, type, ai, aj, ak, al;
2323 real phi, cos_phi, cos_phi2, vid, vtot, dVdphi;
2324 rvec r_ij, r_kl, f_i, f_k = {0, 0, 0};
2325 real st, sth, nrij2, nrkl2, c, cij, ckl;
2328 t2 = 0; /* avoid warning with gcc-3.3. It is never used uninitialized */
2331 ak = al = 0; /* to avoid warnings */
2332 for (i = 0; i < nbonds; )
2334 type = forceatoms[i++];
2335 ai = forceatoms[i++];
2336 aj = forceatoms[i++];
2337 t1 = pbc_rvec_sub(pbc, x[aj], x[ai], r_ij); /* 3 */
2340 ak = forceatoms[i++];
2341 al = forceatoms[i++];
2342 t2 = pbc_rvec_sub(pbc, x[al], x[ak], r_kl); /* 3 */
2351 cos_phi = cos_angle(r_ij, r_kl); /* 25 */
2352 phi = acos(cos_phi); /* 10 */
2354 *dvdlambda += dopdihs_min(forceparams[type].pdihs.cpA,
2355 forceparams[type].pdihs.cpB,
2356 forceparams[type].pdihs.phiA,
2357 forceparams[type].pdihs.phiB,
2358 forceparams[type].pdihs.mult,
2359 phi, lambda, &vid, &dVdphi); /* 40 */
2363 cos_phi2 = sqr(cos_phi); /* 1 */
2366 st = -dVdphi*gmx_invsqrt(1 - cos_phi2); /* 12 */
2367 sth = st*cos_phi; /* 1 */
2368 nrij2 = iprod(r_ij, r_ij); /* 5 */
2369 nrkl2 = iprod(r_kl, r_kl); /* 5 */
2371 c = st*gmx_invsqrt(nrij2*nrkl2); /* 11 */
2372 cij = sth/nrij2; /* 10 */
2373 ckl = sth/nrkl2; /* 10 */
2375 for (m = 0; m < DIM; m++) /* 18+18 */
2377 f_i[m] = (c*r_kl[m]-cij*r_ij[m]);
2382 f_k[m] = (c*r_ij[m]-ckl*r_kl[m]);
2390 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
2393 rvec_inc(fshift[t1], f_i);
2394 rvec_dec(fshift[CENTRAL], f_i);
2399 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, al), dt);
2402 rvec_inc(fshift[t2], f_k);
2403 rvec_dec(fshift[CENTRAL], f_k);
2408 return vtot; /* 184 / 157 (bZAxis) total */
2411 real angres(int nbonds,
2412 const t_iatom forceatoms[], const t_iparams forceparams[],
2413 const rvec x[], rvec f[], rvec fshift[],
2414 const t_pbc *pbc, const t_graph *g,
2415 real lambda, real *dvdlambda,
2416 const t_mdatoms *md, t_fcdata *fcd,
2417 int *global_atom_index)
2419 return low_angres(nbonds, forceatoms, forceparams, x, f, fshift, pbc, g,
2420 lambda, dvdlambda, FALSE);
2423 real angresz(int nbonds,
2424 const t_iatom forceatoms[], const t_iparams forceparams[],
2425 const rvec x[], rvec f[], rvec fshift[],
2426 const t_pbc *pbc, const t_graph *g,
2427 real lambda, real *dvdlambda,
2428 const t_mdatoms *md, t_fcdata *fcd,
2429 int *global_atom_index)
2431 return low_angres(nbonds, forceatoms, forceparams, x, f, fshift, pbc, g,
2432 lambda, dvdlambda, TRUE);
2435 real dihres(int nbonds,
2436 const t_iatom forceatoms[], const t_iparams forceparams[],
2437 const rvec x[], rvec f[], rvec fshift[],
2438 const t_pbc *pbc, const t_graph *g,
2439 real lambda, real *dvdlambda,
2440 const t_mdatoms *md, t_fcdata *fcd,
2441 int *global_atom_index)
2444 int ai, aj, ak, al, i, k, type, t1, t2, t3;
2445 real phi0A, phi0B, dphiA, dphiB, kfacA, kfacB, phi0, dphi, kfac;
2446 real phi, ddphi, ddp, ddp2, dp, sign, d2r, fc, L1;
2447 rvec r_ij, r_kj, r_kl, m, n;
2454 for (i = 0; (i < nbonds); )
2456 type = forceatoms[i++];
2457 ai = forceatoms[i++];
2458 aj = forceatoms[i++];
2459 ak = forceatoms[i++];
2460 al = forceatoms[i++];
2462 phi0A = forceparams[type].dihres.phiA*d2r;
2463 dphiA = forceparams[type].dihres.dphiA*d2r;
2464 kfacA = forceparams[type].dihres.kfacA;
2466 phi0B = forceparams[type].dihres.phiB*d2r;
2467 dphiB = forceparams[type].dihres.dphiB*d2r;
2468 kfacB = forceparams[type].dihres.kfacB;
2470 phi0 = L1*phi0A + lambda*phi0B;
2471 dphi = L1*dphiA + lambda*dphiB;
2472 kfac = L1*kfacA + lambda*kfacB;
2474 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2475 &sign, &t1, &t2, &t3);
2480 fprintf(debug, "dihres[%d]: %d %d %d %d : phi=%f, dphi=%f, kfac=%f\n",
2481 k++, ai, aj, ak, al, phi0, dphi, kfac);
2483 /* phi can jump if phi0 is close to Pi/-Pi, which will cause huge
2484 * force changes if we just apply a normal harmonic.
2485 * Instead, we first calculate phi-phi0 and take it modulo (-Pi,Pi).
2486 * This means we will never have the periodicity problem, unless
2487 * the dihedral is Pi away from phiO, which is very unlikely due to
2491 make_dp_periodic(&dp);
2497 else if (dp < -dphi)
2509 vtot += 0.5*kfac*ddp2;
2512 *dvdlambda += 0.5*(kfacB - kfacA)*ddp2;
2513 /* lambda dependence from changing restraint distances */
2516 *dvdlambda -= kfac*ddp*((dphiB - dphiA)+(phi0B - phi0A));
2520 *dvdlambda += kfac*ddp*((dphiB - dphiA)-(phi0B - phi0A));
2522 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
2523 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2530 real unimplemented(int nbonds,
2531 const t_iatom forceatoms[], const t_iparams forceparams[],
2532 const rvec x[], rvec f[], rvec fshift[],
2533 const t_pbc *pbc, const t_graph *g,
2534 real lambda, real *dvdlambda,
2535 const t_mdatoms *md, t_fcdata *fcd,
2536 int *global_atom_index)
2538 gmx_impl("*** you are using a not implemented function");
2540 return 0.0; /* To make the compiler happy */
2543 real rbdihs(int nbonds,
2544 const t_iatom forceatoms[], const t_iparams forceparams[],
2545 const rvec x[], rvec f[], rvec fshift[],
2546 const t_pbc *pbc, const t_graph *g,
2547 real lambda, real *dvdlambda,
2548 const t_mdatoms *md, t_fcdata *fcd,
2549 int *global_atom_index)
2551 const real c0 = 0.0, c1 = 1.0, c2 = 2.0, c3 = 3.0, c4 = 4.0, c5 = 5.0;
2552 int type, ai, aj, ak, al, i, j;
2554 rvec r_ij, r_kj, r_kl, m, n;
2555 real parmA[NR_RBDIHS];
2556 real parmB[NR_RBDIHS];
2557 real parm[NR_RBDIHS];
2558 real cos_phi, phi, rbp, rbpBA;
2559 real v, sign, ddphi, sin_phi;
2561 real L1 = 1.0-lambda;
2565 for (i = 0; (i < nbonds); )
2567 type = forceatoms[i++];
2568 ai = forceatoms[i++];
2569 aj = forceatoms[i++];
2570 ak = forceatoms[i++];
2571 al = forceatoms[i++];
2573 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2574 &sign, &t1, &t2, &t3); /* 84 */
2576 /* Change to polymer convention */
2583 phi -= M_PI; /* 1 */
2587 /* Beware of accuracy loss, cannot use 1-sqrt(cos^2) ! */
2590 for (j = 0; (j < NR_RBDIHS); j++)
2592 parmA[j] = forceparams[type].rbdihs.rbcA[j];
2593 parmB[j] = forceparams[type].rbdihs.rbcB[j];
2594 parm[j] = L1*parmA[j]+lambda*parmB[j];
2596 /* Calculate cosine powers */
2597 /* Calculate the energy */
2598 /* Calculate the derivative */
2601 dvdl_term += (parmB[0]-parmA[0]);
2606 rbpBA = parmB[1]-parmA[1];
2607 ddphi += rbp*cosfac;
2610 dvdl_term += cosfac*rbpBA;
2612 rbpBA = parmB[2]-parmA[2];
2613 ddphi += c2*rbp*cosfac;
2616 dvdl_term += cosfac*rbpBA;
2618 rbpBA = parmB[3]-parmA[3];
2619 ddphi += c3*rbp*cosfac;
2622 dvdl_term += cosfac*rbpBA;
2624 rbpBA = parmB[4]-parmA[4];
2625 ddphi += c4*rbp*cosfac;
2628 dvdl_term += cosfac*rbpBA;
2630 rbpBA = parmB[5]-parmA[5];
2631 ddphi += c5*rbp*cosfac;
2634 dvdl_term += cosfac*rbpBA;
2636 ddphi = -ddphi*sin_phi; /* 11 */
2638 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
2639 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2642 *dvdlambda += dvdl_term;
2647 int cmap_setup_grid_index(int ip, int grid_spacing, int *ipm1, int *ipp1, int *ipp2)
2653 ip = ip + grid_spacing - 1;
2655 else if (ip > grid_spacing)
2657 ip = ip - grid_spacing - 1;
2666 im1 = grid_spacing - 1;
2668 else if (ip == grid_spacing-2)
2672 else if (ip == grid_spacing-1)
2686 real cmap_dihs(int nbonds,
2687 const t_iatom forceatoms[], const t_iparams forceparams[],
2688 const gmx_cmap_t *cmap_grid,
2689 const rvec x[], rvec f[], rvec fshift[],
2690 const t_pbc *pbc, const t_graph *g,
2691 real lambda, real *dvdlambda,
2692 const t_mdatoms *md, t_fcdata *fcd,
2693 int *global_atom_index)
2695 int i, j, k, n, idx;
2696 int ai, aj, ak, al, am;
2697 int a1i, a1j, a1k, a1l, a2i, a2j, a2k, a2l;
2699 int t11, t21, t31, t12, t22, t32;
2700 int iphi1, ip1m1, ip1p1, ip1p2;
2701 int iphi2, ip2m1, ip2p1, ip2p2;
2703 int pos1, pos2, pos3, pos4, tmp;
2705 real ty[4], ty1[4], ty2[4], ty12[4], tc[16], tx[16];
2706 real phi1, psi1, cos_phi1, sin_phi1, sign1, xphi1;
2707 real phi2, psi2, cos_phi2, sin_phi2, sign2, xphi2;
2708 real dx, xx, tt, tu, e, df1, df2, ddf1, ddf2, ddf12, vtot;
2709 real ra21, rb21, rg21, rg1, rgr1, ra2r1, rb2r1, rabr1;
2710 real ra22, rb22, rg22, rg2, rgr2, ra2r2, rb2r2, rabr2;
2711 real fg1, hg1, fga1, hgb1, gaa1, gbb1;
2712 real fg2, hg2, fga2, hgb2, gaa2, gbb2;
2715 rvec r1_ij, r1_kj, r1_kl, m1, n1;
2716 rvec r2_ij, r2_kj, r2_kl, m2, n2;
2717 rvec f1_i, f1_j, f1_k, f1_l;
2718 rvec f2_i, f2_j, f2_k, f2_l;
2719 rvec a1, b1, a2, b2;
2720 rvec f1, g1, h1, f2, g2, h2;
2721 rvec dtf1, dtg1, dth1, dtf2, dtg2, dth2;
2722 ivec jt1, dt1_ij, dt1_kj, dt1_lj;
2723 ivec jt2, dt2_ij, dt2_kj, dt2_lj;
2727 int loop_index[4][4] = {
2734 /* Total CMAP energy */
2737 for (n = 0; n < nbonds; )
2739 /* Five atoms are involved in the two torsions */
2740 type = forceatoms[n++];
2741 ai = forceatoms[n++];
2742 aj = forceatoms[n++];
2743 ak = forceatoms[n++];
2744 al = forceatoms[n++];
2745 am = forceatoms[n++];
2747 /* Which CMAP type is this */
2748 cmapA = forceparams[type].cmap.cmapA;
2749 cmapd = cmap_grid->cmapdata[cmapA].cmap;
2757 phi1 = dih_angle(x[a1i], x[a1j], x[a1k], x[a1l], pbc, r1_ij, r1_kj, r1_kl, m1, n1,
2758 &sign1, &t11, &t21, &t31); /* 84 */
2760 cos_phi1 = cos(phi1);
2762 a1[0] = r1_ij[1]*r1_kj[2]-r1_ij[2]*r1_kj[1];
2763 a1[1] = r1_ij[2]*r1_kj[0]-r1_ij[0]*r1_kj[2];
2764 a1[2] = r1_ij[0]*r1_kj[1]-r1_ij[1]*r1_kj[0]; /* 9 */
2766 b1[0] = r1_kl[1]*r1_kj[2]-r1_kl[2]*r1_kj[1];
2767 b1[1] = r1_kl[2]*r1_kj[0]-r1_kl[0]*r1_kj[2];
2768 b1[2] = r1_kl[0]*r1_kj[1]-r1_kl[1]*r1_kj[0]; /* 9 */
2770 tmp = pbc_rvec_sub(pbc, x[a1l], x[a1k], h1);
2772 ra21 = iprod(a1, a1); /* 5 */
2773 rb21 = iprod(b1, b1); /* 5 */
2774 rg21 = iprod(r1_kj, r1_kj); /* 5 */
2780 rabr1 = sqrt(ra2r1*rb2r1);
2782 sin_phi1 = rg1 * rabr1 * iprod(a1, h1) * (-1);
2784 if (cos_phi1 < -0.5 || cos_phi1 > 0.5)
2786 phi1 = asin(sin_phi1);
2796 phi1 = -M_PI - phi1;
2802 phi1 = acos(cos_phi1);
2810 xphi1 = phi1 + M_PI; /* 1 */
2812 /* Second torsion */
2818 phi2 = dih_angle(x[a2i], x[a2j], x[a2k], x[a2l], pbc, r2_ij, r2_kj, r2_kl, m2, n2,
2819 &sign2, &t12, &t22, &t32); /* 84 */
2821 cos_phi2 = cos(phi2);
2823 a2[0] = r2_ij[1]*r2_kj[2]-r2_ij[2]*r2_kj[1];
2824 a2[1] = r2_ij[2]*r2_kj[0]-r2_ij[0]*r2_kj[2];
2825 a2[2] = r2_ij[0]*r2_kj[1]-r2_ij[1]*r2_kj[0]; /* 9 */
2827 b2[0] = r2_kl[1]*r2_kj[2]-r2_kl[2]*r2_kj[1];
2828 b2[1] = r2_kl[2]*r2_kj[0]-r2_kl[0]*r2_kj[2];
2829 b2[2] = r2_kl[0]*r2_kj[1]-r2_kl[1]*r2_kj[0]; /* 9 */
2831 tmp = pbc_rvec_sub(pbc, x[a2l], x[a2k], h2);
2833 ra22 = iprod(a2, a2); /* 5 */
2834 rb22 = iprod(b2, b2); /* 5 */
2835 rg22 = iprod(r2_kj, r2_kj); /* 5 */
2841 rabr2 = sqrt(ra2r2*rb2r2);
2843 sin_phi2 = rg2 * rabr2 * iprod(a2, h2) * (-1);
2845 if (cos_phi2 < -0.5 || cos_phi2 > 0.5)
2847 phi2 = asin(sin_phi2);
2857 phi2 = -M_PI - phi2;
2863 phi2 = acos(cos_phi2);
2871 xphi2 = phi2 + M_PI; /* 1 */
2873 /* Range mangling */
2876 xphi1 = xphi1 + 2*M_PI;
2878 else if (xphi1 >= 2*M_PI)
2880 xphi1 = xphi1 - 2*M_PI;
2885 xphi2 = xphi2 + 2*M_PI;
2887 else if (xphi2 >= 2*M_PI)
2889 xphi2 = xphi2 - 2*M_PI;
2892 /* Number of grid points */
2893 dx = 2*M_PI / cmap_grid->grid_spacing;
2895 /* Where on the grid are we */
2896 iphi1 = (int)(xphi1/dx);
2897 iphi2 = (int)(xphi2/dx);
2899 iphi1 = cmap_setup_grid_index(iphi1, cmap_grid->grid_spacing, &ip1m1, &ip1p1, &ip1p2);
2900 iphi2 = cmap_setup_grid_index(iphi2, cmap_grid->grid_spacing, &ip2m1, &ip2p1, &ip2p2);
2902 pos1 = iphi1*cmap_grid->grid_spacing+iphi2;
2903 pos2 = ip1p1*cmap_grid->grid_spacing+iphi2;
2904 pos3 = ip1p1*cmap_grid->grid_spacing+ip2p1;
2905 pos4 = iphi1*cmap_grid->grid_spacing+ip2p1;
2907 ty[0] = cmapd[pos1*4];
2908 ty[1] = cmapd[pos2*4];
2909 ty[2] = cmapd[pos3*4];
2910 ty[3] = cmapd[pos4*4];
2912 ty1[0] = cmapd[pos1*4+1];
2913 ty1[1] = cmapd[pos2*4+1];
2914 ty1[2] = cmapd[pos3*4+1];
2915 ty1[3] = cmapd[pos4*4+1];
2917 ty2[0] = cmapd[pos1*4+2];
2918 ty2[1] = cmapd[pos2*4+2];
2919 ty2[2] = cmapd[pos3*4+2];
2920 ty2[3] = cmapd[pos4*4+2];
2922 ty12[0] = cmapd[pos1*4+3];
2923 ty12[1] = cmapd[pos2*4+3];
2924 ty12[2] = cmapd[pos3*4+3];
2925 ty12[3] = cmapd[pos4*4+3];
2927 /* Switch to degrees */
2928 dx = 360.0 / cmap_grid->grid_spacing;
2929 xphi1 = xphi1 * RAD2DEG;
2930 xphi2 = xphi2 * RAD2DEG;
2932 for (i = 0; i < 4; i++) /* 16 */
2935 tx[i+4] = ty1[i]*dx;
2936 tx[i+8] = ty2[i]*dx;
2937 tx[i+12] = ty12[i]*dx*dx;
2941 for (i = 0; i < 4; i++) /* 1056 */
2943 for (j = 0; j < 4; j++)
2946 for (k = 0; k < 16; k++)
2948 xx = xx + cmap_coeff_matrix[k*16+idx]*tx[k];
2956 tt = (xphi1-iphi1*dx)/dx;
2957 tu = (xphi2-iphi2*dx)/dx;
2966 for (i = 3; i >= 0; i--)
2968 l1 = loop_index[i][3];
2969 l2 = loop_index[i][2];
2970 l3 = loop_index[i][1];
2972 e = tt * e + ((tc[i*4+3]*tu+tc[i*4+2])*tu + tc[i*4+1])*tu+tc[i*4];
2973 df1 = tu * df1 + (3.0*tc[l1]*tt+2.0*tc[l2])*tt+tc[l3];
2974 df2 = tt * df2 + (3.0*tc[i*4+3]*tu+2.0*tc[i*4+2])*tu+tc[i*4+1];
2975 ddf1 = tu * ddf1 + 2.0*3.0*tc[l1]*tt+2.0*tc[l2];
2976 ddf2 = tt * ddf2 + 2.0*3.0*tc[4*i+3]*tu+2.0*tc[4*i+2];
2979 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) +
2980 3.0*tu*tu*(tc[7]+2.0*tc[11]*tt+3.0*tc[15]*tt*tt);
2985 ddf1 = ddf1 * fac * fac;
2986 ddf2 = ddf2 * fac * fac;
2987 ddf12 = ddf12 * fac * fac;
2992 /* Do forces - first torsion */
2993 fg1 = iprod(r1_ij, r1_kj);
2994 hg1 = iprod(r1_kl, r1_kj);
2995 fga1 = fg1*ra2r1*rgr1;
2996 hgb1 = hg1*rb2r1*rgr1;
3000 for (i = 0; i < DIM; i++)
3002 dtf1[i] = gaa1 * a1[i];
3003 dtg1[i] = fga1 * a1[i] - hgb1 * b1[i];
3004 dth1[i] = gbb1 * b1[i];
3006 f1[i] = df1 * dtf1[i];
3007 g1[i] = df1 * dtg1[i];
3008 h1[i] = df1 * dth1[i];
3011 f1_j[i] = -f1[i] - g1[i];
3012 f1_k[i] = h1[i] + g1[i];
3015 f[a1i][i] = f[a1i][i] + f1_i[i];
3016 f[a1j][i] = f[a1j][i] + f1_j[i]; /* - f1[i] - g1[i] */
3017 f[a1k][i] = f[a1k][i] + f1_k[i]; /* h1[i] + g1[i] */
3018 f[a1l][i] = f[a1l][i] + f1_l[i]; /* h1[i] */
3021 /* Do forces - second torsion */
3022 fg2 = iprod(r2_ij, r2_kj);
3023 hg2 = iprod(r2_kl, r2_kj);
3024 fga2 = fg2*ra2r2*rgr2;
3025 hgb2 = hg2*rb2r2*rgr2;
3029 for (i = 0; i < DIM; i++)
3031 dtf2[i] = gaa2 * a2[i];
3032 dtg2[i] = fga2 * a2[i] - hgb2 * b2[i];
3033 dth2[i] = gbb2 * b2[i];
3035 f2[i] = df2 * dtf2[i];
3036 g2[i] = df2 * dtg2[i];
3037 h2[i] = df2 * dth2[i];
3040 f2_j[i] = -f2[i] - g2[i];
3041 f2_k[i] = h2[i] + g2[i];
3044 f[a2i][i] = f[a2i][i] + f2_i[i]; /* f2[i] */
3045 f[a2j][i] = f[a2j][i] + f2_j[i]; /* - f2[i] - g2[i] */
3046 f[a2k][i] = f[a2k][i] + f2_k[i]; /* h2[i] + g2[i] */
3047 f[a2l][i] = f[a2l][i] + f2_l[i]; /* - h2[i] */
3053 copy_ivec(SHIFT_IVEC(g, a1j), jt1);
3054 ivec_sub(SHIFT_IVEC(g, a1i), jt1, dt1_ij);
3055 ivec_sub(SHIFT_IVEC(g, a1k), jt1, dt1_kj);
3056 ivec_sub(SHIFT_IVEC(g, a1l), jt1, dt1_lj);
3057 t11 = IVEC2IS(dt1_ij);
3058 t21 = IVEC2IS(dt1_kj);
3059 t31 = IVEC2IS(dt1_lj);
3061 copy_ivec(SHIFT_IVEC(g, a2j), jt2);
3062 ivec_sub(SHIFT_IVEC(g, a2i), jt2, dt2_ij);
3063 ivec_sub(SHIFT_IVEC(g, a2k), jt2, dt2_kj);
3064 ivec_sub(SHIFT_IVEC(g, a2l), jt2, dt2_lj);
3065 t12 = IVEC2IS(dt2_ij);
3066 t22 = IVEC2IS(dt2_kj);
3067 t32 = IVEC2IS(dt2_lj);
3071 t31 = pbc_rvec_sub(pbc, x[a1l], x[a1j], h1);
3072 t32 = pbc_rvec_sub(pbc, x[a2l], x[a2j], h2);
3080 rvec_inc(fshift[t11], f1_i);
3081 rvec_inc(fshift[CENTRAL], f1_j);
3082 rvec_inc(fshift[t21], f1_k);
3083 rvec_inc(fshift[t31], f1_l);
3085 rvec_inc(fshift[t21], f2_i);
3086 rvec_inc(fshift[CENTRAL], f2_j);
3087 rvec_inc(fshift[t22], f2_k);
3088 rvec_inc(fshift[t32], f2_l);
3095 /***********************************************************
3097 * G R O M O S 9 6 F U N C T I O N S
3099 ***********************************************************/
3100 real g96harmonic(real kA, real kB, real xA, real xB, real x, real lambda,
3103 const real half = 0.5;
3104 real L1, kk, x0, dx, dx2;
3105 real v, f, dvdlambda;
3108 kk = L1*kA+lambda*kB;
3109 x0 = L1*xA+lambda*xB;
3116 dvdlambda = half*(kB-kA)*dx2 + (xA-xB)*kk*dx;
3123 /* That was 21 flops */
3126 real g96bonds(int nbonds,
3127 const t_iatom forceatoms[], const t_iparams forceparams[],
3128 const rvec x[], rvec f[], rvec fshift[],
3129 const t_pbc *pbc, const t_graph *g,
3130 real lambda, real *dvdlambda,
3131 const t_mdatoms *md, t_fcdata *fcd,
3132 int *global_atom_index)
3134 int i, m, ki, ai, aj, type;
3135 real dr2, fbond, vbond, fij, vtot;
3140 for (i = 0; (i < nbonds); )
3142 type = forceatoms[i++];
3143 ai = forceatoms[i++];
3144 aj = forceatoms[i++];
3146 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
3147 dr2 = iprod(dx, dx); /* 5 */
3149 *dvdlambda += g96harmonic(forceparams[type].harmonic.krA,
3150 forceparams[type].harmonic.krB,
3151 forceparams[type].harmonic.rA,
3152 forceparams[type].harmonic.rB,
3153 dr2, lambda, &vbond, &fbond);
3155 vtot += 0.5*vbond; /* 1*/
3159 fprintf(debug, "G96-BONDS: dr = %10g vbond = %10g fbond = %10g\n",
3160 sqrt(dr2), vbond, fbond);
3166 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
3169 for (m = 0; (m < DIM); m++) /* 15 */
3174 fshift[ki][m] += fij;
3175 fshift[CENTRAL][m] -= fij;
3181 real g96bond_angle(const rvec xi, const rvec xj, const rvec xk, const t_pbc *pbc,
3182 rvec r_ij, rvec r_kj,
3184 /* Return value is the angle between the bonds i-j and j-k */
3188 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
3189 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
3191 costh = cos_angle(r_ij, r_kj); /* 25 */
3196 real g96angles(int nbonds,
3197 const t_iatom forceatoms[], const t_iparams forceparams[],
3198 const rvec x[], rvec f[], rvec fshift[],
3199 const t_pbc *pbc, const t_graph *g,
3200 real lambda, real *dvdlambda,
3201 const t_mdatoms *md, t_fcdata *fcd,
3202 int *global_atom_index)
3204 int i, ai, aj, ak, type, m, t1, t2;
3206 real cos_theta, dVdt, va, vtot;
3207 real rij_1, rij_2, rkj_1, rkj_2, rijrkj_1;
3209 ivec jt, dt_ij, dt_kj;
3212 for (i = 0; (i < nbonds); )
3214 type = forceatoms[i++];
3215 ai = forceatoms[i++];
3216 aj = forceatoms[i++];
3217 ak = forceatoms[i++];
3219 cos_theta = g96bond_angle(x[ai], x[aj], x[ak], pbc, r_ij, r_kj, &t1, &t2);
3221 *dvdlambda += g96harmonic(forceparams[type].harmonic.krA,
3222 forceparams[type].harmonic.krB,
3223 forceparams[type].harmonic.rA,
3224 forceparams[type].harmonic.rB,
3225 cos_theta, lambda, &va, &dVdt);
3228 rij_1 = gmx_invsqrt(iprod(r_ij, r_ij));
3229 rkj_1 = gmx_invsqrt(iprod(r_kj, r_kj));
3230 rij_2 = rij_1*rij_1;
3231 rkj_2 = rkj_1*rkj_1;
3232 rijrkj_1 = rij_1*rkj_1; /* 23 */
3237 fprintf(debug, "G96ANGLES: costheta = %10g vth = %10g dV/dct = %10g\n",
3238 cos_theta, va, dVdt);
3241 for (m = 0; (m < DIM); m++) /* 42 */
3243 f_i[m] = dVdt*(r_kj[m]*rijrkj_1 - r_ij[m]*rij_2*cos_theta);
3244 f_k[m] = dVdt*(r_ij[m]*rijrkj_1 - r_kj[m]*rkj_2*cos_theta);
3245 f_j[m] = -f_i[m]-f_k[m];
3253 copy_ivec(SHIFT_IVEC(g, aj), jt);
3255 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3256 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3257 t1 = IVEC2IS(dt_ij);
3258 t2 = IVEC2IS(dt_kj);
3260 rvec_inc(fshift[t1], f_i);
3261 rvec_inc(fshift[CENTRAL], f_j);
3262 rvec_inc(fshift[t2], f_k); /* 9 */
3268 real cross_bond_bond(int nbonds,
3269 const t_iatom forceatoms[], const t_iparams forceparams[],
3270 const rvec x[], rvec f[], rvec fshift[],
3271 const t_pbc *pbc, const t_graph *g,
3272 real lambda, real *dvdlambda,
3273 const t_mdatoms *md, t_fcdata *fcd,
3274 int *global_atom_index)
3276 /* Potential from Lawrence and Skimmer, Chem. Phys. Lett. 372 (2003)
3279 int i, ai, aj, ak, type, m, t1, t2;
3281 real vtot, vrr, s1, s2, r1, r2, r1e, r2e, krr;
3283 ivec jt, dt_ij, dt_kj;
3286 for (i = 0; (i < nbonds); )
3288 type = forceatoms[i++];
3289 ai = forceatoms[i++];
3290 aj = forceatoms[i++];
3291 ak = forceatoms[i++];
3292 r1e = forceparams[type].cross_bb.r1e;
3293 r2e = forceparams[type].cross_bb.r2e;
3294 krr = forceparams[type].cross_bb.krr;
3296 /* Compute distance vectors ... */
3297 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
3298 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
3300 /* ... and their lengths */
3304 /* Deviations from ideality */
3308 /* Energy (can be negative!) */
3313 svmul(-krr*s2/r1, r_ij, f_i);
3314 svmul(-krr*s1/r2, r_kj, f_k);
3316 for (m = 0; (m < DIM); m++) /* 12 */
3318 f_j[m] = -f_i[m] - f_k[m];
3327 copy_ivec(SHIFT_IVEC(g, aj), jt);
3329 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3330 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3331 t1 = IVEC2IS(dt_ij);
3332 t2 = IVEC2IS(dt_kj);
3334 rvec_inc(fshift[t1], f_i);
3335 rvec_inc(fshift[CENTRAL], f_j);
3336 rvec_inc(fshift[t2], f_k); /* 9 */
3342 real cross_bond_angle(int nbonds,
3343 const t_iatom forceatoms[], const t_iparams forceparams[],
3344 const rvec x[], rvec f[], rvec fshift[],
3345 const t_pbc *pbc, const t_graph *g,
3346 real lambda, real *dvdlambda,
3347 const t_mdatoms *md, t_fcdata *fcd,
3348 int *global_atom_index)
3350 /* Potential from Lawrence and Skimmer, Chem. Phys. Lett. 372 (2003)
3353 int i, ai, aj, ak, type, m, t1, t2, t3;
3354 rvec r_ij, r_kj, r_ik;
3355 real vtot, vrt, s1, s2, s3, r1, r2, r3, r1e, r2e, r3e, krt, k1, k2, k3;
3357 ivec jt, dt_ij, dt_kj;
3360 for (i = 0; (i < nbonds); )
3362 type = forceatoms[i++];
3363 ai = forceatoms[i++];
3364 aj = forceatoms[i++];
3365 ak = forceatoms[i++];
3366 r1e = forceparams[type].cross_ba.r1e;
3367 r2e = forceparams[type].cross_ba.r2e;
3368 r3e = forceparams[type].cross_ba.r3e;
3369 krt = forceparams[type].cross_ba.krt;
3371 /* Compute distance vectors ... */
3372 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
3373 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
3374 t3 = pbc_rvec_sub(pbc, x[ai], x[ak], r_ik);
3376 /* ... and their lengths */
3381 /* Deviations from ideality */
3386 /* Energy (can be negative!) */
3387 vrt = krt*s3*(s1+s2);
3393 k3 = -krt*(s1+s2)/r3;
3394 for (m = 0; (m < DIM); m++)
3396 f_i[m] = k1*r_ij[m] + k3*r_ik[m];
3397 f_k[m] = k2*r_kj[m] - k3*r_ik[m];
3398 f_j[m] = -f_i[m] - f_k[m];
3401 for (m = 0; (m < DIM); m++) /* 12 */
3411 copy_ivec(SHIFT_IVEC(g, aj), jt);
3413 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3414 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3415 t1 = IVEC2IS(dt_ij);
3416 t2 = IVEC2IS(dt_kj);
3418 rvec_inc(fshift[t1], f_i);
3419 rvec_inc(fshift[CENTRAL], f_j);
3420 rvec_inc(fshift[t2], f_k); /* 9 */
3426 static real bonded_tab(const char *type, int table_nr,
3427 const bondedtable_t *table, real kA, real kB, real r,
3428 real lambda, real *V, real *F)
3430 real k, tabscale, *VFtab, rt, eps, eps2, Yt, Ft, Geps, Heps2, Fp, VV, FF;
3432 real v, f, dvdlambda;
3434 k = (1.0 - lambda)*kA + lambda*kB;
3436 tabscale = table->scale;
3437 VFtab = table->data;
3443 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",
3444 type, table_nr, r, n0, n0+1, table->n);
3451 Geps = VFtab[nnn+2]*eps;
3452 Heps2 = VFtab[nnn+3]*eps2;
3453 Fp = Ft + Geps + Heps2;
3455 FF = Fp + Geps + 2.0*Heps2;
3457 *F = -k*FF*tabscale;
3459 dvdlambda = (kB - kA)*VV;
3463 /* That was 22 flops */
3466 real tab_bonds(int nbonds,
3467 const t_iatom forceatoms[], const t_iparams forceparams[],
3468 const rvec x[], rvec f[], rvec fshift[],
3469 const t_pbc *pbc, const t_graph *g,
3470 real lambda, real *dvdlambda,
3471 const t_mdatoms *md, t_fcdata *fcd,
3472 int *global_atom_index)
3474 int i, m, ki, ai, aj, type, table;
3475 real dr, dr2, fbond, vbond, fij, vtot;
3480 for (i = 0; (i < nbonds); )
3482 type = forceatoms[i++];
3483 ai = forceatoms[i++];
3484 aj = forceatoms[i++];
3486 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
3487 dr2 = iprod(dx, dx); /* 5 */
3488 dr = dr2*gmx_invsqrt(dr2); /* 10 */
3490 table = forceparams[type].tab.table;
3492 *dvdlambda += bonded_tab("bond", table,
3493 &fcd->bondtab[table],
3494 forceparams[type].tab.kA,
3495 forceparams[type].tab.kB,
3496 dr, lambda, &vbond, &fbond); /* 22 */
3504 vtot += vbond; /* 1*/
3505 fbond *= gmx_invsqrt(dr2); /* 6 */
3509 fprintf(debug, "TABBONDS: dr = %10g vbond = %10g fbond = %10g\n",
3515 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
3518 for (m = 0; (m < DIM); m++) /* 15 */
3523 fshift[ki][m] += fij;
3524 fshift[CENTRAL][m] -= fij;
3530 real tab_angles(int nbonds,
3531 const t_iatom forceatoms[], const t_iparams forceparams[],
3532 const rvec x[], rvec f[], rvec fshift[],
3533 const t_pbc *pbc, const t_graph *g,
3534 real lambda, real *dvdlambda,
3535 const t_mdatoms *md, t_fcdata *fcd,
3536 int *global_atom_index)
3538 int i, ai, aj, ak, t1, t2, type, table;
3540 real cos_theta, cos_theta2, theta, dVdt, va, vtot;
3541 ivec jt, dt_ij, dt_kj;
3544 for (i = 0; (i < nbonds); )
3546 type = forceatoms[i++];
3547 ai = forceatoms[i++];
3548 aj = forceatoms[i++];
3549 ak = forceatoms[i++];
3551 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
3552 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
3554 table = forceparams[type].tab.table;
3556 *dvdlambda += bonded_tab("angle", table,
3557 &fcd->angletab[table],
3558 forceparams[type].tab.kA,
3559 forceparams[type].tab.kB,
3560 theta, lambda, &va, &dVdt); /* 22 */
3563 cos_theta2 = sqr(cos_theta); /* 1 */
3572 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
3573 sth = st*cos_theta; /* 1 */
3577 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
3578 theta*RAD2DEG, va, dVdt);
3581 nrkj2 = iprod(r_kj, r_kj); /* 5 */
3582 nrij2 = iprod(r_ij, r_ij);
3584 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
3585 cii = sth/nrij2; /* 10 */
3586 ckk = sth/nrkj2; /* 10 */
3588 for (m = 0; (m < DIM); m++) /* 39 */
3590 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
3591 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
3592 f_j[m] = -f_i[m]-f_k[m];
3599 copy_ivec(SHIFT_IVEC(g, aj), jt);
3601 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3602 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3603 t1 = IVEC2IS(dt_ij);
3604 t2 = IVEC2IS(dt_kj);
3606 rvec_inc(fshift[t1], f_i);
3607 rvec_inc(fshift[CENTRAL], f_j);
3608 rvec_inc(fshift[t2], f_k);
3614 real tab_dihs(int nbonds,
3615 const t_iatom forceatoms[], const t_iparams forceparams[],
3616 const rvec x[], rvec f[], rvec fshift[],
3617 const t_pbc *pbc, const t_graph *g,
3618 real lambda, real *dvdlambda,
3619 const t_mdatoms *md, t_fcdata *fcd,
3620 int *global_atom_index)
3622 int i, type, ai, aj, ak, al, table;
3624 rvec r_ij, r_kj, r_kl, m, n;
3625 real phi, sign, ddphi, vpd, vtot;
3628 for (i = 0; (i < nbonds); )
3630 type = forceatoms[i++];
3631 ai = forceatoms[i++];
3632 aj = forceatoms[i++];
3633 ak = forceatoms[i++];
3634 al = forceatoms[i++];
3636 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
3637 &sign, &t1, &t2, &t3); /* 84 */
3639 table = forceparams[type].tab.table;
3641 /* Hopefully phi+M_PI never results in values < 0 */
3642 *dvdlambda += bonded_tab("dihedral", table,
3643 &fcd->dihtab[table],
3644 forceparams[type].tab.kA,
3645 forceparams[type].tab.kB,
3646 phi+M_PI, lambda, &vpd, &ddphi);
3649 do_dih_fup(ai, aj, ak, al, -ddphi, r_ij, r_kj, r_kl, m, n,
3650 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
3653 fprintf(debug, "pdih: (%d,%d,%d,%d) phi=%g\n",
3654 ai, aj, ak, al, phi);
3661 /* Return if this is a potential calculated in bondfree.c,
3662 * i.e. an interaction that actually calculates a potential and
3663 * works on multiple atoms (not e.g. a connection or a position restraint).
3665 static gmx_inline gmx_bool ftype_is_bonded_potential(int ftype)
3668 (interaction_function[ftype].flags & IF_BOND) &&
3669 !(ftype == F_CONNBONDS || ftype == F_POSRES) &&
3670 (ftype < F_GB12 || ftype > F_GB14);
3673 static void divide_bondeds_over_threads(t_idef *idef, int nthreads)
3680 idef->nthreads = nthreads;
3682 if (F_NRE*(nthreads+1) > idef->il_thread_division_nalloc)
3684 idef->il_thread_division_nalloc = F_NRE*(nthreads+1);
3685 snew(idef->il_thread_division, idef->il_thread_division_nalloc);
3688 for (ftype = 0; ftype < F_NRE; ftype++)
3690 if (ftype_is_bonded_potential(ftype))
3692 nat1 = interaction_function[ftype].nratoms + 1;
3694 for (t = 0; t <= nthreads; t++)
3696 /* Divide the interactions equally over the threads.
3697 * When the different types of bonded interactions
3698 * are distributed roughly equally over the threads,
3699 * this should lead to well localized output into
3700 * the force buffer on each thread.
3701 * If this is not the case, a more advanced scheme
3702 * (not implemented yet) will do better.
3704 il_nr_thread = (((idef->il[ftype].nr/nat1)*t)/nthreads)*nat1;
3706 /* Ensure that distance restraint pairs with the same label
3707 * end up on the same thread.
3708 * This is slighlty tricky code, since the next for iteration
3709 * may have an initial il_nr_thread lower than the final value
3710 * in the previous iteration, but this will anyhow be increased
3711 * to the approriate value again by this while loop.
3713 while (ftype == F_DISRES &&
3715 il_nr_thread < idef->il[ftype].nr &&
3716 idef->iparams[idef->il[ftype].iatoms[il_nr_thread]].disres.label ==
3717 idef->iparams[idef->il[ftype].iatoms[il_nr_thread-nat1]].disres.label)
3719 il_nr_thread += nat1;
3722 idef->il_thread_division[ftype*(nthreads+1)+t] = il_nr_thread;
3729 calc_bonded_reduction_mask(const t_idef *idef,
3734 int ftype, nb, nat1, nb0, nb1, i, a;
3738 for (ftype = 0; ftype < F_NRE; ftype++)
3740 if (ftype_is_bonded_potential(ftype))
3742 nb = idef->il[ftype].nr;
3745 nat1 = interaction_function[ftype].nratoms + 1;
3747 /* Divide this interaction equally over the threads.
3748 * This is not stored: should match division in calc_bonds.
3750 nb0 = idef->il_thread_division[ftype*(nt+1)+t];
3751 nb1 = idef->il_thread_division[ftype*(nt+1)+t+1];
3753 for (i = nb0; i < nb1; i += nat1)
3755 for (a = 1; a < nat1; a++)
3757 mask |= (1U << (idef->il[ftype].iatoms[i+a]>>shift));
3767 void setup_bonded_threading(t_forcerec *fr, t_idef *idef)
3769 #define MAX_BLOCK_BITS 32
3774 assert(fr->nthreads >= 1);
3777 /* Divide the bonded interaction over the threads */
3778 divide_bondeds_over_threads(idef, fr->nthreads);
3780 if (fr->nthreads == 1)
3787 /* We divide the force array in a maximum of 32 blocks.
3788 * Minimum force block reduction size is 2^6=64.
3791 while (fr->natoms_force > (int)(MAX_BLOCK_BITS*(1U<<fr->red_ashift)))
3797 fprintf(debug, "bonded force buffer block atom shift %d bits\n",
3801 /* Determine to which blocks each thread's bonded force calculation
3802 * contributes. Store this is a mask for each thread.
3804 #pragma omp parallel for num_threads(fr->nthreads) schedule(static)
3805 for (t = 1; t < fr->nthreads; t++)
3807 fr->f_t[t].red_mask =
3808 calc_bonded_reduction_mask(idef, fr->red_ashift, t, fr->nthreads);
3811 /* Determine the maximum number of blocks we need to reduce over */
3814 for (t = 0; t < fr->nthreads; t++)
3817 for (b = 0; b < MAX_BLOCK_BITS; b++)
3819 if (fr->f_t[t].red_mask & (1U<<b))
3821 fr->red_nblock = max(fr->red_nblock, b+1);
3827 fprintf(debug, "thread %d flags %x count %d\n",
3828 t, fr->f_t[t].red_mask, c);
3834 fprintf(debug, "Number of blocks to reduce: %d of size %d\n",
3835 fr->red_nblock, 1<<fr->red_ashift);
3836 fprintf(debug, "Reduction density %.2f density/#thread %.2f\n",
3837 ctot*(1<<fr->red_ashift)/(double)fr->natoms_force,
3838 ctot*(1<<fr->red_ashift)/(double)(fr->natoms_force*fr->nthreads));
3842 static void zero_thread_forces(f_thread_t *f_t, int n,
3843 int nblock, int blocksize)
3845 int b, a0, a1, a, i, j;
3847 if (n > f_t->f_nalloc)
3849 f_t->f_nalloc = over_alloc_large(n);
3850 srenew(f_t->f, f_t->f_nalloc);
3853 if (f_t->red_mask != 0)
3855 for (b = 0; b < nblock; b++)
3857 if (f_t->red_mask && (1U<<b))
3860 a1 = min((b+1)*blocksize, n);
3861 for (a = a0; a < a1; a++)
3863 clear_rvec(f_t->f[a]);
3868 for (i = 0; i < SHIFTS; i++)
3870 clear_rvec(f_t->fshift[i]);
3872 for (i = 0; i < F_NRE; i++)
3876 for (i = 0; i < egNR; i++)
3878 for (j = 0; j < f_t->grpp.nener; j++)
3880 f_t->grpp.ener[i][j] = 0;
3883 for (i = 0; i < efptNR; i++)
3889 static void reduce_thread_force_buffer(int n, rvec *f,
3890 int nthreads, f_thread_t *f_t,
3891 int nblock, int block_size)
3893 /* The max thread number is arbitrary,
3894 * we used a fixed number to avoid memory management.
3895 * Using more than 16 threads is probably never useful performance wise.
3897 #define MAX_BONDED_THREADS 256
3900 if (nthreads > MAX_BONDED_THREADS)
3902 gmx_fatal(FARGS, "Can not reduce bonded forces on more than %d threads",
3903 MAX_BONDED_THREADS);
3906 /* This reduction can run on any number of threads,
3907 * independently of nthreads.
3909 #pragma omp parallel for num_threads(nthreads) schedule(static)
3910 for (b = 0; b < nblock; b++)
3912 rvec *fp[MAX_BONDED_THREADS];
3916 /* Determine which threads contribute to this block */
3918 for (ft = 1; ft < nthreads; ft++)
3920 if (f_t[ft].red_mask & (1U<<b))
3922 fp[nfb++] = f_t[ft].f;
3927 /* Reduce force buffers for threads that contribute */
3929 a1 = (b+1)*block_size;
3931 for (a = a0; a < a1; a++)
3933 for (fb = 0; fb < nfb; fb++)
3935 rvec_inc(f[a], fp[fb][a]);
3942 static void reduce_thread_forces(int n, rvec *f, rvec *fshift,
3943 real *ener, gmx_grppairener_t *grpp, real *dvdl,
3944 int nthreads, f_thread_t *f_t,
3945 int nblock, int block_size,
3946 gmx_bool bCalcEnerVir,
3951 /* Reduce the bonded force buffer */
3952 reduce_thread_force_buffer(n, f, nthreads, f_t, nblock, block_size);
3955 /* When necessary, reduce energy and virial using one thread only */
3960 for (i = 0; i < SHIFTS; i++)
3962 for (t = 1; t < nthreads; t++)
3964 rvec_inc(fshift[i], f_t[t].fshift[i]);
3967 for (i = 0; i < F_NRE; i++)
3969 for (t = 1; t < nthreads; t++)
3971 ener[i] += f_t[t].ener[i];
3974 for (i = 0; i < egNR; i++)
3976 for (j = 0; j < f_t[1].grpp.nener; j++)
3978 for (t = 1; t < nthreads; t++)
3981 grpp->ener[i][j] += f_t[t].grpp.ener[i][j];
3987 for (i = 0; i < efptNR; i++)
3990 for (t = 1; t < nthreads; t++)
3992 dvdl[i] += f_t[t].dvdl[i];
3999 static real calc_one_bond(FILE *fplog, int thread,
4000 int ftype, const t_idef *idef,
4001 rvec x[], rvec f[], rvec fshift[],
4003 const t_pbc *pbc, const t_graph *g,
4004 gmx_grppairener_t *grpp,
4006 real *lambda, real *dvdl,
4007 const t_mdatoms *md, t_fcdata *fcd,
4008 gmx_bool bCalcEnerVir,
4009 int *global_atom_index, gmx_bool bPrintSepPot)
4011 int nat1, nbonds, efptFTYPE;
4016 if (IS_RESTRAINT_TYPE(ftype))
4018 efptFTYPE = efptRESTRAINT;
4022 efptFTYPE = efptBONDED;
4025 nat1 = interaction_function[ftype].nratoms + 1;
4026 nbonds = idef->il[ftype].nr/nat1;
4027 iatoms = idef->il[ftype].iatoms;
4029 nb0 = idef->il_thread_division[ftype*(idef->nthreads+1)+thread];
4030 nbn = idef->il_thread_division[ftype*(idef->nthreads+1)+thread+1] - nb0;
4032 if (!IS_LISTED_LJ_C(ftype))
4034 if (ftype == F_CMAP)
4036 v = cmap_dihs(nbn, iatoms+nb0,
4037 idef->iparams, &idef->cmap_grid,
4038 (const rvec*)x, f, fshift,
4039 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
4040 md, fcd, global_atom_index);
4043 else if (ftype == F_ANGLES &&
4044 !bCalcEnerVir && fr->efep == efepNO)
4046 /* No energies, shift forces, dvdl */
4047 angles_noener_simd(nbn, idef->il[ftype].iatoms+nb0,
4050 pbc, g, lambda[efptFTYPE], md, fcd,
4055 else if (ftype == F_PDIHS &&
4056 !bCalcEnerVir && fr->efep == efepNO)
4058 /* No energies, shift forces, dvdl */
4059 #ifndef SIMD_BONDEDS
4064 (nbn, idef->il[ftype].iatoms+nb0,
4067 pbc, g, lambda[efptFTYPE], md, fcd,
4073 v = interaction_function[ftype].ifunc(nbn, iatoms+nb0,
4075 (const rvec*)x, f, fshift,
4076 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
4077 md, fcd, global_atom_index);
4081 fprintf(fplog, " %-23s #%4d V %12.5e dVdl %12.5e\n",
4082 interaction_function[ftype].longname,
4083 nbonds, v, lambda[efptFTYPE]);
4088 v = do_nonbonded_listed(ftype, nbn, iatoms+nb0, idef->iparams, (const rvec*)x, f, fshift,
4089 pbc, g, lambda, dvdl, md, fr, grpp, global_atom_index);
4093 fprintf(fplog, " %-5s + %-15s #%4d dVdl %12.5e\n",
4094 interaction_function[ftype].longname,
4095 interaction_function[F_LJ14].longname, nbonds, dvdl[efptVDW]);
4096 fprintf(fplog, " %-5s + %-15s #%4d dVdl %12.5e\n",
4097 interaction_function[ftype].longname,
4098 interaction_function[F_COUL14].longname, nbonds, dvdl[efptCOUL]);
4104 inc_nrnb(nrnb, interaction_function[ftype].nrnb_ind, nbonds);
4110 void calc_bonds(FILE *fplog, const gmx_multisim_t *ms,
4112 rvec x[], history_t *hist,
4113 rvec f[], t_forcerec *fr,
4114 const t_pbc *pbc, const t_graph *g,
4115 gmx_enerdata_t *enerd, t_nrnb *nrnb,
4117 const t_mdatoms *md,
4118 t_fcdata *fcd, int *global_atom_index,
4119 t_atomtypes *atype, gmx_genborn_t *born,
4121 gmx_bool bPrintSepPot, gmx_large_int_t step)
4123 gmx_bool bCalcEnerVir;
4125 real v, dvdl[efptNR], dvdl_dum[efptNR]; /* The dummy array is to have a place to store the dhdl at other values
4126 of lambda, which will be thrown away in the end*/
4127 const t_pbc *pbc_null;
4132 assert(fr->nthreads == idef->nthreads);
4135 bCalcEnerVir = (force_flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY));
4137 for (i = 0; i < efptNR; i++)
4151 fprintf(fplog, "Step %s: bonded V and dVdl for this node\n",
4152 gmx_step_str(step, buf));
4158 p_graph(debug, "Bondage is fun", g);
4162 /* Do pre force calculation stuff which might require communication */
4163 if (idef->il[F_ORIRES].nr)
4165 enerd->term[F_ORIRESDEV] =
4166 calc_orires_dev(ms, idef->il[F_ORIRES].nr,
4167 idef->il[F_ORIRES].iatoms,
4168 idef->iparams, md, (const rvec*)x,
4169 pbc_null, fcd, hist);
4171 if (idef->il[F_DISRES].nr)
4173 calc_disres_R_6(ms, idef->il[F_DISRES].nr,
4174 idef->il[F_DISRES].iatoms,
4175 idef->iparams, (const rvec*)x, pbc_null,
4179 #pragma omp parallel for num_threads(fr->nthreads) schedule(static)
4180 for (thread = 0; thread < fr->nthreads; thread++)
4187 gmx_grppairener_t *grpp;
4192 fshift = fr->fshift;
4194 grpp = &enerd->grpp;
4199 zero_thread_forces(&fr->f_t[thread], fr->natoms_force,
4200 fr->red_nblock, 1<<fr->red_ashift);
4202 ft = fr->f_t[thread].f;
4203 fshift = fr->f_t[thread].fshift;
4204 epot = fr->f_t[thread].ener;
4205 grpp = &fr->f_t[thread].grpp;
4206 dvdlt = fr->f_t[thread].dvdl;
4208 /* Loop over all bonded force types to calculate the bonded forces */
4209 for (ftype = 0; (ftype < F_NRE); ftype++)
4211 if (idef->il[ftype].nr > 0 && ftype_is_bonded_potential(ftype))
4213 v = calc_one_bond(fplog, thread, ftype, idef, x,
4214 ft, fshift, fr, pbc_null, g, grpp,
4215 nrnb, lambda, dvdlt,
4216 md, fcd, bCalcEnerVir,
4217 global_atom_index, bPrintSepPot);
4222 if (fr->nthreads > 1)
4224 reduce_thread_forces(fr->natoms_force, f, fr->fshift,
4225 enerd->term, &enerd->grpp, dvdl,
4226 fr->nthreads, fr->f_t,
4227 fr->red_nblock, 1<<fr->red_ashift,
4229 force_flags & GMX_FORCE_DHDL);
4231 if (force_flags & GMX_FORCE_DHDL)
4233 for (i = 0; i < efptNR; i++)
4235 enerd->dvdl_nonlin[i] += dvdl[i];
4239 /* Copy the sum of violations for the distance restraints from fcd */
4242 enerd->term[F_DISRESVIOL] = fcd->disres.sumviol;
4247 void calc_bonds_lambda(FILE *fplog,
4251 const t_pbc *pbc, const t_graph *g,
4252 gmx_grppairener_t *grpp, real *epot, t_nrnb *nrnb,
4254 const t_mdatoms *md,
4256 int *global_atom_index)
4258 int i, ftype, nr_nonperturbed, nr;
4260 real dvdl_dum[efptNR];
4262 const t_pbc *pbc_null;
4274 /* Copy the whole idef, so we can modify the contents locally */
4276 idef_fe.nthreads = 1;
4277 snew(idef_fe.il_thread_division, F_NRE*(idef_fe.nthreads+1));
4279 /* We already have the forces, so we use temp buffers here */
4280 snew(f, fr->natoms_force);
4281 snew(fshift, SHIFTS);
4283 /* Loop over all bonded force types to calculate the bonded energies */
4284 for (ftype = 0; (ftype < F_NRE); ftype++)
4286 if (ftype_is_bonded_potential(ftype))
4288 /* Set the work range of thread 0 to the perturbed bondeds only */
4289 nr_nonperturbed = idef->il[ftype].nr_nonperturbed;
4290 nr = idef->il[ftype].nr;
4291 idef_fe.il_thread_division[ftype*2+0] = nr_nonperturbed;
4292 idef_fe.il_thread_division[ftype*2+1] = nr;
4294 /* This is only to get the flop count correct */
4295 idef_fe.il[ftype].nr = nr - nr_nonperturbed;
4297 if (nr - nr_nonperturbed > 0)
4299 v = calc_one_bond(fplog, 0, ftype, &idef_fe,
4300 x, f, fshift, fr, pbc_null, g,
4301 grpp, nrnb, lambda, dvdl_dum,
4303 global_atom_index, FALSE);
4312 sfree(idef_fe.il_thread_division);