2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team,
6 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
54 #include "gmx_fatal.h"
60 #include "nonbonded.h"
62 /* Include the SIMD macro file and then check for support */
63 #include "gmx_simd_macros.h"
64 #if defined GMX_HAVE_SIMD_MACROS && defined GMX_SIMD_HAVE_TRIGONOMETRIC
66 #include "gmx_simd_vec.h"
69 /* Find a better place for this? */
70 const int cmap_coeff_matrix[] = {
71 1, 0, -3, 2, 0, 0, 0, 0, -3, 0, 9, -6, 2, 0, -6, 4,
72 0, 0, 0, 0, 0, 0, 0, 0, 3, 0, -9, 6, -2, 0, 6, -4,
73 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 9, -6, 0, 0, -6, 4,
74 0, 0, 3, -2, 0, 0, 0, 0, 0, 0, -9, 6, 0, 0, 6, -4,
75 0, 0, 0, 0, 1, 0, -3, 2, -2, 0, 6, -4, 1, 0, -3, 2,
76 0, 0, 0, 0, 0, 0, 0, 0, -1, 0, 3, -2, 1, 0, -3, 2,
77 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -3, 2, 0, 0, 3, -2,
78 0, 0, 0, 0, 0, 0, 3, -2, 0, 0, -6, 4, 0, 0, 3, -2,
79 0, 1, -2, 1, 0, 0, 0, 0, 0, -3, 6, -3, 0, 2, -4, 2,
80 0, 0, 0, 0, 0, 0, 0, 0, 0, 3, -6, 3, 0, -2, 4, -2,
81 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -3, 3, 0, 0, 2, -2,
82 0, 0, -1, 1, 0, 0, 0, 0, 0, 0, 3, -3, 0, 0, -2, 2,
83 0, 0, 0, 0, 0, 1, -2, 1, 0, -2, 4, -2, 0, 1, -2, 1,
84 0, 0, 0, 0, 0, 0, 0, 0, 0, -1, 2, -1, 0, 1, -2, 1,
85 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, -1, 0, 0, -1, 1,
86 0, 0, 0, 0, 0, 0, -1, 1, 0, 0, 2, -2, 0, 0, -1, 1
91 int glatnr(int *global_atom_index, int i)
95 if (global_atom_index == NULL)
101 atnr = global_atom_index[i] + 1;
107 static int pbc_rvec_sub(const t_pbc *pbc, const rvec xi, const rvec xj, rvec dx)
111 return pbc_dx_aiuc(pbc, xi, xj, dx);
115 rvec_sub(xi, xj, dx);
122 /* SIMD PBC data structure, containing 1/boxdiag and the box vectors */
135 /* Set the SIMD pbc data from a normal t_pbc struct */
136 static void set_pbc_simd(const t_pbc *pbc, pbc_simd_t *pbc_simd)
141 /* Setting inv_bdiag to 0 effectively turns off PBC */
142 clear_rvec(inv_bdiag);
145 for (d = 0; d < pbc->ndim_ePBC; d++)
147 inv_bdiag[d] = 1.0/pbc->box[d][d];
151 pbc_simd->inv_bzz = gmx_set1_pr(inv_bdiag[ZZ]);
152 pbc_simd->inv_byy = gmx_set1_pr(inv_bdiag[YY]);
153 pbc_simd->inv_bxx = gmx_set1_pr(inv_bdiag[XX]);
157 pbc_simd->bzx = gmx_set1_pr(pbc->box[ZZ][XX]);
158 pbc_simd->bzy = gmx_set1_pr(pbc->box[ZZ][YY]);
159 pbc_simd->bzz = gmx_set1_pr(pbc->box[ZZ][ZZ]);
160 pbc_simd->byx = gmx_set1_pr(pbc->box[YY][XX]);
161 pbc_simd->byy = gmx_set1_pr(pbc->box[YY][YY]);
162 pbc_simd->bxx = gmx_set1_pr(pbc->box[XX][XX]);
166 pbc_simd->bzx = gmx_setzero_pr();
167 pbc_simd->bzy = gmx_setzero_pr();
168 pbc_simd->bzz = gmx_setzero_pr();
169 pbc_simd->byx = gmx_setzero_pr();
170 pbc_simd->byy = gmx_setzero_pr();
171 pbc_simd->bxx = gmx_setzero_pr();
175 /* Correct distance vector *dx,*dy,*dz for PBC using SIMD */
176 static gmx_inline void
177 pbc_dx_simd(gmx_mm_pr *dx, gmx_mm_pr *dy, gmx_mm_pr *dz,
178 const pbc_simd_t *pbc)
182 sh = gmx_round_pr(gmx_mul_pr(*dz, pbc->inv_bzz));
183 *dx = gmx_nmsub_pr(sh, pbc->bzx, *dx);
184 *dy = gmx_nmsub_pr(sh, pbc->bzy, *dy);
185 *dz = gmx_nmsub_pr(sh, pbc->bzz, *dz);
187 sh = gmx_round_pr(gmx_mul_pr(*dy, pbc->inv_byy));
188 *dx = gmx_nmsub_pr(sh, pbc->byx, *dx);
189 *dy = gmx_nmsub_pr(sh, pbc->byy, *dy);
191 sh = gmx_round_pr(gmx_mul_pr(*dx, pbc->inv_bxx));
192 *dx = gmx_nmsub_pr(sh, pbc->bxx, *dx);
195 #endif /* SIMD_BONDEDS */
198 * Morse potential bond by Frank Everdij
200 * Three parameters needed:
202 * b0 = equilibrium distance in nm
203 * be = beta in nm^-1 (actually, it's nu_e*Sqrt(2*pi*pi*mu/D_e))
204 * cb = well depth in kJ/mol
206 * Note: the potential is referenced to be +cb at infinite separation
207 * and zero at the equilibrium distance!
210 real morse_bonds(int nbonds,
211 const t_iatom forceatoms[], const t_iparams forceparams[],
212 const rvec x[], rvec f[], rvec fshift[],
213 const t_pbc *pbc, const t_graph *g,
214 real lambda, real *dvdlambda,
215 const t_mdatoms *md, t_fcdata *fcd,
216 int *global_atom_index)
218 const real one = 1.0;
219 const real two = 2.0;
220 real dr, dr2, temp, omtemp, cbomtemp, fbond, vbond, fij, vtot;
221 real b0, be, cb, b0A, beA, cbA, b0B, beB, cbB, L1;
223 int i, m, ki, type, ai, aj;
227 for (i = 0; (i < nbonds); )
229 type = forceatoms[i++];
230 ai = forceatoms[i++];
231 aj = forceatoms[i++];
233 b0A = forceparams[type].morse.b0A;
234 beA = forceparams[type].morse.betaA;
235 cbA = forceparams[type].morse.cbA;
237 b0B = forceparams[type].morse.b0B;
238 beB = forceparams[type].morse.betaB;
239 cbB = forceparams[type].morse.cbB;
241 L1 = one-lambda; /* 1 */
242 b0 = L1*b0A + lambda*b0B; /* 3 */
243 be = L1*beA + lambda*beB; /* 3 */
244 cb = L1*cbA + lambda*cbB; /* 3 */
246 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
247 dr2 = iprod(dx, dx); /* 5 */
248 dr = dr2*gmx_invsqrt(dr2); /* 10 */
249 temp = exp(-be*(dr-b0)); /* 12 */
253 /* bonds are constrainted. This may _not_ include bond constraints if they are lambda dependent */
254 *dvdlambda += cbB-cbA;
258 omtemp = one-temp; /* 1 */
259 cbomtemp = cb*omtemp; /* 1 */
260 vbond = cbomtemp*omtemp; /* 1 */
261 fbond = -two*be*temp*cbomtemp*gmx_invsqrt(dr2); /* 9 */
262 vtot += vbond; /* 1 */
264 *dvdlambda += (cbB - cbA) * omtemp * omtemp - (2-2*omtemp)*omtemp * cb * ((b0B-b0A)*be - (beB-beA)*(dr-b0)); /* 15 */
268 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
272 for (m = 0; (m < DIM); m++) /* 15 */
277 fshift[ki][m] += fij;
278 fshift[CENTRAL][m] -= fij;
284 real cubic_bonds(int nbonds,
285 const t_iatom forceatoms[], const t_iparams forceparams[],
286 const rvec x[], rvec f[], rvec fshift[],
287 const t_pbc *pbc, const t_graph *g,
288 real lambda, real *dvdlambda,
289 const t_mdatoms *md, t_fcdata *fcd,
290 int *global_atom_index)
292 const real three = 3.0;
293 const real two = 2.0;
295 real dr, dr2, dist, kdist, kdist2, fbond, vbond, fij, vtot;
297 int i, m, ki, type, ai, aj;
301 for (i = 0; (i < nbonds); )
303 type = forceatoms[i++];
304 ai = forceatoms[i++];
305 aj = forceatoms[i++];
307 b0 = forceparams[type].cubic.b0;
308 kb = forceparams[type].cubic.kb;
309 kcub = forceparams[type].cubic.kcub;
311 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
312 dr2 = iprod(dx, dx); /* 5 */
319 dr = dr2*gmx_invsqrt(dr2); /* 10 */
324 vbond = kdist2 + kcub*kdist2*dist;
325 fbond = -(two*kdist + three*kdist2*kcub)/dr;
327 vtot += vbond; /* 21 */
331 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
334 for (m = 0; (m < DIM); m++) /* 15 */
339 fshift[ki][m] += fij;
340 fshift[CENTRAL][m] -= fij;
346 real FENE_bonds(int nbonds,
347 const t_iatom forceatoms[], const t_iparams forceparams[],
348 const rvec x[], rvec f[], rvec fshift[],
349 const t_pbc *pbc, const t_graph *g,
350 real lambda, real *dvdlambda,
351 const t_mdatoms *md, t_fcdata *fcd,
352 int *global_atom_index)
354 const real half = 0.5;
355 const real one = 1.0;
357 real dr, dr2, bm2, omdr2obm2, fbond, vbond, fij, vtot;
359 int i, m, ki, type, ai, aj;
363 for (i = 0; (i < nbonds); )
365 type = forceatoms[i++];
366 ai = forceatoms[i++];
367 aj = forceatoms[i++];
369 bm = forceparams[type].fene.bm;
370 kb = forceparams[type].fene.kb;
372 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
373 dr2 = iprod(dx, dx); /* 5 */
385 "r^2 (%f) >= bm^2 (%f) in FENE bond between atoms %d and %d",
387 glatnr(global_atom_index, ai),
388 glatnr(global_atom_index, aj));
391 omdr2obm2 = one - dr2/bm2;
393 vbond = -half*kb*bm2*log(omdr2obm2);
394 fbond = -kb/omdr2obm2;
396 vtot += vbond; /* 35 */
400 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
403 for (m = 0; (m < DIM); m++) /* 15 */
408 fshift[ki][m] += fij;
409 fshift[CENTRAL][m] -= fij;
415 real harmonic(real kA, real kB, real xA, real xB, real x, real lambda,
418 const real half = 0.5;
419 real L1, kk, x0, dx, dx2;
420 real v, f, dvdlambda;
423 kk = L1*kA+lambda*kB;
424 x0 = L1*xA+lambda*xB;
431 dvdlambda = half*(kB-kA)*dx2 + (xA-xB)*kk*dx;
438 /* That was 19 flops */
442 real bonds(int nbonds,
443 const t_iatom forceatoms[], const t_iparams forceparams[],
444 const rvec x[], rvec f[], rvec fshift[],
445 const t_pbc *pbc, const t_graph *g,
446 real lambda, real *dvdlambda,
447 const t_mdatoms *md, t_fcdata *fcd,
448 int *global_atom_index)
450 int i, m, ki, ai, aj, type;
451 real dr, dr2, fbond, vbond, fij, vtot;
456 for (i = 0; (i < nbonds); )
458 type = forceatoms[i++];
459 ai = forceatoms[i++];
460 aj = forceatoms[i++];
462 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
463 dr2 = iprod(dx, dx); /* 5 */
464 dr = dr2*gmx_invsqrt(dr2); /* 10 */
466 *dvdlambda += harmonic(forceparams[type].harmonic.krA,
467 forceparams[type].harmonic.krB,
468 forceparams[type].harmonic.rA,
469 forceparams[type].harmonic.rB,
470 dr, lambda, &vbond, &fbond); /* 19 */
478 vtot += vbond; /* 1*/
479 fbond *= gmx_invsqrt(dr2); /* 6 */
483 fprintf(debug, "BONDS: dr = %10g vbond = %10g fbond = %10g\n",
489 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
492 for (m = 0; (m < DIM); m++) /* 15 */
497 fshift[ki][m] += fij;
498 fshift[CENTRAL][m] -= fij;
504 real restraint_bonds(int nbonds,
505 const t_iatom forceatoms[], const t_iparams forceparams[],
506 const rvec x[], rvec f[], rvec fshift[],
507 const t_pbc *pbc, const t_graph *g,
508 real lambda, real *dvdlambda,
509 const t_mdatoms *md, t_fcdata *fcd,
510 int *global_atom_index)
512 int i, m, ki, ai, aj, type;
513 real dr, dr2, fbond, vbond, fij, vtot;
515 real low, dlow, up1, dup1, up2, dup2, k, dk;
523 for (i = 0; (i < nbonds); )
525 type = forceatoms[i++];
526 ai = forceatoms[i++];
527 aj = forceatoms[i++];
529 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
530 dr2 = iprod(dx, dx); /* 5 */
531 dr = dr2*gmx_invsqrt(dr2); /* 10 */
533 low = L1*forceparams[type].restraint.lowA + lambda*forceparams[type].restraint.lowB;
534 dlow = -forceparams[type].restraint.lowA + forceparams[type].restraint.lowB;
535 up1 = L1*forceparams[type].restraint.up1A + lambda*forceparams[type].restraint.up1B;
536 dup1 = -forceparams[type].restraint.up1A + forceparams[type].restraint.up1B;
537 up2 = L1*forceparams[type].restraint.up2A + lambda*forceparams[type].restraint.up2B;
538 dup2 = -forceparams[type].restraint.up2A + forceparams[type].restraint.up2B;
539 k = L1*forceparams[type].restraint.kA + lambda*forceparams[type].restraint.kB;
540 dk = -forceparams[type].restraint.kA + forceparams[type].restraint.kB;
549 *dvdlambda += 0.5*dk*drh2 - k*dlow*drh;
562 *dvdlambda += 0.5*dk*drh2 - k*dup1*drh;
567 vbond = k*(up2 - up1)*(0.5*(up2 - up1) + drh);
568 fbond = -k*(up2 - up1);
569 *dvdlambda += dk*(up2 - up1)*(0.5*(up2 - up1) + drh)
570 + k*(dup2 - dup1)*(up2 - up1 + drh)
571 - k*(up2 - up1)*dup2;
579 vtot += vbond; /* 1*/
580 fbond *= gmx_invsqrt(dr2); /* 6 */
584 fprintf(debug, "BONDS: dr = %10g vbond = %10g fbond = %10g\n",
590 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
593 for (m = 0; (m < DIM); m++) /* 15 */
598 fshift[ki][m] += fij;
599 fshift[CENTRAL][m] -= fij;
606 real polarize(int nbonds,
607 const t_iatom forceatoms[], const t_iparams forceparams[],
608 const rvec x[], rvec f[], rvec fshift[],
609 const t_pbc *pbc, const t_graph *g,
610 real lambda, real *dvdlambda,
611 const t_mdatoms *md, t_fcdata *fcd,
612 int *global_atom_index)
614 int i, m, ki, ai, aj, type;
615 real dr, dr2, fbond, vbond, fij, vtot, ksh;
620 for (i = 0; (i < nbonds); )
622 type = forceatoms[i++];
623 ai = forceatoms[i++];
624 aj = forceatoms[i++];
625 ksh = sqr(md->chargeA[aj])*ONE_4PI_EPS0/forceparams[type].polarize.alpha;
628 fprintf(debug, "POL: local ai = %d aj = %d ksh = %.3f\n", ai, aj, ksh);
631 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
632 dr2 = iprod(dx, dx); /* 5 */
633 dr = dr2*gmx_invsqrt(dr2); /* 10 */
635 *dvdlambda += harmonic(ksh, ksh, 0, 0, dr, lambda, &vbond, &fbond); /* 19 */
642 vtot += vbond; /* 1*/
643 fbond *= gmx_invsqrt(dr2); /* 6 */
647 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
650 for (m = 0; (m < DIM); m++) /* 15 */
655 fshift[ki][m] += fij;
656 fshift[CENTRAL][m] -= fij;
662 real anharm_polarize(int nbonds,
663 const t_iatom forceatoms[], const t_iparams forceparams[],
664 const rvec x[], rvec f[], rvec fshift[],
665 const t_pbc *pbc, const t_graph *g,
666 real lambda, real *dvdlambda,
667 const t_mdatoms *md, t_fcdata *fcd,
668 int *global_atom_index)
670 int i, m, ki, ai, aj, type;
671 real dr, dr2, fbond, vbond, fij, vtot, ksh, khyp, drcut, ddr, ddr3;
676 for (i = 0; (i < nbonds); )
678 type = forceatoms[i++];
679 ai = forceatoms[i++];
680 aj = forceatoms[i++];
681 ksh = sqr(md->chargeA[aj])*ONE_4PI_EPS0/forceparams[type].anharm_polarize.alpha; /* 7*/
682 khyp = forceparams[type].anharm_polarize.khyp;
683 drcut = forceparams[type].anharm_polarize.drcut;
686 fprintf(debug, "POL: local ai = %d aj = %d ksh = %.3f\n", ai, aj, ksh);
689 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
690 dr2 = iprod(dx, dx); /* 5 */
691 dr = dr2*gmx_invsqrt(dr2); /* 10 */
693 *dvdlambda += harmonic(ksh, ksh, 0, 0, dr, lambda, &vbond, &fbond); /* 19 */
704 vbond += khyp*ddr*ddr3;
705 fbond -= 4*khyp*ddr3;
707 fbond *= gmx_invsqrt(dr2); /* 6 */
708 vtot += vbond; /* 1*/
712 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
715 for (m = 0; (m < DIM); m++) /* 15 */
720 fshift[ki][m] += fij;
721 fshift[CENTRAL][m] -= fij;
727 real water_pol(int nbonds,
728 const t_iatom forceatoms[], const t_iparams forceparams[],
729 const rvec x[], rvec f[], rvec fshift[],
730 const t_pbc *pbc, const t_graph *g,
731 real lambda, real *dvdlambda,
732 const t_mdatoms *md, t_fcdata *fcd,
733 int *global_atom_index)
735 /* This routine implements anisotropic polarizibility for water, through
736 * a shell connected to a dummy with spring constant that differ in the
737 * three spatial dimensions in the molecular frame.
739 int i, m, aO, aH1, aH2, aD, aS, type, type0;
740 rvec dOH1, dOH2, dHH, dOD, dDS, nW, kk, dx, kdx, proj;
744 real vtot, fij, r_HH, r_OD, r_nW, tx, ty, tz, qS;
749 type0 = forceatoms[0];
751 qS = md->chargeA[aS];
752 kk[XX] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_x;
753 kk[YY] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_y;
754 kk[ZZ] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_z;
755 r_HH = 1.0/forceparams[type0].wpol.rHH;
756 r_OD = 1.0/forceparams[type0].wpol.rOD;
759 fprintf(debug, "WPOL: qS = %10.5f aS = %5d\n", qS, aS);
760 fprintf(debug, "WPOL: kk = %10.3f %10.3f %10.3f\n",
761 kk[XX], kk[YY], kk[ZZ]);
762 fprintf(debug, "WPOL: rOH = %10.3f rHH = %10.3f rOD = %10.3f\n",
763 forceparams[type0].wpol.rOH,
764 forceparams[type0].wpol.rHH,
765 forceparams[type0].wpol.rOD);
767 for (i = 0; (i < nbonds); i += 6)
769 type = forceatoms[i];
772 gmx_fatal(FARGS, "Sorry, type = %d, type0 = %d, file = %s, line = %d",
773 type, type0, __FILE__, __LINE__);
775 aO = forceatoms[i+1];
776 aH1 = forceatoms[i+2];
777 aH2 = forceatoms[i+3];
778 aD = forceatoms[i+4];
779 aS = forceatoms[i+5];
781 /* Compute vectors describing the water frame */
782 rvec_sub(x[aH1], x[aO], dOH1);
783 rvec_sub(x[aH2], x[aO], dOH2);
784 rvec_sub(x[aH2], x[aH1], dHH);
785 rvec_sub(x[aD], x[aO], dOD);
786 rvec_sub(x[aS], x[aD], dDS);
787 cprod(dOH1, dOH2, nW);
789 /* Compute inverse length of normal vector
790 * (this one could be precomputed, but I'm too lazy now)
792 r_nW = gmx_invsqrt(iprod(nW, nW));
793 /* This is for precision, but does not make a big difference,
796 r_OD = gmx_invsqrt(iprod(dOD, dOD));
798 /* Normalize the vectors in the water frame */
800 svmul(r_HH, dHH, dHH);
801 svmul(r_OD, dOD, dOD);
803 /* Compute displacement of shell along components of the vector */
804 dx[ZZ] = iprod(dDS, dOD);
805 /* Compute projection on the XY plane: dDS - dx[ZZ]*dOD */
806 for (m = 0; (m < DIM); m++)
808 proj[m] = dDS[m]-dx[ZZ]*dOD[m];
811 /*dx[XX] = iprod(dDS,nW);
812 dx[YY] = iprod(dDS,dHH);*/
813 dx[XX] = iprod(proj, nW);
814 for (m = 0; (m < DIM); m++)
816 proj[m] -= dx[XX]*nW[m];
818 dx[YY] = iprod(proj, dHH);
823 fprintf(debug, "WPOL: dx2=%10g dy2=%10g dz2=%10g sum=%10g dDS^2=%10g\n",
824 sqr(dx[XX]), sqr(dx[YY]), sqr(dx[ZZ]), iprod(dx, dx), iprod(dDS, dDS));
825 fprintf(debug, "WPOL: dHH=(%10g,%10g,%10g)\n", dHH[XX], dHH[YY], dHH[ZZ]);
826 fprintf(debug, "WPOL: dOD=(%10g,%10g,%10g), 1/r_OD = %10g\n",
827 dOD[XX], dOD[YY], dOD[ZZ], 1/r_OD);
828 fprintf(debug, "WPOL: nW =(%10g,%10g,%10g), 1/r_nW = %10g\n",
829 nW[XX], nW[YY], nW[ZZ], 1/r_nW);
830 fprintf(debug, "WPOL: dx =%10g, dy =%10g, dz =%10g\n",
831 dx[XX], dx[YY], dx[ZZ]);
832 fprintf(debug, "WPOL: dDSx=%10g, dDSy=%10g, dDSz=%10g\n",
833 dDS[XX], dDS[YY], dDS[ZZ]);
836 /* Now compute the forces and energy */
837 kdx[XX] = kk[XX]*dx[XX];
838 kdx[YY] = kk[YY]*dx[YY];
839 kdx[ZZ] = kk[ZZ]*dx[ZZ];
840 vtot += iprod(dx, kdx);
841 for (m = 0; (m < DIM); m++)
843 /* This is a tensor operation but written out for speed */
857 fprintf(debug, "WPOL: vwpol=%g\n", 0.5*iprod(dx, kdx));
858 fprintf(debug, "WPOL: df = (%10g, %10g, %10g)\n", df[XX], df[YY], df[ZZ]);
866 static real do_1_thole(const rvec xi, const rvec xj, rvec fi, rvec fj,
867 const t_pbc *pbc, real qq,
868 rvec fshift[], real afac)
871 real r12sq, r12_1, r12n, r12bar, v0, v1, fscal, ebar, fff;
874 t = pbc_rvec_sub(pbc, xi, xj, r12); /* 3 */
876 r12sq = iprod(r12, r12); /* 5 */
877 r12_1 = gmx_invsqrt(r12sq); /* 5 */
878 r12bar = afac/r12_1; /* 5 */
879 v0 = qq*ONE_4PI_EPS0*r12_1; /* 2 */
880 ebar = exp(-r12bar); /* 5 */
881 v1 = (1-(1+0.5*r12bar)*ebar); /* 4 */
882 fscal = ((v0*r12_1)*v1 - v0*0.5*afac*ebar*(r12bar+1))*r12_1; /* 9 */
885 fprintf(debug, "THOLE: v0 = %.3f v1 = %.3f r12= % .3f r12bar = %.3f fscal = %.3f ebar = %.3f\n", v0, v1, 1/r12_1, r12bar, fscal, ebar);
888 for (m = 0; (m < DIM); m++)
894 fshift[CENTRAL][m] -= fff;
897 return v0*v1; /* 1 */
901 real thole_pol(int nbonds,
902 const t_iatom forceatoms[], const t_iparams forceparams[],
903 const rvec x[], rvec f[], rvec fshift[],
904 const t_pbc *pbc, const t_graph *g,
905 real lambda, real *dvdlambda,
906 const t_mdatoms *md, t_fcdata *fcd,
907 int *global_atom_index)
909 /* Interaction between two pairs of particles with opposite charge */
910 int i, type, a1, da1, a2, da2;
911 real q1, q2, qq, a, al1, al2, afac;
914 for (i = 0; (i < nbonds); )
916 type = forceatoms[i++];
917 a1 = forceatoms[i++];
918 da1 = forceatoms[i++];
919 a2 = forceatoms[i++];
920 da2 = forceatoms[i++];
921 q1 = md->chargeA[da1];
922 q2 = md->chargeA[da2];
923 a = forceparams[type].thole.a;
924 al1 = forceparams[type].thole.alpha1;
925 al2 = forceparams[type].thole.alpha2;
927 afac = a*pow(al1*al2, -1.0/6.0);
928 V += do_1_thole(x[a1], x[a2], f[a1], f[a2], pbc, qq, fshift, afac);
929 V += do_1_thole(x[da1], x[a2], f[da1], f[a2], pbc, -qq, fshift, afac);
930 V += do_1_thole(x[a1], x[da2], f[a1], f[da2], pbc, -qq, fshift, afac);
931 V += do_1_thole(x[da1], x[da2], f[da1], f[da2], pbc, qq, fshift, afac);
937 real bond_angle(const rvec xi, const rvec xj, const rvec xk, const t_pbc *pbc,
938 rvec r_ij, rvec r_kj, real *costh,
940 /* Return value is the angle between the bonds i-j and j-k */
945 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
946 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
948 *costh = cos_angle(r_ij, r_kj); /* 25 */
949 th = acos(*costh); /* 10 */
954 real angles(int nbonds,
955 const t_iatom forceatoms[], const t_iparams forceparams[],
956 const rvec x[], rvec f[], rvec fshift[],
957 const t_pbc *pbc, const t_graph *g,
958 real lambda, real *dvdlambda,
959 const t_mdatoms *md, t_fcdata *fcd,
960 int *global_atom_index)
962 int i, ai, aj, ak, t1, t2, type;
964 real cos_theta, cos_theta2, theta, dVdt, va, vtot;
965 ivec jt, dt_ij, dt_kj;
968 for (i = 0; i < nbonds; )
970 type = forceatoms[i++];
971 ai = forceatoms[i++];
972 aj = forceatoms[i++];
973 ak = forceatoms[i++];
975 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
976 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
978 *dvdlambda += harmonic(forceparams[type].harmonic.krA,
979 forceparams[type].harmonic.krB,
980 forceparams[type].harmonic.rA*DEG2RAD,
981 forceparams[type].harmonic.rB*DEG2RAD,
982 theta, lambda, &va, &dVdt); /* 21 */
985 cos_theta2 = sqr(cos_theta);
995 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
996 sth = st*cos_theta; /* 1 */
1000 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1001 theta*RAD2DEG, va, dVdt);
1004 nrij2 = iprod(r_ij, r_ij); /* 5 */
1005 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1007 nrij_1 = gmx_invsqrt(nrij2); /* 10 */
1008 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1010 cik = st*nrij_1*nrkj_1; /* 2 */
1011 cii = sth*nrij_1*nrij_1; /* 2 */
1012 ckk = sth*nrkj_1*nrkj_1; /* 2 */
1014 for (m = 0; m < DIM; m++)
1016 f_i[m] = -(cik*r_kj[m] - cii*r_ij[m]);
1017 f_k[m] = -(cik*r_ij[m] - ckk*r_kj[m]);
1018 f_j[m] = -f_i[m] - f_k[m];
1025 copy_ivec(SHIFT_IVEC(g, aj), jt);
1027 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1028 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1029 t1 = IVEC2IS(dt_ij);
1030 t2 = IVEC2IS(dt_kj);
1032 rvec_inc(fshift[t1], f_i);
1033 rvec_inc(fshift[CENTRAL], f_j);
1034 rvec_inc(fshift[t2], f_k);
1043 /* As angles, but using SIMD to calculate many dihedrals at once.
1044 * This routines does not calculate energies and shift forces.
1046 static gmx_inline void
1047 angles_noener_simd(int nbonds,
1048 const t_iatom forceatoms[], const t_iparams forceparams[],
1049 const rvec x[], rvec f[],
1050 const t_pbc *pbc, const t_graph *g,
1052 const t_mdatoms *md, t_fcdata *fcd,
1053 int *global_atom_index)
1055 #define UNROLL GMX_SIMD_WIDTH_HERE
1058 int type, ai[UNROLL], aj[UNROLL], ak[UNROLL];
1059 real coeff_array[2*UNROLL+UNROLL], *coeff;
1060 real dr_array[2*DIM*UNROLL+UNROLL], *dr;
1061 real f_buf_array[6*UNROLL+UNROLL], *f_buf;
1062 gmx_mm_pr k_S, theta0_S;
1063 gmx_mm_pr rijx_S, rijy_S, rijz_S;
1064 gmx_mm_pr rkjx_S, rkjy_S, rkjz_S;
1066 gmx_mm_pr rij_rkj_S;
1067 gmx_mm_pr nrij2_S, nrij_1_S;
1068 gmx_mm_pr nrkj2_S, nrkj_1_S;
1069 gmx_mm_pr cos_S, sin_S;
1071 gmx_mm_pr st_S, sth_S;
1072 gmx_mm_pr cik_S, cii_S, ckk_S;
1073 gmx_mm_pr f_ix_S, f_iy_S, f_iz_S;
1074 gmx_mm_pr f_kx_S, f_ky_S, f_kz_S;
1075 pbc_simd_t pbc_simd;
1077 /* Ensure register memory alignment */
1078 coeff = gmx_simd_align_real(coeff_array);
1079 dr = gmx_simd_align_real(dr_array);
1080 f_buf = gmx_simd_align_real(f_buf_array);
1082 set_pbc_simd(pbc,&pbc_simd);
1084 one_S = gmx_set1_pr(1.0);
1086 /* nbonds is the number of angles times nfa1, here we step UNROLL angles */
1087 for (i = 0; (i < nbonds); i += UNROLL*nfa1)
1089 /* Collect atoms for UNROLL angles.
1090 * iu indexes into forceatoms, we should not let iu go beyond nbonds.
1093 for (s = 0; s < UNROLL; s++)
1095 type = forceatoms[iu];
1096 ai[s] = forceatoms[iu+1];
1097 aj[s] = forceatoms[iu+2];
1098 ak[s] = forceatoms[iu+3];
1100 coeff[s] = forceparams[type].harmonic.krA;
1101 coeff[UNROLL+s] = forceparams[type].harmonic.rA*DEG2RAD;
1103 /* If you can't use pbc_dx_simd below for PBC, e.g. because
1104 * you can't round in SIMD, use pbc_rvec_sub here.
1106 /* Store the non PBC corrected distances packed and aligned */
1107 for (m = 0; m < DIM; m++)
1109 dr[s + m *UNROLL] = x[ai[s]][m] - x[aj[s]][m];
1110 dr[s + (DIM+m)*UNROLL] = x[ak[s]][m] - x[aj[s]][m];
1113 /* At the end fill the arrays with identical entries */
1114 if (iu + nfa1 < nbonds)
1120 k_S = gmx_load_pr(coeff);
1121 theta0_S = gmx_load_pr(coeff+UNROLL);
1123 rijx_S = gmx_load_pr(dr + 0*UNROLL);
1124 rijy_S = gmx_load_pr(dr + 1*UNROLL);
1125 rijz_S = gmx_load_pr(dr + 2*UNROLL);
1126 rkjx_S = gmx_load_pr(dr + 3*UNROLL);
1127 rkjy_S = gmx_load_pr(dr + 4*UNROLL);
1128 rkjz_S = gmx_load_pr(dr + 5*UNROLL);
1130 pbc_dx_simd(&rijx_S, &rijy_S, &rijz_S, &pbc_simd);
1131 pbc_dx_simd(&rkjx_S, &rkjy_S, &rkjz_S, &pbc_simd);
1133 rij_rkj_S = gmx_iprod_pr(rijx_S, rijy_S, rijz_S,
1134 rkjx_S, rkjy_S, rkjz_S);
1136 nrij2_S = gmx_norm2_pr(rijx_S, rijy_S, rijz_S);
1137 nrkj2_S = gmx_norm2_pr(rkjx_S, rkjy_S, rkjz_S);
1139 nrij_1_S = gmx_invsqrt_pr(nrij2_S);
1140 nrkj_1_S = gmx_invsqrt_pr(nrkj2_S);
1142 cos_S = gmx_mul_pr(rij_rkj_S, gmx_mul_pr(nrij_1_S, nrkj_1_S));
1144 theta_S = gmx_acos_pr(cos_S);
1146 sin_S = gmx_invsqrt_pr(gmx_max_pr(gmx_sub_pr(one_S, gmx_mul_pr(cos_S, cos_S)),
1148 st_S = gmx_mul_pr(gmx_mul_pr(k_S, gmx_sub_pr(theta0_S, theta_S)),
1150 sth_S = gmx_mul_pr(st_S, cos_S);
1152 cik_S = gmx_mul_pr(st_S, gmx_mul_pr(nrij_1_S, nrkj_1_S));
1153 cii_S = gmx_mul_pr(sth_S, gmx_mul_pr(nrij_1_S, nrij_1_S));
1154 ckk_S = gmx_mul_pr(sth_S, gmx_mul_pr(nrkj_1_S, nrkj_1_S));
1156 f_ix_S = gmx_mul_pr(cii_S, rijx_S);
1157 f_ix_S = gmx_nmsub_pr(cik_S, rkjx_S, f_ix_S);
1158 f_iy_S = gmx_mul_pr(cii_S, rijy_S);
1159 f_iy_S = gmx_nmsub_pr(cik_S, rkjy_S, f_iy_S);
1160 f_iz_S = gmx_mul_pr(cii_S, rijz_S);
1161 f_iz_S = gmx_nmsub_pr(cik_S, rkjz_S, f_iz_S);
1162 f_kx_S = gmx_mul_pr(ckk_S, rkjx_S);
1163 f_kx_S = gmx_nmsub_pr(cik_S, rijx_S, f_kx_S);
1164 f_ky_S = gmx_mul_pr(ckk_S, rkjy_S);
1165 f_ky_S = gmx_nmsub_pr(cik_S, rijy_S, f_ky_S);
1166 f_kz_S = gmx_mul_pr(ckk_S, rkjz_S);
1167 f_kz_S = gmx_nmsub_pr(cik_S, rijz_S, f_kz_S);
1169 gmx_store_pr(f_buf + 0*UNROLL, f_ix_S);
1170 gmx_store_pr(f_buf + 1*UNROLL, f_iy_S);
1171 gmx_store_pr(f_buf + 2*UNROLL, f_iz_S);
1172 gmx_store_pr(f_buf + 3*UNROLL, f_kx_S);
1173 gmx_store_pr(f_buf + 4*UNROLL, f_ky_S);
1174 gmx_store_pr(f_buf + 5*UNROLL, f_kz_S);
1180 for (m = 0; m < DIM; m++)
1182 f[ai[s]][m] += f_buf[s + m*UNROLL];
1183 f[aj[s]][m] -= f_buf[s + m*UNROLL] + f_buf[s + (DIM+m)*UNROLL];
1184 f[ak[s]][m] += f_buf[s + (DIM+m)*UNROLL];
1189 while (s < UNROLL && iu < nbonds);
1194 #endif /* SIMD_BONDEDS */
1196 real linear_angles(int nbonds,
1197 const t_iatom forceatoms[], const t_iparams forceparams[],
1198 const rvec x[], rvec f[], rvec fshift[],
1199 const t_pbc *pbc, const t_graph *g,
1200 real lambda, real *dvdlambda,
1201 const t_mdatoms *md, t_fcdata *fcd,
1202 int *global_atom_index)
1204 int i, m, ai, aj, ak, t1, t2, type;
1206 real L1, kA, kB, aA, aB, dr, dr2, va, vtot, a, b, klin;
1207 ivec jt, dt_ij, dt_kj;
1208 rvec r_ij, r_kj, r_ik, dx;
1212 for (i = 0; (i < nbonds); )
1214 type = forceatoms[i++];
1215 ai = forceatoms[i++];
1216 aj = forceatoms[i++];
1217 ak = forceatoms[i++];
1219 kA = forceparams[type].linangle.klinA;
1220 kB = forceparams[type].linangle.klinB;
1221 klin = L1*kA + lambda*kB;
1223 aA = forceparams[type].linangle.aA;
1224 aB = forceparams[type].linangle.aB;
1225 a = L1*aA+lambda*aB;
1228 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
1229 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
1230 rvec_sub(r_ij, r_kj, r_ik);
1233 for (m = 0; (m < DIM); m++)
1235 dr = -a * r_ij[m] - b * r_kj[m];
1240 f_j[m] = -(f_i[m]+f_k[m]);
1246 *dvdlambda += 0.5*(kB-kA)*dr2 + klin*(aB-aA)*iprod(dx, r_ik);
1252 copy_ivec(SHIFT_IVEC(g, aj), jt);
1254 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1255 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1256 t1 = IVEC2IS(dt_ij);
1257 t2 = IVEC2IS(dt_kj);
1259 rvec_inc(fshift[t1], f_i);
1260 rvec_inc(fshift[CENTRAL], f_j);
1261 rvec_inc(fshift[t2], f_k);
1266 real urey_bradley(int nbonds,
1267 const t_iatom forceatoms[], const t_iparams forceparams[],
1268 const rvec x[], rvec f[], rvec fshift[],
1269 const t_pbc *pbc, const t_graph *g,
1270 real lambda, real *dvdlambda,
1271 const t_mdatoms *md, t_fcdata *fcd,
1272 int *global_atom_index)
1274 int i, m, ai, aj, ak, t1, t2, type, ki;
1275 rvec r_ij, r_kj, r_ik;
1276 real cos_theta, cos_theta2, theta;
1277 real dVdt, va, vtot, dr, dr2, vbond, fbond, fik;
1278 real kthA, th0A, kUBA, r13A, kthB, th0B, kUBB, r13B;
1279 ivec jt, dt_ij, dt_kj, dt_ik;
1282 for (i = 0; (i < nbonds); )
1284 type = forceatoms[i++];
1285 ai = forceatoms[i++];
1286 aj = forceatoms[i++];
1287 ak = forceatoms[i++];
1288 th0A = forceparams[type].u_b.thetaA*DEG2RAD;
1289 kthA = forceparams[type].u_b.kthetaA;
1290 r13A = forceparams[type].u_b.r13A;
1291 kUBA = forceparams[type].u_b.kUBA;
1292 th0B = forceparams[type].u_b.thetaB*DEG2RAD;
1293 kthB = forceparams[type].u_b.kthetaB;
1294 r13B = forceparams[type].u_b.r13B;
1295 kUBB = forceparams[type].u_b.kUBB;
1297 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
1298 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
1300 *dvdlambda += harmonic(kthA, kthB, th0A, th0B, theta, lambda, &va, &dVdt); /* 21 */
1303 ki = pbc_rvec_sub(pbc, x[ai], x[ak], r_ik); /* 3 */
1304 dr2 = iprod(r_ik, r_ik); /* 5 */
1305 dr = dr2*gmx_invsqrt(dr2); /* 10 */
1307 *dvdlambda += harmonic(kUBA, kUBB, r13A, r13B, dr, lambda, &vbond, &fbond); /* 19 */
1309 cos_theta2 = sqr(cos_theta); /* 1 */
1317 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1318 sth = st*cos_theta; /* 1 */
1322 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1323 theta*RAD2DEG, va, dVdt);
1326 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1327 nrij2 = iprod(r_ij, r_ij);
1329 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
1330 cii = sth/nrij2; /* 10 */
1331 ckk = sth/nrkj2; /* 10 */
1333 for (m = 0; (m < DIM); m++) /* 39 */
1335 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
1336 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
1337 f_j[m] = -f_i[m]-f_k[m];
1344 copy_ivec(SHIFT_IVEC(g, aj), jt);
1346 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1347 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1348 t1 = IVEC2IS(dt_ij);
1349 t2 = IVEC2IS(dt_kj);
1351 rvec_inc(fshift[t1], f_i);
1352 rvec_inc(fshift[CENTRAL], f_j);
1353 rvec_inc(fshift[t2], f_k);
1355 /* Time for the bond calculations */
1361 vtot += vbond; /* 1*/
1362 fbond *= gmx_invsqrt(dr2); /* 6 */
1366 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ak), dt_ik);
1367 ki = IVEC2IS(dt_ik);
1369 for (m = 0; (m < DIM); m++) /* 15 */
1371 fik = fbond*r_ik[m];
1374 fshift[ki][m] += fik;
1375 fshift[CENTRAL][m] -= fik;
1381 real quartic_angles(int nbonds,
1382 const t_iatom forceatoms[], const t_iparams forceparams[],
1383 const rvec x[], rvec f[], rvec fshift[],
1384 const t_pbc *pbc, const t_graph *g,
1385 real lambda, real *dvdlambda,
1386 const t_mdatoms *md, t_fcdata *fcd,
1387 int *global_atom_index)
1389 int i, j, ai, aj, ak, t1, t2, type;
1391 real cos_theta, cos_theta2, theta, dt, dVdt, va, dtp, c, vtot;
1392 ivec jt, dt_ij, dt_kj;
1395 for (i = 0; (i < nbonds); )
1397 type = forceatoms[i++];
1398 ai = forceatoms[i++];
1399 aj = forceatoms[i++];
1400 ak = forceatoms[i++];
1402 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
1403 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
1405 dt = theta - forceparams[type].qangle.theta*DEG2RAD; /* 2 */
1408 va = forceparams[type].qangle.c[0];
1410 for (j = 1; j <= 4; j++)
1412 c = forceparams[type].qangle.c[j];
1421 cos_theta2 = sqr(cos_theta); /* 1 */
1430 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1431 sth = st*cos_theta; /* 1 */
1435 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1436 theta*RAD2DEG, va, dVdt);
1439 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1440 nrij2 = iprod(r_ij, r_ij);
1442 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
1443 cii = sth/nrij2; /* 10 */
1444 ckk = sth/nrkj2; /* 10 */
1446 for (m = 0; (m < DIM); m++) /* 39 */
1448 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
1449 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
1450 f_j[m] = -f_i[m]-f_k[m];
1457 copy_ivec(SHIFT_IVEC(g, aj), jt);
1459 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1460 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1461 t1 = IVEC2IS(dt_ij);
1462 t2 = IVEC2IS(dt_kj);
1464 rvec_inc(fshift[t1], f_i);
1465 rvec_inc(fshift[CENTRAL], f_j);
1466 rvec_inc(fshift[t2], f_k);
1472 real dih_angle(const rvec xi, const rvec xj, const rvec xk, const rvec xl,
1474 rvec r_ij, rvec r_kj, rvec r_kl, rvec m, rvec n,
1475 real *sign, int *t1, int *t2, int *t3)
1479 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
1480 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
1481 *t3 = pbc_rvec_sub(pbc, xk, xl, r_kl); /* 3 */
1483 cprod(r_ij, r_kj, m); /* 9 */
1484 cprod(r_kj, r_kl, n); /* 9 */
1485 phi = gmx_angle(m, n); /* 49 (assuming 25 for atan2) */
1486 ipr = iprod(r_ij, n); /* 5 */
1487 (*sign) = (ipr < 0.0) ? -1.0 : 1.0;
1488 phi = (*sign)*phi; /* 1 */
1496 /* As dih_angle above, but calculates 4 dihedral angles at once using SIMD,
1497 * also calculates the pre-factor required for the dihedral force update.
1498 * Note that bv and buf should be register aligned.
1500 static gmx_inline void
1501 dih_angle_simd(const rvec *x,
1502 const int *ai, const int *aj, const int *ak, const int *al,
1503 const pbc_simd_t *pbc,
1506 gmx_mm_pr *mx_S, gmx_mm_pr *my_S, gmx_mm_pr *mz_S,
1507 gmx_mm_pr *nx_S, gmx_mm_pr *ny_S, gmx_mm_pr *nz_S,
1508 gmx_mm_pr *nrkj_m2_S,
1509 gmx_mm_pr *nrkj_n2_S,
1513 #define UNROLL GMX_SIMD_WIDTH_HERE
1515 gmx_mm_pr rijx_S, rijy_S, rijz_S;
1516 gmx_mm_pr rkjx_S, rkjy_S, rkjz_S;
1517 gmx_mm_pr rklx_S, rkly_S, rklz_S;
1518 gmx_mm_pr cx_S, cy_S, cz_S;
1522 gmx_mm_pr iprm_S, iprn_S;
1523 gmx_mm_pr nrkj2_S, nrkj_1_S, nrkj_2_S, nrkj_S;
1525 gmx_mm_pr fmin_S = gmx_set1_pr(GMX_FLOAT_MIN);
1527 for (s = 0; s < UNROLL; s++)
1529 /* If you can't use pbc_dx_simd below for PBC, e.g. because
1530 * you can't round in SIMD, use pbc_rvec_sub here.
1532 for (m = 0; m < DIM; m++)
1534 dr[s + (0*DIM + m)*UNROLL] = x[ai[s]][m] - x[aj[s]][m];
1535 dr[s + (1*DIM + m)*UNROLL] = x[ak[s]][m] - x[aj[s]][m];
1536 dr[s + (2*DIM + m)*UNROLL] = x[ak[s]][m] - x[al[s]][m];
1540 rijx_S = gmx_load_pr(dr + 0*UNROLL);
1541 rijy_S = gmx_load_pr(dr + 1*UNROLL);
1542 rijz_S = gmx_load_pr(dr + 2*UNROLL);
1543 rkjx_S = gmx_load_pr(dr + 3*UNROLL);
1544 rkjy_S = gmx_load_pr(dr + 4*UNROLL);
1545 rkjz_S = gmx_load_pr(dr + 5*UNROLL);
1546 rklx_S = gmx_load_pr(dr + 6*UNROLL);
1547 rkly_S = gmx_load_pr(dr + 7*UNROLL);
1548 rklz_S = gmx_load_pr(dr + 8*UNROLL);
1550 pbc_dx_simd(&rijx_S, &rijy_S, &rijz_S, pbc);
1551 pbc_dx_simd(&rkjx_S, &rkjy_S, &rkjz_S, pbc);
1552 pbc_dx_simd(&rklx_S, &rkly_S, &rklz_S, pbc);
1554 gmx_cprod_pr(rijx_S, rijy_S, rijz_S,
1555 rkjx_S, rkjy_S, rkjz_S,
1558 gmx_cprod_pr(rkjx_S, rkjy_S, rkjz_S,
1559 rklx_S, rkly_S, rklz_S,
1562 gmx_cprod_pr(*mx_S, *my_S, *mz_S,
1563 *nx_S, *ny_S, *nz_S,
1564 &cx_S, &cy_S, &cz_S);
1566 cn_S = gmx_sqrt_pr(gmx_norm2_pr(cx_S, cy_S, cz_S));
1568 s_S = gmx_iprod_pr(*mx_S, *my_S, *mz_S, *nx_S, *ny_S, *nz_S);
1570 /* Determine the dihedral angle, the sign might need correction */
1571 *phi_S = gmx_atan2_pr(cn_S, s_S);
1573 ipr_S = gmx_iprod_pr(rijx_S, rijy_S, rijz_S,
1574 *nx_S, *ny_S, *nz_S);
1576 iprm_S = gmx_norm2_pr(*mx_S, *my_S, *mz_S);
1577 iprn_S = gmx_norm2_pr(*nx_S, *ny_S, *nz_S);
1579 nrkj2_S = gmx_norm2_pr(rkjx_S, rkjy_S, rkjz_S);
1581 /* Avoid division by zero. When zero, the result is multiplied by 0
1582 * anyhow, so the 3 max below do not affect the final result.
1584 nrkj2_S = gmx_max_pr(nrkj2_S, fmin_S);
1585 nrkj_1_S = gmx_invsqrt_pr(nrkj2_S);
1586 nrkj_2_S = gmx_mul_pr(nrkj_1_S, nrkj_1_S);
1587 nrkj_S = gmx_mul_pr(nrkj2_S, nrkj_1_S);
1589 iprm_S = gmx_max_pr(iprm_S, fmin_S);
1590 iprn_S = gmx_max_pr(iprn_S, fmin_S);
1591 *nrkj_m2_S = gmx_mul_pr(nrkj_S, gmx_inv_pr(iprm_S));
1592 *nrkj_n2_S = gmx_mul_pr(nrkj_S, gmx_inv_pr(iprn_S));
1594 /* Set sign of phi_S with the sign of ipr_S; phi_S is currently positive */
1595 *phi_S = gmx_cpsgn_nonneg_pr(ipr_S, *phi_S);
1597 p_S = gmx_iprod_pr(rijx_S, rijy_S, rijz_S,
1598 rkjx_S, rkjy_S, rkjz_S);
1599 p_S = gmx_mul_pr(p_S, nrkj_2_S);
1601 q_S = gmx_iprod_pr(rklx_S, rkly_S, rklz_S,
1602 rkjx_S, rkjy_S, rkjz_S);
1603 q_S = gmx_mul_pr(q_S, nrkj_2_S);
1605 gmx_store_pr(p, p_S);
1606 gmx_store_pr(q, q_S);
1610 #endif /* SIMD_BONDEDS */
1613 void do_dih_fup(int i, int j, int k, int l, real ddphi,
1614 rvec r_ij, rvec r_kj, rvec r_kl,
1615 rvec m, rvec n, rvec f[], rvec fshift[],
1616 const t_pbc *pbc, const t_graph *g,
1617 const rvec x[], int t1, int t2, int t3)
1620 rvec f_i, f_j, f_k, f_l;
1621 rvec uvec, vvec, svec, dx_jl;
1622 real iprm, iprn, nrkj, nrkj2, nrkj_1, nrkj_2;
1623 real a, b, p, q, toler;
1624 ivec jt, dt_ij, dt_kj, dt_lj;
1626 iprm = iprod(m, m); /* 5 */
1627 iprn = iprod(n, n); /* 5 */
1628 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1629 toler = nrkj2*GMX_REAL_EPS;
1630 if ((iprm > toler) && (iprn > toler))
1632 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1633 nrkj_2 = nrkj_1*nrkj_1; /* 1 */
1634 nrkj = nrkj2*nrkj_1; /* 1 */
1635 a = -ddphi*nrkj/iprm; /* 11 */
1636 svmul(a, m, f_i); /* 3 */
1637 b = ddphi*nrkj/iprn; /* 11 */
1638 svmul(b, n, f_l); /* 3 */
1639 p = iprod(r_ij, r_kj); /* 5 */
1640 p *= nrkj_2; /* 1 */
1641 q = iprod(r_kl, r_kj); /* 5 */
1642 q *= nrkj_2; /* 1 */
1643 svmul(p, f_i, uvec); /* 3 */
1644 svmul(q, f_l, vvec); /* 3 */
1645 rvec_sub(uvec, vvec, svec); /* 3 */
1646 rvec_sub(f_i, svec, f_j); /* 3 */
1647 rvec_add(f_l, svec, f_k); /* 3 */
1648 rvec_inc(f[i], f_i); /* 3 */
1649 rvec_dec(f[j], f_j); /* 3 */
1650 rvec_dec(f[k], f_k); /* 3 */
1651 rvec_inc(f[l], f_l); /* 3 */
1655 copy_ivec(SHIFT_IVEC(g, j), jt);
1656 ivec_sub(SHIFT_IVEC(g, i), jt, dt_ij);
1657 ivec_sub(SHIFT_IVEC(g, k), jt, dt_kj);
1658 ivec_sub(SHIFT_IVEC(g, l), jt, dt_lj);
1659 t1 = IVEC2IS(dt_ij);
1660 t2 = IVEC2IS(dt_kj);
1661 t3 = IVEC2IS(dt_lj);
1665 t3 = pbc_rvec_sub(pbc, x[l], x[j], dx_jl);
1672 rvec_inc(fshift[t1], f_i);
1673 rvec_dec(fshift[CENTRAL], f_j);
1674 rvec_dec(fshift[t2], f_k);
1675 rvec_inc(fshift[t3], f_l);
1680 /* As do_dih_fup above, but without shift forces */
1682 do_dih_fup_noshiftf(int i, int j, int k, int l, real ddphi,
1683 rvec r_ij, rvec r_kj, rvec r_kl,
1684 rvec m, rvec n, rvec f[])
1686 rvec f_i, f_j, f_k, f_l;
1687 rvec uvec, vvec, svec, dx_jl;
1688 real iprm, iprn, nrkj, nrkj2, nrkj_1, nrkj_2;
1689 real a, b, p, q, toler;
1690 ivec jt, dt_ij, dt_kj, dt_lj;
1692 iprm = iprod(m, m); /* 5 */
1693 iprn = iprod(n, n); /* 5 */
1694 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1695 toler = nrkj2*GMX_REAL_EPS;
1696 if ((iprm > toler) && (iprn > toler))
1698 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1699 nrkj_2 = nrkj_1*nrkj_1; /* 1 */
1700 nrkj = nrkj2*nrkj_1; /* 1 */
1701 a = -ddphi*nrkj/iprm; /* 11 */
1702 svmul(a, m, f_i); /* 3 */
1703 b = ddphi*nrkj/iprn; /* 11 */
1704 svmul(b, n, f_l); /* 3 */
1705 p = iprod(r_ij, r_kj); /* 5 */
1706 p *= nrkj_2; /* 1 */
1707 q = iprod(r_kl, r_kj); /* 5 */
1708 q *= nrkj_2; /* 1 */
1709 svmul(p, f_i, uvec); /* 3 */
1710 svmul(q, f_l, vvec); /* 3 */
1711 rvec_sub(uvec, vvec, svec); /* 3 */
1712 rvec_sub(f_i, svec, f_j); /* 3 */
1713 rvec_add(f_l, svec, f_k); /* 3 */
1714 rvec_inc(f[i], f_i); /* 3 */
1715 rvec_dec(f[j], f_j); /* 3 */
1716 rvec_dec(f[k], f_k); /* 3 */
1717 rvec_inc(f[l], f_l); /* 3 */
1721 /* As do_dih_fup_noshiftf above, but with pre-calculated pre-factors */
1722 static gmx_inline void
1723 do_dih_fup_noshiftf_precalc(int i, int j, int k, int l,
1725 real f_i_x, real f_i_y, real f_i_z,
1726 real mf_l_x, real mf_l_y, real mf_l_z,
1729 rvec f_i, f_j, f_k, f_l;
1730 rvec uvec, vvec, svec;
1738 svmul(p, f_i, uvec);
1739 svmul(q, f_l, vvec);
1740 rvec_sub(uvec, vvec, svec);
1741 rvec_sub(f_i, svec, f_j);
1742 rvec_add(f_l, svec, f_k);
1743 rvec_inc(f[i], f_i);
1744 rvec_dec(f[j], f_j);
1745 rvec_dec(f[k], f_k);
1746 rvec_inc(f[l], f_l);
1750 real dopdihs(real cpA, real cpB, real phiA, real phiB, int mult,
1751 real phi, real lambda, real *V, real *F)
1753 real v, dvdlambda, mdphi, v1, sdphi, ddphi;
1754 real L1 = 1.0 - lambda;
1755 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1756 real dph0 = (phiB - phiA)*DEG2RAD;
1757 real cp = L1*cpA + lambda*cpB;
1759 mdphi = mult*phi - ph0;
1761 ddphi = -cp*mult*sdphi;
1762 v1 = 1.0 + cos(mdphi);
1765 dvdlambda = (cpB - cpA)*v1 + cp*dph0*sdphi;
1772 /* That was 40 flops */
1776 dopdihs_noener(real cpA, real cpB, real phiA, real phiB, int mult,
1777 real phi, real lambda, real *F)
1779 real mdphi, sdphi, ddphi;
1780 real L1 = 1.0 - lambda;
1781 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1782 real cp = L1*cpA + lambda*cpB;
1784 mdphi = mult*phi - ph0;
1786 ddphi = -cp*mult*sdphi;
1790 /* That was 20 flops */
1794 dopdihs_mdphi(real cpA, real cpB, real phiA, real phiB, int mult,
1795 real phi, real lambda, real *cp, real *mdphi)
1797 real L1 = 1.0 - lambda;
1798 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1800 *cp = L1*cpA + lambda*cpB;
1802 *mdphi = mult*phi - ph0;
1805 static real dopdihs_min(real cpA, real cpB, real phiA, real phiB, int mult,
1806 real phi, real lambda, real *V, real *F)
1807 /* similar to dopdihs, except for a minus sign *
1808 * and a different treatment of mult/phi0 */
1810 real v, dvdlambda, mdphi, v1, sdphi, ddphi;
1811 real L1 = 1.0 - lambda;
1812 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1813 real dph0 = (phiB - phiA)*DEG2RAD;
1814 real cp = L1*cpA + lambda*cpB;
1816 mdphi = mult*(phi-ph0);
1818 ddphi = cp*mult*sdphi;
1819 v1 = 1.0-cos(mdphi);
1822 dvdlambda = (cpB-cpA)*v1 + cp*dph0*sdphi;
1829 /* That was 40 flops */
1832 real pdihs(int nbonds,
1833 const t_iatom forceatoms[], const t_iparams forceparams[],
1834 const rvec x[], rvec f[], rvec fshift[],
1835 const t_pbc *pbc, const t_graph *g,
1836 real lambda, real *dvdlambda,
1837 const t_mdatoms *md, t_fcdata *fcd,
1838 int *global_atom_index)
1840 int i, type, ai, aj, ak, al;
1842 rvec r_ij, r_kj, r_kl, m, n;
1843 real phi, sign, ddphi, vpd, vtot;
1847 for (i = 0; (i < nbonds); )
1849 type = forceatoms[i++];
1850 ai = forceatoms[i++];
1851 aj = forceatoms[i++];
1852 ak = forceatoms[i++];
1853 al = forceatoms[i++];
1855 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
1856 &sign, &t1, &t2, &t3); /* 84 */
1857 *dvdlambda += dopdihs(forceparams[type].pdihs.cpA,
1858 forceparams[type].pdihs.cpB,
1859 forceparams[type].pdihs.phiA,
1860 forceparams[type].pdihs.phiB,
1861 forceparams[type].pdihs.mult,
1862 phi, lambda, &vpd, &ddphi);
1865 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
1866 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
1869 fprintf(debug, "pdih: (%d,%d,%d,%d) phi=%g\n",
1870 ai, aj, ak, al, phi);
1877 void make_dp_periodic(real *dp) /* 1 flop? */
1879 /* dp cannot be outside (-pi,pi) */
1884 else if (*dp < -M_PI)
1891 /* As pdihs above, but without calculating energies and shift forces */
1893 pdihs_noener(int nbonds,
1894 const t_iatom forceatoms[], const t_iparams forceparams[],
1895 const rvec x[], rvec f[],
1896 const t_pbc *pbc, const t_graph *g,
1898 const t_mdatoms *md, t_fcdata *fcd,
1899 int *global_atom_index)
1901 int i, type, ai, aj, ak, al;
1903 rvec r_ij, r_kj, r_kl, m, n;
1904 real phi, sign, ddphi_tot, ddphi;
1906 for (i = 0; (i < nbonds); )
1908 ai = forceatoms[i+1];
1909 aj = forceatoms[i+2];
1910 ak = forceatoms[i+3];
1911 al = forceatoms[i+4];
1913 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
1914 &sign, &t1, &t2, &t3);
1918 /* Loop over dihedrals working on the same atoms,
1919 * so we avoid recalculating angles and force distributions.
1923 type = forceatoms[i];
1924 dopdihs_noener(forceparams[type].pdihs.cpA,
1925 forceparams[type].pdihs.cpB,
1926 forceparams[type].pdihs.phiA,
1927 forceparams[type].pdihs.phiB,
1928 forceparams[type].pdihs.mult,
1929 phi, lambda, &ddphi);
1934 while (i < nbonds &&
1935 forceatoms[i+1] == ai &&
1936 forceatoms[i+2] == aj &&
1937 forceatoms[i+3] == ak &&
1938 forceatoms[i+4] == al);
1940 do_dih_fup_noshiftf(ai, aj, ak, al, ddphi_tot, r_ij, r_kj, r_kl, m, n, f);
1947 /* As pdihs_noner above, but using SIMD to calculate many dihedrals at once */
1949 pdihs_noener_simd(int nbonds,
1950 const t_iatom forceatoms[], const t_iparams forceparams[],
1951 const rvec x[], rvec f[],
1952 const t_pbc *pbc, const t_graph *g,
1954 const t_mdatoms *md, t_fcdata *fcd,
1955 int *global_atom_index)
1957 #define UNROLL GMX_SIMD_WIDTH_HERE
1960 int type, ai[UNROLL], aj[UNROLL], ak[UNROLL], al[UNROLL];
1961 int t1[UNROLL], t2[UNROLL], t3[UNROLL];
1963 real dr_array[3*DIM*UNROLL+UNROLL], *dr;
1964 real buf_array[7*UNROLL+UNROLL], *buf;
1965 real *cp, *phi0, *mult, *phi, *p, *q, *sf_i, *msf_l;
1966 gmx_mm_pr phi0_S, phi_S;
1967 gmx_mm_pr mx_S, my_S, mz_S;
1968 gmx_mm_pr nx_S, ny_S, nz_S;
1969 gmx_mm_pr nrkj_m2_S, nrkj_n2_S;
1970 gmx_mm_pr cp_S, mdphi_S, mult_S;
1971 gmx_mm_pr sin_S, cos_S;
1973 gmx_mm_pr sf_i_S, msf_l_S;
1974 pbc_simd_t pbc_simd;
1976 /* Ensure SIMD register alignment */
1977 dr = gmx_simd_align_real(dr_array);
1978 buf = gmx_simd_align_real(buf_array);
1980 /* Extract aligned pointer for parameters and variables */
1981 cp = buf + 0*UNROLL;
1982 phi0 = buf + 1*UNROLL;
1983 mult = buf + 2*UNROLL;
1986 sf_i = buf + 5*UNROLL;
1987 msf_l = buf + 6*UNROLL;
1989 set_pbc_simd(pbc, &pbc_simd);
1991 /* nbonds is the number of dihedrals times nfa1, here we step UNROLL dihs */
1992 for (i = 0; (i < nbonds); i += UNROLL*nfa1)
1994 /* Collect atoms quadruplets for UNROLL dihedrals.
1995 * iu indexes into forceatoms, we should not let iu go beyond nbonds.
1998 for (s = 0; s < UNROLL; s++)
2000 type = forceatoms[iu];
2001 ai[s] = forceatoms[iu+1];
2002 aj[s] = forceatoms[iu+2];
2003 ak[s] = forceatoms[iu+3];
2004 al[s] = forceatoms[iu+4];
2006 cp[s] = forceparams[type].pdihs.cpA;
2007 phi0[s] = forceparams[type].pdihs.phiA*DEG2RAD;
2008 mult[s] = forceparams[type].pdihs.mult;
2010 /* At the end fill the arrays with identical entries */
2011 if (iu + nfa1 < nbonds)
2017 /* Caclulate UNROLL dihedral angles at once */
2018 dih_angle_simd(x, ai, aj, ak, al, &pbc_simd,
2021 &mx_S, &my_S, &mz_S,
2022 &nx_S, &ny_S, &nz_S,
2027 cp_S = gmx_load_pr(cp);
2028 phi0_S = gmx_load_pr(phi0);
2029 mult_S = gmx_load_pr(mult);
2031 mdphi_S = gmx_sub_pr(gmx_mul_pr(mult_S, phi_S), phi0_S);
2033 /* Calculate UNROLL sines at once */
2034 gmx_sincos_pr(mdphi_S, &sin_S, &cos_S);
2035 mddphi_S = gmx_mul_pr(gmx_mul_pr(cp_S, mult_S), sin_S);
2036 sf_i_S = gmx_mul_pr(mddphi_S, nrkj_m2_S);
2037 msf_l_S = gmx_mul_pr(mddphi_S, nrkj_n2_S);
2039 /* After this m?_S will contain f[i] */
2040 mx_S = gmx_mul_pr(sf_i_S, mx_S);
2041 my_S = gmx_mul_pr(sf_i_S, my_S);
2042 mz_S = gmx_mul_pr(sf_i_S, mz_S);
2044 /* After this m?_S will contain -f[l] */
2045 nx_S = gmx_mul_pr(msf_l_S, nx_S);
2046 ny_S = gmx_mul_pr(msf_l_S, ny_S);
2047 nz_S = gmx_mul_pr(msf_l_S, nz_S);
2049 gmx_store_pr(dr + 0*UNROLL, mx_S);
2050 gmx_store_pr(dr + 1*UNROLL, my_S);
2051 gmx_store_pr(dr + 2*UNROLL, mz_S);
2052 gmx_store_pr(dr + 3*UNROLL, nx_S);
2053 gmx_store_pr(dr + 4*UNROLL, ny_S);
2054 gmx_store_pr(dr + 5*UNROLL, nz_S);
2060 do_dih_fup_noshiftf_precalc(ai[s], aj[s], ak[s], al[s],
2065 dr[(DIM+XX)*UNROLL+s],
2066 dr[(DIM+YY)*UNROLL+s],
2067 dr[(DIM+ZZ)*UNROLL+s],
2072 while (s < UNROLL && iu < nbonds);
2077 #endif /* SIMD_BONDEDS */
2080 real idihs(int nbonds,
2081 const t_iatom forceatoms[], const t_iparams forceparams[],
2082 const rvec x[], rvec f[], rvec fshift[],
2083 const t_pbc *pbc, const t_graph *g,
2084 real lambda, real *dvdlambda,
2085 const t_mdatoms *md, t_fcdata *fcd,
2086 int *global_atom_index)
2088 int i, type, ai, aj, ak, al;
2090 real phi, phi0, dphi0, ddphi, sign, vtot;
2091 rvec r_ij, r_kj, r_kl, m, n;
2092 real L1, kk, dp, dp2, kA, kB, pA, pB, dvdl_term;
2097 for (i = 0; (i < nbonds); )
2099 type = forceatoms[i++];
2100 ai = forceatoms[i++];
2101 aj = forceatoms[i++];
2102 ak = forceatoms[i++];
2103 al = forceatoms[i++];
2105 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2106 &sign, &t1, &t2, &t3); /* 84 */
2108 /* phi can jump if phi0 is close to Pi/-Pi, which will cause huge
2109 * force changes if we just apply a normal harmonic.
2110 * Instead, we first calculate phi-phi0 and take it modulo (-Pi,Pi).
2111 * This means we will never have the periodicity problem, unless
2112 * the dihedral is Pi away from phiO, which is very unlikely due to
2115 kA = forceparams[type].harmonic.krA;
2116 kB = forceparams[type].harmonic.krB;
2117 pA = forceparams[type].harmonic.rA;
2118 pB = forceparams[type].harmonic.rB;
2120 kk = L1*kA + lambda*kB;
2121 phi0 = (L1*pA + lambda*pB)*DEG2RAD;
2122 dphi0 = (pB - pA)*DEG2RAD;
2126 make_dp_periodic(&dp);
2133 dvdl_term += 0.5*(kB - kA)*dp2 - kk*dphi0*dp;
2135 do_dih_fup(ai, aj, ak, al, (real)(-ddphi), r_ij, r_kj, r_kl, m, n,
2136 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2141 fprintf(debug, "idih: (%d,%d,%d,%d) phi=%g\n",
2142 ai, aj, ak, al, phi);
2147 *dvdlambda += dvdl_term;
2152 real posres(int nbonds,
2153 const t_iatom forceatoms[], const t_iparams forceparams[],
2154 const rvec x[], rvec f[], rvec vir_diag,
2156 real lambda, real *dvdlambda,
2157 int refcoord_scaling, int ePBC, rvec comA, rvec comB)
2159 int i, ai, m, d, type, ki, npbcdim = 0;
2160 const t_iparams *pr;
2163 real posA, posB, ref = 0;
2164 rvec comA_sc, comB_sc, rdist, dpdl, pos, dx;
2165 gmx_bool bForceValid = TRUE;
2167 if ((f == NULL) || (vir_diag == NULL)) /* should both be null together! */
2169 bForceValid = FALSE;
2172 npbcdim = ePBC2npbcdim(ePBC);
2174 if (refcoord_scaling == erscCOM)
2176 clear_rvec(comA_sc);
2177 clear_rvec(comB_sc);
2178 for (m = 0; m < npbcdim; m++)
2180 for (d = m; d < npbcdim; d++)
2182 comA_sc[m] += comA[d]*pbc->box[d][m];
2183 comB_sc[m] += comB[d]*pbc->box[d][m];
2191 for (i = 0; (i < nbonds); )
2193 type = forceatoms[i++];
2194 ai = forceatoms[i++];
2195 pr = &forceparams[type];
2197 for (m = 0; m < DIM; m++)
2199 posA = forceparams[type].posres.pos0A[m];
2200 posB = forceparams[type].posres.pos0B[m];
2203 switch (refcoord_scaling)
2207 rdist[m] = L1*posA + lambda*posB;
2208 dpdl[m] = posB - posA;
2211 /* Box relative coordinates are stored for dimensions with pbc */
2212 posA *= pbc->box[m][m];
2213 posB *= pbc->box[m][m];
2214 for (d = m+1; d < npbcdim; d++)
2216 posA += forceparams[type].posres.pos0A[d]*pbc->box[d][m];
2217 posB += forceparams[type].posres.pos0B[d]*pbc->box[d][m];
2219 ref = L1*posA + lambda*posB;
2221 dpdl[m] = posB - posA;
2224 ref = L1*comA_sc[m] + lambda*comB_sc[m];
2225 rdist[m] = L1*posA + lambda*posB;
2226 dpdl[m] = comB_sc[m] - comA_sc[m] + posB - posA;
2229 gmx_fatal(FARGS, "No such scaling method implemented");
2234 ref = L1*posA + lambda*posB;
2236 dpdl[m] = posB - posA;
2239 /* We do pbc_dx with ref+rdist,
2240 * since with only ref we can be up to half a box vector wrong.
2242 pos[m] = ref + rdist[m];
2247 pbc_dx(pbc, x[ai], pos, dx);
2251 rvec_sub(x[ai], pos, dx);
2254 for (m = 0; (m < DIM); m++)
2256 kk = L1*pr->posres.fcA[m] + lambda*pr->posres.fcB[m];
2258 vtot += 0.5*kk*dx[m]*dx[m];
2260 0.5*(pr->posres.fcB[m] - pr->posres.fcA[m])*dx[m]*dx[m]
2263 /* Here we correct for the pbc_dx which included rdist */
2267 vir_diag[m] -= 0.5*(dx[m] + rdist[m])*fm;
2275 static real low_angres(int nbonds,
2276 const t_iatom forceatoms[], const t_iparams forceparams[],
2277 const rvec x[], rvec f[], rvec fshift[],
2278 const t_pbc *pbc, const t_graph *g,
2279 real lambda, real *dvdlambda,
2282 int i, m, type, ai, aj, ak, al;
2284 real phi, cos_phi, cos_phi2, vid, vtot, dVdphi;
2285 rvec r_ij, r_kl, f_i, f_k = {0, 0, 0};
2286 real st, sth, nrij2, nrkl2, c, cij, ckl;
2289 t2 = 0; /* avoid warning with gcc-3.3. It is never used uninitialized */
2292 ak = al = 0; /* to avoid warnings */
2293 for (i = 0; i < nbonds; )
2295 type = forceatoms[i++];
2296 ai = forceatoms[i++];
2297 aj = forceatoms[i++];
2298 t1 = pbc_rvec_sub(pbc, x[aj], x[ai], r_ij); /* 3 */
2301 ak = forceatoms[i++];
2302 al = forceatoms[i++];
2303 t2 = pbc_rvec_sub(pbc, x[al], x[ak], r_kl); /* 3 */
2312 cos_phi = cos_angle(r_ij, r_kl); /* 25 */
2313 phi = acos(cos_phi); /* 10 */
2315 *dvdlambda += dopdihs_min(forceparams[type].pdihs.cpA,
2316 forceparams[type].pdihs.cpB,
2317 forceparams[type].pdihs.phiA,
2318 forceparams[type].pdihs.phiB,
2319 forceparams[type].pdihs.mult,
2320 phi, lambda, &vid, &dVdphi); /* 40 */
2324 cos_phi2 = sqr(cos_phi); /* 1 */
2327 st = -dVdphi*gmx_invsqrt(1 - cos_phi2); /* 12 */
2328 sth = st*cos_phi; /* 1 */
2329 nrij2 = iprod(r_ij, r_ij); /* 5 */
2330 nrkl2 = iprod(r_kl, r_kl); /* 5 */
2332 c = st*gmx_invsqrt(nrij2*nrkl2); /* 11 */
2333 cij = sth/nrij2; /* 10 */
2334 ckl = sth/nrkl2; /* 10 */
2336 for (m = 0; m < DIM; m++) /* 18+18 */
2338 f_i[m] = (c*r_kl[m]-cij*r_ij[m]);
2343 f_k[m] = (c*r_ij[m]-ckl*r_kl[m]);
2351 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
2354 rvec_inc(fshift[t1], f_i);
2355 rvec_dec(fshift[CENTRAL], f_i);
2360 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, al), dt);
2363 rvec_inc(fshift[t2], f_k);
2364 rvec_dec(fshift[CENTRAL], f_k);
2369 return vtot; /* 184 / 157 (bZAxis) total */
2372 real angres(int nbonds,
2373 const t_iatom forceatoms[], const t_iparams forceparams[],
2374 const rvec x[], rvec f[], rvec fshift[],
2375 const t_pbc *pbc, const t_graph *g,
2376 real lambda, real *dvdlambda,
2377 const t_mdatoms *md, t_fcdata *fcd,
2378 int *global_atom_index)
2380 return low_angres(nbonds, forceatoms, forceparams, x, f, fshift, pbc, g,
2381 lambda, dvdlambda, FALSE);
2384 real angresz(int nbonds,
2385 const t_iatom forceatoms[], const t_iparams forceparams[],
2386 const rvec x[], rvec f[], rvec fshift[],
2387 const t_pbc *pbc, const t_graph *g,
2388 real lambda, real *dvdlambda,
2389 const t_mdatoms *md, t_fcdata *fcd,
2390 int *global_atom_index)
2392 return low_angres(nbonds, forceatoms, forceparams, x, f, fshift, pbc, g,
2393 lambda, dvdlambda, TRUE);
2396 real dihres(int nbonds,
2397 const t_iatom forceatoms[], const t_iparams forceparams[],
2398 const rvec x[], rvec f[], rvec fshift[],
2399 const t_pbc *pbc, const t_graph *g,
2400 real lambda, real *dvdlambda,
2401 const t_mdatoms *md, t_fcdata *fcd,
2402 int *global_atom_index)
2405 int ai, aj, ak, al, i, k, type, t1, t2, t3;
2406 real phi0A, phi0B, dphiA, dphiB, kfacA, kfacB, phi0, dphi, kfac;
2407 real phi, ddphi, ddp, ddp2, dp, sign, d2r, fc, L1;
2408 rvec r_ij, r_kj, r_kl, m, n;
2415 for (i = 0; (i < nbonds); )
2417 type = forceatoms[i++];
2418 ai = forceatoms[i++];
2419 aj = forceatoms[i++];
2420 ak = forceatoms[i++];
2421 al = forceatoms[i++];
2423 phi0A = forceparams[type].dihres.phiA*d2r;
2424 dphiA = forceparams[type].dihres.dphiA*d2r;
2425 kfacA = forceparams[type].dihres.kfacA;
2427 phi0B = forceparams[type].dihres.phiB*d2r;
2428 dphiB = forceparams[type].dihres.dphiB*d2r;
2429 kfacB = forceparams[type].dihres.kfacB;
2431 phi0 = L1*phi0A + lambda*phi0B;
2432 dphi = L1*dphiA + lambda*dphiB;
2433 kfac = L1*kfacA + lambda*kfacB;
2435 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2436 &sign, &t1, &t2, &t3);
2441 fprintf(debug, "dihres[%d]: %d %d %d %d : phi=%f, dphi=%f, kfac=%f\n",
2442 k++, ai, aj, ak, al, phi0, dphi, kfac);
2444 /* phi can jump if phi0 is close to Pi/-Pi, which will cause huge
2445 * force changes if we just apply a normal harmonic.
2446 * Instead, we first calculate phi-phi0 and take it modulo (-Pi,Pi).
2447 * This means we will never have the periodicity problem, unless
2448 * the dihedral is Pi away from phiO, which is very unlikely due to
2452 make_dp_periodic(&dp);
2458 else if (dp < -dphi)
2470 vtot += 0.5*kfac*ddp2;
2473 *dvdlambda += 0.5*(kfacB - kfacA)*ddp2;
2474 /* lambda dependence from changing restraint distances */
2477 *dvdlambda -= kfac*ddp*((dphiB - dphiA)+(phi0B - phi0A));
2481 *dvdlambda += kfac*ddp*((dphiB - dphiA)-(phi0B - phi0A));
2483 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
2484 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2491 real unimplemented(int nbonds,
2492 const t_iatom forceatoms[], const t_iparams forceparams[],
2493 const rvec x[], rvec f[], rvec fshift[],
2494 const t_pbc *pbc, const t_graph *g,
2495 real lambda, real *dvdlambda,
2496 const t_mdatoms *md, t_fcdata *fcd,
2497 int *global_atom_index)
2499 gmx_impl("*** you are using a not implemented function");
2501 return 0.0; /* To make the compiler happy */
2504 real rbdihs(int nbonds,
2505 const t_iatom forceatoms[], const t_iparams forceparams[],
2506 const rvec x[], rvec f[], rvec fshift[],
2507 const t_pbc *pbc, const t_graph *g,
2508 real lambda, real *dvdlambda,
2509 const t_mdatoms *md, t_fcdata *fcd,
2510 int *global_atom_index)
2512 const real c0 = 0.0, c1 = 1.0, c2 = 2.0, c3 = 3.0, c4 = 4.0, c5 = 5.0;
2513 int type, ai, aj, ak, al, i, j;
2515 rvec r_ij, r_kj, r_kl, m, n;
2516 real parmA[NR_RBDIHS];
2517 real parmB[NR_RBDIHS];
2518 real parm[NR_RBDIHS];
2519 real cos_phi, phi, rbp, rbpBA;
2520 real v, sign, ddphi, sin_phi;
2522 real L1 = 1.0-lambda;
2526 for (i = 0; (i < nbonds); )
2528 type = forceatoms[i++];
2529 ai = forceatoms[i++];
2530 aj = forceatoms[i++];
2531 ak = forceatoms[i++];
2532 al = forceatoms[i++];
2534 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2535 &sign, &t1, &t2, &t3); /* 84 */
2537 /* Change to polymer convention */
2544 phi -= M_PI; /* 1 */
2548 /* Beware of accuracy loss, cannot use 1-sqrt(cos^2) ! */
2551 for (j = 0; (j < NR_RBDIHS); j++)
2553 parmA[j] = forceparams[type].rbdihs.rbcA[j];
2554 parmB[j] = forceparams[type].rbdihs.rbcB[j];
2555 parm[j] = L1*parmA[j]+lambda*parmB[j];
2557 /* Calculate cosine powers */
2558 /* Calculate the energy */
2559 /* Calculate the derivative */
2562 dvdl_term += (parmB[0]-parmA[0]);
2567 rbpBA = parmB[1]-parmA[1];
2568 ddphi += rbp*cosfac;
2571 dvdl_term += cosfac*rbpBA;
2573 rbpBA = parmB[2]-parmA[2];
2574 ddphi += c2*rbp*cosfac;
2577 dvdl_term += cosfac*rbpBA;
2579 rbpBA = parmB[3]-parmA[3];
2580 ddphi += c3*rbp*cosfac;
2583 dvdl_term += cosfac*rbpBA;
2585 rbpBA = parmB[4]-parmA[4];
2586 ddphi += c4*rbp*cosfac;
2589 dvdl_term += cosfac*rbpBA;
2591 rbpBA = parmB[5]-parmA[5];
2592 ddphi += c5*rbp*cosfac;
2595 dvdl_term += cosfac*rbpBA;
2597 ddphi = -ddphi*sin_phi; /* 11 */
2599 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
2600 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2603 *dvdlambda += dvdl_term;
2608 int cmap_setup_grid_index(int ip, int grid_spacing, int *ipm1, int *ipp1, int *ipp2)
2614 ip = ip + grid_spacing - 1;
2616 else if (ip > grid_spacing)
2618 ip = ip - grid_spacing - 1;
2627 im1 = grid_spacing - 1;
2629 else if (ip == grid_spacing-2)
2633 else if (ip == grid_spacing-1)
2647 real cmap_dihs(int nbonds,
2648 const t_iatom forceatoms[], const t_iparams forceparams[],
2649 const gmx_cmap_t *cmap_grid,
2650 const rvec x[], rvec f[], rvec fshift[],
2651 const t_pbc *pbc, const t_graph *g,
2652 real lambda, real *dvdlambda,
2653 const t_mdatoms *md, t_fcdata *fcd,
2654 int *global_atom_index)
2656 int i, j, k, n, idx;
2657 int ai, aj, ak, al, am;
2658 int a1i, a1j, a1k, a1l, a2i, a2j, a2k, a2l;
2660 int t11, t21, t31, t12, t22, t32;
2661 int iphi1, ip1m1, ip1p1, ip1p2;
2662 int iphi2, ip2m1, ip2p1, ip2p2;
2664 int pos1, pos2, pos3, pos4, tmp;
2666 real ty[4], ty1[4], ty2[4], ty12[4], tc[16], tx[16];
2667 real phi1, psi1, cos_phi1, sin_phi1, sign1, xphi1;
2668 real phi2, psi2, cos_phi2, sin_phi2, sign2, xphi2;
2669 real dx, xx, tt, tu, e, df1, df2, ddf1, ddf2, ddf12, vtot;
2670 real ra21, rb21, rg21, rg1, rgr1, ra2r1, rb2r1, rabr1;
2671 real ra22, rb22, rg22, rg2, rgr2, ra2r2, rb2r2, rabr2;
2672 real fg1, hg1, fga1, hgb1, gaa1, gbb1;
2673 real fg2, hg2, fga2, hgb2, gaa2, gbb2;
2676 rvec r1_ij, r1_kj, r1_kl, m1, n1;
2677 rvec r2_ij, r2_kj, r2_kl, m2, n2;
2678 rvec f1_i, f1_j, f1_k, f1_l;
2679 rvec f2_i, f2_j, f2_k, f2_l;
2680 rvec a1, b1, a2, b2;
2681 rvec f1, g1, h1, f2, g2, h2;
2682 rvec dtf1, dtg1, dth1, dtf2, dtg2, dth2;
2683 ivec jt1, dt1_ij, dt1_kj, dt1_lj;
2684 ivec jt2, dt2_ij, dt2_kj, dt2_lj;
2688 int loop_index[4][4] = {
2695 /* Total CMAP energy */
2698 for (n = 0; n < nbonds; )
2700 /* Five atoms are involved in the two torsions */
2701 type = forceatoms[n++];
2702 ai = forceatoms[n++];
2703 aj = forceatoms[n++];
2704 ak = forceatoms[n++];
2705 al = forceatoms[n++];
2706 am = forceatoms[n++];
2708 /* Which CMAP type is this */
2709 cmapA = forceparams[type].cmap.cmapA;
2710 cmapd = cmap_grid->cmapdata[cmapA].cmap;
2718 phi1 = dih_angle(x[a1i], x[a1j], x[a1k], x[a1l], pbc, r1_ij, r1_kj, r1_kl, m1, n1,
2719 &sign1, &t11, &t21, &t31); /* 84 */
2721 cos_phi1 = cos(phi1);
2723 a1[0] = r1_ij[1]*r1_kj[2]-r1_ij[2]*r1_kj[1];
2724 a1[1] = r1_ij[2]*r1_kj[0]-r1_ij[0]*r1_kj[2];
2725 a1[2] = r1_ij[0]*r1_kj[1]-r1_ij[1]*r1_kj[0]; /* 9 */
2727 b1[0] = r1_kl[1]*r1_kj[2]-r1_kl[2]*r1_kj[1];
2728 b1[1] = r1_kl[2]*r1_kj[0]-r1_kl[0]*r1_kj[2];
2729 b1[2] = r1_kl[0]*r1_kj[1]-r1_kl[1]*r1_kj[0]; /* 9 */
2731 tmp = pbc_rvec_sub(pbc, x[a1l], x[a1k], h1);
2733 ra21 = iprod(a1, a1); /* 5 */
2734 rb21 = iprod(b1, b1); /* 5 */
2735 rg21 = iprod(r1_kj, r1_kj); /* 5 */
2741 rabr1 = sqrt(ra2r1*rb2r1);
2743 sin_phi1 = rg1 * rabr1 * iprod(a1, h1) * (-1);
2745 if (cos_phi1 < -0.5 || cos_phi1 > 0.5)
2747 phi1 = asin(sin_phi1);
2757 phi1 = -M_PI - phi1;
2763 phi1 = acos(cos_phi1);
2771 xphi1 = phi1 + M_PI; /* 1 */
2773 /* Second torsion */
2779 phi2 = dih_angle(x[a2i], x[a2j], x[a2k], x[a2l], pbc, r2_ij, r2_kj, r2_kl, m2, n2,
2780 &sign2, &t12, &t22, &t32); /* 84 */
2782 cos_phi2 = cos(phi2);
2784 a2[0] = r2_ij[1]*r2_kj[2]-r2_ij[2]*r2_kj[1];
2785 a2[1] = r2_ij[2]*r2_kj[0]-r2_ij[0]*r2_kj[2];
2786 a2[2] = r2_ij[0]*r2_kj[1]-r2_ij[1]*r2_kj[0]; /* 9 */
2788 b2[0] = r2_kl[1]*r2_kj[2]-r2_kl[2]*r2_kj[1];
2789 b2[1] = r2_kl[2]*r2_kj[0]-r2_kl[0]*r2_kj[2];
2790 b2[2] = r2_kl[0]*r2_kj[1]-r2_kl[1]*r2_kj[0]; /* 9 */
2792 tmp = pbc_rvec_sub(pbc, x[a2l], x[a2k], h2);
2794 ra22 = iprod(a2, a2); /* 5 */
2795 rb22 = iprod(b2, b2); /* 5 */
2796 rg22 = iprod(r2_kj, r2_kj); /* 5 */
2802 rabr2 = sqrt(ra2r2*rb2r2);
2804 sin_phi2 = rg2 * rabr2 * iprod(a2, h2) * (-1);
2806 if (cos_phi2 < -0.5 || cos_phi2 > 0.5)
2808 phi2 = asin(sin_phi2);
2818 phi2 = -M_PI - phi2;
2824 phi2 = acos(cos_phi2);
2832 xphi2 = phi2 + M_PI; /* 1 */
2834 /* Range mangling */
2837 xphi1 = xphi1 + 2*M_PI;
2839 else if (xphi1 >= 2*M_PI)
2841 xphi1 = xphi1 - 2*M_PI;
2846 xphi2 = xphi2 + 2*M_PI;
2848 else if (xphi2 >= 2*M_PI)
2850 xphi2 = xphi2 - 2*M_PI;
2853 /* Number of grid points */
2854 dx = 2*M_PI / cmap_grid->grid_spacing;
2856 /* Where on the grid are we */
2857 iphi1 = (int)(xphi1/dx);
2858 iphi2 = (int)(xphi2/dx);
2860 iphi1 = cmap_setup_grid_index(iphi1, cmap_grid->grid_spacing, &ip1m1, &ip1p1, &ip1p2);
2861 iphi2 = cmap_setup_grid_index(iphi2, cmap_grid->grid_spacing, &ip2m1, &ip2p1, &ip2p2);
2863 pos1 = iphi1*cmap_grid->grid_spacing+iphi2;
2864 pos2 = ip1p1*cmap_grid->grid_spacing+iphi2;
2865 pos3 = ip1p1*cmap_grid->grid_spacing+ip2p1;
2866 pos4 = iphi1*cmap_grid->grid_spacing+ip2p1;
2868 ty[0] = cmapd[pos1*4];
2869 ty[1] = cmapd[pos2*4];
2870 ty[2] = cmapd[pos3*4];
2871 ty[3] = cmapd[pos4*4];
2873 ty1[0] = cmapd[pos1*4+1];
2874 ty1[1] = cmapd[pos2*4+1];
2875 ty1[2] = cmapd[pos3*4+1];
2876 ty1[3] = cmapd[pos4*4+1];
2878 ty2[0] = cmapd[pos1*4+2];
2879 ty2[1] = cmapd[pos2*4+2];
2880 ty2[2] = cmapd[pos3*4+2];
2881 ty2[3] = cmapd[pos4*4+2];
2883 ty12[0] = cmapd[pos1*4+3];
2884 ty12[1] = cmapd[pos2*4+3];
2885 ty12[2] = cmapd[pos3*4+3];
2886 ty12[3] = cmapd[pos4*4+3];
2888 /* Switch to degrees */
2889 dx = 360.0 / cmap_grid->grid_spacing;
2890 xphi1 = xphi1 * RAD2DEG;
2891 xphi2 = xphi2 * RAD2DEG;
2893 for (i = 0; i < 4; i++) /* 16 */
2896 tx[i+4] = ty1[i]*dx;
2897 tx[i+8] = ty2[i]*dx;
2898 tx[i+12] = ty12[i]*dx*dx;
2902 for (i = 0; i < 4; i++) /* 1056 */
2904 for (j = 0; j < 4; j++)
2907 for (k = 0; k < 16; k++)
2909 xx = xx + cmap_coeff_matrix[k*16+idx]*tx[k];
2917 tt = (xphi1-iphi1*dx)/dx;
2918 tu = (xphi2-iphi2*dx)/dx;
2927 for (i = 3; i >= 0; i--)
2929 l1 = loop_index[i][3];
2930 l2 = loop_index[i][2];
2931 l3 = loop_index[i][1];
2933 e = tt * e + ((tc[i*4+3]*tu+tc[i*4+2])*tu + tc[i*4+1])*tu+tc[i*4];
2934 df1 = tu * df1 + (3.0*tc[l1]*tt+2.0*tc[l2])*tt+tc[l3];
2935 df2 = tt * df2 + (3.0*tc[i*4+3]*tu+2.0*tc[i*4+2])*tu+tc[i*4+1];
2936 ddf1 = tu * ddf1 + 2.0*3.0*tc[l1]*tt+2.0*tc[l2];
2937 ddf2 = tt * ddf2 + 2.0*3.0*tc[4*i+3]*tu+2.0*tc[4*i+2];
2940 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) +
2941 3.0*tu*tu*(tc[7]+2.0*tc[11]*tt+3.0*tc[15]*tt*tt);
2946 ddf1 = ddf1 * fac * fac;
2947 ddf2 = ddf2 * fac * fac;
2948 ddf12 = ddf12 * fac * fac;
2953 /* Do forces - first torsion */
2954 fg1 = iprod(r1_ij, r1_kj);
2955 hg1 = iprod(r1_kl, r1_kj);
2956 fga1 = fg1*ra2r1*rgr1;
2957 hgb1 = hg1*rb2r1*rgr1;
2961 for (i = 0; i < DIM; i++)
2963 dtf1[i] = gaa1 * a1[i];
2964 dtg1[i] = fga1 * a1[i] - hgb1 * b1[i];
2965 dth1[i] = gbb1 * b1[i];
2967 f1[i] = df1 * dtf1[i];
2968 g1[i] = df1 * dtg1[i];
2969 h1[i] = df1 * dth1[i];
2972 f1_j[i] = -f1[i] - g1[i];
2973 f1_k[i] = h1[i] + g1[i];
2976 f[a1i][i] = f[a1i][i] + f1_i[i];
2977 f[a1j][i] = f[a1j][i] + f1_j[i]; /* - f1[i] - g1[i] */
2978 f[a1k][i] = f[a1k][i] + f1_k[i]; /* h1[i] + g1[i] */
2979 f[a1l][i] = f[a1l][i] + f1_l[i]; /* h1[i] */
2982 /* Do forces - second torsion */
2983 fg2 = iprod(r2_ij, r2_kj);
2984 hg2 = iprod(r2_kl, r2_kj);
2985 fga2 = fg2*ra2r2*rgr2;
2986 hgb2 = hg2*rb2r2*rgr2;
2990 for (i = 0; i < DIM; i++)
2992 dtf2[i] = gaa2 * a2[i];
2993 dtg2[i] = fga2 * a2[i] - hgb2 * b2[i];
2994 dth2[i] = gbb2 * b2[i];
2996 f2[i] = df2 * dtf2[i];
2997 g2[i] = df2 * dtg2[i];
2998 h2[i] = df2 * dth2[i];
3001 f2_j[i] = -f2[i] - g2[i];
3002 f2_k[i] = h2[i] + g2[i];
3005 f[a2i][i] = f[a2i][i] + f2_i[i]; /* f2[i] */
3006 f[a2j][i] = f[a2j][i] + f2_j[i]; /* - f2[i] - g2[i] */
3007 f[a2k][i] = f[a2k][i] + f2_k[i]; /* h2[i] + g2[i] */
3008 f[a2l][i] = f[a2l][i] + f2_l[i]; /* - h2[i] */
3014 copy_ivec(SHIFT_IVEC(g, a1j), jt1);
3015 ivec_sub(SHIFT_IVEC(g, a1i), jt1, dt1_ij);
3016 ivec_sub(SHIFT_IVEC(g, a1k), jt1, dt1_kj);
3017 ivec_sub(SHIFT_IVEC(g, a1l), jt1, dt1_lj);
3018 t11 = IVEC2IS(dt1_ij);
3019 t21 = IVEC2IS(dt1_kj);
3020 t31 = IVEC2IS(dt1_lj);
3022 copy_ivec(SHIFT_IVEC(g, a2j), jt2);
3023 ivec_sub(SHIFT_IVEC(g, a2i), jt2, dt2_ij);
3024 ivec_sub(SHIFT_IVEC(g, a2k), jt2, dt2_kj);
3025 ivec_sub(SHIFT_IVEC(g, a2l), jt2, dt2_lj);
3026 t12 = IVEC2IS(dt2_ij);
3027 t22 = IVEC2IS(dt2_kj);
3028 t32 = IVEC2IS(dt2_lj);
3032 t31 = pbc_rvec_sub(pbc, x[a1l], x[a1j], h1);
3033 t32 = pbc_rvec_sub(pbc, x[a2l], x[a2j], h2);
3041 rvec_inc(fshift[t11], f1_i);
3042 rvec_inc(fshift[CENTRAL], f1_j);
3043 rvec_inc(fshift[t21], f1_k);
3044 rvec_inc(fshift[t31], f1_l);
3046 rvec_inc(fshift[t21], f2_i);
3047 rvec_inc(fshift[CENTRAL], f2_j);
3048 rvec_inc(fshift[t22], f2_k);
3049 rvec_inc(fshift[t32], f2_l);
3056 /***********************************************************
3058 * G R O M O S 9 6 F U N C T I O N S
3060 ***********************************************************/
3061 real g96harmonic(real kA, real kB, real xA, real xB, real x, real lambda,
3064 const real half = 0.5;
3065 real L1, kk, x0, dx, dx2;
3066 real v, f, dvdlambda;
3069 kk = L1*kA+lambda*kB;
3070 x0 = L1*xA+lambda*xB;
3077 dvdlambda = half*(kB-kA)*dx2 + (xA-xB)*kk*dx;
3084 /* That was 21 flops */
3087 real g96bonds(int nbonds,
3088 const t_iatom forceatoms[], const t_iparams forceparams[],
3089 const rvec x[], rvec f[], rvec fshift[],
3090 const t_pbc *pbc, const t_graph *g,
3091 real lambda, real *dvdlambda,
3092 const t_mdatoms *md, t_fcdata *fcd,
3093 int *global_atom_index)
3095 int i, m, ki, ai, aj, type;
3096 real dr2, fbond, vbond, fij, vtot;
3101 for (i = 0; (i < nbonds); )
3103 type = forceatoms[i++];
3104 ai = forceatoms[i++];
3105 aj = forceatoms[i++];
3107 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
3108 dr2 = iprod(dx, dx); /* 5 */
3110 *dvdlambda += g96harmonic(forceparams[type].harmonic.krA,
3111 forceparams[type].harmonic.krB,
3112 forceparams[type].harmonic.rA,
3113 forceparams[type].harmonic.rB,
3114 dr2, lambda, &vbond, &fbond);
3116 vtot += 0.5*vbond; /* 1*/
3120 fprintf(debug, "G96-BONDS: dr = %10g vbond = %10g fbond = %10g\n",
3121 sqrt(dr2), vbond, fbond);
3127 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
3130 for (m = 0; (m < DIM); m++) /* 15 */
3135 fshift[ki][m] += fij;
3136 fshift[CENTRAL][m] -= fij;
3142 real g96bond_angle(const rvec xi, const rvec xj, const rvec xk, const t_pbc *pbc,
3143 rvec r_ij, rvec r_kj,
3145 /* Return value is the angle between the bonds i-j and j-k */
3149 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
3150 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
3152 costh = cos_angle(r_ij, r_kj); /* 25 */
3157 real g96angles(int nbonds,
3158 const t_iatom forceatoms[], const t_iparams forceparams[],
3159 const rvec x[], rvec f[], rvec fshift[],
3160 const t_pbc *pbc, const t_graph *g,
3161 real lambda, real *dvdlambda,
3162 const t_mdatoms *md, t_fcdata *fcd,
3163 int *global_atom_index)
3165 int i, ai, aj, ak, type, m, t1, t2;
3167 real cos_theta, dVdt, va, vtot;
3168 real rij_1, rij_2, rkj_1, rkj_2, rijrkj_1;
3170 ivec jt, dt_ij, dt_kj;
3173 for (i = 0; (i < nbonds); )
3175 type = forceatoms[i++];
3176 ai = forceatoms[i++];
3177 aj = forceatoms[i++];
3178 ak = forceatoms[i++];
3180 cos_theta = g96bond_angle(x[ai], x[aj], x[ak], pbc, r_ij, r_kj, &t1, &t2);
3182 *dvdlambda += g96harmonic(forceparams[type].harmonic.krA,
3183 forceparams[type].harmonic.krB,
3184 forceparams[type].harmonic.rA,
3185 forceparams[type].harmonic.rB,
3186 cos_theta, lambda, &va, &dVdt);
3189 rij_1 = gmx_invsqrt(iprod(r_ij, r_ij));
3190 rkj_1 = gmx_invsqrt(iprod(r_kj, r_kj));
3191 rij_2 = rij_1*rij_1;
3192 rkj_2 = rkj_1*rkj_1;
3193 rijrkj_1 = rij_1*rkj_1; /* 23 */
3198 fprintf(debug, "G96ANGLES: costheta = %10g vth = %10g dV/dct = %10g\n",
3199 cos_theta, va, dVdt);
3202 for (m = 0; (m < DIM); m++) /* 42 */
3204 f_i[m] = dVdt*(r_kj[m]*rijrkj_1 - r_ij[m]*rij_2*cos_theta);
3205 f_k[m] = dVdt*(r_ij[m]*rijrkj_1 - r_kj[m]*rkj_2*cos_theta);
3206 f_j[m] = -f_i[m]-f_k[m];
3214 copy_ivec(SHIFT_IVEC(g, aj), jt);
3216 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3217 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3218 t1 = IVEC2IS(dt_ij);
3219 t2 = IVEC2IS(dt_kj);
3221 rvec_inc(fshift[t1], f_i);
3222 rvec_inc(fshift[CENTRAL], f_j);
3223 rvec_inc(fshift[t2], f_k); /* 9 */
3229 real cross_bond_bond(int nbonds,
3230 const t_iatom forceatoms[], const t_iparams forceparams[],
3231 const rvec x[], rvec f[], rvec fshift[],
3232 const t_pbc *pbc, const t_graph *g,
3233 real lambda, real *dvdlambda,
3234 const t_mdatoms *md, t_fcdata *fcd,
3235 int *global_atom_index)
3237 /* Potential from Lawrence and Skimmer, Chem. Phys. Lett. 372 (2003)
3240 int i, ai, aj, ak, type, m, t1, t2;
3242 real vtot, vrr, s1, s2, r1, r2, r1e, r2e, krr;
3244 ivec jt, dt_ij, dt_kj;
3247 for (i = 0; (i < nbonds); )
3249 type = forceatoms[i++];
3250 ai = forceatoms[i++];
3251 aj = forceatoms[i++];
3252 ak = forceatoms[i++];
3253 r1e = forceparams[type].cross_bb.r1e;
3254 r2e = forceparams[type].cross_bb.r2e;
3255 krr = forceparams[type].cross_bb.krr;
3257 /* Compute distance vectors ... */
3258 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
3259 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
3261 /* ... and their lengths */
3265 /* Deviations from ideality */
3269 /* Energy (can be negative!) */
3274 svmul(-krr*s2/r1, r_ij, f_i);
3275 svmul(-krr*s1/r2, r_kj, f_k);
3277 for (m = 0; (m < DIM); m++) /* 12 */
3279 f_j[m] = -f_i[m] - f_k[m];
3288 copy_ivec(SHIFT_IVEC(g, aj), jt);
3290 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3291 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3292 t1 = IVEC2IS(dt_ij);
3293 t2 = IVEC2IS(dt_kj);
3295 rvec_inc(fshift[t1], f_i);
3296 rvec_inc(fshift[CENTRAL], f_j);
3297 rvec_inc(fshift[t2], f_k); /* 9 */
3303 real cross_bond_angle(int nbonds,
3304 const t_iatom forceatoms[], const t_iparams forceparams[],
3305 const rvec x[], rvec f[], rvec fshift[],
3306 const t_pbc *pbc, const t_graph *g,
3307 real lambda, real *dvdlambda,
3308 const t_mdatoms *md, t_fcdata *fcd,
3309 int *global_atom_index)
3311 /* Potential from Lawrence and Skimmer, Chem. Phys. Lett. 372 (2003)
3314 int i, ai, aj, ak, type, m, t1, t2, t3;
3315 rvec r_ij, r_kj, r_ik;
3316 real vtot, vrt, s1, s2, s3, r1, r2, r3, r1e, r2e, r3e, krt, k1, k2, k3;
3318 ivec jt, dt_ij, dt_kj;
3321 for (i = 0; (i < nbonds); )
3323 type = forceatoms[i++];
3324 ai = forceatoms[i++];
3325 aj = forceatoms[i++];
3326 ak = forceatoms[i++];
3327 r1e = forceparams[type].cross_ba.r1e;
3328 r2e = forceparams[type].cross_ba.r2e;
3329 r3e = forceparams[type].cross_ba.r3e;
3330 krt = forceparams[type].cross_ba.krt;
3332 /* Compute distance vectors ... */
3333 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
3334 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
3335 t3 = pbc_rvec_sub(pbc, x[ai], x[ak], r_ik);
3337 /* ... and their lengths */
3342 /* Deviations from ideality */
3347 /* Energy (can be negative!) */
3348 vrt = krt*s3*(s1+s2);
3354 k3 = -krt*(s1+s2)/r3;
3355 for (m = 0; (m < DIM); m++)
3357 f_i[m] = k1*r_ij[m] + k3*r_ik[m];
3358 f_k[m] = k2*r_kj[m] - k3*r_ik[m];
3359 f_j[m] = -f_i[m] - f_k[m];
3362 for (m = 0; (m < DIM); m++) /* 12 */
3372 copy_ivec(SHIFT_IVEC(g, aj), jt);
3374 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3375 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3376 t1 = IVEC2IS(dt_ij);
3377 t2 = IVEC2IS(dt_kj);
3379 rvec_inc(fshift[t1], f_i);
3380 rvec_inc(fshift[CENTRAL], f_j);
3381 rvec_inc(fshift[t2], f_k); /* 9 */
3387 static real bonded_tab(const char *type, int table_nr,
3388 const bondedtable_t *table, real kA, real kB, real r,
3389 real lambda, real *V, real *F)
3391 real k, tabscale, *VFtab, rt, eps, eps2, Yt, Ft, Geps, Heps2, Fp, VV, FF;
3393 real v, f, dvdlambda;
3395 k = (1.0 - lambda)*kA + lambda*kB;
3397 tabscale = table->scale;
3398 VFtab = table->data;
3404 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",
3405 type, table_nr, r, n0, n0+1, table->n);
3412 Geps = VFtab[nnn+2]*eps;
3413 Heps2 = VFtab[nnn+3]*eps2;
3414 Fp = Ft + Geps + Heps2;
3416 FF = Fp + Geps + 2.0*Heps2;
3418 *F = -k*FF*tabscale;
3420 dvdlambda = (kB - kA)*VV;
3424 /* That was 22 flops */
3427 real tab_bonds(int nbonds,
3428 const t_iatom forceatoms[], const t_iparams forceparams[],
3429 const rvec x[], rvec f[], rvec fshift[],
3430 const t_pbc *pbc, const t_graph *g,
3431 real lambda, real *dvdlambda,
3432 const t_mdatoms *md, t_fcdata *fcd,
3433 int *global_atom_index)
3435 int i, m, ki, ai, aj, type, table;
3436 real dr, dr2, fbond, vbond, fij, vtot;
3441 for (i = 0; (i < nbonds); )
3443 type = forceatoms[i++];
3444 ai = forceatoms[i++];
3445 aj = forceatoms[i++];
3447 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
3448 dr2 = iprod(dx, dx); /* 5 */
3449 dr = dr2*gmx_invsqrt(dr2); /* 10 */
3451 table = forceparams[type].tab.table;
3453 *dvdlambda += bonded_tab("bond", table,
3454 &fcd->bondtab[table],
3455 forceparams[type].tab.kA,
3456 forceparams[type].tab.kB,
3457 dr, lambda, &vbond, &fbond); /* 22 */
3465 vtot += vbond; /* 1*/
3466 fbond *= gmx_invsqrt(dr2); /* 6 */
3470 fprintf(debug, "TABBONDS: dr = %10g vbond = %10g fbond = %10g\n",
3476 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
3479 for (m = 0; (m < DIM); m++) /* 15 */
3484 fshift[ki][m] += fij;
3485 fshift[CENTRAL][m] -= fij;
3491 real tab_angles(int nbonds,
3492 const t_iatom forceatoms[], const t_iparams forceparams[],
3493 const rvec x[], rvec f[], rvec fshift[],
3494 const t_pbc *pbc, const t_graph *g,
3495 real lambda, real *dvdlambda,
3496 const t_mdatoms *md, t_fcdata *fcd,
3497 int *global_atom_index)
3499 int i, ai, aj, ak, t1, t2, type, table;
3501 real cos_theta, cos_theta2, theta, dVdt, va, vtot;
3502 ivec jt, dt_ij, dt_kj;
3505 for (i = 0; (i < nbonds); )
3507 type = forceatoms[i++];
3508 ai = forceatoms[i++];
3509 aj = forceatoms[i++];
3510 ak = forceatoms[i++];
3512 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
3513 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
3515 table = forceparams[type].tab.table;
3517 *dvdlambda += bonded_tab("angle", table,
3518 &fcd->angletab[table],
3519 forceparams[type].tab.kA,
3520 forceparams[type].tab.kB,
3521 theta, lambda, &va, &dVdt); /* 22 */
3524 cos_theta2 = sqr(cos_theta); /* 1 */
3533 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
3534 sth = st*cos_theta; /* 1 */
3538 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
3539 theta*RAD2DEG, va, dVdt);
3542 nrkj2 = iprod(r_kj, r_kj); /* 5 */
3543 nrij2 = iprod(r_ij, r_ij);
3545 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
3546 cii = sth/nrij2; /* 10 */
3547 ckk = sth/nrkj2; /* 10 */
3549 for (m = 0; (m < DIM); m++) /* 39 */
3551 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
3552 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
3553 f_j[m] = -f_i[m]-f_k[m];
3560 copy_ivec(SHIFT_IVEC(g, aj), jt);
3562 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3563 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3564 t1 = IVEC2IS(dt_ij);
3565 t2 = IVEC2IS(dt_kj);
3567 rvec_inc(fshift[t1], f_i);
3568 rvec_inc(fshift[CENTRAL], f_j);
3569 rvec_inc(fshift[t2], f_k);
3575 real tab_dihs(int nbonds,
3576 const t_iatom forceatoms[], const t_iparams forceparams[],
3577 const rvec x[], rvec f[], rvec fshift[],
3578 const t_pbc *pbc, const t_graph *g,
3579 real lambda, real *dvdlambda,
3580 const t_mdatoms *md, t_fcdata *fcd,
3581 int *global_atom_index)
3583 int i, type, ai, aj, ak, al, table;
3585 rvec r_ij, r_kj, r_kl, m, n;
3586 real phi, sign, ddphi, vpd, vtot;
3589 for (i = 0; (i < nbonds); )
3591 type = forceatoms[i++];
3592 ai = forceatoms[i++];
3593 aj = forceatoms[i++];
3594 ak = forceatoms[i++];
3595 al = forceatoms[i++];
3597 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
3598 &sign, &t1, &t2, &t3); /* 84 */
3600 table = forceparams[type].tab.table;
3602 /* Hopefully phi+M_PI never results in values < 0 */
3603 *dvdlambda += bonded_tab("dihedral", table,
3604 &fcd->dihtab[table],
3605 forceparams[type].tab.kA,
3606 forceparams[type].tab.kB,
3607 phi+M_PI, lambda, &vpd, &ddphi);
3610 do_dih_fup(ai, aj, ak, al, -ddphi, r_ij, r_kj, r_kl, m, n,
3611 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
3614 fprintf(debug, "pdih: (%d,%d,%d,%d) phi=%g\n",
3615 ai, aj, ak, al, phi);
3622 /* Return if this is a potential calculated in bondfree.c,
3623 * i.e. an interaction that actually calculates a potential and
3624 * works on multiple atoms (not e.g. a connection or a position restraint).
3626 static gmx_inline gmx_bool ftype_is_bonded_potential(int ftype)
3629 (interaction_function[ftype].flags & IF_BOND) &&
3630 !(ftype == F_CONNBONDS || ftype == F_POSRES) &&
3631 (ftype < F_GB12 || ftype > F_GB14);
3634 static void divide_bondeds_over_threads(t_idef *idef, int nthreads)
3641 idef->nthreads = nthreads;
3643 if (F_NRE*(nthreads+1) > idef->il_thread_division_nalloc)
3645 idef->il_thread_division_nalloc = F_NRE*(nthreads+1);
3646 snew(idef->il_thread_division, idef->il_thread_division_nalloc);
3649 for (ftype = 0; ftype < F_NRE; ftype++)
3651 if (ftype_is_bonded_potential(ftype))
3653 nat1 = interaction_function[ftype].nratoms + 1;
3655 for (t = 0; t <= nthreads; t++)
3657 /* Divide the interactions equally over the threads.
3658 * When the different types of bonded interactions
3659 * are distributed roughly equally over the threads,
3660 * this should lead to well localized output into
3661 * the force buffer on each thread.
3662 * If this is not the case, a more advanced scheme
3663 * (not implemented yet) will do better.
3665 il_nr_thread = (((idef->il[ftype].nr/nat1)*t)/nthreads)*nat1;
3667 /* Ensure that distance restraint pairs with the same label
3668 * end up on the same thread.
3669 * This is slighlty tricky code, since the next for iteration
3670 * may have an initial il_nr_thread lower than the final value
3671 * in the previous iteration, but this will anyhow be increased
3672 * to the approriate value again by this while loop.
3674 while (ftype == F_DISRES &&
3676 il_nr_thread < idef->il[ftype].nr &&
3677 idef->iparams[idef->il[ftype].iatoms[il_nr_thread]].disres.label ==
3678 idef->iparams[idef->il[ftype].iatoms[il_nr_thread-nat1]].disres.label)
3680 il_nr_thread += nat1;
3683 idef->il_thread_division[ftype*(nthreads+1)+t] = il_nr_thread;
3690 calc_bonded_reduction_mask(const t_idef *idef,
3695 int ftype, nb, nat1, nb0, nb1, i, a;
3699 for (ftype = 0; ftype < F_NRE; ftype++)
3701 if (ftype_is_bonded_potential(ftype))
3703 nb = idef->il[ftype].nr;
3706 nat1 = interaction_function[ftype].nratoms + 1;
3708 /* Divide this interaction equally over the threads.
3709 * This is not stored: should match division in calc_bonds.
3711 nb0 = idef->il_thread_division[ftype*(nt+1)+t];
3712 nb1 = idef->il_thread_division[ftype*(nt+1)+t+1];
3714 for (i = nb0; i < nb1; i += nat1)
3716 for (a = 1; a < nat1; a++)
3718 mask |= (1U << (idef->il[ftype].iatoms[i+a]>>shift));
3728 void setup_bonded_threading(t_forcerec *fr, t_idef *idef)
3730 #define MAX_BLOCK_BITS 32
3735 assert(fr->nthreads >= 1);
3738 /* Divide the bonded interaction over the threads */
3739 divide_bondeds_over_threads(idef, fr->nthreads);
3741 if (fr->nthreads == 1)
3748 /* We divide the force array in a maximum of 32 blocks.
3749 * Minimum force block reduction size is 2^6=64.
3752 while (fr->natoms_force > (int)(MAX_BLOCK_BITS*(1U<<fr->red_ashift)))
3758 fprintf(debug, "bonded force buffer block atom shift %d bits\n",
3762 /* Determine to which blocks each thread's bonded force calculation
3763 * contributes. Store this is a mask for each thread.
3765 #pragma omp parallel for num_threads(fr->nthreads) schedule(static)
3766 for (t = 1; t < fr->nthreads; t++)
3768 fr->f_t[t].red_mask =
3769 calc_bonded_reduction_mask(idef, fr->red_ashift, t, fr->nthreads);
3772 /* Determine the maximum number of blocks we need to reduce over */
3775 for (t = 0; t < fr->nthreads; t++)
3778 for (b = 0; b < MAX_BLOCK_BITS; b++)
3780 if (fr->f_t[t].red_mask & (1U<<b))
3782 fr->red_nblock = max(fr->red_nblock, b+1);
3788 fprintf(debug, "thread %d flags %x count %d\n",
3789 t, fr->f_t[t].red_mask, c);
3795 fprintf(debug, "Number of blocks to reduce: %d of size %d\n",
3796 fr->red_nblock, 1<<fr->red_ashift);
3797 fprintf(debug, "Reduction density %.2f density/#thread %.2f\n",
3798 ctot*(1<<fr->red_ashift)/(double)fr->natoms_force,
3799 ctot*(1<<fr->red_ashift)/(double)(fr->natoms_force*fr->nthreads));
3803 static void zero_thread_forces(f_thread_t *f_t, int n,
3804 int nblock, int blocksize)
3806 int b, a0, a1, a, i, j;
3808 if (n > f_t->f_nalloc)
3810 f_t->f_nalloc = over_alloc_large(n);
3811 srenew(f_t->f, f_t->f_nalloc);
3814 if (f_t->red_mask != 0)
3816 for (b = 0; b < nblock; b++)
3818 if (f_t->red_mask && (1U<<b))
3821 a1 = min((b+1)*blocksize, n);
3822 for (a = a0; a < a1; a++)
3824 clear_rvec(f_t->f[a]);
3829 for (i = 0; i < SHIFTS; i++)
3831 clear_rvec(f_t->fshift[i]);
3833 for (i = 0; i < F_NRE; i++)
3837 for (i = 0; i < egNR; i++)
3839 for (j = 0; j < f_t->grpp.nener; j++)
3841 f_t->grpp.ener[i][j] = 0;
3844 for (i = 0; i < efptNR; i++)
3850 static void reduce_thread_force_buffer(int n, rvec *f,
3851 int nthreads, f_thread_t *f_t,
3852 int nblock, int block_size)
3854 /* The max thread number is arbitrary,
3855 * we used a fixed number to avoid memory management.
3856 * Using more than 16 threads is probably never useful performance wise.
3858 #define MAX_BONDED_THREADS 256
3861 if (nthreads > MAX_BONDED_THREADS)
3863 gmx_fatal(FARGS, "Can not reduce bonded forces on more than %d threads",
3864 MAX_BONDED_THREADS);
3867 /* This reduction can run on any number of threads,
3868 * independently of nthreads.
3870 #pragma omp parallel for num_threads(nthreads) schedule(static)
3871 for (b = 0; b < nblock; b++)
3873 rvec *fp[MAX_BONDED_THREADS];
3877 /* Determine which threads contribute to this block */
3879 for (ft = 1; ft < nthreads; ft++)
3881 if (f_t[ft].red_mask & (1U<<b))
3883 fp[nfb++] = f_t[ft].f;
3888 /* Reduce force buffers for threads that contribute */
3890 a1 = (b+1)*block_size;
3892 for (a = a0; a < a1; a++)
3894 for (fb = 0; fb < nfb; fb++)
3896 rvec_inc(f[a], fp[fb][a]);
3903 static void reduce_thread_forces(int n, rvec *f, rvec *fshift,
3904 real *ener, gmx_grppairener_t *grpp, real *dvdl,
3905 int nthreads, f_thread_t *f_t,
3906 int nblock, int block_size,
3907 gmx_bool bCalcEnerVir,
3912 /* Reduce the bonded force buffer */
3913 reduce_thread_force_buffer(n, f, nthreads, f_t, nblock, block_size);
3916 /* When necessary, reduce energy and virial using one thread only */
3921 for (i = 0; i < SHIFTS; i++)
3923 for (t = 1; t < nthreads; t++)
3925 rvec_inc(fshift[i], f_t[t].fshift[i]);
3928 for (i = 0; i < F_NRE; i++)
3930 for (t = 1; t < nthreads; t++)
3932 ener[i] += f_t[t].ener[i];
3935 for (i = 0; i < egNR; i++)
3937 for (j = 0; j < f_t[1].grpp.nener; j++)
3939 for (t = 1; t < nthreads; t++)
3942 grpp->ener[i][j] += f_t[t].grpp.ener[i][j];
3948 for (i = 0; i < efptNR; i++)
3951 for (t = 1; t < nthreads; t++)
3953 dvdl[i] += f_t[t].dvdl[i];
3960 static real calc_one_bond(FILE *fplog, int thread,
3961 int ftype, const t_idef *idef,
3962 rvec x[], rvec f[], rvec fshift[],
3964 const t_pbc *pbc, const t_graph *g,
3965 gmx_grppairener_t *grpp,
3967 real *lambda, real *dvdl,
3968 const t_mdatoms *md, t_fcdata *fcd,
3969 gmx_bool bCalcEnerVir,
3970 int *global_atom_index, gmx_bool bPrintSepPot)
3972 int nat1, nbonds, efptFTYPE;
3977 if (IS_RESTRAINT_TYPE(ftype))
3979 efptFTYPE = efptRESTRAINT;
3983 efptFTYPE = efptBONDED;
3986 nat1 = interaction_function[ftype].nratoms + 1;
3987 nbonds = idef->il[ftype].nr/nat1;
3988 iatoms = idef->il[ftype].iatoms;
3990 nb0 = idef->il_thread_division[ftype*(idef->nthreads+1)+thread];
3991 nbn = idef->il_thread_division[ftype*(idef->nthreads+1)+thread+1] - nb0;
3993 if (!IS_LISTED_LJ_C(ftype))
3995 if (ftype == F_CMAP)
3997 v = cmap_dihs(nbn, iatoms+nb0,
3998 idef->iparams, &idef->cmap_grid,
3999 (const rvec*)x, f, fshift,
4000 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
4001 md, fcd, global_atom_index);
4004 else if (ftype == F_ANGLES &&
4005 !bCalcEnerVir && fr->efep == efepNO)
4007 /* No energies, shift forces, dvdl */
4008 angles_noener_simd(nbn, idef->il[ftype].iatoms+nb0,
4011 pbc, g, lambda[efptFTYPE], md, fcd,
4016 else if (ftype == F_PDIHS &&
4017 !bCalcEnerVir && fr->efep == efepNO)
4019 /* No energies, shift forces, dvdl */
4020 #ifndef SIMD_BONDEDS
4025 (nbn, idef->il[ftype].iatoms+nb0,
4028 pbc, g, lambda[efptFTYPE], md, fcd,
4034 v = interaction_function[ftype].ifunc(nbn, iatoms+nb0,
4036 (const rvec*)x, f, fshift,
4037 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
4038 md, fcd, global_atom_index);
4042 fprintf(fplog, " %-23s #%4d V %12.5e dVdl %12.5e\n",
4043 interaction_function[ftype].longname,
4044 nbonds, v, lambda[efptFTYPE]);
4049 v = do_nonbonded_listed(ftype, nbn, iatoms+nb0, idef->iparams, (const rvec*)x, f, fshift,
4050 pbc, g, lambda, dvdl, md, fr, grpp, global_atom_index);
4054 fprintf(fplog, " %-5s + %-15s #%4d dVdl %12.5e\n",
4055 interaction_function[ftype].longname,
4056 interaction_function[F_LJ14].longname, nbonds, dvdl[efptVDW]);
4057 fprintf(fplog, " %-5s + %-15s #%4d dVdl %12.5e\n",
4058 interaction_function[ftype].longname,
4059 interaction_function[F_COUL14].longname, nbonds, dvdl[efptCOUL]);
4065 inc_nrnb(nrnb, interaction_function[ftype].nrnb_ind, nbonds);
4071 void calc_bonds(FILE *fplog, const gmx_multisim_t *ms,
4073 rvec x[], history_t *hist,
4074 rvec f[], t_forcerec *fr,
4075 const t_pbc *pbc, const t_graph *g,
4076 gmx_enerdata_t *enerd, t_nrnb *nrnb,
4078 const t_mdatoms *md,
4079 t_fcdata *fcd, int *global_atom_index,
4080 t_atomtypes *atype, gmx_genborn_t *born,
4082 gmx_bool bPrintSepPot, gmx_large_int_t step)
4084 gmx_bool bCalcEnerVir;
4086 real v, dvdl[efptNR], dvdl_dum[efptNR]; /* The dummy array is to have a place to store the dhdl at other values
4087 of lambda, which will be thrown away in the end*/
4088 const t_pbc *pbc_null;
4093 assert(fr->nthreads == idef->nthreads);
4096 bCalcEnerVir = (force_flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY));
4098 for (i = 0; i < efptNR; i++)
4112 fprintf(fplog, "Step %s: bonded V and dVdl for this node\n",
4113 gmx_step_str(step, buf));
4119 p_graph(debug, "Bondage is fun", g);
4123 /* Do pre force calculation stuff which might require communication */
4124 if (idef->il[F_ORIRES].nr)
4126 enerd->term[F_ORIRESDEV] =
4127 calc_orires_dev(ms, idef->il[F_ORIRES].nr,
4128 idef->il[F_ORIRES].iatoms,
4129 idef->iparams, md, (const rvec*)x,
4130 pbc_null, fcd, hist);
4132 if (idef->il[F_DISRES].nr)
4134 calc_disres_R_6(ms, idef->il[F_DISRES].nr,
4135 idef->il[F_DISRES].iatoms,
4136 idef->iparams, (const rvec*)x, pbc_null,
4140 #pragma omp parallel for num_threads(fr->nthreads) schedule(static)
4141 for (thread = 0; thread < fr->nthreads; thread++)
4148 gmx_grppairener_t *grpp;
4153 fshift = fr->fshift;
4155 grpp = &enerd->grpp;
4160 zero_thread_forces(&fr->f_t[thread], fr->natoms_force,
4161 fr->red_nblock, 1<<fr->red_ashift);
4163 ft = fr->f_t[thread].f;
4164 fshift = fr->f_t[thread].fshift;
4165 epot = fr->f_t[thread].ener;
4166 grpp = &fr->f_t[thread].grpp;
4167 dvdlt = fr->f_t[thread].dvdl;
4169 /* Loop over all bonded force types to calculate the bonded forces */
4170 for (ftype = 0; (ftype < F_NRE); ftype++)
4172 if (idef->il[ftype].nr > 0 && ftype_is_bonded_potential(ftype))
4174 v = calc_one_bond(fplog, thread, ftype, idef, x,
4175 ft, fshift, fr, pbc_null, g, grpp,
4176 nrnb, lambda, dvdlt,
4177 md, fcd, bCalcEnerVir,
4178 global_atom_index, bPrintSepPot);
4183 if (fr->nthreads > 1)
4185 reduce_thread_forces(fr->natoms_force, f, fr->fshift,
4186 enerd->term, &enerd->grpp, dvdl,
4187 fr->nthreads, fr->f_t,
4188 fr->red_nblock, 1<<fr->red_ashift,
4190 force_flags & GMX_FORCE_DHDL);
4192 if (force_flags & GMX_FORCE_DHDL)
4194 for (i = 0; i < efptNR; i++)
4196 enerd->dvdl_nonlin[i] += dvdl[i];
4200 /* Copy the sum of violations for the distance restraints from fcd */
4203 enerd->term[F_DISRESVIOL] = fcd->disres.sumviol;
4208 void calc_bonds_lambda(FILE *fplog,
4212 const t_pbc *pbc, const t_graph *g,
4213 gmx_grppairener_t *grpp, real *epot, t_nrnb *nrnb,
4215 const t_mdatoms *md,
4217 int *global_atom_index)
4219 int i, ftype, nr_nonperturbed, nr;
4221 real dvdl_dum[efptNR];
4223 const t_pbc *pbc_null;
4235 /* Copy the whole idef, so we can modify the contents locally */
4237 idef_fe.nthreads = 1;
4238 snew(idef_fe.il_thread_division, F_NRE*(idef_fe.nthreads+1));
4240 /* We already have the forces, so we use temp buffers here */
4241 snew(f, fr->natoms_force);
4242 snew(fshift, SHIFTS);
4244 /* Loop over all bonded force types to calculate the bonded energies */
4245 for (ftype = 0; (ftype < F_NRE); ftype++)
4247 if (ftype_is_bonded_potential(ftype))
4249 /* Set the work range of thread 0 to the perturbed bondeds only */
4250 nr_nonperturbed = idef->il[ftype].nr_nonperturbed;
4251 nr = idef->il[ftype].nr;
4252 idef_fe.il_thread_division[ftype*2+0] = nr_nonperturbed;
4253 idef_fe.il_thread_division[ftype*2+1] = nr;
4255 /* This is only to get the flop count correct */
4256 idef_fe.il[ftype].nr = nr - nr_nonperturbed;
4258 if (nr - nr_nonperturbed > 0)
4260 v = calc_one_bond(fplog, 0, ftype, &idef_fe,
4261 x, f, fshift, fr, pbc_null, g,
4262 grpp, nrnb, lambda, dvdl_dum,
4264 global_atom_index, FALSE);
4273 sfree(idef_fe.il_thread_division);