2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
45 #include "gromacs/math/utilities.h"
53 #include "gmx_fatal.h"
59 #include "nonbonded.h"
61 /* Include the SIMD macro file and then check for support */
62 #include "gromacs/simd/macros.h"
63 #if defined GMX_HAVE_SIMD_MACROS && defined GMX_SIMD_HAVE_TRIGONOMETRIC
65 #include "gromacs/simd/vector_operations.h"
68 /* Find a better place for this? */
69 const int cmap_coeff_matrix[] = {
70 1, 0, -3, 2, 0, 0, 0, 0, -3, 0, 9, -6, 2, 0, -6, 4,
71 0, 0, 0, 0, 0, 0, 0, 0, 3, 0, -9, 6, -2, 0, 6, -4,
72 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 9, -6, 0, 0, -6, 4,
73 0, 0, 3, -2, 0, 0, 0, 0, 0, 0, -9, 6, 0, 0, 6, -4,
74 0, 0, 0, 0, 1, 0, -3, 2, -2, 0, 6, -4, 1, 0, -3, 2,
75 0, 0, 0, 0, 0, 0, 0, 0, -1, 0, 3, -2, 1, 0, -3, 2,
76 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -3, 2, 0, 0, 3, -2,
77 0, 0, 0, 0, 0, 0, 3, -2, 0, 0, -6, 4, 0, 0, 3, -2,
78 0, 1, -2, 1, 0, 0, 0, 0, 0, -3, 6, -3, 0, 2, -4, 2,
79 0, 0, 0, 0, 0, 0, 0, 0, 0, 3, -6, 3, 0, -2, 4, -2,
80 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -3, 3, 0, 0, 2, -2,
81 0, 0, -1, 1, 0, 0, 0, 0, 0, 0, 3, -3, 0, 0, -2, 2,
82 0, 0, 0, 0, 0, 1, -2, 1, 0, -2, 4, -2, 0, 1, -2, 1,
83 0, 0, 0, 0, 0, 0, 0, 0, 0, -1, 2, -1, 0, 1, -2, 1,
84 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, -1, 0, 0, -1, 1,
85 0, 0, 0, 0, 0, 0, -1, 1, 0, 0, 2, -2, 0, 0, -1, 1
90 int glatnr(int *global_atom_index, int i)
94 if (global_atom_index == NULL)
100 atnr = global_atom_index[i] + 1;
106 static int pbc_rvec_sub(const t_pbc *pbc, const rvec xi, const rvec xj, rvec dx)
110 return pbc_dx_aiuc(pbc, xi, xj, dx);
114 rvec_sub(xi, xj, dx);
121 /* SIMD PBC data structure, containing 1/boxdiag and the box vectors */
123 gmx_simd_real_t inv_bzz;
124 gmx_simd_real_t inv_byy;
125 gmx_simd_real_t inv_bxx;
134 /* Set the SIMD pbc data from a normal t_pbc struct */
135 static void set_pbc_simd(const t_pbc *pbc, pbc_simd_t *pbc_simd)
140 /* Setting inv_bdiag to 0 effectively turns off PBC */
141 clear_rvec(inv_bdiag);
144 for (d = 0; d < pbc->ndim_ePBC; d++)
146 inv_bdiag[d] = 1.0/pbc->box[d][d];
150 pbc_simd->inv_bzz = gmx_simd_set1_r(inv_bdiag[ZZ]);
151 pbc_simd->inv_byy = gmx_simd_set1_r(inv_bdiag[YY]);
152 pbc_simd->inv_bxx = gmx_simd_set1_r(inv_bdiag[XX]);
156 pbc_simd->bzx = gmx_simd_set1_r(pbc->box[ZZ][XX]);
157 pbc_simd->bzy = gmx_simd_set1_r(pbc->box[ZZ][YY]);
158 pbc_simd->bzz = gmx_simd_set1_r(pbc->box[ZZ][ZZ]);
159 pbc_simd->byx = gmx_simd_set1_r(pbc->box[YY][XX]);
160 pbc_simd->byy = gmx_simd_set1_r(pbc->box[YY][YY]);
161 pbc_simd->bxx = gmx_simd_set1_r(pbc->box[XX][XX]);
165 pbc_simd->bzx = gmx_simd_setzero_r();
166 pbc_simd->bzy = gmx_simd_setzero_r();
167 pbc_simd->bzz = gmx_simd_setzero_r();
168 pbc_simd->byx = gmx_simd_setzero_r();
169 pbc_simd->byy = gmx_simd_setzero_r();
170 pbc_simd->bxx = gmx_simd_setzero_r();
174 /* Correct distance vector *dx,*dy,*dz for PBC using SIMD */
175 static gmx_inline void
176 pbc_dx_simd(gmx_simd_real_t *dx, gmx_simd_real_t *dy, gmx_simd_real_t *dz,
177 const pbc_simd_t *pbc)
181 sh = gmx_simd_round_r(gmx_simd_mul_r(*dz, pbc->inv_bzz));
182 *dx = gmx_simd_fnmadd_r(sh, pbc->bzx, *dx);
183 *dy = gmx_simd_fnmadd_r(sh, pbc->bzy, *dy);
184 *dz = gmx_simd_fnmadd_r(sh, pbc->bzz, *dz);
186 sh = gmx_simd_round_r(gmx_simd_mul_r(*dy, pbc->inv_byy));
187 *dx = gmx_simd_fnmadd_r(sh, pbc->byx, *dx);
188 *dy = gmx_simd_fnmadd_r(sh, pbc->byy, *dy);
190 sh = gmx_simd_round_r(gmx_simd_mul_r(*dx, pbc->inv_bxx));
191 *dx = gmx_simd_fnmadd_r(sh, pbc->bxx, *dx);
194 #endif /* SIMD_BONDEDS */
197 * Morse potential bond by Frank Everdij
199 * Three parameters needed:
201 * b0 = equilibrium distance in nm
202 * be = beta in nm^-1 (actually, it's nu_e*Sqrt(2*pi*pi*mu/D_e))
203 * cb = well depth in kJ/mol
205 * Note: the potential is referenced to be +cb at infinite separation
206 * and zero at the equilibrium distance!
209 real morse_bonds(int nbonds,
210 const t_iatom forceatoms[], const t_iparams forceparams[],
211 const rvec x[], rvec f[], rvec fshift[],
212 const t_pbc *pbc, const t_graph *g,
213 real lambda, real *dvdlambda,
214 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
215 int gmx_unused *global_atom_index)
217 const real one = 1.0;
218 const real two = 2.0;
219 real dr, dr2, temp, omtemp, cbomtemp, fbond, vbond, fij, vtot;
220 real b0, be, cb, b0A, beA, cbA, b0B, beB, cbB, L1;
222 int i, m, ki, type, ai, aj;
226 for (i = 0; (i < nbonds); )
228 type = forceatoms[i++];
229 ai = forceatoms[i++];
230 aj = forceatoms[i++];
232 b0A = forceparams[type].morse.b0A;
233 beA = forceparams[type].morse.betaA;
234 cbA = forceparams[type].morse.cbA;
236 b0B = forceparams[type].morse.b0B;
237 beB = forceparams[type].morse.betaB;
238 cbB = forceparams[type].morse.cbB;
240 L1 = one-lambda; /* 1 */
241 b0 = L1*b0A + lambda*b0B; /* 3 */
242 be = L1*beA + lambda*beB; /* 3 */
243 cb = L1*cbA + lambda*cbB; /* 3 */
245 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
246 dr2 = iprod(dx, dx); /* 5 */
247 dr = dr2*gmx_invsqrt(dr2); /* 10 */
248 temp = exp(-be*(dr-b0)); /* 12 */
252 /* bonds are constrainted. This may _not_ include bond constraints if they are lambda dependent */
253 *dvdlambda += cbB-cbA;
257 omtemp = one-temp; /* 1 */
258 cbomtemp = cb*omtemp; /* 1 */
259 vbond = cbomtemp*omtemp; /* 1 */
260 fbond = -two*be*temp*cbomtemp*gmx_invsqrt(dr2); /* 9 */
261 vtot += vbond; /* 1 */
263 *dvdlambda += (cbB - cbA) * omtemp * omtemp - (2-2*omtemp)*omtemp * cb * ((b0B-b0A)*be - (beB-beA)*(dr-b0)); /* 15 */
267 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
271 for (m = 0; (m < DIM); m++) /* 15 */
276 fshift[ki][m] += fij;
277 fshift[CENTRAL][m] -= fij;
283 real cubic_bonds(int nbonds,
284 const t_iatom forceatoms[], const t_iparams forceparams[],
285 const rvec x[], rvec f[], rvec fshift[],
286 const t_pbc *pbc, const t_graph *g,
287 real gmx_unused lambda, real gmx_unused *dvdlambda,
288 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
289 int gmx_unused *global_atom_index)
291 const real three = 3.0;
292 const real two = 2.0;
294 real dr, dr2, dist, kdist, kdist2, fbond, vbond, fij, vtot;
296 int i, m, ki, type, ai, aj;
300 for (i = 0; (i < nbonds); )
302 type = forceatoms[i++];
303 ai = forceatoms[i++];
304 aj = forceatoms[i++];
306 b0 = forceparams[type].cubic.b0;
307 kb = forceparams[type].cubic.kb;
308 kcub = forceparams[type].cubic.kcub;
310 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
311 dr2 = iprod(dx, dx); /* 5 */
318 dr = dr2*gmx_invsqrt(dr2); /* 10 */
323 vbond = kdist2 + kcub*kdist2*dist;
324 fbond = -(two*kdist + three*kdist2*kcub)/dr;
326 vtot += vbond; /* 21 */
330 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
333 for (m = 0; (m < DIM); m++) /* 15 */
338 fshift[ki][m] += fij;
339 fshift[CENTRAL][m] -= fij;
345 real FENE_bonds(int nbonds,
346 const t_iatom forceatoms[], const t_iparams forceparams[],
347 const rvec x[], rvec f[], rvec fshift[],
348 const t_pbc *pbc, const t_graph *g,
349 real gmx_unused lambda, real gmx_unused *dvdlambda,
350 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
351 int *global_atom_index)
353 const real half = 0.5;
354 const real one = 1.0;
356 real dr, dr2, bm2, omdr2obm2, fbond, vbond, fij, vtot;
358 int i, m, ki, type, ai, aj;
362 for (i = 0; (i < nbonds); )
364 type = forceatoms[i++];
365 ai = forceatoms[i++];
366 aj = forceatoms[i++];
368 bm = forceparams[type].fene.bm;
369 kb = forceparams[type].fene.kb;
371 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
372 dr2 = iprod(dx, dx); /* 5 */
384 "r^2 (%f) >= bm^2 (%f) in FENE bond between atoms %d and %d",
386 glatnr(global_atom_index, ai),
387 glatnr(global_atom_index, aj));
390 omdr2obm2 = one - dr2/bm2;
392 vbond = -half*kb*bm2*log(omdr2obm2);
393 fbond = -kb/omdr2obm2;
395 vtot += vbond; /* 35 */
399 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
402 for (m = 0; (m < DIM); m++) /* 15 */
407 fshift[ki][m] += fij;
408 fshift[CENTRAL][m] -= fij;
414 real harmonic(real kA, real kB, real xA, real xB, real x, real lambda,
417 const real half = 0.5;
418 real L1, kk, x0, dx, dx2;
419 real v, f, dvdlambda;
422 kk = L1*kA+lambda*kB;
423 x0 = L1*xA+lambda*xB;
430 dvdlambda = half*(kB-kA)*dx2 + (xA-xB)*kk*dx;
437 /* That was 19 flops */
441 real bonds(int nbonds,
442 const t_iatom forceatoms[], const t_iparams forceparams[],
443 const rvec x[], rvec f[], rvec fshift[],
444 const t_pbc *pbc, const t_graph *g,
445 real lambda, real *dvdlambda,
446 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
447 int gmx_unused *global_atom_index)
449 int i, m, ki, ai, aj, type;
450 real dr, dr2, fbond, vbond, fij, vtot;
455 for (i = 0; (i < nbonds); )
457 type = forceatoms[i++];
458 ai = forceatoms[i++];
459 aj = forceatoms[i++];
461 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
462 dr2 = iprod(dx, dx); /* 5 */
463 dr = dr2*gmx_invsqrt(dr2); /* 10 */
465 *dvdlambda += harmonic(forceparams[type].harmonic.krA,
466 forceparams[type].harmonic.krB,
467 forceparams[type].harmonic.rA,
468 forceparams[type].harmonic.rB,
469 dr, lambda, &vbond, &fbond); /* 19 */
477 vtot += vbond; /* 1*/
478 fbond *= gmx_invsqrt(dr2); /* 6 */
482 fprintf(debug, "BONDS: dr = %10g vbond = %10g fbond = %10g\n",
488 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
491 for (m = 0; (m < DIM); m++) /* 15 */
496 fshift[ki][m] += fij;
497 fshift[CENTRAL][m] -= fij;
503 real restraint_bonds(int nbonds,
504 const t_iatom forceatoms[], const t_iparams forceparams[],
505 const rvec x[], rvec f[], rvec fshift[],
506 const t_pbc *pbc, const t_graph *g,
507 real lambda, real *dvdlambda,
508 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
509 int gmx_unused *global_atom_index)
511 int i, m, ki, ai, aj, type;
512 real dr, dr2, fbond, vbond, fij, vtot;
514 real low, dlow, up1, dup1, up2, dup2, k, dk;
522 for (i = 0; (i < nbonds); )
524 type = forceatoms[i++];
525 ai = forceatoms[i++];
526 aj = forceatoms[i++];
528 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
529 dr2 = iprod(dx, dx); /* 5 */
530 dr = dr2*gmx_invsqrt(dr2); /* 10 */
532 low = L1*forceparams[type].restraint.lowA + lambda*forceparams[type].restraint.lowB;
533 dlow = -forceparams[type].restraint.lowA + forceparams[type].restraint.lowB;
534 up1 = L1*forceparams[type].restraint.up1A + lambda*forceparams[type].restraint.up1B;
535 dup1 = -forceparams[type].restraint.up1A + forceparams[type].restraint.up1B;
536 up2 = L1*forceparams[type].restraint.up2A + lambda*forceparams[type].restraint.up2B;
537 dup2 = -forceparams[type].restraint.up2A + forceparams[type].restraint.up2B;
538 k = L1*forceparams[type].restraint.kA + lambda*forceparams[type].restraint.kB;
539 dk = -forceparams[type].restraint.kA + forceparams[type].restraint.kB;
548 *dvdlambda += 0.5*dk*drh2 - k*dlow*drh;
561 *dvdlambda += 0.5*dk*drh2 - k*dup1*drh;
566 vbond = k*(up2 - up1)*(0.5*(up2 - up1) + drh);
567 fbond = -k*(up2 - up1);
568 *dvdlambda += dk*(up2 - up1)*(0.5*(up2 - up1) + drh)
569 + k*(dup2 - dup1)*(up2 - up1 + drh)
570 - k*(up2 - up1)*dup2;
578 vtot += vbond; /* 1*/
579 fbond *= gmx_invsqrt(dr2); /* 6 */
583 fprintf(debug, "BONDS: dr = %10g vbond = %10g fbond = %10g\n",
589 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
592 for (m = 0; (m < DIM); m++) /* 15 */
597 fshift[ki][m] += fij;
598 fshift[CENTRAL][m] -= fij;
605 real polarize(int nbonds,
606 const t_iatom forceatoms[], const t_iparams forceparams[],
607 const rvec x[], rvec f[], rvec fshift[],
608 const t_pbc *pbc, const t_graph *g,
609 real lambda, real *dvdlambda,
610 const t_mdatoms *md, t_fcdata gmx_unused *fcd,
611 int gmx_unused *global_atom_index)
613 int i, m, ki, ai, aj, type;
614 real dr, dr2, fbond, vbond, fij, vtot, ksh;
619 for (i = 0; (i < nbonds); )
621 type = forceatoms[i++];
622 ai = forceatoms[i++];
623 aj = forceatoms[i++];
624 ksh = sqr(md->chargeA[aj])*ONE_4PI_EPS0/forceparams[type].polarize.alpha;
627 fprintf(debug, "POL: local ai = %d aj = %d ksh = %.3f\n", ai, aj, ksh);
630 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
631 dr2 = iprod(dx, dx); /* 5 */
632 dr = dr2*gmx_invsqrt(dr2); /* 10 */
634 *dvdlambda += harmonic(ksh, ksh, 0, 0, dr, lambda, &vbond, &fbond); /* 19 */
641 vtot += vbond; /* 1*/
642 fbond *= gmx_invsqrt(dr2); /* 6 */
646 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
649 for (m = 0; (m < DIM); m++) /* 15 */
654 fshift[ki][m] += fij;
655 fshift[CENTRAL][m] -= fij;
661 real anharm_polarize(int nbonds,
662 const t_iatom forceatoms[], const t_iparams forceparams[],
663 const rvec x[], rvec f[], rvec fshift[],
664 const t_pbc *pbc, const t_graph *g,
665 real lambda, real *dvdlambda,
666 const t_mdatoms *md, t_fcdata gmx_unused *fcd,
667 int gmx_unused *global_atom_index)
669 int i, m, ki, ai, aj, type;
670 real dr, dr2, fbond, vbond, fij, vtot, ksh, khyp, drcut, ddr, ddr3;
675 for (i = 0; (i < nbonds); )
677 type = forceatoms[i++];
678 ai = forceatoms[i++];
679 aj = forceatoms[i++];
680 ksh = sqr(md->chargeA[aj])*ONE_4PI_EPS0/forceparams[type].anharm_polarize.alpha; /* 7*/
681 khyp = forceparams[type].anharm_polarize.khyp;
682 drcut = forceparams[type].anharm_polarize.drcut;
685 fprintf(debug, "POL: local ai = %d aj = %d ksh = %.3f\n", ai, aj, ksh);
688 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
689 dr2 = iprod(dx, dx); /* 5 */
690 dr = dr2*gmx_invsqrt(dr2); /* 10 */
692 *dvdlambda += harmonic(ksh, ksh, 0, 0, dr, lambda, &vbond, &fbond); /* 19 */
703 vbond += khyp*ddr*ddr3;
704 fbond -= 4*khyp*ddr3;
706 fbond *= gmx_invsqrt(dr2); /* 6 */
707 vtot += vbond; /* 1*/
711 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
714 for (m = 0; (m < DIM); m++) /* 15 */
719 fshift[ki][m] += fij;
720 fshift[CENTRAL][m] -= fij;
726 real water_pol(int nbonds,
727 const t_iatom forceatoms[], const t_iparams forceparams[],
728 const rvec x[], rvec f[], rvec gmx_unused fshift[],
729 const t_pbc gmx_unused *pbc, const t_graph gmx_unused *g,
730 real gmx_unused lambda, real gmx_unused *dvdlambda,
731 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
732 int gmx_unused *global_atom_index)
734 /* This routine implements anisotropic polarizibility for water, through
735 * a shell connected to a dummy with spring constant that differ in the
736 * three spatial dimensions in the molecular frame.
738 int i, m, aO, aH1, aH2, aD, aS, type, type0;
739 rvec dOH1, dOH2, dHH, dOD, dDS, nW, kk, dx, kdx, proj;
743 real vtot, fij, r_HH, r_OD, r_nW, tx, ty, tz, qS;
748 type0 = forceatoms[0];
750 qS = md->chargeA[aS];
751 kk[XX] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_x;
752 kk[YY] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_y;
753 kk[ZZ] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_z;
754 r_HH = 1.0/forceparams[type0].wpol.rHH;
755 r_OD = 1.0/forceparams[type0].wpol.rOD;
758 fprintf(debug, "WPOL: qS = %10.5f aS = %5d\n", qS, aS);
759 fprintf(debug, "WPOL: kk = %10.3f %10.3f %10.3f\n",
760 kk[XX], kk[YY], kk[ZZ]);
761 fprintf(debug, "WPOL: rOH = %10.3f rHH = %10.3f rOD = %10.3f\n",
762 forceparams[type0].wpol.rOH,
763 forceparams[type0].wpol.rHH,
764 forceparams[type0].wpol.rOD);
766 for (i = 0; (i < nbonds); i += 6)
768 type = forceatoms[i];
771 gmx_fatal(FARGS, "Sorry, type = %d, type0 = %d, file = %s, line = %d",
772 type, type0, __FILE__, __LINE__);
774 aO = forceatoms[i+1];
775 aH1 = forceatoms[i+2];
776 aH2 = forceatoms[i+3];
777 aD = forceatoms[i+4];
778 aS = forceatoms[i+5];
780 /* Compute vectors describing the water frame */
781 rvec_sub(x[aH1], x[aO], dOH1);
782 rvec_sub(x[aH2], x[aO], dOH2);
783 rvec_sub(x[aH2], x[aH1], dHH);
784 rvec_sub(x[aD], x[aO], dOD);
785 rvec_sub(x[aS], x[aD], dDS);
786 cprod(dOH1, dOH2, nW);
788 /* Compute inverse length of normal vector
789 * (this one could be precomputed, but I'm too lazy now)
791 r_nW = gmx_invsqrt(iprod(nW, nW));
792 /* This is for precision, but does not make a big difference,
795 r_OD = gmx_invsqrt(iprod(dOD, dOD));
797 /* Normalize the vectors in the water frame */
799 svmul(r_HH, dHH, dHH);
800 svmul(r_OD, dOD, dOD);
802 /* Compute displacement of shell along components of the vector */
803 dx[ZZ] = iprod(dDS, dOD);
804 /* Compute projection on the XY plane: dDS - dx[ZZ]*dOD */
805 for (m = 0; (m < DIM); m++)
807 proj[m] = dDS[m]-dx[ZZ]*dOD[m];
810 /*dx[XX] = iprod(dDS,nW);
811 dx[YY] = iprod(dDS,dHH);*/
812 dx[XX] = iprod(proj, nW);
813 for (m = 0; (m < DIM); m++)
815 proj[m] -= dx[XX]*nW[m];
817 dx[YY] = iprod(proj, dHH);
822 fprintf(debug, "WPOL: dx2=%10g dy2=%10g dz2=%10g sum=%10g dDS^2=%10g\n",
823 sqr(dx[XX]), sqr(dx[YY]), sqr(dx[ZZ]), iprod(dx, dx), iprod(dDS, dDS));
824 fprintf(debug, "WPOL: dHH=(%10g,%10g,%10g)\n", dHH[XX], dHH[YY], dHH[ZZ]);
825 fprintf(debug, "WPOL: dOD=(%10g,%10g,%10g), 1/r_OD = %10g\n",
826 dOD[XX], dOD[YY], dOD[ZZ], 1/r_OD);
827 fprintf(debug, "WPOL: nW =(%10g,%10g,%10g), 1/r_nW = %10g\n",
828 nW[XX], nW[YY], nW[ZZ], 1/r_nW);
829 fprintf(debug, "WPOL: dx =%10g, dy =%10g, dz =%10g\n",
830 dx[XX], dx[YY], dx[ZZ]);
831 fprintf(debug, "WPOL: dDSx=%10g, dDSy=%10g, dDSz=%10g\n",
832 dDS[XX], dDS[YY], dDS[ZZ]);
835 /* Now compute the forces and energy */
836 kdx[XX] = kk[XX]*dx[XX];
837 kdx[YY] = kk[YY]*dx[YY];
838 kdx[ZZ] = kk[ZZ]*dx[ZZ];
839 vtot += iprod(dx, kdx);
840 for (m = 0; (m < DIM); m++)
842 /* This is a tensor operation but written out for speed */
856 fprintf(debug, "WPOL: vwpol=%g\n", 0.5*iprod(dx, kdx));
857 fprintf(debug, "WPOL: df = (%10g, %10g, %10g)\n", df[XX], df[YY], df[ZZ]);
865 static real do_1_thole(const rvec xi, const rvec xj, rvec fi, rvec fj,
866 const t_pbc *pbc, real qq,
867 rvec fshift[], real afac)
870 real r12sq, r12_1, r12n, r12bar, v0, v1, fscal, ebar, fff;
873 t = pbc_rvec_sub(pbc, xi, xj, r12); /* 3 */
875 r12sq = iprod(r12, r12); /* 5 */
876 r12_1 = gmx_invsqrt(r12sq); /* 5 */
877 r12bar = afac/r12_1; /* 5 */
878 v0 = qq*ONE_4PI_EPS0*r12_1; /* 2 */
879 ebar = exp(-r12bar); /* 5 */
880 v1 = (1-(1+0.5*r12bar)*ebar); /* 4 */
881 fscal = ((v0*r12_1)*v1 - v0*0.5*afac*ebar*(r12bar+1))*r12_1; /* 9 */
884 fprintf(debug, "THOLE: v0 = %.3f v1 = %.3f r12= % .3f r12bar = %.3f fscal = %.3f ebar = %.3f\n", v0, v1, 1/r12_1, r12bar, fscal, ebar);
887 for (m = 0; (m < DIM); m++)
893 fshift[CENTRAL][m] -= fff;
896 return v0*v1; /* 1 */
900 real thole_pol(int nbonds,
901 const t_iatom forceatoms[], const t_iparams forceparams[],
902 const rvec x[], rvec f[], rvec fshift[],
903 const t_pbc *pbc, const t_graph gmx_unused *g,
904 real gmx_unused lambda, real gmx_unused *dvdlambda,
905 const t_mdatoms *md, t_fcdata gmx_unused *fcd,
906 int gmx_unused *global_atom_index)
908 /* Interaction between two pairs of particles with opposite charge */
909 int i, type, a1, da1, a2, da2;
910 real q1, q2, qq, a, al1, al2, afac;
913 for (i = 0; (i < nbonds); )
915 type = forceatoms[i++];
916 a1 = forceatoms[i++];
917 da1 = forceatoms[i++];
918 a2 = forceatoms[i++];
919 da2 = forceatoms[i++];
920 q1 = md->chargeA[da1];
921 q2 = md->chargeA[da2];
922 a = forceparams[type].thole.a;
923 al1 = forceparams[type].thole.alpha1;
924 al2 = forceparams[type].thole.alpha2;
926 afac = a*pow(al1*al2, -1.0/6.0);
927 V += do_1_thole(x[a1], x[a2], f[a1], f[a2], pbc, qq, fshift, afac);
928 V += do_1_thole(x[da1], x[a2], f[da1], f[a2], pbc, -qq, fshift, afac);
929 V += do_1_thole(x[a1], x[da2], f[a1], f[da2], pbc, -qq, fshift, afac);
930 V += do_1_thole(x[da1], x[da2], f[da1], f[da2], pbc, qq, fshift, afac);
936 real bond_angle(const rvec xi, const rvec xj, const rvec xk, const t_pbc *pbc,
937 rvec r_ij, rvec r_kj, real *costh,
939 /* Return value is the angle between the bonds i-j and j-k */
944 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
945 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
947 *costh = cos_angle(r_ij, r_kj); /* 25 */
948 th = acos(*costh); /* 10 */
953 real angles(int nbonds,
954 const t_iatom forceatoms[], const t_iparams forceparams[],
955 const rvec x[], rvec f[], rvec fshift[],
956 const t_pbc *pbc, const t_graph *g,
957 real lambda, real *dvdlambda,
958 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
959 int gmx_unused *global_atom_index)
961 int i, ai, aj, ak, t1, t2, type;
963 real cos_theta, cos_theta2, theta, dVdt, va, vtot;
964 ivec jt, dt_ij, dt_kj;
967 for (i = 0; i < nbonds; )
969 type = forceatoms[i++];
970 ai = forceatoms[i++];
971 aj = forceatoms[i++];
972 ak = forceatoms[i++];
974 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
975 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
977 *dvdlambda += harmonic(forceparams[type].harmonic.krA,
978 forceparams[type].harmonic.krB,
979 forceparams[type].harmonic.rA*DEG2RAD,
980 forceparams[type].harmonic.rB*DEG2RAD,
981 theta, lambda, &va, &dVdt); /* 21 */
984 cos_theta2 = sqr(cos_theta);
994 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
995 sth = st*cos_theta; /* 1 */
999 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1000 theta*RAD2DEG, va, dVdt);
1003 nrij2 = iprod(r_ij, r_ij); /* 5 */
1004 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1006 nrij_1 = gmx_invsqrt(nrij2); /* 10 */
1007 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1009 cik = st*nrij_1*nrkj_1; /* 2 */
1010 cii = sth*nrij_1*nrij_1; /* 2 */
1011 ckk = sth*nrkj_1*nrkj_1; /* 2 */
1013 for (m = 0; m < DIM; m++)
1015 f_i[m] = -(cik*r_kj[m] - cii*r_ij[m]);
1016 f_k[m] = -(cik*r_ij[m] - ckk*r_kj[m]);
1017 f_j[m] = -f_i[m] - f_k[m];
1024 copy_ivec(SHIFT_IVEC(g, aj), jt);
1026 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1027 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1028 t1 = IVEC2IS(dt_ij);
1029 t2 = IVEC2IS(dt_kj);
1031 rvec_inc(fshift[t1], f_i);
1032 rvec_inc(fshift[CENTRAL], f_j);
1033 rvec_inc(fshift[t2], f_k);
1042 /* As angles, but using SIMD to calculate many dihedrals at once.
1043 * This routines does not calculate energies and shift forces.
1045 static gmx_inline void
1046 angles_noener_simd(int nbonds,
1047 const t_iatom forceatoms[], const t_iparams forceparams[],
1048 const rvec x[], rvec f[],
1049 const t_pbc *pbc, const t_graph gmx_unused *g,
1050 real gmx_unused lambda,
1051 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1052 int gmx_unused *global_atom_index)
1056 int type, ai[GMX_SIMD_REAL_WIDTH], aj[GMX_SIMD_REAL_WIDTH];
1057 int ak[GMX_SIMD_REAL_WIDTH];
1058 real coeff_array[2*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *coeff;
1059 real dr_array[2*DIM*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *dr;
1060 real f_buf_array[6*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *f_buf;
1061 gmx_simd_real_t k_S, theta0_S;
1062 gmx_simd_real_t rijx_S, rijy_S, rijz_S;
1063 gmx_simd_real_t rkjx_S, rkjy_S, rkjz_S;
1064 gmx_simd_real_t one_S;
1065 gmx_simd_real_t min_one_plus_eps_S;
1066 gmx_simd_real_t rij_rkj_S;
1067 gmx_simd_real_t nrij2_S, nrij_1_S;
1068 gmx_simd_real_t nrkj2_S, nrkj_1_S;
1069 gmx_simd_real_t cos_S, invsin_S;
1070 gmx_simd_real_t theta_S;
1071 gmx_simd_real_t st_S, sth_S;
1072 gmx_simd_real_t cik_S, cii_S, ckk_S;
1073 gmx_simd_real_t f_ix_S, f_iy_S, f_iz_S;
1074 gmx_simd_real_t 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_r(coeff_array);
1079 dr = gmx_simd_align_r(dr_array);
1080 f_buf = gmx_simd_align_r(f_buf_array);
1082 set_pbc_simd(pbc, &pbc_simd);
1084 one_S = gmx_simd_set1_r(1.0);
1086 /* The smallest number > -1 */
1087 min_one_plus_eps_S = gmx_simd_set1_r(-1.0 + 2*GMX_REAL_EPS);
1089 /* nbonds is the number of angles times nfa1, here we step GMX_SIMD_REAL_WIDTH angles */
1090 for (i = 0; (i < nbonds); i += GMX_SIMD_REAL_WIDTH*nfa1)
1092 /* Collect atoms for GMX_SIMD_REAL_WIDTH angles.
1093 * iu indexes into forceatoms, we should not let iu go beyond nbonds.
1096 for (s = 0; s < GMX_SIMD_REAL_WIDTH; s++)
1098 type = forceatoms[iu];
1099 ai[s] = forceatoms[iu+1];
1100 aj[s] = forceatoms[iu+2];
1101 ak[s] = forceatoms[iu+3];
1103 coeff[s] = forceparams[type].harmonic.krA;
1104 coeff[GMX_SIMD_REAL_WIDTH+s] = forceparams[type].harmonic.rA*DEG2RAD;
1106 /* If you can't use pbc_dx_simd below for PBC, e.g. because
1107 * you can't round in SIMD, use pbc_rvec_sub here.
1109 /* Store the non PBC corrected distances packed and aligned */
1110 for (m = 0; m < DIM; m++)
1112 dr[s + m *GMX_SIMD_REAL_WIDTH] = x[ai[s]][m] - x[aj[s]][m];
1113 dr[s + (DIM+m)*GMX_SIMD_REAL_WIDTH] = x[ak[s]][m] - x[aj[s]][m];
1116 /* At the end fill the arrays with identical entries */
1117 if (iu + nfa1 < nbonds)
1123 k_S = gmx_simd_load_r(coeff);
1124 theta0_S = gmx_simd_load_r(coeff+GMX_SIMD_REAL_WIDTH);
1126 rijx_S = gmx_simd_load_r(dr + 0*GMX_SIMD_REAL_WIDTH);
1127 rijy_S = gmx_simd_load_r(dr + 1*GMX_SIMD_REAL_WIDTH);
1128 rijz_S = gmx_simd_load_r(dr + 2*GMX_SIMD_REAL_WIDTH);
1129 rkjx_S = gmx_simd_load_r(dr + 3*GMX_SIMD_REAL_WIDTH);
1130 rkjy_S = gmx_simd_load_r(dr + 4*GMX_SIMD_REAL_WIDTH);
1131 rkjz_S = gmx_simd_load_r(dr + 5*GMX_SIMD_REAL_WIDTH);
1133 pbc_dx_simd(&rijx_S, &rijy_S, &rijz_S, &pbc_simd);
1134 pbc_dx_simd(&rkjx_S, &rkjy_S, &rkjz_S, &pbc_simd);
1136 rij_rkj_S = gmx_simd_iprod_r(rijx_S, rijy_S, rijz_S,
1137 rkjx_S, rkjy_S, rkjz_S);
1139 nrij2_S = gmx_simd_norm2_r(rijx_S, rijy_S, rijz_S);
1140 nrkj2_S = gmx_simd_norm2_r(rkjx_S, rkjy_S, rkjz_S);
1142 nrij_1_S = gmx_simd_invsqrt_r(nrij2_S);
1143 nrkj_1_S = gmx_simd_invsqrt_r(nrkj2_S);
1145 cos_S = gmx_simd_mul_r(rij_rkj_S, gmx_simd_mul_r(nrij_1_S, nrkj_1_S));
1147 /* To allow for 180 degrees, we take the max of cos and -1 + 1bit,
1148 * so we can safely get the 1/sin from 1/sqrt(1 - cos^2).
1149 * This also ensures that rounding errors would cause the argument
1150 * of gmx_simd_acos_r to be < -1.
1151 * Note that we do not take precautions for cos(0)=1, so the outer
1152 * atoms in an angle should not be on top of each other.
1154 cos_S = gmx_simd_max_r(cos_S, min_one_plus_eps_S);
1156 theta_S = gmx_simd_acos_r(cos_S);
1158 invsin_S = gmx_simd_invsqrt_r(gmx_simd_sub_r(one_S, gmx_simd_mul_r(cos_S, cos_S)));
1160 st_S = gmx_simd_mul_r(gmx_simd_mul_r(k_S, gmx_simd_sub_r(theta0_S, theta_S)),
1162 sth_S = gmx_simd_mul_r(st_S, cos_S);
1164 cik_S = gmx_simd_mul_r(st_S, gmx_simd_mul_r(nrij_1_S, nrkj_1_S));
1165 cii_S = gmx_simd_mul_r(sth_S, gmx_simd_mul_r(nrij_1_S, nrij_1_S));
1166 ckk_S = gmx_simd_mul_r(sth_S, gmx_simd_mul_r(nrkj_1_S, nrkj_1_S));
1168 f_ix_S = gmx_simd_mul_r(cii_S, rijx_S);
1169 f_ix_S = gmx_simd_fnmadd_r(cik_S, rkjx_S, f_ix_S);
1170 f_iy_S = gmx_simd_mul_r(cii_S, rijy_S);
1171 f_iy_S = gmx_simd_fnmadd_r(cik_S, rkjy_S, f_iy_S);
1172 f_iz_S = gmx_simd_mul_r(cii_S, rijz_S);
1173 f_iz_S = gmx_simd_fnmadd_r(cik_S, rkjz_S, f_iz_S);
1174 f_kx_S = gmx_simd_mul_r(ckk_S, rkjx_S);
1175 f_kx_S = gmx_simd_fnmadd_r(cik_S, rijx_S, f_kx_S);
1176 f_ky_S = gmx_simd_mul_r(ckk_S, rkjy_S);
1177 f_ky_S = gmx_simd_fnmadd_r(cik_S, rijy_S, f_ky_S);
1178 f_kz_S = gmx_simd_mul_r(ckk_S, rkjz_S);
1179 f_kz_S = gmx_simd_fnmadd_r(cik_S, rijz_S, f_kz_S);
1181 gmx_simd_store_r(f_buf + 0*GMX_SIMD_REAL_WIDTH, f_ix_S);
1182 gmx_simd_store_r(f_buf + 1*GMX_SIMD_REAL_WIDTH, f_iy_S);
1183 gmx_simd_store_r(f_buf + 2*GMX_SIMD_REAL_WIDTH, f_iz_S);
1184 gmx_simd_store_r(f_buf + 3*GMX_SIMD_REAL_WIDTH, f_kx_S);
1185 gmx_simd_store_r(f_buf + 4*GMX_SIMD_REAL_WIDTH, f_ky_S);
1186 gmx_simd_store_r(f_buf + 5*GMX_SIMD_REAL_WIDTH, f_kz_S);
1192 for (m = 0; m < DIM; m++)
1194 f[ai[s]][m] += f_buf[s + m*GMX_SIMD_REAL_WIDTH];
1195 f[aj[s]][m] -= f_buf[s + m*GMX_SIMD_REAL_WIDTH] + f_buf[s + (DIM+m)*GMX_SIMD_REAL_WIDTH];
1196 f[ak[s]][m] += f_buf[s + (DIM+m)*GMX_SIMD_REAL_WIDTH];
1201 while (s < GMX_SIMD_REAL_WIDTH && iu < nbonds);
1205 #endif /* SIMD_BONDEDS */
1207 real linear_angles(int nbonds,
1208 const t_iatom forceatoms[], const t_iparams forceparams[],
1209 const rvec x[], rvec f[], rvec fshift[],
1210 const t_pbc *pbc, const t_graph *g,
1211 real lambda, real *dvdlambda,
1212 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1213 int gmx_unused *global_atom_index)
1215 int i, m, ai, aj, ak, t1, t2, type;
1217 real L1, kA, kB, aA, aB, dr, dr2, va, vtot, a, b, klin;
1218 ivec jt, dt_ij, dt_kj;
1219 rvec r_ij, r_kj, r_ik, dx;
1223 for (i = 0; (i < nbonds); )
1225 type = forceatoms[i++];
1226 ai = forceatoms[i++];
1227 aj = forceatoms[i++];
1228 ak = forceatoms[i++];
1230 kA = forceparams[type].linangle.klinA;
1231 kB = forceparams[type].linangle.klinB;
1232 klin = L1*kA + lambda*kB;
1234 aA = forceparams[type].linangle.aA;
1235 aB = forceparams[type].linangle.aB;
1236 a = L1*aA+lambda*aB;
1239 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
1240 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
1241 rvec_sub(r_ij, r_kj, r_ik);
1244 for (m = 0; (m < DIM); m++)
1246 dr = -a * r_ij[m] - b * r_kj[m];
1251 f_j[m] = -(f_i[m]+f_k[m]);
1257 *dvdlambda += 0.5*(kB-kA)*dr2 + klin*(aB-aA)*iprod(dx, r_ik);
1263 copy_ivec(SHIFT_IVEC(g, aj), jt);
1265 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1266 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1267 t1 = IVEC2IS(dt_ij);
1268 t2 = IVEC2IS(dt_kj);
1270 rvec_inc(fshift[t1], f_i);
1271 rvec_inc(fshift[CENTRAL], f_j);
1272 rvec_inc(fshift[t2], f_k);
1277 real urey_bradley(int nbonds,
1278 const t_iatom forceatoms[], const t_iparams forceparams[],
1279 const rvec x[], rvec f[], rvec fshift[],
1280 const t_pbc *pbc, const t_graph *g,
1281 real lambda, real *dvdlambda,
1282 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1283 int gmx_unused *global_atom_index)
1285 int i, m, ai, aj, ak, t1, t2, type, ki;
1286 rvec r_ij, r_kj, r_ik;
1287 real cos_theta, cos_theta2, theta;
1288 real dVdt, va, vtot, dr, dr2, vbond, fbond, fik;
1289 real kthA, th0A, kUBA, r13A, kthB, th0B, kUBB, r13B;
1290 ivec jt, dt_ij, dt_kj, dt_ik;
1293 for (i = 0; (i < nbonds); )
1295 type = forceatoms[i++];
1296 ai = forceatoms[i++];
1297 aj = forceatoms[i++];
1298 ak = forceatoms[i++];
1299 th0A = forceparams[type].u_b.thetaA*DEG2RAD;
1300 kthA = forceparams[type].u_b.kthetaA;
1301 r13A = forceparams[type].u_b.r13A;
1302 kUBA = forceparams[type].u_b.kUBA;
1303 th0B = forceparams[type].u_b.thetaB*DEG2RAD;
1304 kthB = forceparams[type].u_b.kthetaB;
1305 r13B = forceparams[type].u_b.r13B;
1306 kUBB = forceparams[type].u_b.kUBB;
1308 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
1309 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
1311 *dvdlambda += harmonic(kthA, kthB, th0A, th0B, theta, lambda, &va, &dVdt); /* 21 */
1314 ki = pbc_rvec_sub(pbc, x[ai], x[ak], r_ik); /* 3 */
1315 dr2 = iprod(r_ik, r_ik); /* 5 */
1316 dr = dr2*gmx_invsqrt(dr2); /* 10 */
1318 *dvdlambda += harmonic(kUBA, kUBB, r13A, r13B, dr, lambda, &vbond, &fbond); /* 19 */
1320 cos_theta2 = sqr(cos_theta); /* 1 */
1328 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1329 sth = st*cos_theta; /* 1 */
1333 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1334 theta*RAD2DEG, va, dVdt);
1337 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1338 nrij2 = iprod(r_ij, r_ij);
1340 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
1341 cii = sth/nrij2; /* 10 */
1342 ckk = sth/nrkj2; /* 10 */
1344 for (m = 0; (m < DIM); m++) /* 39 */
1346 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
1347 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
1348 f_j[m] = -f_i[m]-f_k[m];
1355 copy_ivec(SHIFT_IVEC(g, aj), jt);
1357 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1358 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1359 t1 = IVEC2IS(dt_ij);
1360 t2 = IVEC2IS(dt_kj);
1362 rvec_inc(fshift[t1], f_i);
1363 rvec_inc(fshift[CENTRAL], f_j);
1364 rvec_inc(fshift[t2], f_k);
1366 /* Time for the bond calculations */
1372 vtot += vbond; /* 1*/
1373 fbond *= gmx_invsqrt(dr2); /* 6 */
1377 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ak), dt_ik);
1378 ki = IVEC2IS(dt_ik);
1380 for (m = 0; (m < DIM); m++) /* 15 */
1382 fik = fbond*r_ik[m];
1385 fshift[ki][m] += fik;
1386 fshift[CENTRAL][m] -= fik;
1392 real quartic_angles(int nbonds,
1393 const t_iatom forceatoms[], const t_iparams forceparams[],
1394 const rvec x[], rvec f[], rvec fshift[],
1395 const t_pbc *pbc, const t_graph *g,
1396 real gmx_unused lambda, real gmx_unused *dvdlambda,
1397 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1398 int gmx_unused *global_atom_index)
1400 int i, j, ai, aj, ak, t1, t2, type;
1402 real cos_theta, cos_theta2, theta, dt, dVdt, va, dtp, c, vtot;
1403 ivec jt, dt_ij, dt_kj;
1406 for (i = 0; (i < nbonds); )
1408 type = forceatoms[i++];
1409 ai = forceatoms[i++];
1410 aj = forceatoms[i++];
1411 ak = forceatoms[i++];
1413 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
1414 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
1416 dt = theta - forceparams[type].qangle.theta*DEG2RAD; /* 2 */
1419 va = forceparams[type].qangle.c[0];
1421 for (j = 1; j <= 4; j++)
1423 c = forceparams[type].qangle.c[j];
1432 cos_theta2 = sqr(cos_theta); /* 1 */
1441 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1442 sth = st*cos_theta; /* 1 */
1446 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1447 theta*RAD2DEG, va, dVdt);
1450 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1451 nrij2 = iprod(r_ij, r_ij);
1453 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
1454 cii = sth/nrij2; /* 10 */
1455 ckk = sth/nrkj2; /* 10 */
1457 for (m = 0; (m < DIM); m++) /* 39 */
1459 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
1460 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
1461 f_j[m] = -f_i[m]-f_k[m];
1468 copy_ivec(SHIFT_IVEC(g, aj), jt);
1470 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1471 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1472 t1 = IVEC2IS(dt_ij);
1473 t2 = IVEC2IS(dt_kj);
1475 rvec_inc(fshift[t1], f_i);
1476 rvec_inc(fshift[CENTRAL], f_j);
1477 rvec_inc(fshift[t2], f_k);
1483 real dih_angle(const rvec xi, const rvec xj, const rvec xk, const rvec xl,
1485 rvec r_ij, rvec r_kj, rvec r_kl, rvec m, rvec n,
1486 real *sign, int *t1, int *t2, int *t3)
1490 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
1491 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
1492 *t3 = pbc_rvec_sub(pbc, xk, xl, r_kl); /* 3 */
1494 cprod(r_ij, r_kj, m); /* 9 */
1495 cprod(r_kj, r_kl, n); /* 9 */
1496 phi = gmx_angle(m, n); /* 49 (assuming 25 for atan2) */
1497 ipr = iprod(r_ij, n); /* 5 */
1498 (*sign) = (ipr < 0.0) ? -1.0 : 1.0;
1499 phi = (*sign)*phi; /* 1 */
1507 /* As dih_angle above, but calculates 4 dihedral angles at once using SIMD,
1508 * also calculates the pre-factor required for the dihedral force update.
1509 * Note that bv and buf should be register aligned.
1511 static gmx_inline void
1512 dih_angle_simd(const rvec *x,
1513 const int *ai, const int *aj, const int *ak, const int *al,
1514 const pbc_simd_t *pbc,
1516 gmx_simd_real_t *phi_S,
1517 gmx_simd_real_t *mx_S, gmx_simd_real_t *my_S, gmx_simd_real_t *mz_S,
1518 gmx_simd_real_t *nx_S, gmx_simd_real_t *ny_S, gmx_simd_real_t *nz_S,
1519 gmx_simd_real_t *nrkj_m2_S,
1520 gmx_simd_real_t *nrkj_n2_S,
1525 gmx_simd_real_t rijx_S, rijy_S, rijz_S;
1526 gmx_simd_real_t rkjx_S, rkjy_S, rkjz_S;
1527 gmx_simd_real_t rklx_S, rkly_S, rklz_S;
1528 gmx_simd_real_t cx_S, cy_S, cz_S;
1529 gmx_simd_real_t cn_S;
1530 gmx_simd_real_t s_S;
1531 gmx_simd_real_t ipr_S;
1532 gmx_simd_real_t iprm_S, iprn_S;
1533 gmx_simd_real_t nrkj2_S, nrkj_1_S, nrkj_2_S, nrkj_S;
1534 gmx_simd_real_t toler_S;
1535 gmx_simd_real_t p_S, q_S;
1536 gmx_simd_real_t nrkj2_min_S;
1537 gmx_simd_real_t real_eps_S;
1539 /* Used to avoid division by zero.
1540 * We take into acount that we multiply the result by real_eps_S.
1542 nrkj2_min_S = gmx_simd_set1_r(GMX_REAL_MIN/(2*GMX_REAL_EPS));
1544 /* The value of the last significant bit (GMX_REAL_EPS is half of that) */
1545 real_eps_S = gmx_simd_set1_r(2*GMX_REAL_EPS);
1547 for (s = 0; s < GMX_SIMD_REAL_WIDTH; s++)
1549 /* If you can't use pbc_dx_simd below for PBC, e.g. because
1550 * you can't round in SIMD, use pbc_rvec_sub here.
1552 for (m = 0; m < DIM; m++)
1554 dr[s + (0*DIM + m)*GMX_SIMD_REAL_WIDTH] = x[ai[s]][m] - x[aj[s]][m];
1555 dr[s + (1*DIM + m)*GMX_SIMD_REAL_WIDTH] = x[ak[s]][m] - x[aj[s]][m];
1556 dr[s + (2*DIM + m)*GMX_SIMD_REAL_WIDTH] = x[ak[s]][m] - x[al[s]][m];
1560 rijx_S = gmx_simd_load_r(dr + 0*GMX_SIMD_REAL_WIDTH);
1561 rijy_S = gmx_simd_load_r(dr + 1*GMX_SIMD_REAL_WIDTH);
1562 rijz_S = gmx_simd_load_r(dr + 2*GMX_SIMD_REAL_WIDTH);
1563 rkjx_S = gmx_simd_load_r(dr + 3*GMX_SIMD_REAL_WIDTH);
1564 rkjy_S = gmx_simd_load_r(dr + 4*GMX_SIMD_REAL_WIDTH);
1565 rkjz_S = gmx_simd_load_r(dr + 5*GMX_SIMD_REAL_WIDTH);
1566 rklx_S = gmx_simd_load_r(dr + 6*GMX_SIMD_REAL_WIDTH);
1567 rkly_S = gmx_simd_load_r(dr + 7*GMX_SIMD_REAL_WIDTH);
1568 rklz_S = gmx_simd_load_r(dr + 8*GMX_SIMD_REAL_WIDTH);
1570 pbc_dx_simd(&rijx_S, &rijy_S, &rijz_S, pbc);
1571 pbc_dx_simd(&rkjx_S, &rkjy_S, &rkjz_S, pbc);
1572 pbc_dx_simd(&rklx_S, &rkly_S, &rklz_S, pbc);
1574 gmx_simd_cprod_r(rijx_S, rijy_S, rijz_S,
1575 rkjx_S, rkjy_S, rkjz_S,
1578 gmx_simd_cprod_r(rkjx_S, rkjy_S, rkjz_S,
1579 rklx_S, rkly_S, rklz_S,
1582 gmx_simd_cprod_r(*mx_S, *my_S, *mz_S,
1583 *nx_S, *ny_S, *nz_S,
1584 &cx_S, &cy_S, &cz_S);
1586 cn_S = gmx_simd_sqrt_r(gmx_simd_norm2_r(cx_S, cy_S, cz_S));
1588 s_S = gmx_simd_iprod_r(*mx_S, *my_S, *mz_S, *nx_S, *ny_S, *nz_S);
1590 /* Determine the dihedral angle, the sign might need correction */
1591 *phi_S = gmx_simd_atan2_r(cn_S, s_S);
1593 ipr_S = gmx_simd_iprod_r(rijx_S, rijy_S, rijz_S,
1594 *nx_S, *ny_S, *nz_S);
1596 iprm_S = gmx_simd_norm2_r(*mx_S, *my_S, *mz_S);
1597 iprn_S = gmx_simd_norm2_r(*nx_S, *ny_S, *nz_S);
1599 nrkj2_S = gmx_simd_norm2_r(rkjx_S, rkjy_S, rkjz_S);
1601 /* Avoid division by zero. When zero, the result is multiplied by 0
1602 * anyhow, so the 3 max below do not affect the final result.
1604 nrkj2_S = gmx_simd_max_r(nrkj2_S, nrkj2_min_S);
1605 nrkj_1_S = gmx_simd_invsqrt_r(nrkj2_S);
1606 nrkj_2_S = gmx_simd_mul_r(nrkj_1_S, nrkj_1_S);
1607 nrkj_S = gmx_simd_mul_r(nrkj2_S, nrkj_1_S);
1609 toler_S = gmx_simd_mul_r(nrkj2_S, real_eps_S);
1611 /* Here the plain-C code uses a conditional, but we can't do that in SIMD.
1612 * So we take a max with the tolerance instead. Since we multiply with
1613 * m or n later, the max does not affect the results.
1615 iprm_S = gmx_simd_max_r(iprm_S, toler_S);
1616 iprn_S = gmx_simd_max_r(iprn_S, toler_S);
1617 *nrkj_m2_S = gmx_simd_mul_r(nrkj_S, gmx_simd_inv_r(iprm_S));
1618 *nrkj_n2_S = gmx_simd_mul_r(nrkj_S, gmx_simd_inv_r(iprn_S));
1620 /* Set sign of phi_S with the sign of ipr_S; phi_S is currently positive */
1621 *phi_S = gmx_cpsgn_nonneg_pr(ipr_S, *phi_S);
1623 p_S = gmx_simd_iprod_r(rijx_S, rijy_S, rijz_S,
1624 rkjx_S, rkjy_S, rkjz_S);
1625 p_S = gmx_simd_mul_r(p_S, nrkj_2_S);
1627 q_S = gmx_simd_iprod_r(rklx_S, rkly_S, rklz_S,
1628 rkjx_S, rkjy_S, rkjz_S);
1629 q_S = gmx_simd_mul_r(q_S, nrkj_2_S);
1631 gmx_simd_store_r(p, p_S);
1632 gmx_simd_store_r(q, q_S);
1635 #endif /* SIMD_BONDEDS */
1638 void do_dih_fup(int i, int j, int k, int l, real ddphi,
1639 rvec r_ij, rvec r_kj, rvec r_kl,
1640 rvec m, rvec n, rvec f[], rvec fshift[],
1641 const t_pbc *pbc, const t_graph *g,
1642 const rvec x[], int t1, int t2, int t3)
1645 rvec f_i, f_j, f_k, f_l;
1646 rvec uvec, vvec, svec, dx_jl;
1647 real iprm, iprn, nrkj, nrkj2, nrkj_1, nrkj_2;
1648 real a, b, p, q, toler;
1649 ivec jt, dt_ij, dt_kj, dt_lj;
1651 iprm = iprod(m, m); /* 5 */
1652 iprn = iprod(n, n); /* 5 */
1653 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1654 toler = nrkj2*GMX_REAL_EPS;
1655 if ((iprm > toler) && (iprn > toler))
1657 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1658 nrkj_2 = nrkj_1*nrkj_1; /* 1 */
1659 nrkj = nrkj2*nrkj_1; /* 1 */
1660 a = -ddphi*nrkj/iprm; /* 11 */
1661 svmul(a, m, f_i); /* 3 */
1662 b = ddphi*nrkj/iprn; /* 11 */
1663 svmul(b, n, f_l); /* 3 */
1664 p = iprod(r_ij, r_kj); /* 5 */
1665 p *= nrkj_2; /* 1 */
1666 q = iprod(r_kl, r_kj); /* 5 */
1667 q *= nrkj_2; /* 1 */
1668 svmul(p, f_i, uvec); /* 3 */
1669 svmul(q, f_l, vvec); /* 3 */
1670 rvec_sub(uvec, vvec, svec); /* 3 */
1671 rvec_sub(f_i, svec, f_j); /* 3 */
1672 rvec_add(f_l, svec, f_k); /* 3 */
1673 rvec_inc(f[i], f_i); /* 3 */
1674 rvec_dec(f[j], f_j); /* 3 */
1675 rvec_dec(f[k], f_k); /* 3 */
1676 rvec_inc(f[l], f_l); /* 3 */
1680 copy_ivec(SHIFT_IVEC(g, j), jt);
1681 ivec_sub(SHIFT_IVEC(g, i), jt, dt_ij);
1682 ivec_sub(SHIFT_IVEC(g, k), jt, dt_kj);
1683 ivec_sub(SHIFT_IVEC(g, l), jt, dt_lj);
1684 t1 = IVEC2IS(dt_ij);
1685 t2 = IVEC2IS(dt_kj);
1686 t3 = IVEC2IS(dt_lj);
1690 t3 = pbc_rvec_sub(pbc, x[l], x[j], dx_jl);
1697 rvec_inc(fshift[t1], f_i);
1698 rvec_dec(fshift[CENTRAL], f_j);
1699 rvec_dec(fshift[t2], f_k);
1700 rvec_inc(fshift[t3], f_l);
1705 /* As do_dih_fup above, but without shift forces */
1707 do_dih_fup_noshiftf(int i, int j, int k, int l, real ddphi,
1708 rvec r_ij, rvec r_kj, rvec r_kl,
1709 rvec m, rvec n, rvec f[])
1711 rvec f_i, f_j, f_k, f_l;
1712 rvec uvec, vvec, svec, dx_jl;
1713 real iprm, iprn, nrkj, nrkj2, nrkj_1, nrkj_2;
1714 real a, b, p, q, toler;
1715 ivec jt, dt_ij, dt_kj, dt_lj;
1717 iprm = iprod(m, m); /* 5 */
1718 iprn = iprod(n, n); /* 5 */
1719 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1720 toler = nrkj2*GMX_REAL_EPS;
1721 if ((iprm > toler) && (iprn > toler))
1723 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1724 nrkj_2 = nrkj_1*nrkj_1; /* 1 */
1725 nrkj = nrkj2*nrkj_1; /* 1 */
1726 a = -ddphi*nrkj/iprm; /* 11 */
1727 svmul(a, m, f_i); /* 3 */
1728 b = ddphi*nrkj/iprn; /* 11 */
1729 svmul(b, n, f_l); /* 3 */
1730 p = iprod(r_ij, r_kj); /* 5 */
1731 p *= nrkj_2; /* 1 */
1732 q = iprod(r_kl, r_kj); /* 5 */
1733 q *= nrkj_2; /* 1 */
1734 svmul(p, f_i, uvec); /* 3 */
1735 svmul(q, f_l, vvec); /* 3 */
1736 rvec_sub(uvec, vvec, svec); /* 3 */
1737 rvec_sub(f_i, svec, f_j); /* 3 */
1738 rvec_add(f_l, svec, f_k); /* 3 */
1739 rvec_inc(f[i], f_i); /* 3 */
1740 rvec_dec(f[j], f_j); /* 3 */
1741 rvec_dec(f[k], f_k); /* 3 */
1742 rvec_inc(f[l], f_l); /* 3 */
1746 /* As do_dih_fup_noshiftf above, but with pre-calculated pre-factors */
1747 static gmx_inline void
1748 do_dih_fup_noshiftf_precalc(int i, int j, int k, int l,
1750 real f_i_x, real f_i_y, real f_i_z,
1751 real mf_l_x, real mf_l_y, real mf_l_z,
1754 rvec f_i, f_j, f_k, f_l;
1755 rvec uvec, vvec, svec;
1763 svmul(p, f_i, uvec);
1764 svmul(q, f_l, vvec);
1765 rvec_sub(uvec, vvec, svec);
1766 rvec_sub(f_i, svec, f_j);
1767 rvec_add(f_l, svec, f_k);
1768 rvec_inc(f[i], f_i);
1769 rvec_dec(f[j], f_j);
1770 rvec_dec(f[k], f_k);
1771 rvec_inc(f[l], f_l);
1775 real dopdihs(real cpA, real cpB, real phiA, real phiB, int mult,
1776 real phi, real lambda, real *V, real *F)
1778 real v, dvdlambda, mdphi, v1, sdphi, ddphi;
1779 real L1 = 1.0 - lambda;
1780 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1781 real dph0 = (phiB - phiA)*DEG2RAD;
1782 real cp = L1*cpA + lambda*cpB;
1784 mdphi = mult*phi - ph0;
1786 ddphi = -cp*mult*sdphi;
1787 v1 = 1.0 + cos(mdphi);
1790 dvdlambda = (cpB - cpA)*v1 + cp*dph0*sdphi;
1797 /* That was 40 flops */
1801 dopdihs_noener(real cpA, real cpB, real phiA, real phiB, int mult,
1802 real phi, real lambda, real *F)
1804 real mdphi, sdphi, ddphi;
1805 real L1 = 1.0 - lambda;
1806 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1807 real cp = L1*cpA + lambda*cpB;
1809 mdphi = mult*phi - ph0;
1811 ddphi = -cp*mult*sdphi;
1815 /* That was 20 flops */
1819 dopdihs_mdphi(real cpA, real cpB, real phiA, real phiB, int mult,
1820 real phi, real lambda, real *cp, real *mdphi)
1822 real L1 = 1.0 - lambda;
1823 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1825 *cp = L1*cpA + lambda*cpB;
1827 *mdphi = mult*phi - ph0;
1830 static real dopdihs_min(real cpA, real cpB, real phiA, real phiB, int mult,
1831 real phi, real lambda, real *V, real *F)
1832 /* similar to dopdihs, except for a minus sign *
1833 * and a different treatment of mult/phi0 */
1835 real v, dvdlambda, mdphi, v1, sdphi, ddphi;
1836 real L1 = 1.0 - lambda;
1837 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1838 real dph0 = (phiB - phiA)*DEG2RAD;
1839 real cp = L1*cpA + lambda*cpB;
1841 mdphi = mult*(phi-ph0);
1843 ddphi = cp*mult*sdphi;
1844 v1 = 1.0-cos(mdphi);
1847 dvdlambda = (cpB-cpA)*v1 + cp*dph0*sdphi;
1854 /* That was 40 flops */
1857 real pdihs(int nbonds,
1858 const t_iatom forceatoms[], const t_iparams forceparams[],
1859 const rvec x[], rvec f[], rvec fshift[],
1860 const t_pbc *pbc, const t_graph *g,
1861 real lambda, real *dvdlambda,
1862 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1863 int gmx_unused *global_atom_index)
1865 int i, type, ai, aj, ak, al;
1867 rvec r_ij, r_kj, r_kl, m, n;
1868 real phi, sign, ddphi, vpd, vtot;
1872 for (i = 0; (i < nbonds); )
1874 type = forceatoms[i++];
1875 ai = forceatoms[i++];
1876 aj = forceatoms[i++];
1877 ak = forceatoms[i++];
1878 al = forceatoms[i++];
1880 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
1881 &sign, &t1, &t2, &t3); /* 84 */
1882 *dvdlambda += dopdihs(forceparams[type].pdihs.cpA,
1883 forceparams[type].pdihs.cpB,
1884 forceparams[type].pdihs.phiA,
1885 forceparams[type].pdihs.phiB,
1886 forceparams[type].pdihs.mult,
1887 phi, lambda, &vpd, &ddphi);
1890 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
1891 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
1894 fprintf(debug, "pdih: (%d,%d,%d,%d) phi=%g\n",
1895 ai, aj, ak, al, phi);
1902 void make_dp_periodic(real *dp) /* 1 flop? */
1904 /* dp cannot be outside (-pi,pi) */
1909 else if (*dp < -M_PI)
1916 /* As pdihs above, but without calculating energies and shift forces */
1918 pdihs_noener(int nbonds,
1919 const t_iatom forceatoms[], const t_iparams forceparams[],
1920 const rvec x[], rvec f[],
1921 const t_pbc gmx_unused *pbc, const t_graph gmx_unused *g,
1923 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1924 int gmx_unused *global_atom_index)
1926 int i, type, ai, aj, ak, al;
1928 rvec r_ij, r_kj, r_kl, m, n;
1929 real phi, sign, ddphi_tot, ddphi;
1931 for (i = 0; (i < nbonds); )
1933 ai = forceatoms[i+1];
1934 aj = forceatoms[i+2];
1935 ak = forceatoms[i+3];
1936 al = forceatoms[i+4];
1938 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
1939 &sign, &t1, &t2, &t3);
1943 /* Loop over dihedrals working on the same atoms,
1944 * so we avoid recalculating angles and force distributions.
1948 type = forceatoms[i];
1949 dopdihs_noener(forceparams[type].pdihs.cpA,
1950 forceparams[type].pdihs.cpB,
1951 forceparams[type].pdihs.phiA,
1952 forceparams[type].pdihs.phiB,
1953 forceparams[type].pdihs.mult,
1954 phi, lambda, &ddphi);
1959 while (i < nbonds &&
1960 forceatoms[i+1] == ai &&
1961 forceatoms[i+2] == aj &&
1962 forceatoms[i+3] == ak &&
1963 forceatoms[i+4] == al);
1965 do_dih_fup_noshiftf(ai, aj, ak, al, ddphi_tot, r_ij, r_kj, r_kl, m, n, f);
1972 /* As pdihs_noner above, but using SIMD to calculate many dihedrals at once */
1974 pdihs_noener_simd(int nbonds,
1975 const t_iatom forceatoms[], const t_iparams forceparams[],
1976 const rvec x[], rvec f[],
1977 const t_pbc *pbc, const t_graph gmx_unused *g,
1978 real gmx_unused lambda,
1979 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1980 int gmx_unused *global_atom_index)
1984 int type, ai[GMX_SIMD_REAL_WIDTH], aj[GMX_SIMD_REAL_WIDTH], ak[GMX_SIMD_REAL_WIDTH], al[GMX_SIMD_REAL_WIDTH];
1985 int t1[GMX_SIMD_REAL_WIDTH], t2[GMX_SIMD_REAL_WIDTH], t3[GMX_SIMD_REAL_WIDTH];
1987 real dr_array[3*DIM*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *dr;
1988 real buf_array[7*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *buf;
1989 real *cp, *phi0, *mult, *phi, *p, *q, *sf_i, *msf_l;
1990 gmx_simd_real_t phi0_S, phi_S;
1991 gmx_simd_real_t mx_S, my_S, mz_S;
1992 gmx_simd_real_t nx_S, ny_S, nz_S;
1993 gmx_simd_real_t nrkj_m2_S, nrkj_n2_S;
1994 gmx_simd_real_t cp_S, mdphi_S, mult_S;
1995 gmx_simd_real_t sin_S, cos_S;
1996 gmx_simd_real_t mddphi_S;
1997 gmx_simd_real_t sf_i_S, msf_l_S;
1998 pbc_simd_t pbc_simd;
2000 /* Ensure SIMD register alignment */
2001 dr = gmx_simd_align_r(dr_array);
2002 buf = gmx_simd_align_r(buf_array);
2004 /* Extract aligned pointer for parameters and variables */
2005 cp = buf + 0*GMX_SIMD_REAL_WIDTH;
2006 phi0 = buf + 1*GMX_SIMD_REAL_WIDTH;
2007 mult = buf + 2*GMX_SIMD_REAL_WIDTH;
2008 p = buf + 3*GMX_SIMD_REAL_WIDTH;
2009 q = buf + 4*GMX_SIMD_REAL_WIDTH;
2010 sf_i = buf + 5*GMX_SIMD_REAL_WIDTH;
2011 msf_l = buf + 6*GMX_SIMD_REAL_WIDTH;
2013 set_pbc_simd(pbc, &pbc_simd);
2015 /* nbonds is the number of dihedrals times nfa1, here we step GMX_SIMD_REAL_WIDTH dihs */
2016 for (i = 0; (i < nbonds); i += GMX_SIMD_REAL_WIDTH*nfa1)
2018 /* Collect atoms quadruplets for GMX_SIMD_REAL_WIDTH dihedrals.
2019 * iu indexes into forceatoms, we should not let iu go beyond nbonds.
2022 for (s = 0; s < GMX_SIMD_REAL_WIDTH; s++)
2024 type = forceatoms[iu];
2025 ai[s] = forceatoms[iu+1];
2026 aj[s] = forceatoms[iu+2];
2027 ak[s] = forceatoms[iu+3];
2028 al[s] = forceatoms[iu+4];
2030 cp[s] = forceparams[type].pdihs.cpA;
2031 phi0[s] = forceparams[type].pdihs.phiA*DEG2RAD;
2032 mult[s] = forceparams[type].pdihs.mult;
2034 /* At the end fill the arrays with identical entries */
2035 if (iu + nfa1 < nbonds)
2041 /* Caclulate GMX_SIMD_REAL_WIDTH dihedral angles at once */
2042 dih_angle_simd(x, ai, aj, ak, al, &pbc_simd,
2045 &mx_S, &my_S, &mz_S,
2046 &nx_S, &ny_S, &nz_S,
2051 cp_S = gmx_simd_load_r(cp);
2052 phi0_S = gmx_simd_load_r(phi0);
2053 mult_S = gmx_simd_load_r(mult);
2055 mdphi_S = gmx_simd_sub_r(gmx_simd_mul_r(mult_S, phi_S), phi0_S);
2057 /* Calculate GMX_SIMD_REAL_WIDTH sines at once */
2058 gmx_simd_sincos_r(mdphi_S, &sin_S, &cos_S);
2059 mddphi_S = gmx_simd_mul_r(gmx_simd_mul_r(cp_S, mult_S), sin_S);
2060 sf_i_S = gmx_simd_mul_r(mddphi_S, nrkj_m2_S);
2061 msf_l_S = gmx_simd_mul_r(mddphi_S, nrkj_n2_S);
2063 /* After this m?_S will contain f[i] */
2064 mx_S = gmx_simd_mul_r(sf_i_S, mx_S);
2065 my_S = gmx_simd_mul_r(sf_i_S, my_S);
2066 mz_S = gmx_simd_mul_r(sf_i_S, mz_S);
2068 /* After this m?_S will contain -f[l] */
2069 nx_S = gmx_simd_mul_r(msf_l_S, nx_S);
2070 ny_S = gmx_simd_mul_r(msf_l_S, ny_S);
2071 nz_S = gmx_simd_mul_r(msf_l_S, nz_S);
2073 gmx_simd_store_r(dr + 0*GMX_SIMD_REAL_WIDTH, mx_S);
2074 gmx_simd_store_r(dr + 1*GMX_SIMD_REAL_WIDTH, my_S);
2075 gmx_simd_store_r(dr + 2*GMX_SIMD_REAL_WIDTH, mz_S);
2076 gmx_simd_store_r(dr + 3*GMX_SIMD_REAL_WIDTH, nx_S);
2077 gmx_simd_store_r(dr + 4*GMX_SIMD_REAL_WIDTH, ny_S);
2078 gmx_simd_store_r(dr + 5*GMX_SIMD_REAL_WIDTH, nz_S);
2084 do_dih_fup_noshiftf_precalc(ai[s], aj[s], ak[s], al[s],
2086 dr[ XX *GMX_SIMD_REAL_WIDTH+s],
2087 dr[ YY *GMX_SIMD_REAL_WIDTH+s],
2088 dr[ ZZ *GMX_SIMD_REAL_WIDTH+s],
2089 dr[(DIM+XX)*GMX_SIMD_REAL_WIDTH+s],
2090 dr[(DIM+YY)*GMX_SIMD_REAL_WIDTH+s],
2091 dr[(DIM+ZZ)*GMX_SIMD_REAL_WIDTH+s],
2096 while (s < GMX_SIMD_REAL_WIDTH && iu < nbonds);
2100 #endif /* SIMD_BONDEDS */
2103 real idihs(int nbonds,
2104 const t_iatom forceatoms[], const t_iparams forceparams[],
2105 const rvec x[], rvec f[], rvec fshift[],
2106 const t_pbc *pbc, const t_graph *g,
2107 real lambda, real *dvdlambda,
2108 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2109 int gmx_unused *global_atom_index)
2111 int i, type, ai, aj, ak, al;
2113 real phi, phi0, dphi0, ddphi, sign, vtot;
2114 rvec r_ij, r_kj, r_kl, m, n;
2115 real L1, kk, dp, dp2, kA, kB, pA, pB, dvdl_term;
2120 for (i = 0; (i < nbonds); )
2122 type = forceatoms[i++];
2123 ai = forceatoms[i++];
2124 aj = forceatoms[i++];
2125 ak = forceatoms[i++];
2126 al = forceatoms[i++];
2128 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2129 &sign, &t1, &t2, &t3); /* 84 */
2131 /* phi can jump if phi0 is close to Pi/-Pi, which will cause huge
2132 * force changes if we just apply a normal harmonic.
2133 * Instead, we first calculate phi-phi0 and take it modulo (-Pi,Pi).
2134 * This means we will never have the periodicity problem, unless
2135 * the dihedral is Pi away from phiO, which is very unlikely due to
2138 kA = forceparams[type].harmonic.krA;
2139 kB = forceparams[type].harmonic.krB;
2140 pA = forceparams[type].harmonic.rA;
2141 pB = forceparams[type].harmonic.rB;
2143 kk = L1*kA + lambda*kB;
2144 phi0 = (L1*pA + lambda*pB)*DEG2RAD;
2145 dphi0 = (pB - pA)*DEG2RAD;
2149 make_dp_periodic(&dp);
2156 dvdl_term += 0.5*(kB - kA)*dp2 - kk*dphi0*dp;
2158 do_dih_fup(ai, aj, ak, al, (real)(-ddphi), r_ij, r_kj, r_kl, m, n,
2159 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2164 fprintf(debug, "idih: (%d,%d,%d,%d) phi=%g\n",
2165 ai, aj, ak, al, phi);
2170 *dvdlambda += dvdl_term;
2175 /*! \brief returns dx, rdist, and dpdl for functions posres() and fbposres()
2177 static void posres_dx(const rvec x, const rvec pos0A, const rvec pos0B,
2178 const rvec comA_sc, const rvec comB_sc,
2180 t_pbc *pbc, int refcoord_scaling, int npbcdim,
2181 rvec dx, rvec rdist, rvec dpdl)
2184 real posA, posB, L1, ref = 0.;
2189 for (m = 0; m < DIM; m++)
2195 switch (refcoord_scaling)
2199 rdist[m] = L1*posA + lambda*posB;
2200 dpdl[m] = posB - posA;
2203 /* Box relative coordinates are stored for dimensions with pbc */
2204 posA *= pbc->box[m][m];
2205 posB *= pbc->box[m][m];
2206 for (d = m+1; d < npbcdim; d++)
2208 posA += pos0A[d]*pbc->box[d][m];
2209 posB += pos0B[d]*pbc->box[d][m];
2211 ref = L1*posA + lambda*posB;
2213 dpdl[m] = posB - posA;
2216 ref = L1*comA_sc[m] + lambda*comB_sc[m];
2217 rdist[m] = L1*posA + lambda*posB;
2218 dpdl[m] = comB_sc[m] - comA_sc[m] + posB - posA;
2221 gmx_fatal(FARGS, "No such scaling method implemented");
2226 ref = L1*posA + lambda*posB;
2228 dpdl[m] = posB - posA;
2231 /* We do pbc_dx with ref+rdist,
2232 * since with only ref we can be up to half a box vector wrong.
2234 pos[m] = ref + rdist[m];
2239 pbc_dx(pbc, x, pos, dx);
2243 rvec_sub(x, pos, dx);
2247 /*! \brief Adds forces of flat-bottomed positions restraints to f[]
2248 * and fixes vir_diag. Returns the flat-bottomed potential. */
2249 real fbposres(int nbonds,
2250 const t_iatom forceatoms[], const t_iparams forceparams[],
2251 const rvec x[], rvec f[], rvec vir_diag,
2253 int refcoord_scaling, int ePBC, rvec com)
2254 /* compute flat-bottomed positions restraints */
2256 int i, ai, m, d, type, npbcdim = 0, fbdim;
2257 const t_iparams *pr;
2259 real ref = 0, dr, dr2, rpot, rfb, rfb2, fact, invdr;
2260 rvec com_sc, rdist, pos, dx, dpdl, fm;
2263 npbcdim = ePBC2npbcdim(ePBC);
2265 if (refcoord_scaling == erscCOM)
2268 for (m = 0; m < npbcdim; m++)
2270 for (d = m; d < npbcdim; d++)
2272 com_sc[m] += com[d]*pbc->box[d][m];
2278 for (i = 0; (i < nbonds); )
2280 type = forceatoms[i++];
2281 ai = forceatoms[i++];
2282 pr = &forceparams[type];
2284 /* same calculation as for normal posres, but with identical A and B states, and lambda==0 */
2285 posres_dx(x[ai], forceparams[type].fbposres.pos0, forceparams[type].fbposres.pos0,
2286 com_sc, com_sc, 0.0,
2287 pbc, refcoord_scaling, npbcdim,
2293 kk = pr->fbposres.k;
2294 rfb = pr->fbposres.r;
2297 /* with rfb<0, push particle out of the sphere/cylinder/layer */
2305 switch (pr->fbposres.geom)
2307 case efbposresSPHERE:
2308 /* spherical flat-bottom posres */
2311 ( (dr2 > rfb2 && bInvert == FALSE ) || (dr2 < rfb2 && bInvert == TRUE ) )
2315 v = 0.5*kk*sqr(dr - rfb);
2316 fact = -kk*(dr-rfb)/dr; /* Force pointing to the center pos0 */
2317 svmul(fact, dx, fm);
2320 case efbposresCYLINDER:
2321 /* cylidrical flat-bottom posres in x-y plane. fm[ZZ] = 0. */
2322 dr2 = sqr(dx[XX])+sqr(dx[YY]);
2324 ( (dr2 > rfb2 && bInvert == FALSE ) || (dr2 < rfb2 && bInvert == TRUE ) )
2329 v = 0.5*kk*sqr(dr - rfb);
2330 fm[XX] = -kk*(dr-rfb)*dx[XX]*invdr; /* Force pointing to the center */
2331 fm[YY] = -kk*(dr-rfb)*dx[YY]*invdr;
2334 case efbposresX: /* fbdim=XX */
2335 case efbposresY: /* fbdim=YY */
2336 case efbposresZ: /* fbdim=ZZ */
2337 /* 1D flat-bottom potential */
2338 fbdim = pr->fbposres.geom - efbposresX;
2340 if ( ( dr > rfb && bInvert == FALSE ) || ( 0 < dr && dr < rfb && bInvert == TRUE ) )
2342 v = 0.5*kk*sqr(dr - rfb);
2343 fm[fbdim] = -kk*(dr - rfb);
2345 else if ( (dr < (-rfb) && bInvert == FALSE ) || ( (-rfb) < dr && dr < 0 && bInvert == TRUE ))
2347 v = 0.5*kk*sqr(dr + rfb);
2348 fm[fbdim] = -kk*(dr + rfb);
2355 for (m = 0; (m < DIM); m++)
2358 /* Here we correct for the pbc_dx which included rdist */
2359 vir_diag[m] -= 0.5*(dx[m] + rdist[m])*fm[m];
2367 real posres(int nbonds,
2368 const t_iatom forceatoms[], const t_iparams forceparams[],
2369 const rvec x[], rvec f[], rvec vir_diag,
2371 real lambda, real *dvdlambda,
2372 int refcoord_scaling, int ePBC, rvec comA, rvec comB)
2374 int i, ai, m, d, type, ki, npbcdim = 0;
2375 const t_iparams *pr;
2378 real posA, posB, ref = 0;
2379 rvec comA_sc, comB_sc, rdist, dpdl, pos, dx;
2380 gmx_bool bForceValid = TRUE;
2382 if ((f == NULL) || (vir_diag == NULL)) /* should both be null together! */
2384 bForceValid = FALSE;
2387 npbcdim = ePBC2npbcdim(ePBC);
2389 if (refcoord_scaling == erscCOM)
2391 clear_rvec(comA_sc);
2392 clear_rvec(comB_sc);
2393 for (m = 0; m < npbcdim; m++)
2395 for (d = m; d < npbcdim; d++)
2397 comA_sc[m] += comA[d]*pbc->box[d][m];
2398 comB_sc[m] += comB[d]*pbc->box[d][m];
2406 for (i = 0; (i < nbonds); )
2408 type = forceatoms[i++];
2409 ai = forceatoms[i++];
2410 pr = &forceparams[type];
2412 /* return dx, rdist, and dpdl */
2413 posres_dx(x[ai], forceparams[type].posres.pos0A, forceparams[type].posres.pos0B,
2414 comA_sc, comB_sc, lambda,
2415 pbc, refcoord_scaling, npbcdim,
2418 for (m = 0; (m < DIM); m++)
2420 kk = L1*pr->posres.fcA[m] + lambda*pr->posres.fcB[m];
2422 vtot += 0.5*kk*dx[m]*dx[m];
2424 0.5*(pr->posres.fcB[m] - pr->posres.fcA[m])*dx[m]*dx[m]
2427 /* Here we correct for the pbc_dx which included rdist */
2431 vir_diag[m] -= 0.5*(dx[m] + rdist[m])*fm;
2439 static real low_angres(int nbonds,
2440 const t_iatom forceatoms[], const t_iparams forceparams[],
2441 const rvec x[], rvec f[], rvec fshift[],
2442 const t_pbc *pbc, const t_graph *g,
2443 real lambda, real *dvdlambda,
2446 int i, m, type, ai, aj, ak, al;
2448 real phi, cos_phi, cos_phi2, vid, vtot, dVdphi;
2449 rvec r_ij, r_kl, f_i, f_k = {0, 0, 0};
2450 real st, sth, nrij2, nrkl2, c, cij, ckl;
2453 t2 = 0; /* avoid warning with gcc-3.3. It is never used uninitialized */
2456 ak = al = 0; /* to avoid warnings */
2457 for (i = 0; i < nbonds; )
2459 type = forceatoms[i++];
2460 ai = forceatoms[i++];
2461 aj = forceatoms[i++];
2462 t1 = pbc_rvec_sub(pbc, x[aj], x[ai], r_ij); /* 3 */
2465 ak = forceatoms[i++];
2466 al = forceatoms[i++];
2467 t2 = pbc_rvec_sub(pbc, x[al], x[ak], r_kl); /* 3 */
2476 cos_phi = cos_angle(r_ij, r_kl); /* 25 */
2477 phi = acos(cos_phi); /* 10 */
2479 *dvdlambda += dopdihs_min(forceparams[type].pdihs.cpA,
2480 forceparams[type].pdihs.cpB,
2481 forceparams[type].pdihs.phiA,
2482 forceparams[type].pdihs.phiB,
2483 forceparams[type].pdihs.mult,
2484 phi, lambda, &vid, &dVdphi); /* 40 */
2488 cos_phi2 = sqr(cos_phi); /* 1 */
2491 st = -dVdphi*gmx_invsqrt(1 - cos_phi2); /* 12 */
2492 sth = st*cos_phi; /* 1 */
2493 nrij2 = iprod(r_ij, r_ij); /* 5 */
2494 nrkl2 = iprod(r_kl, r_kl); /* 5 */
2496 c = st*gmx_invsqrt(nrij2*nrkl2); /* 11 */
2497 cij = sth/nrij2; /* 10 */
2498 ckl = sth/nrkl2; /* 10 */
2500 for (m = 0; m < DIM; m++) /* 18+18 */
2502 f_i[m] = (c*r_kl[m]-cij*r_ij[m]);
2507 f_k[m] = (c*r_ij[m]-ckl*r_kl[m]);
2515 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
2518 rvec_inc(fshift[t1], f_i);
2519 rvec_dec(fshift[CENTRAL], f_i);
2524 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, al), dt);
2527 rvec_inc(fshift[t2], f_k);
2528 rvec_dec(fshift[CENTRAL], f_k);
2533 return vtot; /* 184 / 157 (bZAxis) total */
2536 real angres(int nbonds,
2537 const t_iatom forceatoms[], const t_iparams forceparams[],
2538 const rvec x[], rvec f[], rvec fshift[],
2539 const t_pbc *pbc, const t_graph *g,
2540 real lambda, real *dvdlambda,
2541 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2542 int gmx_unused *global_atom_index)
2544 return low_angres(nbonds, forceatoms, forceparams, x, f, fshift, pbc, g,
2545 lambda, dvdlambda, FALSE);
2548 real angresz(int nbonds,
2549 const t_iatom forceatoms[], const t_iparams forceparams[],
2550 const rvec x[], rvec f[], rvec fshift[],
2551 const t_pbc *pbc, const t_graph *g,
2552 real lambda, real *dvdlambda,
2553 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2554 int gmx_unused *global_atom_index)
2556 return low_angres(nbonds, forceatoms, forceparams, x, f, fshift, pbc, g,
2557 lambda, dvdlambda, TRUE);
2560 real dihres(int nbonds,
2561 const t_iatom forceatoms[], const t_iparams forceparams[],
2562 const rvec x[], rvec f[], rvec fshift[],
2563 const t_pbc *pbc, const t_graph *g,
2564 real lambda, real *dvdlambda,
2565 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2566 int gmx_unused *global_atom_index)
2569 int ai, aj, ak, al, i, k, type, t1, t2, t3;
2570 real phi0A, phi0B, dphiA, dphiB, kfacA, kfacB, phi0, dphi, kfac;
2571 real phi, ddphi, ddp, ddp2, dp, sign, d2r, fc, L1;
2572 rvec r_ij, r_kj, r_kl, m, n;
2579 for (i = 0; (i < nbonds); )
2581 type = forceatoms[i++];
2582 ai = forceatoms[i++];
2583 aj = forceatoms[i++];
2584 ak = forceatoms[i++];
2585 al = forceatoms[i++];
2587 phi0A = forceparams[type].dihres.phiA*d2r;
2588 dphiA = forceparams[type].dihres.dphiA*d2r;
2589 kfacA = forceparams[type].dihres.kfacA;
2591 phi0B = forceparams[type].dihres.phiB*d2r;
2592 dphiB = forceparams[type].dihres.dphiB*d2r;
2593 kfacB = forceparams[type].dihres.kfacB;
2595 phi0 = L1*phi0A + lambda*phi0B;
2596 dphi = L1*dphiA + lambda*dphiB;
2597 kfac = L1*kfacA + lambda*kfacB;
2599 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2600 &sign, &t1, &t2, &t3);
2605 fprintf(debug, "dihres[%d]: %d %d %d %d : phi=%f, dphi=%f, kfac=%f\n",
2606 k++, ai, aj, ak, al, phi0, dphi, kfac);
2608 /* phi can jump if phi0 is close to Pi/-Pi, which will cause huge
2609 * force changes if we just apply a normal harmonic.
2610 * Instead, we first calculate phi-phi0 and take it modulo (-Pi,Pi).
2611 * This means we will never have the periodicity problem, unless
2612 * the dihedral is Pi away from phiO, which is very unlikely due to
2616 make_dp_periodic(&dp);
2622 else if (dp < -dphi)
2634 vtot += 0.5*kfac*ddp2;
2637 *dvdlambda += 0.5*(kfacB - kfacA)*ddp2;
2638 /* lambda dependence from changing restraint distances */
2641 *dvdlambda -= kfac*ddp*((dphiB - dphiA)+(phi0B - phi0A));
2645 *dvdlambda += kfac*ddp*((dphiB - dphiA)-(phi0B - phi0A));
2647 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
2648 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2655 real unimplemented(int gmx_unused nbonds,
2656 const t_iatom gmx_unused forceatoms[], const t_iparams gmx_unused forceparams[],
2657 const rvec gmx_unused x[], rvec gmx_unused f[], rvec gmx_unused fshift[],
2658 const t_pbc gmx_unused *pbc, const t_graph gmx_unused *g,
2659 real gmx_unused lambda, real gmx_unused *dvdlambda,
2660 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2661 int gmx_unused *global_atom_index)
2663 gmx_impl("*** you are using a not implemented function");
2665 return 0.0; /* To make the compiler happy */
2668 real rbdihs(int nbonds,
2669 const t_iatom forceatoms[], const t_iparams forceparams[],
2670 const rvec x[], rvec f[], rvec fshift[],
2671 const t_pbc *pbc, const t_graph *g,
2672 real lambda, real *dvdlambda,
2673 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2674 int gmx_unused *global_atom_index)
2676 const real c0 = 0.0, c1 = 1.0, c2 = 2.0, c3 = 3.0, c4 = 4.0, c5 = 5.0;
2677 int type, ai, aj, ak, al, i, j;
2679 rvec r_ij, r_kj, r_kl, m, n;
2680 real parmA[NR_RBDIHS];
2681 real parmB[NR_RBDIHS];
2682 real parm[NR_RBDIHS];
2683 real cos_phi, phi, rbp, rbpBA;
2684 real v, sign, ddphi, sin_phi;
2686 real L1 = 1.0-lambda;
2690 for (i = 0; (i < nbonds); )
2692 type = forceatoms[i++];
2693 ai = forceatoms[i++];
2694 aj = forceatoms[i++];
2695 ak = forceatoms[i++];
2696 al = forceatoms[i++];
2698 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2699 &sign, &t1, &t2, &t3); /* 84 */
2701 /* Change to polymer convention */
2708 phi -= M_PI; /* 1 */
2712 /* Beware of accuracy loss, cannot use 1-sqrt(cos^2) ! */
2715 for (j = 0; (j < NR_RBDIHS); j++)
2717 parmA[j] = forceparams[type].rbdihs.rbcA[j];
2718 parmB[j] = forceparams[type].rbdihs.rbcB[j];
2719 parm[j] = L1*parmA[j]+lambda*parmB[j];
2721 /* Calculate cosine powers */
2722 /* Calculate the energy */
2723 /* Calculate the derivative */
2726 dvdl_term += (parmB[0]-parmA[0]);
2731 rbpBA = parmB[1]-parmA[1];
2732 ddphi += rbp*cosfac;
2735 dvdl_term += cosfac*rbpBA;
2737 rbpBA = parmB[2]-parmA[2];
2738 ddphi += c2*rbp*cosfac;
2741 dvdl_term += cosfac*rbpBA;
2743 rbpBA = parmB[3]-parmA[3];
2744 ddphi += c3*rbp*cosfac;
2747 dvdl_term += cosfac*rbpBA;
2749 rbpBA = parmB[4]-parmA[4];
2750 ddphi += c4*rbp*cosfac;
2753 dvdl_term += cosfac*rbpBA;
2755 rbpBA = parmB[5]-parmA[5];
2756 ddphi += c5*rbp*cosfac;
2759 dvdl_term += cosfac*rbpBA;
2761 ddphi = -ddphi*sin_phi; /* 11 */
2763 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
2764 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2767 *dvdlambda += dvdl_term;
2772 int cmap_setup_grid_index(int ip, int grid_spacing, int *ipm1, int *ipp1, int *ipp2)
2778 ip = ip + grid_spacing - 1;
2780 else if (ip > grid_spacing)
2782 ip = ip - grid_spacing - 1;
2791 im1 = grid_spacing - 1;
2793 else if (ip == grid_spacing-2)
2797 else if (ip == grid_spacing-1)
2811 real cmap_dihs(int nbonds,
2812 const t_iatom forceatoms[], const t_iparams forceparams[],
2813 const gmx_cmap_t *cmap_grid,
2814 const rvec x[], rvec f[], rvec fshift[],
2815 const t_pbc *pbc, const t_graph *g,
2816 real gmx_unused lambda, real gmx_unused *dvdlambda,
2817 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2818 int gmx_unused *global_atom_index)
2820 int i, j, k, n, idx;
2821 int ai, aj, ak, al, am;
2822 int a1i, a1j, a1k, a1l, a2i, a2j, a2k, a2l;
2824 int t11, t21, t31, t12, t22, t32;
2825 int iphi1, ip1m1, ip1p1, ip1p2;
2826 int iphi2, ip2m1, ip2p1, ip2p2;
2828 int pos1, pos2, pos3, pos4, tmp;
2830 real ty[4], ty1[4], ty2[4], ty12[4], tc[16], tx[16];
2831 real phi1, psi1, cos_phi1, sin_phi1, sign1, xphi1;
2832 real phi2, psi2, cos_phi2, sin_phi2, sign2, xphi2;
2833 real dx, xx, tt, tu, e, df1, df2, ddf1, ddf2, ddf12, vtot;
2834 real ra21, rb21, rg21, rg1, rgr1, ra2r1, rb2r1, rabr1;
2835 real ra22, rb22, rg22, rg2, rgr2, ra2r2, rb2r2, rabr2;
2836 real fg1, hg1, fga1, hgb1, gaa1, gbb1;
2837 real fg2, hg2, fga2, hgb2, gaa2, gbb2;
2840 rvec r1_ij, r1_kj, r1_kl, m1, n1;
2841 rvec r2_ij, r2_kj, r2_kl, m2, n2;
2842 rvec f1_i, f1_j, f1_k, f1_l;
2843 rvec f2_i, f2_j, f2_k, f2_l;
2844 rvec a1, b1, a2, b2;
2845 rvec f1, g1, h1, f2, g2, h2;
2846 rvec dtf1, dtg1, dth1, dtf2, dtg2, dth2;
2847 ivec jt1, dt1_ij, dt1_kj, dt1_lj;
2848 ivec jt2, dt2_ij, dt2_kj, dt2_lj;
2852 int loop_index[4][4] = {
2859 /* Total CMAP energy */
2862 for (n = 0; n < nbonds; )
2864 /* Five atoms are involved in the two torsions */
2865 type = forceatoms[n++];
2866 ai = forceatoms[n++];
2867 aj = forceatoms[n++];
2868 ak = forceatoms[n++];
2869 al = forceatoms[n++];
2870 am = forceatoms[n++];
2872 /* Which CMAP type is this */
2873 cmapA = forceparams[type].cmap.cmapA;
2874 cmapd = cmap_grid->cmapdata[cmapA].cmap;
2882 phi1 = dih_angle(x[a1i], x[a1j], x[a1k], x[a1l], pbc, r1_ij, r1_kj, r1_kl, m1, n1,
2883 &sign1, &t11, &t21, &t31); /* 84 */
2885 cos_phi1 = cos(phi1);
2887 a1[0] = r1_ij[1]*r1_kj[2]-r1_ij[2]*r1_kj[1];
2888 a1[1] = r1_ij[2]*r1_kj[0]-r1_ij[0]*r1_kj[2];
2889 a1[2] = r1_ij[0]*r1_kj[1]-r1_ij[1]*r1_kj[0]; /* 9 */
2891 b1[0] = r1_kl[1]*r1_kj[2]-r1_kl[2]*r1_kj[1];
2892 b1[1] = r1_kl[2]*r1_kj[0]-r1_kl[0]*r1_kj[2];
2893 b1[2] = r1_kl[0]*r1_kj[1]-r1_kl[1]*r1_kj[0]; /* 9 */
2895 tmp = pbc_rvec_sub(pbc, x[a1l], x[a1k], h1);
2897 ra21 = iprod(a1, a1); /* 5 */
2898 rb21 = iprod(b1, b1); /* 5 */
2899 rg21 = iprod(r1_kj, r1_kj); /* 5 */
2905 rabr1 = sqrt(ra2r1*rb2r1);
2907 sin_phi1 = rg1 * rabr1 * iprod(a1, h1) * (-1);
2909 if (cos_phi1 < -0.5 || cos_phi1 > 0.5)
2911 phi1 = asin(sin_phi1);
2921 phi1 = -M_PI - phi1;
2927 phi1 = acos(cos_phi1);
2935 xphi1 = phi1 + M_PI; /* 1 */
2937 /* Second torsion */
2943 phi2 = dih_angle(x[a2i], x[a2j], x[a2k], x[a2l], pbc, r2_ij, r2_kj, r2_kl, m2, n2,
2944 &sign2, &t12, &t22, &t32); /* 84 */
2946 cos_phi2 = cos(phi2);
2948 a2[0] = r2_ij[1]*r2_kj[2]-r2_ij[2]*r2_kj[1];
2949 a2[1] = r2_ij[2]*r2_kj[0]-r2_ij[0]*r2_kj[2];
2950 a2[2] = r2_ij[0]*r2_kj[1]-r2_ij[1]*r2_kj[0]; /* 9 */
2952 b2[0] = r2_kl[1]*r2_kj[2]-r2_kl[2]*r2_kj[1];
2953 b2[1] = r2_kl[2]*r2_kj[0]-r2_kl[0]*r2_kj[2];
2954 b2[2] = r2_kl[0]*r2_kj[1]-r2_kl[1]*r2_kj[0]; /* 9 */
2956 tmp = pbc_rvec_sub(pbc, x[a2l], x[a2k], h2);
2958 ra22 = iprod(a2, a2); /* 5 */
2959 rb22 = iprod(b2, b2); /* 5 */
2960 rg22 = iprod(r2_kj, r2_kj); /* 5 */
2966 rabr2 = sqrt(ra2r2*rb2r2);
2968 sin_phi2 = rg2 * rabr2 * iprod(a2, h2) * (-1);
2970 if (cos_phi2 < -0.5 || cos_phi2 > 0.5)
2972 phi2 = asin(sin_phi2);
2982 phi2 = -M_PI - phi2;
2988 phi2 = acos(cos_phi2);
2996 xphi2 = phi2 + M_PI; /* 1 */
2998 /* Range mangling */
3001 xphi1 = xphi1 + 2*M_PI;
3003 else if (xphi1 >= 2*M_PI)
3005 xphi1 = xphi1 - 2*M_PI;
3010 xphi2 = xphi2 + 2*M_PI;
3012 else if (xphi2 >= 2*M_PI)
3014 xphi2 = xphi2 - 2*M_PI;
3017 /* Number of grid points */
3018 dx = 2*M_PI / cmap_grid->grid_spacing;
3020 /* Where on the grid are we */
3021 iphi1 = (int)(xphi1/dx);
3022 iphi2 = (int)(xphi2/dx);
3024 iphi1 = cmap_setup_grid_index(iphi1, cmap_grid->grid_spacing, &ip1m1, &ip1p1, &ip1p2);
3025 iphi2 = cmap_setup_grid_index(iphi2, cmap_grid->grid_spacing, &ip2m1, &ip2p1, &ip2p2);
3027 pos1 = iphi1*cmap_grid->grid_spacing+iphi2;
3028 pos2 = ip1p1*cmap_grid->grid_spacing+iphi2;
3029 pos3 = ip1p1*cmap_grid->grid_spacing+ip2p1;
3030 pos4 = iphi1*cmap_grid->grid_spacing+ip2p1;
3032 ty[0] = cmapd[pos1*4];
3033 ty[1] = cmapd[pos2*4];
3034 ty[2] = cmapd[pos3*4];
3035 ty[3] = cmapd[pos4*4];
3037 ty1[0] = cmapd[pos1*4+1];
3038 ty1[1] = cmapd[pos2*4+1];
3039 ty1[2] = cmapd[pos3*4+1];
3040 ty1[3] = cmapd[pos4*4+1];
3042 ty2[0] = cmapd[pos1*4+2];
3043 ty2[1] = cmapd[pos2*4+2];
3044 ty2[2] = cmapd[pos3*4+2];
3045 ty2[3] = cmapd[pos4*4+2];
3047 ty12[0] = cmapd[pos1*4+3];
3048 ty12[1] = cmapd[pos2*4+3];
3049 ty12[2] = cmapd[pos3*4+3];
3050 ty12[3] = cmapd[pos4*4+3];
3052 /* Switch to degrees */
3053 dx = 360.0 / cmap_grid->grid_spacing;
3054 xphi1 = xphi1 * RAD2DEG;
3055 xphi2 = xphi2 * RAD2DEG;
3057 for (i = 0; i < 4; i++) /* 16 */
3060 tx[i+4] = ty1[i]*dx;
3061 tx[i+8] = ty2[i]*dx;
3062 tx[i+12] = ty12[i]*dx*dx;
3066 for (i = 0; i < 4; i++) /* 1056 */
3068 for (j = 0; j < 4; j++)
3071 for (k = 0; k < 16; k++)
3073 xx = xx + cmap_coeff_matrix[k*16+idx]*tx[k];
3081 tt = (xphi1-iphi1*dx)/dx;
3082 tu = (xphi2-iphi2*dx)/dx;
3091 for (i = 3; i >= 0; i--)
3093 l1 = loop_index[i][3];
3094 l2 = loop_index[i][2];
3095 l3 = loop_index[i][1];
3097 e = tt * e + ((tc[i*4+3]*tu+tc[i*4+2])*tu + tc[i*4+1])*tu+tc[i*4];
3098 df1 = tu * df1 + (3.0*tc[l1]*tt+2.0*tc[l2])*tt+tc[l3];
3099 df2 = tt * df2 + (3.0*tc[i*4+3]*tu+2.0*tc[i*4+2])*tu+tc[i*4+1];
3100 ddf1 = tu * ddf1 + 2.0*3.0*tc[l1]*tt+2.0*tc[l2];
3101 ddf2 = tt * ddf2 + 2.0*3.0*tc[4*i+3]*tu+2.0*tc[4*i+2];
3104 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) +
3105 3.0*tu*tu*(tc[7]+2.0*tc[11]*tt+3.0*tc[15]*tt*tt);
3110 ddf1 = ddf1 * fac * fac;
3111 ddf2 = ddf2 * fac * fac;
3112 ddf12 = ddf12 * fac * fac;
3117 /* Do forces - first torsion */
3118 fg1 = iprod(r1_ij, r1_kj);
3119 hg1 = iprod(r1_kl, r1_kj);
3120 fga1 = fg1*ra2r1*rgr1;
3121 hgb1 = hg1*rb2r1*rgr1;
3125 for (i = 0; i < DIM; i++)
3127 dtf1[i] = gaa1 * a1[i];
3128 dtg1[i] = fga1 * a1[i] - hgb1 * b1[i];
3129 dth1[i] = gbb1 * b1[i];
3131 f1[i] = df1 * dtf1[i];
3132 g1[i] = df1 * dtg1[i];
3133 h1[i] = df1 * dth1[i];
3136 f1_j[i] = -f1[i] - g1[i];
3137 f1_k[i] = h1[i] + g1[i];
3140 f[a1i][i] = f[a1i][i] + f1_i[i];
3141 f[a1j][i] = f[a1j][i] + f1_j[i]; /* - f1[i] - g1[i] */
3142 f[a1k][i] = f[a1k][i] + f1_k[i]; /* h1[i] + g1[i] */
3143 f[a1l][i] = f[a1l][i] + f1_l[i]; /* h1[i] */
3146 /* Do forces - second torsion */
3147 fg2 = iprod(r2_ij, r2_kj);
3148 hg2 = iprod(r2_kl, r2_kj);
3149 fga2 = fg2*ra2r2*rgr2;
3150 hgb2 = hg2*rb2r2*rgr2;
3154 for (i = 0; i < DIM; i++)
3156 dtf2[i] = gaa2 * a2[i];
3157 dtg2[i] = fga2 * a2[i] - hgb2 * b2[i];
3158 dth2[i] = gbb2 * b2[i];
3160 f2[i] = df2 * dtf2[i];
3161 g2[i] = df2 * dtg2[i];
3162 h2[i] = df2 * dth2[i];
3165 f2_j[i] = -f2[i] - g2[i];
3166 f2_k[i] = h2[i] + g2[i];
3169 f[a2i][i] = f[a2i][i] + f2_i[i]; /* f2[i] */
3170 f[a2j][i] = f[a2j][i] + f2_j[i]; /* - f2[i] - g2[i] */
3171 f[a2k][i] = f[a2k][i] + f2_k[i]; /* h2[i] + g2[i] */
3172 f[a2l][i] = f[a2l][i] + f2_l[i]; /* - h2[i] */
3178 copy_ivec(SHIFT_IVEC(g, a1j), jt1);
3179 ivec_sub(SHIFT_IVEC(g, a1i), jt1, dt1_ij);
3180 ivec_sub(SHIFT_IVEC(g, a1k), jt1, dt1_kj);
3181 ivec_sub(SHIFT_IVEC(g, a1l), jt1, dt1_lj);
3182 t11 = IVEC2IS(dt1_ij);
3183 t21 = IVEC2IS(dt1_kj);
3184 t31 = IVEC2IS(dt1_lj);
3186 copy_ivec(SHIFT_IVEC(g, a2j), jt2);
3187 ivec_sub(SHIFT_IVEC(g, a2i), jt2, dt2_ij);
3188 ivec_sub(SHIFT_IVEC(g, a2k), jt2, dt2_kj);
3189 ivec_sub(SHIFT_IVEC(g, a2l), jt2, dt2_lj);
3190 t12 = IVEC2IS(dt2_ij);
3191 t22 = IVEC2IS(dt2_kj);
3192 t32 = IVEC2IS(dt2_lj);
3196 t31 = pbc_rvec_sub(pbc, x[a1l], x[a1j], h1);
3197 t32 = pbc_rvec_sub(pbc, x[a2l], x[a2j], h2);
3205 rvec_inc(fshift[t11], f1_i);
3206 rvec_inc(fshift[CENTRAL], f1_j);
3207 rvec_inc(fshift[t21], f1_k);
3208 rvec_inc(fshift[t31], f1_l);
3210 rvec_inc(fshift[t21], f2_i);
3211 rvec_inc(fshift[CENTRAL], f2_j);
3212 rvec_inc(fshift[t22], f2_k);
3213 rvec_inc(fshift[t32], f2_l);
3220 /***********************************************************
3222 * G R O M O S 9 6 F U N C T I O N S
3224 ***********************************************************/
3225 real g96harmonic(real kA, real kB, real xA, real xB, real x, real lambda,
3228 const real half = 0.5;
3229 real L1, kk, x0, dx, dx2;
3230 real v, f, dvdlambda;
3233 kk = L1*kA+lambda*kB;
3234 x0 = L1*xA+lambda*xB;
3241 dvdlambda = half*(kB-kA)*dx2 + (xA-xB)*kk*dx;
3248 /* That was 21 flops */
3251 real g96bonds(int nbonds,
3252 const t_iatom forceatoms[], const t_iparams forceparams[],
3253 const rvec x[], rvec f[], rvec fshift[],
3254 const t_pbc *pbc, const t_graph *g,
3255 real lambda, real *dvdlambda,
3256 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3257 int gmx_unused *global_atom_index)
3259 int i, m, ki, ai, aj, type;
3260 real dr2, fbond, vbond, fij, vtot;
3265 for (i = 0; (i < nbonds); )
3267 type = forceatoms[i++];
3268 ai = forceatoms[i++];
3269 aj = forceatoms[i++];
3271 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
3272 dr2 = iprod(dx, dx); /* 5 */
3274 *dvdlambda += g96harmonic(forceparams[type].harmonic.krA,
3275 forceparams[type].harmonic.krB,
3276 forceparams[type].harmonic.rA,
3277 forceparams[type].harmonic.rB,
3278 dr2, lambda, &vbond, &fbond);
3280 vtot += 0.5*vbond; /* 1*/
3284 fprintf(debug, "G96-BONDS: dr = %10g vbond = %10g fbond = %10g\n",
3285 sqrt(dr2), vbond, fbond);
3291 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
3294 for (m = 0; (m < DIM); m++) /* 15 */
3299 fshift[ki][m] += fij;
3300 fshift[CENTRAL][m] -= fij;
3306 real g96bond_angle(const rvec xi, const rvec xj, const rvec xk, const t_pbc *pbc,
3307 rvec r_ij, rvec r_kj,
3309 /* Return value is the angle between the bonds i-j and j-k */
3313 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
3314 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
3316 costh = cos_angle(r_ij, r_kj); /* 25 */
3321 real g96angles(int nbonds,
3322 const t_iatom forceatoms[], const t_iparams forceparams[],
3323 const rvec x[], rvec f[], rvec fshift[],
3324 const t_pbc *pbc, const t_graph *g,
3325 real lambda, real *dvdlambda,
3326 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3327 int gmx_unused *global_atom_index)
3329 int i, ai, aj, ak, type, m, t1, t2;
3331 real cos_theta, dVdt, va, vtot;
3332 real rij_1, rij_2, rkj_1, rkj_2, rijrkj_1;
3334 ivec jt, dt_ij, dt_kj;
3337 for (i = 0; (i < nbonds); )
3339 type = forceatoms[i++];
3340 ai = forceatoms[i++];
3341 aj = forceatoms[i++];
3342 ak = forceatoms[i++];
3344 cos_theta = g96bond_angle(x[ai], x[aj], x[ak], pbc, r_ij, r_kj, &t1, &t2);
3346 *dvdlambda += g96harmonic(forceparams[type].harmonic.krA,
3347 forceparams[type].harmonic.krB,
3348 forceparams[type].harmonic.rA,
3349 forceparams[type].harmonic.rB,
3350 cos_theta, lambda, &va, &dVdt);
3353 rij_1 = gmx_invsqrt(iprod(r_ij, r_ij));
3354 rkj_1 = gmx_invsqrt(iprod(r_kj, r_kj));
3355 rij_2 = rij_1*rij_1;
3356 rkj_2 = rkj_1*rkj_1;
3357 rijrkj_1 = rij_1*rkj_1; /* 23 */
3362 fprintf(debug, "G96ANGLES: costheta = %10g vth = %10g dV/dct = %10g\n",
3363 cos_theta, va, dVdt);
3366 for (m = 0; (m < DIM); m++) /* 42 */
3368 f_i[m] = dVdt*(r_kj[m]*rijrkj_1 - r_ij[m]*rij_2*cos_theta);
3369 f_k[m] = dVdt*(r_ij[m]*rijrkj_1 - r_kj[m]*rkj_2*cos_theta);
3370 f_j[m] = -f_i[m]-f_k[m];
3378 copy_ivec(SHIFT_IVEC(g, aj), jt);
3380 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3381 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3382 t1 = IVEC2IS(dt_ij);
3383 t2 = IVEC2IS(dt_kj);
3385 rvec_inc(fshift[t1], f_i);
3386 rvec_inc(fshift[CENTRAL], f_j);
3387 rvec_inc(fshift[t2], f_k); /* 9 */
3393 real cross_bond_bond(int nbonds,
3394 const t_iatom forceatoms[], const t_iparams forceparams[],
3395 const rvec x[], rvec f[], rvec fshift[],
3396 const t_pbc *pbc, const t_graph *g,
3397 real gmx_unused lambda, real gmx_unused *dvdlambda,
3398 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3399 int gmx_unused *global_atom_index)
3401 /* Potential from Lawrence and Skimmer, Chem. Phys. Lett. 372 (2003)
3404 int i, ai, aj, ak, type, m, t1, t2;
3406 real vtot, vrr, s1, s2, r1, r2, r1e, r2e, krr;
3408 ivec jt, dt_ij, dt_kj;
3411 for (i = 0; (i < nbonds); )
3413 type = forceatoms[i++];
3414 ai = forceatoms[i++];
3415 aj = forceatoms[i++];
3416 ak = forceatoms[i++];
3417 r1e = forceparams[type].cross_bb.r1e;
3418 r2e = forceparams[type].cross_bb.r2e;
3419 krr = forceparams[type].cross_bb.krr;
3421 /* Compute distance vectors ... */
3422 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
3423 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
3425 /* ... and their lengths */
3429 /* Deviations from ideality */
3433 /* Energy (can be negative!) */
3438 svmul(-krr*s2/r1, r_ij, f_i);
3439 svmul(-krr*s1/r2, r_kj, f_k);
3441 for (m = 0; (m < DIM); m++) /* 12 */
3443 f_j[m] = -f_i[m] - f_k[m];
3452 copy_ivec(SHIFT_IVEC(g, aj), jt);
3454 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3455 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3456 t1 = IVEC2IS(dt_ij);
3457 t2 = IVEC2IS(dt_kj);
3459 rvec_inc(fshift[t1], f_i);
3460 rvec_inc(fshift[CENTRAL], f_j);
3461 rvec_inc(fshift[t2], f_k); /* 9 */
3467 real cross_bond_angle(int nbonds,
3468 const t_iatom forceatoms[], const t_iparams forceparams[],
3469 const rvec x[], rvec f[], rvec fshift[],
3470 const t_pbc *pbc, const t_graph *g,
3471 real gmx_unused lambda, real gmx_unused *dvdlambda,
3472 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3473 int gmx_unused *global_atom_index)
3475 /* Potential from Lawrence and Skimmer, Chem. Phys. Lett. 372 (2003)
3478 int i, ai, aj, ak, type, m, t1, t2, t3;
3479 rvec r_ij, r_kj, r_ik;
3480 real vtot, vrt, s1, s2, s3, r1, r2, r3, r1e, r2e, r3e, krt, k1, k2, k3;
3482 ivec jt, dt_ij, dt_kj;
3485 for (i = 0; (i < nbonds); )
3487 type = forceatoms[i++];
3488 ai = forceatoms[i++];
3489 aj = forceatoms[i++];
3490 ak = forceatoms[i++];
3491 r1e = forceparams[type].cross_ba.r1e;
3492 r2e = forceparams[type].cross_ba.r2e;
3493 r3e = forceparams[type].cross_ba.r3e;
3494 krt = forceparams[type].cross_ba.krt;
3496 /* Compute distance vectors ... */
3497 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
3498 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
3499 t3 = pbc_rvec_sub(pbc, x[ai], x[ak], r_ik);
3501 /* ... and their lengths */
3506 /* Deviations from ideality */
3511 /* Energy (can be negative!) */
3512 vrt = krt*s3*(s1+s2);
3518 k3 = -krt*(s1+s2)/r3;
3519 for (m = 0; (m < DIM); m++)
3521 f_i[m] = k1*r_ij[m] + k3*r_ik[m];
3522 f_k[m] = k2*r_kj[m] - k3*r_ik[m];
3523 f_j[m] = -f_i[m] - f_k[m];
3526 for (m = 0; (m < DIM); m++) /* 12 */
3536 copy_ivec(SHIFT_IVEC(g, aj), jt);
3538 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3539 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3540 t1 = IVEC2IS(dt_ij);
3541 t2 = IVEC2IS(dt_kj);
3543 rvec_inc(fshift[t1], f_i);
3544 rvec_inc(fshift[CENTRAL], f_j);
3545 rvec_inc(fshift[t2], f_k); /* 9 */
3551 static real bonded_tab(const char *type, int table_nr,
3552 const bondedtable_t *table, real kA, real kB, real r,
3553 real lambda, real *V, real *F)
3555 real k, tabscale, *VFtab, rt, eps, eps2, Yt, Ft, Geps, Heps2, Fp, VV, FF;
3557 real v, f, dvdlambda;
3559 k = (1.0 - lambda)*kA + lambda*kB;
3561 tabscale = table->scale;
3562 VFtab = table->data;
3568 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",
3569 type, table_nr, r, n0, n0+1, table->n);
3576 Geps = VFtab[nnn+2]*eps;
3577 Heps2 = VFtab[nnn+3]*eps2;
3578 Fp = Ft + Geps + Heps2;
3580 FF = Fp + Geps + 2.0*Heps2;
3582 *F = -k*FF*tabscale;
3584 dvdlambda = (kB - kA)*VV;
3588 /* That was 22 flops */
3591 real tab_bonds(int nbonds,
3592 const t_iatom forceatoms[], const t_iparams forceparams[],
3593 const rvec x[], rvec f[], rvec fshift[],
3594 const t_pbc *pbc, const t_graph *g,
3595 real lambda, real *dvdlambda,
3596 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
3597 int gmx_unused *global_atom_index)
3599 int i, m, ki, ai, aj, type, table;
3600 real dr, dr2, fbond, vbond, fij, vtot;
3605 for (i = 0; (i < nbonds); )
3607 type = forceatoms[i++];
3608 ai = forceatoms[i++];
3609 aj = forceatoms[i++];
3611 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
3612 dr2 = iprod(dx, dx); /* 5 */
3613 dr = dr2*gmx_invsqrt(dr2); /* 10 */
3615 table = forceparams[type].tab.table;
3617 *dvdlambda += bonded_tab("bond", table,
3618 &fcd->bondtab[table],
3619 forceparams[type].tab.kA,
3620 forceparams[type].tab.kB,
3621 dr, lambda, &vbond, &fbond); /* 22 */
3629 vtot += vbond; /* 1*/
3630 fbond *= gmx_invsqrt(dr2); /* 6 */
3634 fprintf(debug, "TABBONDS: dr = %10g vbond = %10g fbond = %10g\n",
3640 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
3643 for (m = 0; (m < DIM); m++) /* 15 */
3648 fshift[ki][m] += fij;
3649 fshift[CENTRAL][m] -= fij;
3655 real tab_angles(int nbonds,
3656 const t_iatom forceatoms[], const t_iparams forceparams[],
3657 const rvec x[], rvec f[], rvec fshift[],
3658 const t_pbc *pbc, const t_graph *g,
3659 real lambda, real *dvdlambda,
3660 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
3661 int gmx_unused *global_atom_index)
3663 int i, ai, aj, ak, t1, t2, type, table;
3665 real cos_theta, cos_theta2, theta, dVdt, va, vtot;
3666 ivec jt, dt_ij, dt_kj;
3669 for (i = 0; (i < nbonds); )
3671 type = forceatoms[i++];
3672 ai = forceatoms[i++];
3673 aj = forceatoms[i++];
3674 ak = forceatoms[i++];
3676 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
3677 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
3679 table = forceparams[type].tab.table;
3681 *dvdlambda += bonded_tab("angle", table,
3682 &fcd->angletab[table],
3683 forceparams[type].tab.kA,
3684 forceparams[type].tab.kB,
3685 theta, lambda, &va, &dVdt); /* 22 */
3688 cos_theta2 = sqr(cos_theta); /* 1 */
3697 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
3698 sth = st*cos_theta; /* 1 */
3702 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
3703 theta*RAD2DEG, va, dVdt);
3706 nrkj2 = iprod(r_kj, r_kj); /* 5 */
3707 nrij2 = iprod(r_ij, r_ij);
3709 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
3710 cii = sth/nrij2; /* 10 */
3711 ckk = sth/nrkj2; /* 10 */
3713 for (m = 0; (m < DIM); m++) /* 39 */
3715 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
3716 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
3717 f_j[m] = -f_i[m]-f_k[m];
3724 copy_ivec(SHIFT_IVEC(g, aj), jt);
3726 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3727 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3728 t1 = IVEC2IS(dt_ij);
3729 t2 = IVEC2IS(dt_kj);
3731 rvec_inc(fshift[t1], f_i);
3732 rvec_inc(fshift[CENTRAL], f_j);
3733 rvec_inc(fshift[t2], f_k);
3739 real tab_dihs(int nbonds,
3740 const t_iatom forceatoms[], const t_iparams forceparams[],
3741 const rvec x[], rvec f[], rvec fshift[],
3742 const t_pbc *pbc, const t_graph *g,
3743 real lambda, real *dvdlambda,
3744 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
3745 int gmx_unused *global_atom_index)
3747 int i, type, ai, aj, ak, al, table;
3749 rvec r_ij, r_kj, r_kl, m, n;
3750 real phi, sign, ddphi, vpd, vtot;
3753 for (i = 0; (i < nbonds); )
3755 type = forceatoms[i++];
3756 ai = forceatoms[i++];
3757 aj = forceatoms[i++];
3758 ak = forceatoms[i++];
3759 al = forceatoms[i++];
3761 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
3762 &sign, &t1, &t2, &t3); /* 84 */
3764 table = forceparams[type].tab.table;
3766 /* Hopefully phi+M_PI never results in values < 0 */
3767 *dvdlambda += bonded_tab("dihedral", table,
3768 &fcd->dihtab[table],
3769 forceparams[type].tab.kA,
3770 forceparams[type].tab.kB,
3771 phi+M_PI, lambda, &vpd, &ddphi);
3774 do_dih_fup(ai, aj, ak, al, -ddphi, r_ij, r_kj, r_kl, m, n,
3775 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
3778 fprintf(debug, "pdih: (%d,%d,%d,%d) phi=%g\n",
3779 ai, aj, ak, al, phi);
3786 /* Return if this is a potential calculated in bondfree.c,
3787 * i.e. an interaction that actually calculates a potential and
3788 * works on multiple atoms (not e.g. a connection or a position restraint).
3790 static gmx_inline gmx_bool ftype_is_bonded_potential(int ftype)
3793 (interaction_function[ftype].flags & IF_BOND) &&
3794 !(ftype == F_CONNBONDS || ftype == F_POSRES || ftype == F_FBPOSRES) &&
3795 (ftype < F_GB12 || ftype > F_GB14);
3798 static void divide_bondeds_over_threads(t_idef *idef, int nthreads)
3805 idef->nthreads = nthreads;
3807 if (F_NRE*(nthreads+1) > idef->il_thread_division_nalloc)
3809 idef->il_thread_division_nalloc = F_NRE*(nthreads+1);
3810 snew(idef->il_thread_division, idef->il_thread_division_nalloc);
3813 for (ftype = 0; ftype < F_NRE; ftype++)
3815 if (ftype_is_bonded_potential(ftype))
3817 nat1 = interaction_function[ftype].nratoms + 1;
3819 for (t = 0; t <= nthreads; t++)
3821 /* Divide the interactions equally over the threads.
3822 * When the different types of bonded interactions
3823 * are distributed roughly equally over the threads,
3824 * this should lead to well localized output into
3825 * the force buffer on each thread.
3826 * If this is not the case, a more advanced scheme
3827 * (not implemented yet) will do better.
3829 il_nr_thread = (((idef->il[ftype].nr/nat1)*t)/nthreads)*nat1;
3831 /* Ensure that distance restraint pairs with the same label
3832 * end up on the same thread.
3833 * This is slighlty tricky code, since the next for iteration
3834 * may have an initial il_nr_thread lower than the final value
3835 * in the previous iteration, but this will anyhow be increased
3836 * to the approriate value again by this while loop.
3838 while (ftype == F_DISRES &&
3840 il_nr_thread < idef->il[ftype].nr &&
3841 idef->iparams[idef->il[ftype].iatoms[il_nr_thread]].disres.label ==
3842 idef->iparams[idef->il[ftype].iatoms[il_nr_thread-nat1]].disres.label)
3844 il_nr_thread += nat1;
3847 idef->il_thread_division[ftype*(nthreads+1)+t] = il_nr_thread;
3854 calc_bonded_reduction_mask(const t_idef *idef,
3859 int ftype, nb, nat1, nb0, nb1, i, a;
3863 for (ftype = 0; ftype < F_NRE; ftype++)
3865 if (ftype_is_bonded_potential(ftype))
3867 nb = idef->il[ftype].nr;
3870 nat1 = interaction_function[ftype].nratoms + 1;
3872 /* Divide this interaction equally over the threads.
3873 * This is not stored: should match division in calc_bonds.
3875 nb0 = idef->il_thread_division[ftype*(nt+1)+t];
3876 nb1 = idef->il_thread_division[ftype*(nt+1)+t+1];
3878 for (i = nb0; i < nb1; i += nat1)
3880 for (a = 1; a < nat1; a++)
3882 mask |= (1U << (idef->il[ftype].iatoms[i+a]>>shift));
3892 void setup_bonded_threading(t_forcerec *fr, t_idef *idef)
3894 #define MAX_BLOCK_BITS 32
3898 assert(fr->nthreads >= 1);
3900 /* Divide the bonded interaction over the threads */
3901 divide_bondeds_over_threads(idef, fr->nthreads);
3903 if (fr->nthreads == 1)
3910 /* We divide the force array in a maximum of 32 blocks.
3911 * Minimum force block reduction size is 2^6=64.
3914 while (fr->natoms_force > (int)(MAX_BLOCK_BITS*(1U<<fr->red_ashift)))
3920 fprintf(debug, "bonded force buffer block atom shift %d bits\n",
3924 /* Determine to which blocks each thread's bonded force calculation
3925 * contributes. Store this is a mask for each thread.
3927 #pragma omp parallel for num_threads(fr->nthreads) schedule(static)
3928 for (t = 1; t < fr->nthreads; t++)
3930 fr->f_t[t].red_mask =
3931 calc_bonded_reduction_mask(idef, fr->red_ashift, t, fr->nthreads);
3934 /* Determine the maximum number of blocks we need to reduce over */
3937 for (t = 0; t < fr->nthreads; t++)
3940 for (b = 0; b < MAX_BLOCK_BITS; b++)
3942 if (fr->f_t[t].red_mask & (1U<<b))
3944 fr->red_nblock = max(fr->red_nblock, b+1);
3950 fprintf(debug, "thread %d flags %x count %d\n",
3951 t, fr->f_t[t].red_mask, c);
3957 fprintf(debug, "Number of blocks to reduce: %d of size %d\n",
3958 fr->red_nblock, 1<<fr->red_ashift);
3959 fprintf(debug, "Reduction density %.2f density/#thread %.2f\n",
3960 ctot*(1<<fr->red_ashift)/(double)fr->natoms_force,
3961 ctot*(1<<fr->red_ashift)/(double)(fr->natoms_force*fr->nthreads));
3965 static void zero_thread_forces(f_thread_t *f_t, int n,
3966 int nblock, int blocksize)
3968 int b, a0, a1, a, i, j;
3970 if (n > f_t->f_nalloc)
3972 f_t->f_nalloc = over_alloc_large(n);
3973 srenew(f_t->f, f_t->f_nalloc);
3976 if (f_t->red_mask != 0)
3978 for (b = 0; b < nblock; b++)
3980 if (f_t->red_mask && (1U<<b))
3983 a1 = min((b+1)*blocksize, n);
3984 for (a = a0; a < a1; a++)
3986 clear_rvec(f_t->f[a]);
3991 for (i = 0; i < SHIFTS; i++)
3993 clear_rvec(f_t->fshift[i]);
3995 for (i = 0; i < F_NRE; i++)
3999 for (i = 0; i < egNR; i++)
4001 for (j = 0; j < f_t->grpp.nener; j++)
4003 f_t->grpp.ener[i][j] = 0;
4006 for (i = 0; i < efptNR; i++)
4012 static void reduce_thread_force_buffer(int n, rvec *f,
4013 int nthreads, f_thread_t *f_t,
4014 int nblock, int block_size)
4016 /* The max thread number is arbitrary,
4017 * we used a fixed number to avoid memory management.
4018 * Using more than 16 threads is probably never useful performance wise.
4020 #define MAX_BONDED_THREADS 256
4023 if (nthreads > MAX_BONDED_THREADS)
4025 gmx_fatal(FARGS, "Can not reduce bonded forces on more than %d threads",
4026 MAX_BONDED_THREADS);
4029 /* This reduction can run on any number of threads,
4030 * independently of nthreads.
4032 #pragma omp parallel for num_threads(nthreads) schedule(static)
4033 for (b = 0; b < nblock; b++)
4035 rvec *fp[MAX_BONDED_THREADS];
4039 /* Determine which threads contribute to this block */
4041 for (ft = 1; ft < nthreads; ft++)
4043 if (f_t[ft].red_mask & (1U<<b))
4045 fp[nfb++] = f_t[ft].f;
4050 /* Reduce force buffers for threads that contribute */
4052 a1 = (b+1)*block_size;
4054 for (a = a0; a < a1; a++)
4056 for (fb = 0; fb < nfb; fb++)
4058 rvec_inc(f[a], fp[fb][a]);
4065 static void reduce_thread_forces(int n, rvec *f, rvec *fshift,
4066 real *ener, gmx_grppairener_t *grpp, real *dvdl,
4067 int nthreads, f_thread_t *f_t,
4068 int nblock, int block_size,
4069 gmx_bool bCalcEnerVir,
4074 /* Reduce the bonded force buffer */
4075 reduce_thread_force_buffer(n, f, nthreads, f_t, nblock, block_size);
4078 /* When necessary, reduce energy and virial using one thread only */
4083 for (i = 0; i < SHIFTS; i++)
4085 for (t = 1; t < nthreads; t++)
4087 rvec_inc(fshift[i], f_t[t].fshift[i]);
4090 for (i = 0; i < F_NRE; i++)
4092 for (t = 1; t < nthreads; t++)
4094 ener[i] += f_t[t].ener[i];
4097 for (i = 0; i < egNR; i++)
4099 for (j = 0; j < f_t[1].grpp.nener; j++)
4101 for (t = 1; t < nthreads; t++)
4104 grpp->ener[i][j] += f_t[t].grpp.ener[i][j];
4110 for (i = 0; i < efptNR; i++)
4113 for (t = 1; t < nthreads; t++)
4115 dvdl[i] += f_t[t].dvdl[i];
4122 static real calc_one_bond(FILE *fplog, int thread,
4123 int ftype, const t_idef *idef,
4124 rvec x[], rvec f[], rvec fshift[],
4126 const t_pbc *pbc, const t_graph *g,
4127 gmx_grppairener_t *grpp,
4129 real *lambda, real *dvdl,
4130 const t_mdatoms *md, t_fcdata *fcd,
4131 gmx_bool bCalcEnerVir,
4132 int *global_atom_index, gmx_bool bPrintSepPot)
4134 int nat1, nbonds, efptFTYPE;
4139 if (IS_RESTRAINT_TYPE(ftype))
4141 efptFTYPE = efptRESTRAINT;
4145 efptFTYPE = efptBONDED;
4148 nat1 = interaction_function[ftype].nratoms + 1;
4149 nbonds = idef->il[ftype].nr/nat1;
4150 iatoms = idef->il[ftype].iatoms;
4152 nb0 = idef->il_thread_division[ftype*(idef->nthreads+1)+thread];
4153 nbn = idef->il_thread_division[ftype*(idef->nthreads+1)+thread+1] - nb0;
4155 if (!IS_LISTED_LJ_C(ftype))
4157 if (ftype == F_CMAP)
4159 v = cmap_dihs(nbn, iatoms+nb0,
4160 idef->iparams, &idef->cmap_grid,
4161 (const rvec*)x, f, fshift,
4162 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
4163 md, fcd, global_atom_index);
4166 else if (ftype == F_ANGLES &&
4167 !bCalcEnerVir && fr->efep == efepNO)
4169 /* No energies, shift forces, dvdl */
4170 angles_noener_simd(nbn, idef->il[ftype].iatoms+nb0,
4173 pbc, g, lambda[efptFTYPE], md, fcd,
4178 else if (ftype == F_PDIHS &&
4179 !bCalcEnerVir && fr->efep == efepNO)
4181 /* No energies, shift forces, dvdl */
4182 #ifndef SIMD_BONDEDS
4187 (nbn, idef->il[ftype].iatoms+nb0,
4190 pbc, g, lambda[efptFTYPE], md, fcd,
4196 v = interaction_function[ftype].ifunc(nbn, iatoms+nb0,
4198 (const rvec*)x, f, fshift,
4199 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
4200 md, fcd, global_atom_index);
4204 fprintf(fplog, " %-23s #%4d V %12.5e dVdl %12.5e\n",
4205 interaction_function[ftype].longname,
4206 nbonds, v, lambda[efptFTYPE]);
4211 v = do_nonbonded_listed(ftype, nbn, iatoms+nb0, idef->iparams, (const rvec*)x, f, fshift,
4212 pbc, g, lambda, dvdl, md, fr, grpp, global_atom_index);
4216 fprintf(fplog, " %-5s + %-15s #%4d dVdl %12.5e\n",
4217 interaction_function[ftype].longname,
4218 interaction_function[F_LJ14].longname, nbonds, dvdl[efptVDW]);
4219 fprintf(fplog, " %-5s + %-15s #%4d dVdl %12.5e\n",
4220 interaction_function[ftype].longname,
4221 interaction_function[F_COUL14].longname, nbonds, dvdl[efptCOUL]);
4227 inc_nrnb(nrnb, interaction_function[ftype].nrnb_ind, nbonds);
4233 void calc_bonds(FILE *fplog, const gmx_multisim_t *ms,
4235 rvec x[], history_t *hist,
4236 rvec f[], t_forcerec *fr,
4237 const t_pbc *pbc, const t_graph *g,
4238 gmx_enerdata_t *enerd, t_nrnb *nrnb,
4240 const t_mdatoms *md,
4241 t_fcdata *fcd, int *global_atom_index,
4242 t_atomtypes gmx_unused *atype, gmx_genborn_t gmx_unused *born,
4244 gmx_bool bPrintSepPot, gmx_int64_t step)
4246 gmx_bool bCalcEnerVir;
4248 real v, dvdl[efptNR], dvdl_dum[efptNR]; /* The dummy array is to have a place to store the dhdl at other values
4249 of lambda, which will be thrown away in the end*/
4250 const t_pbc *pbc_null;
4254 assert(fr->nthreads == idef->nthreads);
4256 bCalcEnerVir = (force_flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY));
4258 for (i = 0; i < efptNR; i++)
4272 fprintf(fplog, "Step %s: bonded V and dVdl for this node\n",
4273 gmx_step_str(step, buf));
4279 p_graph(debug, "Bondage is fun", g);
4283 /* Do pre force calculation stuff which might require communication */
4284 if (idef->il[F_ORIRES].nr)
4286 enerd->term[F_ORIRESDEV] =
4287 calc_orires_dev(ms, idef->il[F_ORIRES].nr,
4288 idef->il[F_ORIRES].iatoms,
4289 idef->iparams, md, (const rvec*)x,
4290 pbc_null, fcd, hist);
4292 if (idef->il[F_DISRES].nr)
4294 calc_disres_R_6(idef->il[F_DISRES].nr,
4295 idef->il[F_DISRES].iatoms,
4296 idef->iparams, (const rvec*)x, pbc_null,
4299 if (fcd->disres.nsystems > 1)
4301 gmx_sum_sim(2*fcd->disres.nres, fcd->disres.Rt_6, ms);
4306 #pragma omp parallel for num_threads(fr->nthreads) schedule(static)
4307 for (thread = 0; thread < fr->nthreads; thread++)
4314 gmx_grppairener_t *grpp;
4319 fshift = fr->fshift;
4321 grpp = &enerd->grpp;
4326 zero_thread_forces(&fr->f_t[thread], fr->natoms_force,
4327 fr->red_nblock, 1<<fr->red_ashift);
4329 ft = fr->f_t[thread].f;
4330 fshift = fr->f_t[thread].fshift;
4331 epot = fr->f_t[thread].ener;
4332 grpp = &fr->f_t[thread].grpp;
4333 dvdlt = fr->f_t[thread].dvdl;
4335 /* Loop over all bonded force types to calculate the bonded forces */
4336 for (ftype = 0; (ftype < F_NRE); ftype++)
4338 if (idef->il[ftype].nr > 0 && ftype_is_bonded_potential(ftype))
4340 v = calc_one_bond(fplog, thread, ftype, idef, x,
4341 ft, fshift, fr, pbc_null, g, grpp,
4342 nrnb, lambda, dvdlt,
4343 md, fcd, bCalcEnerVir,
4344 global_atom_index, bPrintSepPot);
4349 if (fr->nthreads > 1)
4351 reduce_thread_forces(fr->natoms_force, f, fr->fshift,
4352 enerd->term, &enerd->grpp, dvdl,
4353 fr->nthreads, fr->f_t,
4354 fr->red_nblock, 1<<fr->red_ashift,
4356 force_flags & GMX_FORCE_DHDL);
4358 if (force_flags & GMX_FORCE_DHDL)
4360 for (i = 0; i < efptNR; i++)
4362 enerd->dvdl_nonlin[i] += dvdl[i];
4366 /* Copy the sum of violations for the distance restraints from fcd */
4369 enerd->term[F_DISRESVIOL] = fcd->disres.sumviol;
4374 void calc_bonds_lambda(FILE *fplog,
4378 const t_pbc *pbc, const t_graph *g,
4379 gmx_grppairener_t *grpp, real *epot, t_nrnb *nrnb,
4381 const t_mdatoms *md,
4383 int *global_atom_index)
4385 int i, ftype, nr_nonperturbed, nr;
4387 real dvdl_dum[efptNR];
4389 const t_pbc *pbc_null;
4401 /* Copy the whole idef, so we can modify the contents locally */
4403 idef_fe.nthreads = 1;
4404 snew(idef_fe.il_thread_division, F_NRE*(idef_fe.nthreads+1));
4406 /* We already have the forces, so we use temp buffers here */
4407 snew(f, fr->natoms_force);
4408 snew(fshift, SHIFTS);
4410 /* Loop over all bonded force types to calculate the bonded energies */
4411 for (ftype = 0; (ftype < F_NRE); ftype++)
4413 if (ftype_is_bonded_potential(ftype))
4415 /* Set the work range of thread 0 to the perturbed bondeds only */
4416 nr_nonperturbed = idef->il[ftype].nr_nonperturbed;
4417 nr = idef->il[ftype].nr;
4418 idef_fe.il_thread_division[ftype*2+0] = nr_nonperturbed;
4419 idef_fe.il_thread_division[ftype*2+1] = nr;
4421 /* This is only to get the flop count correct */
4422 idef_fe.il[ftype].nr = nr - nr_nonperturbed;
4424 if (nr - nr_nonperturbed > 0)
4426 v = calc_one_bond(fplog, 0, ftype, &idef_fe,
4427 x, f, fshift, fr, pbc_null, g,
4428 grpp, nrnb, lambda, dvdl_dum,
4430 global_atom_index, FALSE);
4439 sfree(idef_fe.il_thread_division);