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"
62 #include "gromacs/simd/simd.h"
63 #include "gromacs/simd/simd_math.h"
64 #include "gromacs/simd/vector_operations.h"
66 /* Find a better place for this? */
67 const int cmap_coeff_matrix[] = {
68 1, 0, -3, 2, 0, 0, 0, 0, -3, 0, 9, -6, 2, 0, -6, 4,
69 0, 0, 0, 0, 0, 0, 0, 0, 3, 0, -9, 6, -2, 0, 6, -4,
70 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 9, -6, 0, 0, -6, 4,
71 0, 0, 3, -2, 0, 0, 0, 0, 0, 0, -9, 6, 0, 0, 6, -4,
72 0, 0, 0, 0, 1, 0, -3, 2, -2, 0, 6, -4, 1, 0, -3, 2,
73 0, 0, 0, 0, 0, 0, 0, 0, -1, 0, 3, -2, 1, 0, -3, 2,
74 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -3, 2, 0, 0, 3, -2,
75 0, 0, 0, 0, 0, 0, 3, -2, 0, 0, -6, 4, 0, 0, 3, -2,
76 0, 1, -2, 1, 0, 0, 0, 0, 0, -3, 6, -3, 0, 2, -4, 2,
77 0, 0, 0, 0, 0, 0, 0, 0, 0, 3, -6, 3, 0, -2, 4, -2,
78 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -3, 3, 0, 0, 2, -2,
79 0, 0, -1, 1, 0, 0, 0, 0, 0, 0, 3, -3, 0, 0, -2, 2,
80 0, 0, 0, 0, 0, 1, -2, 1, 0, -2, 4, -2, 0, 1, -2, 1,
81 0, 0, 0, 0, 0, 0, 0, 0, 0, -1, 2, -1, 0, 1, -2, 1,
82 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, -1, 0, 0, -1, 1,
83 0, 0, 0, 0, 0, 0, -1, 1, 0, 0, 2, -2, 0, 0, -1, 1
88 int glatnr(int *global_atom_index, int i)
92 if (global_atom_index == NULL)
98 atnr = global_atom_index[i] + 1;
104 static int pbc_rvec_sub(const t_pbc *pbc, const rvec xi, const rvec xj, rvec dx)
108 return pbc_dx_aiuc(pbc, xi, xj, dx);
112 rvec_sub(xi, xj, dx);
117 #ifdef GMX_SIMD_HAVE_REAL
119 /* SIMD PBC data structure, containing 1/boxdiag and the box vectors */
121 gmx_simd_real_t inv_bzz;
122 gmx_simd_real_t inv_byy;
123 gmx_simd_real_t inv_bxx;
132 /* Set the SIMD pbc data from a normal t_pbc struct */
133 static void set_pbc_simd(const t_pbc *pbc, pbc_simd_t *pbc_simd)
138 /* Setting inv_bdiag to 0 effectively turns off PBC */
139 clear_rvec(inv_bdiag);
142 for (d = 0; d < pbc->ndim_ePBC; d++)
144 inv_bdiag[d] = 1.0/pbc->box[d][d];
148 pbc_simd->inv_bzz = gmx_simd_set1_r(inv_bdiag[ZZ]);
149 pbc_simd->inv_byy = gmx_simd_set1_r(inv_bdiag[YY]);
150 pbc_simd->inv_bxx = gmx_simd_set1_r(inv_bdiag[XX]);
154 pbc_simd->bzx = gmx_simd_set1_r(pbc->box[ZZ][XX]);
155 pbc_simd->bzy = gmx_simd_set1_r(pbc->box[ZZ][YY]);
156 pbc_simd->bzz = gmx_simd_set1_r(pbc->box[ZZ][ZZ]);
157 pbc_simd->byx = gmx_simd_set1_r(pbc->box[YY][XX]);
158 pbc_simd->byy = gmx_simd_set1_r(pbc->box[YY][YY]);
159 pbc_simd->bxx = gmx_simd_set1_r(pbc->box[XX][XX]);
163 pbc_simd->bzx = gmx_simd_setzero_r();
164 pbc_simd->bzy = gmx_simd_setzero_r();
165 pbc_simd->bzz = gmx_simd_setzero_r();
166 pbc_simd->byx = gmx_simd_setzero_r();
167 pbc_simd->byy = gmx_simd_setzero_r();
168 pbc_simd->bxx = gmx_simd_setzero_r();
172 /* Correct distance vector *dx,*dy,*dz for PBC using SIMD */
173 static gmx_inline void
174 pbc_dx_simd(gmx_simd_real_t *dx, gmx_simd_real_t *dy, gmx_simd_real_t *dz,
175 const pbc_simd_t *pbc)
179 sh = gmx_simd_round_r(gmx_simd_mul_r(*dz, pbc->inv_bzz));
180 *dx = gmx_simd_fnmadd_r(sh, pbc->bzx, *dx);
181 *dy = gmx_simd_fnmadd_r(sh, pbc->bzy, *dy);
182 *dz = gmx_simd_fnmadd_r(sh, pbc->bzz, *dz);
184 sh = gmx_simd_round_r(gmx_simd_mul_r(*dy, pbc->inv_byy));
185 *dx = gmx_simd_fnmadd_r(sh, pbc->byx, *dx);
186 *dy = gmx_simd_fnmadd_r(sh, pbc->byy, *dy);
188 sh = gmx_simd_round_r(gmx_simd_mul_r(*dx, pbc->inv_bxx));
189 *dx = gmx_simd_fnmadd_r(sh, pbc->bxx, *dx);
192 #endif /* GMX_SIMD_HAVE_REAL */
195 * Morse potential bond by Frank Everdij
197 * Three parameters needed:
199 * b0 = equilibrium distance in nm
200 * be = beta in nm^-1 (actually, it's nu_e*Sqrt(2*pi*pi*mu/D_e))
201 * cb = well depth in kJ/mol
203 * Note: the potential is referenced to be +cb at infinite separation
204 * and zero at the equilibrium distance!
207 real morse_bonds(int nbonds,
208 const t_iatom forceatoms[], const t_iparams forceparams[],
209 const rvec x[], rvec f[], rvec fshift[],
210 const t_pbc *pbc, const t_graph *g,
211 real lambda, real *dvdlambda,
212 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
213 int gmx_unused *global_atom_index)
215 const real one = 1.0;
216 const real two = 2.0;
217 real dr, dr2, temp, omtemp, cbomtemp, fbond, vbond, fij, vtot;
218 real b0, be, cb, b0A, beA, cbA, b0B, beB, cbB, L1;
220 int i, m, ki, type, ai, aj;
224 for (i = 0; (i < nbonds); )
226 type = forceatoms[i++];
227 ai = forceatoms[i++];
228 aj = forceatoms[i++];
230 b0A = forceparams[type].morse.b0A;
231 beA = forceparams[type].morse.betaA;
232 cbA = forceparams[type].morse.cbA;
234 b0B = forceparams[type].morse.b0B;
235 beB = forceparams[type].morse.betaB;
236 cbB = forceparams[type].morse.cbB;
238 L1 = one-lambda; /* 1 */
239 b0 = L1*b0A + lambda*b0B; /* 3 */
240 be = L1*beA + lambda*beB; /* 3 */
241 cb = L1*cbA + lambda*cbB; /* 3 */
243 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
244 dr2 = iprod(dx, dx); /* 5 */
245 dr = dr2*gmx_invsqrt(dr2); /* 10 */
246 temp = exp(-be*(dr-b0)); /* 12 */
250 /* bonds are constrainted. This may _not_ include bond constraints if they are lambda dependent */
251 *dvdlambda += cbB-cbA;
255 omtemp = one-temp; /* 1 */
256 cbomtemp = cb*omtemp; /* 1 */
257 vbond = cbomtemp*omtemp; /* 1 */
258 fbond = -two*be*temp*cbomtemp*gmx_invsqrt(dr2); /* 9 */
259 vtot += vbond; /* 1 */
261 *dvdlambda += (cbB - cbA) * omtemp * omtemp - (2-2*omtemp)*omtemp * cb * ((b0B-b0A)*be - (beB-beA)*(dr-b0)); /* 15 */
265 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
269 for (m = 0; (m < DIM); m++) /* 15 */
274 fshift[ki][m] += fij;
275 fshift[CENTRAL][m] -= fij;
281 real cubic_bonds(int nbonds,
282 const t_iatom forceatoms[], const t_iparams forceparams[],
283 const rvec x[], rvec f[], rvec fshift[],
284 const t_pbc *pbc, const t_graph *g,
285 real gmx_unused lambda, real gmx_unused *dvdlambda,
286 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
287 int gmx_unused *global_atom_index)
289 const real three = 3.0;
290 const real two = 2.0;
292 real dr, dr2, dist, kdist, kdist2, fbond, vbond, fij, vtot;
294 int i, m, ki, type, ai, aj;
298 for (i = 0; (i < nbonds); )
300 type = forceatoms[i++];
301 ai = forceatoms[i++];
302 aj = forceatoms[i++];
304 b0 = forceparams[type].cubic.b0;
305 kb = forceparams[type].cubic.kb;
306 kcub = forceparams[type].cubic.kcub;
308 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
309 dr2 = iprod(dx, dx); /* 5 */
316 dr = dr2*gmx_invsqrt(dr2); /* 10 */
321 vbond = kdist2 + kcub*kdist2*dist;
322 fbond = -(two*kdist + three*kdist2*kcub)/dr;
324 vtot += vbond; /* 21 */
328 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
331 for (m = 0; (m < DIM); m++) /* 15 */
336 fshift[ki][m] += fij;
337 fshift[CENTRAL][m] -= fij;
343 real FENE_bonds(int nbonds,
344 const t_iatom forceatoms[], const t_iparams forceparams[],
345 const rvec x[], rvec f[], rvec fshift[],
346 const t_pbc *pbc, const t_graph *g,
347 real gmx_unused lambda, real gmx_unused *dvdlambda,
348 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
349 int *global_atom_index)
351 const real half = 0.5;
352 const real one = 1.0;
354 real dr, dr2, bm2, omdr2obm2, fbond, vbond, fij, vtot;
356 int i, m, ki, type, ai, aj;
360 for (i = 0; (i < nbonds); )
362 type = forceatoms[i++];
363 ai = forceatoms[i++];
364 aj = forceatoms[i++];
366 bm = forceparams[type].fene.bm;
367 kb = forceparams[type].fene.kb;
369 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
370 dr2 = iprod(dx, dx); /* 5 */
382 "r^2 (%f) >= bm^2 (%f) in FENE bond between atoms %d and %d",
384 glatnr(global_atom_index, ai),
385 glatnr(global_atom_index, aj));
388 omdr2obm2 = one - dr2/bm2;
390 vbond = -half*kb*bm2*log(omdr2obm2);
391 fbond = -kb/omdr2obm2;
393 vtot += vbond; /* 35 */
397 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
400 for (m = 0; (m < DIM); m++) /* 15 */
405 fshift[ki][m] += fij;
406 fshift[CENTRAL][m] -= fij;
412 real harmonic(real kA, real kB, real xA, real xB, real x, real lambda,
415 const real half = 0.5;
416 real L1, kk, x0, dx, dx2;
417 real v, f, dvdlambda;
420 kk = L1*kA+lambda*kB;
421 x0 = L1*xA+lambda*xB;
428 dvdlambda = half*(kB-kA)*dx2 + (xA-xB)*kk*dx;
435 /* That was 19 flops */
439 real bonds(int nbonds,
440 const t_iatom forceatoms[], const t_iparams forceparams[],
441 const rvec x[], rvec f[], rvec fshift[],
442 const t_pbc *pbc, const t_graph *g,
443 real lambda, real *dvdlambda,
444 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
445 int gmx_unused *global_atom_index)
447 int i, m, ki, ai, aj, type;
448 real dr, dr2, fbond, vbond, fij, vtot;
453 for (i = 0; (i < nbonds); )
455 type = forceatoms[i++];
456 ai = forceatoms[i++];
457 aj = forceatoms[i++];
459 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
460 dr2 = iprod(dx, dx); /* 5 */
461 dr = dr2*gmx_invsqrt(dr2); /* 10 */
463 *dvdlambda += harmonic(forceparams[type].harmonic.krA,
464 forceparams[type].harmonic.krB,
465 forceparams[type].harmonic.rA,
466 forceparams[type].harmonic.rB,
467 dr, lambda, &vbond, &fbond); /* 19 */
475 vtot += vbond; /* 1*/
476 fbond *= gmx_invsqrt(dr2); /* 6 */
480 fprintf(debug, "BONDS: dr = %10g vbond = %10g fbond = %10g\n",
486 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
489 for (m = 0; (m < DIM); m++) /* 15 */
494 fshift[ki][m] += fij;
495 fshift[CENTRAL][m] -= fij;
501 real restraint_bonds(int nbonds,
502 const t_iatom forceatoms[], const t_iparams forceparams[],
503 const rvec x[], rvec f[], rvec fshift[],
504 const t_pbc *pbc, const t_graph *g,
505 real lambda, real *dvdlambda,
506 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
507 int gmx_unused *global_atom_index)
509 int i, m, ki, ai, aj, type;
510 real dr, dr2, fbond, vbond, fij, vtot;
512 real low, dlow, up1, dup1, up2, dup2, k, dk;
520 for (i = 0; (i < nbonds); )
522 type = forceatoms[i++];
523 ai = forceatoms[i++];
524 aj = forceatoms[i++];
526 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
527 dr2 = iprod(dx, dx); /* 5 */
528 dr = dr2*gmx_invsqrt(dr2); /* 10 */
530 low = L1*forceparams[type].restraint.lowA + lambda*forceparams[type].restraint.lowB;
531 dlow = -forceparams[type].restraint.lowA + forceparams[type].restraint.lowB;
532 up1 = L1*forceparams[type].restraint.up1A + lambda*forceparams[type].restraint.up1B;
533 dup1 = -forceparams[type].restraint.up1A + forceparams[type].restraint.up1B;
534 up2 = L1*forceparams[type].restraint.up2A + lambda*forceparams[type].restraint.up2B;
535 dup2 = -forceparams[type].restraint.up2A + forceparams[type].restraint.up2B;
536 k = L1*forceparams[type].restraint.kA + lambda*forceparams[type].restraint.kB;
537 dk = -forceparams[type].restraint.kA + forceparams[type].restraint.kB;
546 *dvdlambda += 0.5*dk*drh2 - k*dlow*drh;
559 *dvdlambda += 0.5*dk*drh2 - k*dup1*drh;
564 vbond = k*(up2 - up1)*(0.5*(up2 - up1) + drh);
565 fbond = -k*(up2 - up1);
566 *dvdlambda += dk*(up2 - up1)*(0.5*(up2 - up1) + drh)
567 + k*(dup2 - dup1)*(up2 - up1 + drh)
568 - k*(up2 - up1)*dup2;
576 vtot += vbond; /* 1*/
577 fbond *= gmx_invsqrt(dr2); /* 6 */
581 fprintf(debug, "BONDS: dr = %10g vbond = %10g fbond = %10g\n",
587 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
590 for (m = 0; (m < DIM); m++) /* 15 */
595 fshift[ki][m] += fij;
596 fshift[CENTRAL][m] -= fij;
603 real polarize(int nbonds,
604 const t_iatom forceatoms[], const t_iparams forceparams[],
605 const rvec x[], rvec f[], rvec fshift[],
606 const t_pbc *pbc, const t_graph *g,
607 real lambda, real *dvdlambda,
608 const t_mdatoms *md, t_fcdata gmx_unused *fcd,
609 int gmx_unused *global_atom_index)
611 int i, m, ki, ai, aj, type;
612 real dr, dr2, fbond, vbond, fij, vtot, ksh;
617 for (i = 0; (i < nbonds); )
619 type = forceatoms[i++];
620 ai = forceatoms[i++];
621 aj = forceatoms[i++];
622 ksh = sqr(md->chargeA[aj])*ONE_4PI_EPS0/forceparams[type].polarize.alpha;
625 fprintf(debug, "POL: local ai = %d aj = %d ksh = %.3f\n", ai, aj, ksh);
628 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
629 dr2 = iprod(dx, dx); /* 5 */
630 dr = dr2*gmx_invsqrt(dr2); /* 10 */
632 *dvdlambda += harmonic(ksh, ksh, 0, 0, dr, lambda, &vbond, &fbond); /* 19 */
639 vtot += vbond; /* 1*/
640 fbond *= gmx_invsqrt(dr2); /* 6 */
644 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
647 for (m = 0; (m < DIM); m++) /* 15 */
652 fshift[ki][m] += fij;
653 fshift[CENTRAL][m] -= fij;
659 real anharm_polarize(int nbonds,
660 const t_iatom forceatoms[], const t_iparams forceparams[],
661 const rvec x[], rvec f[], rvec fshift[],
662 const t_pbc *pbc, const t_graph *g,
663 real lambda, real *dvdlambda,
664 const t_mdatoms *md, t_fcdata gmx_unused *fcd,
665 int gmx_unused *global_atom_index)
667 int i, m, ki, ai, aj, type;
668 real dr, dr2, fbond, vbond, fij, vtot, ksh, khyp, drcut, ddr, ddr3;
673 for (i = 0; (i < nbonds); )
675 type = forceatoms[i++];
676 ai = forceatoms[i++];
677 aj = forceatoms[i++];
678 ksh = sqr(md->chargeA[aj])*ONE_4PI_EPS0/forceparams[type].anharm_polarize.alpha; /* 7*/
679 khyp = forceparams[type].anharm_polarize.khyp;
680 drcut = forceparams[type].anharm_polarize.drcut;
683 fprintf(debug, "POL: local ai = %d aj = %d ksh = %.3f\n", ai, aj, ksh);
686 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
687 dr2 = iprod(dx, dx); /* 5 */
688 dr = dr2*gmx_invsqrt(dr2); /* 10 */
690 *dvdlambda += harmonic(ksh, ksh, 0, 0, dr, lambda, &vbond, &fbond); /* 19 */
701 vbond += khyp*ddr*ddr3;
702 fbond -= 4*khyp*ddr3;
704 fbond *= gmx_invsqrt(dr2); /* 6 */
705 vtot += vbond; /* 1*/
709 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
712 for (m = 0; (m < DIM); m++) /* 15 */
717 fshift[ki][m] += fij;
718 fshift[CENTRAL][m] -= fij;
724 real water_pol(int nbonds,
725 const t_iatom forceatoms[], const t_iparams forceparams[],
726 const rvec x[], rvec f[], rvec gmx_unused fshift[],
727 const t_pbc gmx_unused *pbc, const t_graph gmx_unused *g,
728 real gmx_unused lambda, real gmx_unused *dvdlambda,
729 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
730 int gmx_unused *global_atom_index)
732 /* This routine implements anisotropic polarizibility for water, through
733 * a shell connected to a dummy with spring constant that differ in the
734 * three spatial dimensions in the molecular frame.
736 int i, m, aO, aH1, aH2, aD, aS, type, type0;
737 rvec dOH1, dOH2, dHH, dOD, dDS, nW, kk, dx, kdx, proj;
741 real vtot, fij, r_HH, r_OD, r_nW, tx, ty, tz, qS;
746 type0 = forceatoms[0];
748 qS = md->chargeA[aS];
749 kk[XX] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_x;
750 kk[YY] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_y;
751 kk[ZZ] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_z;
752 r_HH = 1.0/forceparams[type0].wpol.rHH;
753 r_OD = 1.0/forceparams[type0].wpol.rOD;
756 fprintf(debug, "WPOL: qS = %10.5f aS = %5d\n", qS, aS);
757 fprintf(debug, "WPOL: kk = %10.3f %10.3f %10.3f\n",
758 kk[XX], kk[YY], kk[ZZ]);
759 fprintf(debug, "WPOL: rOH = %10.3f rHH = %10.3f rOD = %10.3f\n",
760 forceparams[type0].wpol.rOH,
761 forceparams[type0].wpol.rHH,
762 forceparams[type0].wpol.rOD);
764 for (i = 0; (i < nbonds); i += 6)
766 type = forceatoms[i];
769 gmx_fatal(FARGS, "Sorry, type = %d, type0 = %d, file = %s, line = %d",
770 type, type0, __FILE__, __LINE__);
772 aO = forceatoms[i+1];
773 aH1 = forceatoms[i+2];
774 aH2 = forceatoms[i+3];
775 aD = forceatoms[i+4];
776 aS = forceatoms[i+5];
778 /* Compute vectors describing the water frame */
779 rvec_sub(x[aH1], x[aO], dOH1);
780 rvec_sub(x[aH2], x[aO], dOH2);
781 rvec_sub(x[aH2], x[aH1], dHH);
782 rvec_sub(x[aD], x[aO], dOD);
783 rvec_sub(x[aS], x[aD], dDS);
784 cprod(dOH1, dOH2, nW);
786 /* Compute inverse length of normal vector
787 * (this one could be precomputed, but I'm too lazy now)
789 r_nW = gmx_invsqrt(iprod(nW, nW));
790 /* This is for precision, but does not make a big difference,
793 r_OD = gmx_invsqrt(iprod(dOD, dOD));
795 /* Normalize the vectors in the water frame */
797 svmul(r_HH, dHH, dHH);
798 svmul(r_OD, dOD, dOD);
800 /* Compute displacement of shell along components of the vector */
801 dx[ZZ] = iprod(dDS, dOD);
802 /* Compute projection on the XY plane: dDS - dx[ZZ]*dOD */
803 for (m = 0; (m < DIM); m++)
805 proj[m] = dDS[m]-dx[ZZ]*dOD[m];
808 /*dx[XX] = iprod(dDS,nW);
809 dx[YY] = iprod(dDS,dHH);*/
810 dx[XX] = iprod(proj, nW);
811 for (m = 0; (m < DIM); m++)
813 proj[m] -= dx[XX]*nW[m];
815 dx[YY] = iprod(proj, dHH);
820 fprintf(debug, "WPOL: dx2=%10g dy2=%10g dz2=%10g sum=%10g dDS^2=%10g\n",
821 sqr(dx[XX]), sqr(dx[YY]), sqr(dx[ZZ]), iprod(dx, dx), iprod(dDS, dDS));
822 fprintf(debug, "WPOL: dHH=(%10g,%10g,%10g)\n", dHH[XX], dHH[YY], dHH[ZZ]);
823 fprintf(debug, "WPOL: dOD=(%10g,%10g,%10g), 1/r_OD = %10g\n",
824 dOD[XX], dOD[YY], dOD[ZZ], 1/r_OD);
825 fprintf(debug, "WPOL: nW =(%10g,%10g,%10g), 1/r_nW = %10g\n",
826 nW[XX], nW[YY], nW[ZZ], 1/r_nW);
827 fprintf(debug, "WPOL: dx =%10g, dy =%10g, dz =%10g\n",
828 dx[XX], dx[YY], dx[ZZ]);
829 fprintf(debug, "WPOL: dDSx=%10g, dDSy=%10g, dDSz=%10g\n",
830 dDS[XX], dDS[YY], dDS[ZZ]);
833 /* Now compute the forces and energy */
834 kdx[XX] = kk[XX]*dx[XX];
835 kdx[YY] = kk[YY]*dx[YY];
836 kdx[ZZ] = kk[ZZ]*dx[ZZ];
837 vtot += iprod(dx, kdx);
838 for (m = 0; (m < DIM); m++)
840 /* This is a tensor operation but written out for speed */
854 fprintf(debug, "WPOL: vwpol=%g\n", 0.5*iprod(dx, kdx));
855 fprintf(debug, "WPOL: df = (%10g, %10g, %10g)\n", df[XX], df[YY], df[ZZ]);
863 static real do_1_thole(const rvec xi, const rvec xj, rvec fi, rvec fj,
864 const t_pbc *pbc, real qq,
865 rvec fshift[], real afac)
868 real r12sq, r12_1, r12n, r12bar, v0, v1, fscal, ebar, fff;
871 t = pbc_rvec_sub(pbc, xi, xj, r12); /* 3 */
873 r12sq = iprod(r12, r12); /* 5 */
874 r12_1 = gmx_invsqrt(r12sq); /* 5 */
875 r12bar = afac/r12_1; /* 5 */
876 v0 = qq*ONE_4PI_EPS0*r12_1; /* 2 */
877 ebar = exp(-r12bar); /* 5 */
878 v1 = (1-(1+0.5*r12bar)*ebar); /* 4 */
879 fscal = ((v0*r12_1)*v1 - v0*0.5*afac*ebar*(r12bar+1))*r12_1; /* 9 */
882 fprintf(debug, "THOLE: v0 = %.3f v1 = %.3f r12= % .3f r12bar = %.3f fscal = %.3f ebar = %.3f\n", v0, v1, 1/r12_1, r12bar, fscal, ebar);
885 for (m = 0; (m < DIM); m++)
891 fshift[CENTRAL][m] -= fff;
894 return v0*v1; /* 1 */
898 real thole_pol(int nbonds,
899 const t_iatom forceatoms[], const t_iparams forceparams[],
900 const rvec x[], rvec f[], rvec fshift[],
901 const t_pbc *pbc, const t_graph gmx_unused *g,
902 real gmx_unused lambda, real gmx_unused *dvdlambda,
903 const t_mdatoms *md, t_fcdata gmx_unused *fcd,
904 int gmx_unused *global_atom_index)
906 /* Interaction between two pairs of particles with opposite charge */
907 int i, type, a1, da1, a2, da2;
908 real q1, q2, qq, a, al1, al2, afac;
911 for (i = 0; (i < nbonds); )
913 type = forceatoms[i++];
914 a1 = forceatoms[i++];
915 da1 = forceatoms[i++];
916 a2 = forceatoms[i++];
917 da2 = forceatoms[i++];
918 q1 = md->chargeA[da1];
919 q2 = md->chargeA[da2];
920 a = forceparams[type].thole.a;
921 al1 = forceparams[type].thole.alpha1;
922 al2 = forceparams[type].thole.alpha2;
924 afac = a*pow(al1*al2, -1.0/6.0);
925 V += do_1_thole(x[a1], x[a2], f[a1], f[a2], pbc, qq, fshift, afac);
926 V += do_1_thole(x[da1], x[a2], f[da1], f[a2], pbc, -qq, fshift, afac);
927 V += do_1_thole(x[a1], x[da2], f[a1], f[da2], pbc, -qq, fshift, afac);
928 V += do_1_thole(x[da1], x[da2], f[da1], f[da2], pbc, qq, fshift, afac);
934 real bond_angle(const rvec xi, const rvec xj, const rvec xk, const t_pbc *pbc,
935 rvec r_ij, rvec r_kj, real *costh,
937 /* Return value is the angle between the bonds i-j and j-k */
942 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
943 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
945 *costh = cos_angle(r_ij, r_kj); /* 25 */
946 th = acos(*costh); /* 10 */
951 real angles(int nbonds,
952 const t_iatom forceatoms[], const t_iparams forceparams[],
953 const rvec x[], rvec f[], rvec fshift[],
954 const t_pbc *pbc, const t_graph *g,
955 real lambda, real *dvdlambda,
956 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
957 int gmx_unused *global_atom_index)
959 int i, ai, aj, ak, t1, t2, type;
961 real cos_theta, cos_theta2, theta, dVdt, va, vtot;
962 ivec jt, dt_ij, dt_kj;
965 for (i = 0; i < nbonds; )
967 type = forceatoms[i++];
968 ai = forceatoms[i++];
969 aj = forceatoms[i++];
970 ak = forceatoms[i++];
972 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
973 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
975 *dvdlambda += harmonic(forceparams[type].harmonic.krA,
976 forceparams[type].harmonic.krB,
977 forceparams[type].harmonic.rA*DEG2RAD,
978 forceparams[type].harmonic.rB*DEG2RAD,
979 theta, lambda, &va, &dVdt); /* 21 */
982 cos_theta2 = sqr(cos_theta);
992 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
993 sth = st*cos_theta; /* 1 */
997 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
998 theta*RAD2DEG, va, dVdt);
1001 nrij2 = iprod(r_ij, r_ij); /* 5 */
1002 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1004 nrij_1 = gmx_invsqrt(nrij2); /* 10 */
1005 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1007 cik = st*nrij_1*nrkj_1; /* 2 */
1008 cii = sth*nrij_1*nrij_1; /* 2 */
1009 ckk = sth*nrkj_1*nrkj_1; /* 2 */
1011 for (m = 0; m < DIM; m++)
1013 f_i[m] = -(cik*r_kj[m] - cii*r_ij[m]);
1014 f_k[m] = -(cik*r_ij[m] - ckk*r_kj[m]);
1015 f_j[m] = -f_i[m] - f_k[m];
1022 copy_ivec(SHIFT_IVEC(g, aj), jt);
1024 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1025 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1026 t1 = IVEC2IS(dt_ij);
1027 t2 = IVEC2IS(dt_kj);
1029 rvec_inc(fshift[t1], f_i);
1030 rvec_inc(fshift[CENTRAL], f_j);
1031 rvec_inc(fshift[t2], f_k);
1038 #ifdef GMX_SIMD_HAVE_REAL
1040 /* As angles, but using SIMD to calculate many dihedrals at once.
1041 * This routines does not calculate energies and shift forces.
1043 static gmx_inline void
1044 angles_noener_simd(int nbonds,
1045 const t_iatom forceatoms[], const t_iparams forceparams[],
1046 const rvec x[], rvec f[],
1047 const t_pbc *pbc, const t_graph gmx_unused *g,
1048 real gmx_unused lambda,
1049 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1050 int gmx_unused *global_atom_index)
1054 int type, ai[GMX_SIMD_REAL_WIDTH], aj[GMX_SIMD_REAL_WIDTH];
1055 int ak[GMX_SIMD_REAL_WIDTH];
1056 real coeff_array[2*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *coeff;
1057 real dr_array[2*DIM*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *dr;
1058 real f_buf_array[6*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *f_buf;
1059 gmx_simd_real_t k_S, theta0_S;
1060 gmx_simd_real_t rijx_S, rijy_S, rijz_S;
1061 gmx_simd_real_t rkjx_S, rkjy_S, rkjz_S;
1062 gmx_simd_real_t one_S;
1063 gmx_simd_real_t min_one_plus_eps_S;
1064 gmx_simd_real_t rij_rkj_S;
1065 gmx_simd_real_t nrij2_S, nrij_1_S;
1066 gmx_simd_real_t nrkj2_S, nrkj_1_S;
1067 gmx_simd_real_t cos_S, invsin_S;
1068 gmx_simd_real_t theta_S;
1069 gmx_simd_real_t st_S, sth_S;
1070 gmx_simd_real_t cik_S, cii_S, ckk_S;
1071 gmx_simd_real_t f_ix_S, f_iy_S, f_iz_S;
1072 gmx_simd_real_t f_kx_S, f_ky_S, f_kz_S;
1073 pbc_simd_t pbc_simd;
1075 /* Ensure register memory alignment */
1076 coeff = gmx_simd_align_r(coeff_array);
1077 dr = gmx_simd_align_r(dr_array);
1078 f_buf = gmx_simd_align_r(f_buf_array);
1080 set_pbc_simd(pbc, &pbc_simd);
1082 one_S = gmx_simd_set1_r(1.0);
1084 /* The smallest number > -1 */
1085 min_one_plus_eps_S = gmx_simd_set1_r(-1.0 + 2*GMX_REAL_EPS);
1087 /* nbonds is the number of angles times nfa1, here we step GMX_SIMD_REAL_WIDTH angles */
1088 for (i = 0; (i < nbonds); i += GMX_SIMD_REAL_WIDTH*nfa1)
1090 /* Collect atoms for GMX_SIMD_REAL_WIDTH angles.
1091 * iu indexes into forceatoms, we should not let iu go beyond nbonds.
1094 for (s = 0; s < GMX_SIMD_REAL_WIDTH; s++)
1096 type = forceatoms[iu];
1097 ai[s] = forceatoms[iu+1];
1098 aj[s] = forceatoms[iu+2];
1099 ak[s] = forceatoms[iu+3];
1101 coeff[s] = forceparams[type].harmonic.krA;
1102 coeff[GMX_SIMD_REAL_WIDTH+s] = forceparams[type].harmonic.rA*DEG2RAD;
1104 /* If you can't use pbc_dx_simd below for PBC, e.g. because
1105 * you can't round in SIMD, use pbc_rvec_sub here.
1107 /* Store the non PBC corrected distances packed and aligned */
1108 for (m = 0; m < DIM; m++)
1110 dr[s + m *GMX_SIMD_REAL_WIDTH] = x[ai[s]][m] - x[aj[s]][m];
1111 dr[s + (DIM+m)*GMX_SIMD_REAL_WIDTH] = x[ak[s]][m] - x[aj[s]][m];
1114 /* At the end fill the arrays with identical entries */
1115 if (iu + nfa1 < nbonds)
1121 k_S = gmx_simd_load_r(coeff);
1122 theta0_S = gmx_simd_load_r(coeff+GMX_SIMD_REAL_WIDTH);
1124 rijx_S = gmx_simd_load_r(dr + 0*GMX_SIMD_REAL_WIDTH);
1125 rijy_S = gmx_simd_load_r(dr + 1*GMX_SIMD_REAL_WIDTH);
1126 rijz_S = gmx_simd_load_r(dr + 2*GMX_SIMD_REAL_WIDTH);
1127 rkjx_S = gmx_simd_load_r(dr + 3*GMX_SIMD_REAL_WIDTH);
1128 rkjy_S = gmx_simd_load_r(dr + 4*GMX_SIMD_REAL_WIDTH);
1129 rkjz_S = gmx_simd_load_r(dr + 5*GMX_SIMD_REAL_WIDTH);
1131 pbc_dx_simd(&rijx_S, &rijy_S, &rijz_S, &pbc_simd);
1132 pbc_dx_simd(&rkjx_S, &rkjy_S, &rkjz_S, &pbc_simd);
1134 rij_rkj_S = gmx_simd_iprod_r(rijx_S, rijy_S, rijz_S,
1135 rkjx_S, rkjy_S, rkjz_S);
1137 nrij2_S = gmx_simd_norm2_r(rijx_S, rijy_S, rijz_S);
1138 nrkj2_S = gmx_simd_norm2_r(rkjx_S, rkjy_S, rkjz_S);
1140 nrij_1_S = gmx_simd_invsqrt_r(nrij2_S);
1141 nrkj_1_S = gmx_simd_invsqrt_r(nrkj2_S);
1143 cos_S = gmx_simd_mul_r(rij_rkj_S, gmx_simd_mul_r(nrij_1_S, nrkj_1_S));
1145 /* To allow for 180 degrees, we take the max of cos and -1 + 1bit,
1146 * so we can safely get the 1/sin from 1/sqrt(1 - cos^2).
1147 * This also ensures that rounding errors would cause the argument
1148 * of gmx_simd_acos_r to be < -1.
1149 * Note that we do not take precautions for cos(0)=1, so the outer
1150 * atoms in an angle should not be on top of each other.
1152 cos_S = gmx_simd_max_r(cos_S, min_one_plus_eps_S);
1154 theta_S = gmx_simd_acos_r(cos_S);
1156 invsin_S = gmx_simd_invsqrt_r(gmx_simd_sub_r(one_S, gmx_simd_mul_r(cos_S, cos_S)));
1158 st_S = gmx_simd_mul_r(gmx_simd_mul_r(k_S, gmx_simd_sub_r(theta0_S, theta_S)),
1160 sth_S = gmx_simd_mul_r(st_S, cos_S);
1162 cik_S = gmx_simd_mul_r(st_S, gmx_simd_mul_r(nrij_1_S, nrkj_1_S));
1163 cii_S = gmx_simd_mul_r(sth_S, gmx_simd_mul_r(nrij_1_S, nrij_1_S));
1164 ckk_S = gmx_simd_mul_r(sth_S, gmx_simd_mul_r(nrkj_1_S, nrkj_1_S));
1166 f_ix_S = gmx_simd_mul_r(cii_S, rijx_S);
1167 f_ix_S = gmx_simd_fnmadd_r(cik_S, rkjx_S, f_ix_S);
1168 f_iy_S = gmx_simd_mul_r(cii_S, rijy_S);
1169 f_iy_S = gmx_simd_fnmadd_r(cik_S, rkjy_S, f_iy_S);
1170 f_iz_S = gmx_simd_mul_r(cii_S, rijz_S);
1171 f_iz_S = gmx_simd_fnmadd_r(cik_S, rkjz_S, f_iz_S);
1172 f_kx_S = gmx_simd_mul_r(ckk_S, rkjx_S);
1173 f_kx_S = gmx_simd_fnmadd_r(cik_S, rijx_S, f_kx_S);
1174 f_ky_S = gmx_simd_mul_r(ckk_S, rkjy_S);
1175 f_ky_S = gmx_simd_fnmadd_r(cik_S, rijy_S, f_ky_S);
1176 f_kz_S = gmx_simd_mul_r(ckk_S, rkjz_S);
1177 f_kz_S = gmx_simd_fnmadd_r(cik_S, rijz_S, f_kz_S);
1179 gmx_simd_store_r(f_buf + 0*GMX_SIMD_REAL_WIDTH, f_ix_S);
1180 gmx_simd_store_r(f_buf + 1*GMX_SIMD_REAL_WIDTH, f_iy_S);
1181 gmx_simd_store_r(f_buf + 2*GMX_SIMD_REAL_WIDTH, f_iz_S);
1182 gmx_simd_store_r(f_buf + 3*GMX_SIMD_REAL_WIDTH, f_kx_S);
1183 gmx_simd_store_r(f_buf + 4*GMX_SIMD_REAL_WIDTH, f_ky_S);
1184 gmx_simd_store_r(f_buf + 5*GMX_SIMD_REAL_WIDTH, f_kz_S);
1190 for (m = 0; m < DIM; m++)
1192 f[ai[s]][m] += f_buf[s + m*GMX_SIMD_REAL_WIDTH];
1193 f[aj[s]][m] -= f_buf[s + m*GMX_SIMD_REAL_WIDTH] + f_buf[s + (DIM+m)*GMX_SIMD_REAL_WIDTH];
1194 f[ak[s]][m] += f_buf[s + (DIM+m)*GMX_SIMD_REAL_WIDTH];
1199 while (s < GMX_SIMD_REAL_WIDTH && iu < nbonds);
1203 #endif /* GMX_SIMD_HAVE_REAL */
1205 real linear_angles(int nbonds,
1206 const t_iatom forceatoms[], const t_iparams forceparams[],
1207 const rvec x[], rvec f[], rvec fshift[],
1208 const t_pbc *pbc, const t_graph *g,
1209 real lambda, real *dvdlambda,
1210 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1211 int gmx_unused *global_atom_index)
1213 int i, m, ai, aj, ak, t1, t2, type;
1215 real L1, kA, kB, aA, aB, dr, dr2, va, vtot, a, b, klin;
1216 ivec jt, dt_ij, dt_kj;
1217 rvec r_ij, r_kj, r_ik, dx;
1221 for (i = 0; (i < nbonds); )
1223 type = forceatoms[i++];
1224 ai = forceatoms[i++];
1225 aj = forceatoms[i++];
1226 ak = forceatoms[i++];
1228 kA = forceparams[type].linangle.klinA;
1229 kB = forceparams[type].linangle.klinB;
1230 klin = L1*kA + lambda*kB;
1232 aA = forceparams[type].linangle.aA;
1233 aB = forceparams[type].linangle.aB;
1234 a = L1*aA+lambda*aB;
1237 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
1238 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
1239 rvec_sub(r_ij, r_kj, r_ik);
1242 for (m = 0; (m < DIM); m++)
1244 dr = -a * r_ij[m] - b * r_kj[m];
1249 f_j[m] = -(f_i[m]+f_k[m]);
1255 *dvdlambda += 0.5*(kB-kA)*dr2 + klin*(aB-aA)*iprod(dx, r_ik);
1261 copy_ivec(SHIFT_IVEC(g, aj), jt);
1263 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1264 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1265 t1 = IVEC2IS(dt_ij);
1266 t2 = IVEC2IS(dt_kj);
1268 rvec_inc(fshift[t1], f_i);
1269 rvec_inc(fshift[CENTRAL], f_j);
1270 rvec_inc(fshift[t2], f_k);
1275 real urey_bradley(int nbonds,
1276 const t_iatom forceatoms[], const t_iparams forceparams[],
1277 const rvec x[], rvec f[], rvec fshift[],
1278 const t_pbc *pbc, const t_graph *g,
1279 real lambda, real *dvdlambda,
1280 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1281 int gmx_unused *global_atom_index)
1283 int i, m, ai, aj, ak, t1, t2, type, ki;
1284 rvec r_ij, r_kj, r_ik;
1285 real cos_theta, cos_theta2, theta;
1286 real dVdt, va, vtot, dr, dr2, vbond, fbond, fik;
1287 real kthA, th0A, kUBA, r13A, kthB, th0B, kUBB, r13B;
1288 ivec jt, dt_ij, dt_kj, dt_ik;
1291 for (i = 0; (i < nbonds); )
1293 type = forceatoms[i++];
1294 ai = forceatoms[i++];
1295 aj = forceatoms[i++];
1296 ak = forceatoms[i++];
1297 th0A = forceparams[type].u_b.thetaA*DEG2RAD;
1298 kthA = forceparams[type].u_b.kthetaA;
1299 r13A = forceparams[type].u_b.r13A;
1300 kUBA = forceparams[type].u_b.kUBA;
1301 th0B = forceparams[type].u_b.thetaB*DEG2RAD;
1302 kthB = forceparams[type].u_b.kthetaB;
1303 r13B = forceparams[type].u_b.r13B;
1304 kUBB = forceparams[type].u_b.kUBB;
1306 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
1307 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
1309 *dvdlambda += harmonic(kthA, kthB, th0A, th0B, theta, lambda, &va, &dVdt); /* 21 */
1312 ki = pbc_rvec_sub(pbc, x[ai], x[ak], r_ik); /* 3 */
1313 dr2 = iprod(r_ik, r_ik); /* 5 */
1314 dr = dr2*gmx_invsqrt(dr2); /* 10 */
1316 *dvdlambda += harmonic(kUBA, kUBB, r13A, r13B, dr, lambda, &vbond, &fbond); /* 19 */
1318 cos_theta2 = sqr(cos_theta); /* 1 */
1326 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1327 sth = st*cos_theta; /* 1 */
1331 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1332 theta*RAD2DEG, va, dVdt);
1335 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1336 nrij2 = iprod(r_ij, r_ij);
1338 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
1339 cii = sth/nrij2; /* 10 */
1340 ckk = sth/nrkj2; /* 10 */
1342 for (m = 0; (m < DIM); m++) /* 39 */
1344 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
1345 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
1346 f_j[m] = -f_i[m]-f_k[m];
1353 copy_ivec(SHIFT_IVEC(g, aj), jt);
1355 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1356 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1357 t1 = IVEC2IS(dt_ij);
1358 t2 = IVEC2IS(dt_kj);
1360 rvec_inc(fshift[t1], f_i);
1361 rvec_inc(fshift[CENTRAL], f_j);
1362 rvec_inc(fshift[t2], f_k);
1364 /* Time for the bond calculations */
1370 vtot += vbond; /* 1*/
1371 fbond *= gmx_invsqrt(dr2); /* 6 */
1375 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ak), dt_ik);
1376 ki = IVEC2IS(dt_ik);
1378 for (m = 0; (m < DIM); m++) /* 15 */
1380 fik = fbond*r_ik[m];
1383 fshift[ki][m] += fik;
1384 fshift[CENTRAL][m] -= fik;
1390 real quartic_angles(int nbonds,
1391 const t_iatom forceatoms[], const t_iparams forceparams[],
1392 const rvec x[], rvec f[], rvec fshift[],
1393 const t_pbc *pbc, const t_graph *g,
1394 real gmx_unused lambda, real gmx_unused *dvdlambda,
1395 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1396 int gmx_unused *global_atom_index)
1398 int i, j, ai, aj, ak, t1, t2, type;
1400 real cos_theta, cos_theta2, theta, dt, dVdt, va, dtp, c, vtot;
1401 ivec jt, dt_ij, dt_kj;
1404 for (i = 0; (i < nbonds); )
1406 type = forceatoms[i++];
1407 ai = forceatoms[i++];
1408 aj = forceatoms[i++];
1409 ak = forceatoms[i++];
1411 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
1412 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
1414 dt = theta - forceparams[type].qangle.theta*DEG2RAD; /* 2 */
1417 va = forceparams[type].qangle.c[0];
1419 for (j = 1; j <= 4; j++)
1421 c = forceparams[type].qangle.c[j];
1430 cos_theta2 = sqr(cos_theta); /* 1 */
1439 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1440 sth = st*cos_theta; /* 1 */
1444 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1445 theta*RAD2DEG, va, dVdt);
1448 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1449 nrij2 = iprod(r_ij, r_ij);
1451 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
1452 cii = sth/nrij2; /* 10 */
1453 ckk = sth/nrkj2; /* 10 */
1455 for (m = 0; (m < DIM); m++) /* 39 */
1457 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
1458 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
1459 f_j[m] = -f_i[m]-f_k[m];
1466 copy_ivec(SHIFT_IVEC(g, aj), jt);
1468 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1469 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1470 t1 = IVEC2IS(dt_ij);
1471 t2 = IVEC2IS(dt_kj);
1473 rvec_inc(fshift[t1], f_i);
1474 rvec_inc(fshift[CENTRAL], f_j);
1475 rvec_inc(fshift[t2], f_k);
1481 real dih_angle(const rvec xi, const rvec xj, const rvec xk, const rvec xl,
1483 rvec r_ij, rvec r_kj, rvec r_kl, rvec m, rvec n,
1484 real *sign, int *t1, int *t2, int *t3)
1488 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
1489 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
1490 *t3 = pbc_rvec_sub(pbc, xk, xl, r_kl); /* 3 */
1492 cprod(r_ij, r_kj, m); /* 9 */
1493 cprod(r_kj, r_kl, n); /* 9 */
1494 phi = gmx_angle(m, n); /* 49 (assuming 25 for atan2) */
1495 ipr = iprod(r_ij, n); /* 5 */
1496 (*sign) = (ipr < 0.0) ? -1.0 : 1.0;
1497 phi = (*sign)*phi; /* 1 */
1503 #ifdef GMX_SIMD_HAVE_REAL
1505 /* As dih_angle above, but calculates 4 dihedral angles at once using SIMD,
1506 * also calculates the pre-factor required for the dihedral force update.
1507 * Note that bv and buf should be register aligned.
1509 static gmx_inline void
1510 dih_angle_simd(const rvec *x,
1511 const int *ai, const int *aj, const int *ak, const int *al,
1512 const pbc_simd_t *pbc,
1514 gmx_simd_real_t *phi_S,
1515 gmx_simd_real_t *mx_S, gmx_simd_real_t *my_S, gmx_simd_real_t *mz_S,
1516 gmx_simd_real_t *nx_S, gmx_simd_real_t *ny_S, gmx_simd_real_t *nz_S,
1517 gmx_simd_real_t *nrkj_m2_S,
1518 gmx_simd_real_t *nrkj_n2_S,
1523 gmx_simd_real_t rijx_S, rijy_S, rijz_S;
1524 gmx_simd_real_t rkjx_S, rkjy_S, rkjz_S;
1525 gmx_simd_real_t rklx_S, rkly_S, rklz_S;
1526 gmx_simd_real_t cx_S, cy_S, cz_S;
1527 gmx_simd_real_t cn_S;
1528 gmx_simd_real_t s_S;
1529 gmx_simd_real_t ipr_S;
1530 gmx_simd_real_t iprm_S, iprn_S;
1531 gmx_simd_real_t nrkj2_S, nrkj_1_S, nrkj_2_S, nrkj_S;
1532 gmx_simd_real_t toler_S;
1533 gmx_simd_real_t p_S, q_S;
1534 gmx_simd_real_t nrkj2_min_S;
1535 gmx_simd_real_t real_eps_S;
1537 /* Used to avoid division by zero.
1538 * We take into acount that we multiply the result by real_eps_S.
1540 nrkj2_min_S = gmx_simd_set1_r(GMX_REAL_MIN/(2*GMX_REAL_EPS));
1542 /* The value of the last significant bit (GMX_REAL_EPS is half of that) */
1543 real_eps_S = gmx_simd_set1_r(2*GMX_REAL_EPS);
1545 for (s = 0; s < GMX_SIMD_REAL_WIDTH; s++)
1547 /* If you can't use pbc_dx_simd below for PBC, e.g. because
1548 * you can't round in SIMD, use pbc_rvec_sub here.
1550 for (m = 0; m < DIM; m++)
1552 dr[s + (0*DIM + m)*GMX_SIMD_REAL_WIDTH] = x[ai[s]][m] - x[aj[s]][m];
1553 dr[s + (1*DIM + m)*GMX_SIMD_REAL_WIDTH] = x[ak[s]][m] - x[aj[s]][m];
1554 dr[s + (2*DIM + m)*GMX_SIMD_REAL_WIDTH] = x[ak[s]][m] - x[al[s]][m];
1558 rijx_S = gmx_simd_load_r(dr + 0*GMX_SIMD_REAL_WIDTH);
1559 rijy_S = gmx_simd_load_r(dr + 1*GMX_SIMD_REAL_WIDTH);
1560 rijz_S = gmx_simd_load_r(dr + 2*GMX_SIMD_REAL_WIDTH);
1561 rkjx_S = gmx_simd_load_r(dr + 3*GMX_SIMD_REAL_WIDTH);
1562 rkjy_S = gmx_simd_load_r(dr + 4*GMX_SIMD_REAL_WIDTH);
1563 rkjz_S = gmx_simd_load_r(dr + 5*GMX_SIMD_REAL_WIDTH);
1564 rklx_S = gmx_simd_load_r(dr + 6*GMX_SIMD_REAL_WIDTH);
1565 rkly_S = gmx_simd_load_r(dr + 7*GMX_SIMD_REAL_WIDTH);
1566 rklz_S = gmx_simd_load_r(dr + 8*GMX_SIMD_REAL_WIDTH);
1568 pbc_dx_simd(&rijx_S, &rijy_S, &rijz_S, pbc);
1569 pbc_dx_simd(&rkjx_S, &rkjy_S, &rkjz_S, pbc);
1570 pbc_dx_simd(&rklx_S, &rkly_S, &rklz_S, pbc);
1572 gmx_simd_cprod_r(rijx_S, rijy_S, rijz_S,
1573 rkjx_S, rkjy_S, rkjz_S,
1576 gmx_simd_cprod_r(rkjx_S, rkjy_S, rkjz_S,
1577 rklx_S, rkly_S, rklz_S,
1580 gmx_simd_cprod_r(*mx_S, *my_S, *mz_S,
1581 *nx_S, *ny_S, *nz_S,
1582 &cx_S, &cy_S, &cz_S);
1584 cn_S = gmx_simd_sqrt_r(gmx_simd_norm2_r(cx_S, cy_S, cz_S));
1586 s_S = gmx_simd_iprod_r(*mx_S, *my_S, *mz_S, *nx_S, *ny_S, *nz_S);
1588 /* Determine the dihedral angle, the sign might need correction */
1589 *phi_S = gmx_simd_atan2_r(cn_S, s_S);
1591 ipr_S = gmx_simd_iprod_r(rijx_S, rijy_S, rijz_S,
1592 *nx_S, *ny_S, *nz_S);
1594 iprm_S = gmx_simd_norm2_r(*mx_S, *my_S, *mz_S);
1595 iprn_S = gmx_simd_norm2_r(*nx_S, *ny_S, *nz_S);
1597 nrkj2_S = gmx_simd_norm2_r(rkjx_S, rkjy_S, rkjz_S);
1599 /* Avoid division by zero. When zero, the result is multiplied by 0
1600 * anyhow, so the 3 max below do not affect the final result.
1602 nrkj2_S = gmx_simd_max_r(nrkj2_S, nrkj2_min_S);
1603 nrkj_1_S = gmx_simd_invsqrt_r(nrkj2_S);
1604 nrkj_2_S = gmx_simd_mul_r(nrkj_1_S, nrkj_1_S);
1605 nrkj_S = gmx_simd_mul_r(nrkj2_S, nrkj_1_S);
1607 toler_S = gmx_simd_mul_r(nrkj2_S, real_eps_S);
1609 /* Here the plain-C code uses a conditional, but we can't do that in SIMD.
1610 * So we take a max with the tolerance instead. Since we multiply with
1611 * m or n later, the max does not affect the results.
1613 iprm_S = gmx_simd_max_r(iprm_S, toler_S);
1614 iprn_S = gmx_simd_max_r(iprn_S, toler_S);
1615 *nrkj_m2_S = gmx_simd_mul_r(nrkj_S, gmx_simd_inv_r(iprm_S));
1616 *nrkj_n2_S = gmx_simd_mul_r(nrkj_S, gmx_simd_inv_r(iprn_S));
1618 /* Set sign of phi_S with the sign of ipr_S; phi_S is currently positive */
1619 *phi_S = gmx_simd_xor_sign_r(*phi_S, ipr_S);
1620 p_S = gmx_simd_iprod_r(rijx_S, rijy_S, rijz_S,
1621 rkjx_S, rkjy_S, rkjz_S);
1622 p_S = gmx_simd_mul_r(p_S, nrkj_2_S);
1624 q_S = gmx_simd_iprod_r(rklx_S, rkly_S, rklz_S,
1625 rkjx_S, rkjy_S, rkjz_S);
1626 q_S = gmx_simd_mul_r(q_S, nrkj_2_S);
1628 gmx_simd_store_r(p, p_S);
1629 gmx_simd_store_r(q, q_S);
1632 #endif /* GMX_SIMD_HAVE_REAL */
1635 void do_dih_fup(int i, int j, int k, int l, real ddphi,
1636 rvec r_ij, rvec r_kj, rvec r_kl,
1637 rvec m, rvec n, rvec f[], rvec fshift[],
1638 const t_pbc *pbc, const t_graph *g,
1639 const rvec x[], int t1, int t2, int t3)
1642 rvec f_i, f_j, f_k, f_l;
1643 rvec uvec, vvec, svec, dx_jl;
1644 real iprm, iprn, nrkj, nrkj2, nrkj_1, nrkj_2;
1645 real a, b, p, q, toler;
1646 ivec jt, dt_ij, dt_kj, dt_lj;
1648 iprm = iprod(m, m); /* 5 */
1649 iprn = iprod(n, n); /* 5 */
1650 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1651 toler = nrkj2*GMX_REAL_EPS;
1652 if ((iprm > toler) && (iprn > toler))
1654 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1655 nrkj_2 = nrkj_1*nrkj_1; /* 1 */
1656 nrkj = nrkj2*nrkj_1; /* 1 */
1657 a = -ddphi*nrkj/iprm; /* 11 */
1658 svmul(a, m, f_i); /* 3 */
1659 b = ddphi*nrkj/iprn; /* 11 */
1660 svmul(b, n, f_l); /* 3 */
1661 p = iprod(r_ij, r_kj); /* 5 */
1662 p *= nrkj_2; /* 1 */
1663 q = iprod(r_kl, r_kj); /* 5 */
1664 q *= nrkj_2; /* 1 */
1665 svmul(p, f_i, uvec); /* 3 */
1666 svmul(q, f_l, vvec); /* 3 */
1667 rvec_sub(uvec, vvec, svec); /* 3 */
1668 rvec_sub(f_i, svec, f_j); /* 3 */
1669 rvec_add(f_l, svec, f_k); /* 3 */
1670 rvec_inc(f[i], f_i); /* 3 */
1671 rvec_dec(f[j], f_j); /* 3 */
1672 rvec_dec(f[k], f_k); /* 3 */
1673 rvec_inc(f[l], f_l); /* 3 */
1677 copy_ivec(SHIFT_IVEC(g, j), jt);
1678 ivec_sub(SHIFT_IVEC(g, i), jt, dt_ij);
1679 ivec_sub(SHIFT_IVEC(g, k), jt, dt_kj);
1680 ivec_sub(SHIFT_IVEC(g, l), jt, dt_lj);
1681 t1 = IVEC2IS(dt_ij);
1682 t2 = IVEC2IS(dt_kj);
1683 t3 = IVEC2IS(dt_lj);
1687 t3 = pbc_rvec_sub(pbc, x[l], x[j], dx_jl);
1694 rvec_inc(fshift[t1], f_i);
1695 rvec_dec(fshift[CENTRAL], f_j);
1696 rvec_dec(fshift[t2], f_k);
1697 rvec_inc(fshift[t3], f_l);
1702 /* As do_dih_fup above, but without shift forces */
1704 do_dih_fup_noshiftf(int i, int j, int k, int l, real ddphi,
1705 rvec r_ij, rvec r_kj, rvec r_kl,
1706 rvec m, rvec n, rvec f[])
1708 rvec f_i, f_j, f_k, f_l;
1709 rvec uvec, vvec, svec, dx_jl;
1710 real iprm, iprn, nrkj, nrkj2, nrkj_1, nrkj_2;
1711 real a, b, p, q, toler;
1712 ivec jt, dt_ij, dt_kj, dt_lj;
1714 iprm = iprod(m, m); /* 5 */
1715 iprn = iprod(n, n); /* 5 */
1716 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1717 toler = nrkj2*GMX_REAL_EPS;
1718 if ((iprm > toler) && (iprn > toler))
1720 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1721 nrkj_2 = nrkj_1*nrkj_1; /* 1 */
1722 nrkj = nrkj2*nrkj_1; /* 1 */
1723 a = -ddphi*nrkj/iprm; /* 11 */
1724 svmul(a, m, f_i); /* 3 */
1725 b = ddphi*nrkj/iprn; /* 11 */
1726 svmul(b, n, f_l); /* 3 */
1727 p = iprod(r_ij, r_kj); /* 5 */
1728 p *= nrkj_2; /* 1 */
1729 q = iprod(r_kl, r_kj); /* 5 */
1730 q *= nrkj_2; /* 1 */
1731 svmul(p, f_i, uvec); /* 3 */
1732 svmul(q, f_l, vvec); /* 3 */
1733 rvec_sub(uvec, vvec, svec); /* 3 */
1734 rvec_sub(f_i, svec, f_j); /* 3 */
1735 rvec_add(f_l, svec, f_k); /* 3 */
1736 rvec_inc(f[i], f_i); /* 3 */
1737 rvec_dec(f[j], f_j); /* 3 */
1738 rvec_dec(f[k], f_k); /* 3 */
1739 rvec_inc(f[l], f_l); /* 3 */
1743 /* As do_dih_fup_noshiftf above, but with pre-calculated pre-factors */
1744 static gmx_inline void
1745 do_dih_fup_noshiftf_precalc(int i, int j, int k, int l,
1747 real f_i_x, real f_i_y, real f_i_z,
1748 real mf_l_x, real mf_l_y, real mf_l_z,
1751 rvec f_i, f_j, f_k, f_l;
1752 rvec uvec, vvec, svec;
1760 svmul(p, f_i, uvec);
1761 svmul(q, f_l, vvec);
1762 rvec_sub(uvec, vvec, svec);
1763 rvec_sub(f_i, svec, f_j);
1764 rvec_add(f_l, svec, f_k);
1765 rvec_inc(f[i], f_i);
1766 rvec_dec(f[j], f_j);
1767 rvec_dec(f[k], f_k);
1768 rvec_inc(f[l], f_l);
1772 real dopdihs(real cpA, real cpB, real phiA, real phiB, int mult,
1773 real phi, real lambda, real *V, real *F)
1775 real v, dvdlambda, mdphi, v1, sdphi, ddphi;
1776 real L1 = 1.0 - lambda;
1777 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1778 real dph0 = (phiB - phiA)*DEG2RAD;
1779 real cp = L1*cpA + lambda*cpB;
1781 mdphi = mult*phi - ph0;
1783 ddphi = -cp*mult*sdphi;
1784 v1 = 1.0 + cos(mdphi);
1787 dvdlambda = (cpB - cpA)*v1 + cp*dph0*sdphi;
1794 /* That was 40 flops */
1798 dopdihs_noener(real cpA, real cpB, real phiA, real phiB, int mult,
1799 real phi, real lambda, real *F)
1801 real mdphi, sdphi, ddphi;
1802 real L1 = 1.0 - lambda;
1803 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1804 real cp = L1*cpA + lambda*cpB;
1806 mdphi = mult*phi - ph0;
1808 ddphi = -cp*mult*sdphi;
1812 /* That was 20 flops */
1816 dopdihs_mdphi(real cpA, real cpB, real phiA, real phiB, int mult,
1817 real phi, real lambda, real *cp, real *mdphi)
1819 real L1 = 1.0 - lambda;
1820 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1822 *cp = L1*cpA + lambda*cpB;
1824 *mdphi = mult*phi - ph0;
1827 static real dopdihs_min(real cpA, real cpB, real phiA, real phiB, int mult,
1828 real phi, real lambda, real *V, real *F)
1829 /* similar to dopdihs, except for a minus sign *
1830 * and a different treatment of mult/phi0 */
1832 real v, dvdlambda, mdphi, v1, sdphi, ddphi;
1833 real L1 = 1.0 - lambda;
1834 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1835 real dph0 = (phiB - phiA)*DEG2RAD;
1836 real cp = L1*cpA + lambda*cpB;
1838 mdphi = mult*(phi-ph0);
1840 ddphi = cp*mult*sdphi;
1841 v1 = 1.0-cos(mdphi);
1844 dvdlambda = (cpB-cpA)*v1 + cp*dph0*sdphi;
1851 /* That was 40 flops */
1854 real pdihs(int nbonds,
1855 const t_iatom forceatoms[], const t_iparams forceparams[],
1856 const rvec x[], rvec f[], rvec fshift[],
1857 const t_pbc *pbc, const t_graph *g,
1858 real lambda, real *dvdlambda,
1859 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1860 int gmx_unused *global_atom_index)
1862 int i, type, ai, aj, ak, al;
1864 rvec r_ij, r_kj, r_kl, m, n;
1865 real phi, sign, ddphi, vpd, vtot;
1869 for (i = 0; (i < nbonds); )
1871 type = forceatoms[i++];
1872 ai = forceatoms[i++];
1873 aj = forceatoms[i++];
1874 ak = forceatoms[i++];
1875 al = forceatoms[i++];
1877 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
1878 &sign, &t1, &t2, &t3); /* 84 */
1879 *dvdlambda += dopdihs(forceparams[type].pdihs.cpA,
1880 forceparams[type].pdihs.cpB,
1881 forceparams[type].pdihs.phiA,
1882 forceparams[type].pdihs.phiB,
1883 forceparams[type].pdihs.mult,
1884 phi, lambda, &vpd, &ddphi);
1887 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
1888 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
1891 fprintf(debug, "pdih: (%d,%d,%d,%d) phi=%g\n",
1892 ai, aj, ak, al, phi);
1899 void make_dp_periodic(real *dp) /* 1 flop? */
1901 /* dp cannot be outside (-pi,pi) */
1906 else if (*dp < -M_PI)
1913 /* As pdihs above, but without calculating energies and shift forces */
1915 pdihs_noener(int nbonds,
1916 const t_iatom forceatoms[], const t_iparams forceparams[],
1917 const rvec x[], rvec f[],
1918 const t_pbc gmx_unused *pbc, const t_graph gmx_unused *g,
1920 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1921 int gmx_unused *global_atom_index)
1923 int i, type, ai, aj, ak, al;
1925 rvec r_ij, r_kj, r_kl, m, n;
1926 real phi, sign, ddphi_tot, ddphi;
1928 for (i = 0; (i < nbonds); )
1930 ai = forceatoms[i+1];
1931 aj = forceatoms[i+2];
1932 ak = forceatoms[i+3];
1933 al = forceatoms[i+4];
1935 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
1936 &sign, &t1, &t2, &t3);
1940 /* Loop over dihedrals working on the same atoms,
1941 * so we avoid recalculating angles and force distributions.
1945 type = forceatoms[i];
1946 dopdihs_noener(forceparams[type].pdihs.cpA,
1947 forceparams[type].pdihs.cpB,
1948 forceparams[type].pdihs.phiA,
1949 forceparams[type].pdihs.phiB,
1950 forceparams[type].pdihs.mult,
1951 phi, lambda, &ddphi);
1956 while (i < nbonds &&
1957 forceatoms[i+1] == ai &&
1958 forceatoms[i+2] == aj &&
1959 forceatoms[i+3] == ak &&
1960 forceatoms[i+4] == al);
1962 do_dih_fup_noshiftf(ai, aj, ak, al, ddphi_tot, r_ij, r_kj, r_kl, m, n, f);
1967 #ifdef GMX_SIMD_HAVE_REAL
1969 /* As pdihs_noner above, but using SIMD to calculate many dihedrals at once */
1971 pdihs_noener_simd(int nbonds,
1972 const t_iatom forceatoms[], const t_iparams forceparams[],
1973 const rvec x[], rvec f[],
1974 const t_pbc *pbc, const t_graph gmx_unused *g,
1975 real gmx_unused lambda,
1976 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1977 int gmx_unused *global_atom_index)
1981 int type, ai[GMX_SIMD_REAL_WIDTH], aj[GMX_SIMD_REAL_WIDTH], ak[GMX_SIMD_REAL_WIDTH], al[GMX_SIMD_REAL_WIDTH];
1982 int t1[GMX_SIMD_REAL_WIDTH], t2[GMX_SIMD_REAL_WIDTH], t3[GMX_SIMD_REAL_WIDTH];
1984 real dr_array[3*DIM*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *dr;
1985 real buf_array[7*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *buf;
1986 real *cp, *phi0, *mult, *phi, *p, *q, *sf_i, *msf_l;
1987 gmx_simd_real_t phi0_S, phi_S;
1988 gmx_simd_real_t mx_S, my_S, mz_S;
1989 gmx_simd_real_t nx_S, ny_S, nz_S;
1990 gmx_simd_real_t nrkj_m2_S, nrkj_n2_S;
1991 gmx_simd_real_t cp_S, mdphi_S, mult_S;
1992 gmx_simd_real_t sin_S, cos_S;
1993 gmx_simd_real_t mddphi_S;
1994 gmx_simd_real_t sf_i_S, msf_l_S;
1995 pbc_simd_t pbc_simd;
1997 /* Ensure SIMD register alignment */
1998 dr = gmx_simd_align_r(dr_array);
1999 buf = gmx_simd_align_r(buf_array);
2001 /* Extract aligned pointer for parameters and variables */
2002 cp = buf + 0*GMX_SIMD_REAL_WIDTH;
2003 phi0 = buf + 1*GMX_SIMD_REAL_WIDTH;
2004 mult = buf + 2*GMX_SIMD_REAL_WIDTH;
2005 p = buf + 3*GMX_SIMD_REAL_WIDTH;
2006 q = buf + 4*GMX_SIMD_REAL_WIDTH;
2007 sf_i = buf + 5*GMX_SIMD_REAL_WIDTH;
2008 msf_l = buf + 6*GMX_SIMD_REAL_WIDTH;
2010 set_pbc_simd(pbc, &pbc_simd);
2012 /* nbonds is the number of dihedrals times nfa1, here we step GMX_SIMD_REAL_WIDTH dihs */
2013 for (i = 0; (i < nbonds); i += GMX_SIMD_REAL_WIDTH*nfa1)
2015 /* Collect atoms quadruplets for GMX_SIMD_REAL_WIDTH dihedrals.
2016 * iu indexes into forceatoms, we should not let iu go beyond nbonds.
2019 for (s = 0; s < GMX_SIMD_REAL_WIDTH; s++)
2021 type = forceatoms[iu];
2022 ai[s] = forceatoms[iu+1];
2023 aj[s] = forceatoms[iu+2];
2024 ak[s] = forceatoms[iu+3];
2025 al[s] = forceatoms[iu+4];
2027 cp[s] = forceparams[type].pdihs.cpA;
2028 phi0[s] = forceparams[type].pdihs.phiA*DEG2RAD;
2029 mult[s] = forceparams[type].pdihs.mult;
2031 /* At the end fill the arrays with identical entries */
2032 if (iu + nfa1 < nbonds)
2038 /* Caclulate GMX_SIMD_REAL_WIDTH dihedral angles at once */
2039 dih_angle_simd(x, ai, aj, ak, al, &pbc_simd,
2042 &mx_S, &my_S, &mz_S,
2043 &nx_S, &ny_S, &nz_S,
2048 cp_S = gmx_simd_load_r(cp);
2049 phi0_S = gmx_simd_load_r(phi0);
2050 mult_S = gmx_simd_load_r(mult);
2052 mdphi_S = gmx_simd_sub_r(gmx_simd_mul_r(mult_S, phi_S), phi0_S);
2054 /* Calculate GMX_SIMD_REAL_WIDTH sines at once */
2055 gmx_simd_sincos_r(mdphi_S, &sin_S, &cos_S);
2056 mddphi_S = gmx_simd_mul_r(gmx_simd_mul_r(cp_S, mult_S), sin_S);
2057 sf_i_S = gmx_simd_mul_r(mddphi_S, nrkj_m2_S);
2058 msf_l_S = gmx_simd_mul_r(mddphi_S, nrkj_n2_S);
2060 /* After this m?_S will contain f[i] */
2061 mx_S = gmx_simd_mul_r(sf_i_S, mx_S);
2062 my_S = gmx_simd_mul_r(sf_i_S, my_S);
2063 mz_S = gmx_simd_mul_r(sf_i_S, mz_S);
2065 /* After this m?_S will contain -f[l] */
2066 nx_S = gmx_simd_mul_r(msf_l_S, nx_S);
2067 ny_S = gmx_simd_mul_r(msf_l_S, ny_S);
2068 nz_S = gmx_simd_mul_r(msf_l_S, nz_S);
2070 gmx_simd_store_r(dr + 0*GMX_SIMD_REAL_WIDTH, mx_S);
2071 gmx_simd_store_r(dr + 1*GMX_SIMD_REAL_WIDTH, my_S);
2072 gmx_simd_store_r(dr + 2*GMX_SIMD_REAL_WIDTH, mz_S);
2073 gmx_simd_store_r(dr + 3*GMX_SIMD_REAL_WIDTH, nx_S);
2074 gmx_simd_store_r(dr + 4*GMX_SIMD_REAL_WIDTH, ny_S);
2075 gmx_simd_store_r(dr + 5*GMX_SIMD_REAL_WIDTH, nz_S);
2081 do_dih_fup_noshiftf_precalc(ai[s], aj[s], ak[s], al[s],
2083 dr[ XX *GMX_SIMD_REAL_WIDTH+s],
2084 dr[ YY *GMX_SIMD_REAL_WIDTH+s],
2085 dr[ ZZ *GMX_SIMD_REAL_WIDTH+s],
2086 dr[(DIM+XX)*GMX_SIMD_REAL_WIDTH+s],
2087 dr[(DIM+YY)*GMX_SIMD_REAL_WIDTH+s],
2088 dr[(DIM+ZZ)*GMX_SIMD_REAL_WIDTH+s],
2093 while (s < GMX_SIMD_REAL_WIDTH && iu < nbonds);
2097 #endif /* GMX_SIMD_HAVE_REAL */
2100 real idihs(int nbonds,
2101 const t_iatom forceatoms[], const t_iparams forceparams[],
2102 const rvec x[], rvec f[], rvec fshift[],
2103 const t_pbc *pbc, const t_graph *g,
2104 real lambda, real *dvdlambda,
2105 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2106 int gmx_unused *global_atom_index)
2108 int i, type, ai, aj, ak, al;
2110 real phi, phi0, dphi0, ddphi, sign, vtot;
2111 rvec r_ij, r_kj, r_kl, m, n;
2112 real L1, kk, dp, dp2, kA, kB, pA, pB, dvdl_term;
2117 for (i = 0; (i < nbonds); )
2119 type = forceatoms[i++];
2120 ai = forceatoms[i++];
2121 aj = forceatoms[i++];
2122 ak = forceatoms[i++];
2123 al = forceatoms[i++];
2125 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2126 &sign, &t1, &t2, &t3); /* 84 */
2128 /* phi can jump if phi0 is close to Pi/-Pi, which will cause huge
2129 * force changes if we just apply a normal harmonic.
2130 * Instead, we first calculate phi-phi0 and take it modulo (-Pi,Pi).
2131 * This means we will never have the periodicity problem, unless
2132 * the dihedral is Pi away from phiO, which is very unlikely due to
2135 kA = forceparams[type].harmonic.krA;
2136 kB = forceparams[type].harmonic.krB;
2137 pA = forceparams[type].harmonic.rA;
2138 pB = forceparams[type].harmonic.rB;
2140 kk = L1*kA + lambda*kB;
2141 phi0 = (L1*pA + lambda*pB)*DEG2RAD;
2142 dphi0 = (pB - pA)*DEG2RAD;
2146 make_dp_periodic(&dp);
2153 dvdl_term += 0.5*(kB - kA)*dp2 - kk*dphi0*dp;
2155 do_dih_fup(ai, aj, ak, al, (real)(-ddphi), r_ij, r_kj, r_kl, m, n,
2156 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2161 fprintf(debug, "idih: (%d,%d,%d,%d) phi=%g\n",
2162 ai, aj, ak, al, phi);
2167 *dvdlambda += dvdl_term;
2172 /*! \brief returns dx, rdist, and dpdl for functions posres() and fbposres()
2174 static void posres_dx(const rvec x, const rvec pos0A, const rvec pos0B,
2175 const rvec comA_sc, const rvec comB_sc,
2177 t_pbc *pbc, int refcoord_scaling, int npbcdim,
2178 rvec dx, rvec rdist, rvec dpdl)
2181 real posA, posB, L1, ref = 0.;
2186 for (m = 0; m < DIM; m++)
2192 switch (refcoord_scaling)
2196 rdist[m] = L1*posA + lambda*posB;
2197 dpdl[m] = posB - posA;
2200 /* Box relative coordinates are stored for dimensions with pbc */
2201 posA *= pbc->box[m][m];
2202 posB *= pbc->box[m][m];
2203 assert(npbcdim <= DIM);
2204 for (d = m+1; d < npbcdim; d++)
2206 posA += pos0A[d]*pbc->box[d][m];
2207 posB += pos0B[d]*pbc->box[d][m];
2209 ref = L1*posA + lambda*posB;
2211 dpdl[m] = posB - posA;
2214 ref = L1*comA_sc[m] + lambda*comB_sc[m];
2215 rdist[m] = L1*posA + lambda*posB;
2216 dpdl[m] = comB_sc[m] - comA_sc[m] + posB - posA;
2219 gmx_fatal(FARGS, "No such scaling method implemented");
2224 ref = L1*posA + lambda*posB;
2226 dpdl[m] = posB - posA;
2229 /* We do pbc_dx with ref+rdist,
2230 * since with only ref we can be up to half a box vector wrong.
2232 pos[m] = ref + rdist[m];
2237 pbc_dx(pbc, x, pos, dx);
2241 rvec_sub(x, pos, dx);
2245 /*! \brief Adds forces of flat-bottomed positions restraints to f[]
2246 * and fixes vir_diag. Returns the flat-bottomed potential. */
2247 real fbposres(int nbonds,
2248 const t_iatom forceatoms[], const t_iparams forceparams[],
2249 const rvec x[], rvec f[], rvec vir_diag,
2251 int refcoord_scaling, int ePBC, rvec com)
2252 /* compute flat-bottomed positions restraints */
2254 int i, ai, m, d, type, npbcdim = 0, fbdim;
2255 const t_iparams *pr;
2257 real ref = 0, dr, dr2, rpot, rfb, rfb2, fact, invdr;
2258 rvec com_sc, rdist, pos, dx, dpdl, fm;
2261 npbcdim = ePBC2npbcdim(ePBC);
2263 if (refcoord_scaling == erscCOM)
2266 for (m = 0; m < npbcdim; m++)
2268 assert(npbcdim <= DIM);
2269 for (d = m; d < npbcdim; d++)
2271 com_sc[m] += com[d]*pbc->box[d][m];
2277 for (i = 0; (i < nbonds); )
2279 type = forceatoms[i++];
2280 ai = forceatoms[i++];
2281 pr = &forceparams[type];
2283 /* same calculation as for normal posres, but with identical A and B states, and lambda==0 */
2284 posres_dx(x[ai], forceparams[type].fbposres.pos0, forceparams[type].fbposres.pos0,
2285 com_sc, com_sc, 0.0,
2286 pbc, refcoord_scaling, npbcdim,
2292 kk = pr->fbposres.k;
2293 rfb = pr->fbposres.r;
2296 /* with rfb<0, push particle out of the sphere/cylinder/layer */
2304 switch (pr->fbposres.geom)
2306 case efbposresSPHERE:
2307 /* spherical flat-bottom posres */
2310 ( (dr2 > rfb2 && bInvert == FALSE ) || (dr2 < rfb2 && bInvert == TRUE ) )
2314 v = 0.5*kk*sqr(dr - rfb);
2315 fact = -kk*(dr-rfb)/dr; /* Force pointing to the center pos0 */
2316 svmul(fact, dx, fm);
2319 case efbposresCYLINDER:
2320 /* cylidrical flat-bottom posres in x-y plane. fm[ZZ] = 0. */
2321 dr2 = sqr(dx[XX])+sqr(dx[YY]);
2323 ( (dr2 > rfb2 && bInvert == FALSE ) || (dr2 < rfb2 && bInvert == TRUE ) )
2328 v = 0.5*kk*sqr(dr - rfb);
2329 fm[XX] = -kk*(dr-rfb)*dx[XX]*invdr; /* Force pointing to the center */
2330 fm[YY] = -kk*(dr-rfb)*dx[YY]*invdr;
2333 case efbposresX: /* fbdim=XX */
2334 case efbposresY: /* fbdim=YY */
2335 case efbposresZ: /* fbdim=ZZ */
2336 /* 1D flat-bottom potential */
2337 fbdim = pr->fbposres.geom - efbposresX;
2339 if ( ( dr > rfb && bInvert == FALSE ) || ( 0 < dr && dr < rfb && bInvert == TRUE ) )
2341 v = 0.5*kk*sqr(dr - rfb);
2342 fm[fbdim] = -kk*(dr - rfb);
2344 else if ( (dr < (-rfb) && bInvert == FALSE ) || ( (-rfb) < dr && dr < 0 && bInvert == TRUE ))
2346 v = 0.5*kk*sqr(dr + rfb);
2347 fm[fbdim] = -kk*(dr + rfb);
2354 for (m = 0; (m < DIM); m++)
2357 /* Here we correct for the pbc_dx which included rdist */
2358 vir_diag[m] -= 0.5*(dx[m] + rdist[m])*fm[m];
2366 real posres(int nbonds,
2367 const t_iatom forceatoms[], const t_iparams forceparams[],
2368 const rvec x[], rvec f[], rvec vir_diag,
2370 real lambda, real *dvdlambda,
2371 int refcoord_scaling, int ePBC, rvec comA, rvec comB)
2373 int i, ai, m, d, type, ki, npbcdim = 0;
2374 const t_iparams *pr;
2377 real posA, posB, ref = 0;
2378 rvec comA_sc, comB_sc, rdist, dpdl, pos, dx;
2379 gmx_bool bForceValid = TRUE;
2381 if ((f == NULL) || (vir_diag == NULL)) /* should both be null together! */
2383 bForceValid = FALSE;
2386 npbcdim = ePBC2npbcdim(ePBC);
2388 if (refcoord_scaling == erscCOM)
2390 clear_rvec(comA_sc);
2391 clear_rvec(comB_sc);
2392 for (m = 0; m < npbcdim; m++)
2394 assert(npbcdim <= DIM);
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 restrangles(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 gmx_unused lambda, real gmx_unused *dvdlambda,
2673 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2674 int gmx_unused *global_atom_index)
2676 int i, d, ai, aj, ak, type, m;
2680 ivec jt, dt_ij, dt_kj;
2682 real prefactor, ratio_ante, ratio_post;
2683 rvec delta_ante, delta_post, vec_temp;
2686 for (i = 0; (i < nbonds); )
2688 type = forceatoms[i++];
2689 ai = forceatoms[i++];
2690 aj = forceatoms[i++];
2691 ak = forceatoms[i++];
2693 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], vec_temp);
2694 pbc_rvec_sub(pbc, x[aj], x[ai], delta_ante);
2695 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], delta_post);
2698 /* This function computes factors needed for restricted angle potential.
2699 * The restricted angle potential is used in coarse-grained simulations to avoid singularities
2700 * when three particles align and the dihedral angle and dihedral potential
2701 * cannot be calculated. This potential is calculated using the formula:
2702 real restrangles(int nbonds,
2703 const t_iatom forceatoms[],const t_iparams forceparams[],
2704 const rvec x[],rvec f[],rvec fshift[],
2705 const t_pbc *pbc,const t_graph *g,
2706 real gmx_unused lambda,real gmx_unused *dvdlambda,
2707 const t_mdatoms gmx_unused *md,t_fcdata gmx_unused *fcd,
2708 int gmx_unused *global_atom_index)
2710 int i, d, ai, aj, ak, type, m;
2714 ivec jt,dt_ij,dt_kj;
2716 real prefactor, ratio_ante, ratio_post;
2717 rvec delta_ante, delta_post, vec_temp;
2720 for(i=0; (i<nbonds); )
2722 type = forceatoms[i++];
2723 ai = forceatoms[i++];
2724 aj = forceatoms[i++];
2725 ak = forceatoms[i++];
2727 * \f[V_{\rm ReB}(\theta_i) = \frac{1}{2} k_{\theta} \frac{(\cos\theta_i - \cos\theta_0)^2}
2728 * {\sin^2\theta_i}\f] ({eq:ReB} and ref \cite{MonicaGoga2013} from the manual).
2729 * For more explanations see comments file "restcbt.h". */
2731 compute_factors_restangles(type, forceparams, delta_ante, delta_post,
2732 &prefactor, &ratio_ante, &ratio_post, &v);
2734 /* Forces are computed per component */
2735 for (d = 0; d < DIM; d++)
2737 f_i[d] = prefactor * (ratio_ante * delta_ante[d] - delta_post[d]);
2738 f_j[d] = prefactor * ((ratio_post + 1.0) * delta_post[d] - (ratio_ante + 1.0) * delta_ante[d]);
2739 f_k[d] = prefactor * (delta_ante[d] - ratio_post * delta_post[d]);
2742 /* Computation of potential energy */
2748 for (m = 0; (m < DIM); m++)
2757 copy_ivec(SHIFT_IVEC(g, aj), jt);
2758 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
2759 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
2760 t1 = IVEC2IS(dt_ij);
2761 t2 = IVEC2IS(dt_kj);
2764 rvec_inc(fshift[t1], f_i);
2765 rvec_inc(fshift[CENTRAL], f_j);
2766 rvec_inc(fshift[t2], f_k);
2772 real restrdihs(int nbonds,
2773 const t_iatom forceatoms[], const t_iparams forceparams[],
2774 const rvec x[], rvec f[], rvec fshift[],
2775 const t_pbc *pbc, const t_graph *g,
2776 real gmx_unused lambda, real gmx_unused *dvlambda,
2777 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2778 int gmx_unused *global_atom_index)
2780 int i, d, type, ai, aj, ak, al;
2781 rvec f_i, f_j, f_k, f_l;
2783 ivec jt, dt_ij, dt_kj, dt_lj;
2786 rvec delta_ante, delta_crnt, delta_post, vec_temp;
2787 real factor_phi_ai_ante, factor_phi_ai_crnt, factor_phi_ai_post;
2788 real factor_phi_aj_ante, factor_phi_aj_crnt, factor_phi_aj_post;
2789 real factor_phi_ak_ante, factor_phi_ak_crnt, factor_phi_ak_post;
2790 real factor_phi_al_ante, factor_phi_al_crnt, factor_phi_al_post;
2795 for (i = 0; (i < nbonds); )
2797 type = forceatoms[i++];
2798 ai = forceatoms[i++];
2799 aj = forceatoms[i++];
2800 ak = forceatoms[i++];
2801 al = forceatoms[i++];
2803 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], vec_temp);
2804 pbc_rvec_sub(pbc, x[aj], x[ai], delta_ante);
2805 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], delta_crnt);
2806 t3 = pbc_rvec_sub(pbc, x[ak], x[al], vec_temp);
2807 pbc_rvec_sub(pbc, x[al], x[ak], delta_post);
2809 /* This function computes factors needed for restricted angle potential.
2810 * The restricted angle potential is used in coarse-grained simulations to avoid singularities
2811 * when three particles align and the dihedral angle and dihedral potential cannot be calculated.
2812 * This potential is calculated using the formula:
2813 * \f[V_{\rm ReB}(\theta_i) = \frac{1}{2} k_{\theta}
2814 * \frac{(\cos\theta_i - \cos\theta_0)^2}{\sin^2\theta_i}\f]
2815 * ({eq:ReB} and ref \cite{MonicaGoga2013} from the manual).
2816 * For more explanations see comments file "restcbt.h" */
2818 compute_factors_restrdihs(type, forceparams,
2819 delta_ante, delta_crnt, delta_post,
2820 &factor_phi_ai_ante, &factor_phi_ai_crnt, &factor_phi_ai_post,
2821 &factor_phi_aj_ante, &factor_phi_aj_crnt, &factor_phi_aj_post,
2822 &factor_phi_ak_ante, &factor_phi_ak_crnt, &factor_phi_ak_post,
2823 &factor_phi_al_ante, &factor_phi_al_crnt, &factor_phi_al_post,
2824 &prefactor_phi, &v);
2827 /* Computation of forces per component */
2828 for (d = 0; d < DIM; d++)
2830 f_i[d] = prefactor_phi * (factor_phi_ai_ante * delta_ante[d] + factor_phi_ai_crnt * delta_crnt[d] + factor_phi_ai_post * delta_post[d]);
2831 f_j[d] = prefactor_phi * (factor_phi_aj_ante * delta_ante[d] + factor_phi_aj_crnt * delta_crnt[d] + factor_phi_aj_post * delta_post[d]);
2832 f_k[d] = prefactor_phi * (factor_phi_ak_ante * delta_ante[d] + factor_phi_ak_crnt * delta_crnt[d] + factor_phi_ak_post * delta_post[d]);
2833 f_l[d] = prefactor_phi * (factor_phi_al_ante * delta_ante[d] + factor_phi_al_crnt * delta_crnt[d] + factor_phi_al_post * delta_post[d]);
2835 /* Computation of the energy */
2841 /* Updating the forces */
2843 rvec_inc(f[ai], f_i);
2844 rvec_inc(f[aj], f_j);
2845 rvec_inc(f[ak], f_k);
2846 rvec_inc(f[al], f_l);
2849 /* Updating the fshift forces for the pressure coupling */
2852 copy_ivec(SHIFT_IVEC(g, aj), jt);
2853 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
2854 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
2855 ivec_sub(SHIFT_IVEC(g, al), jt, dt_lj);
2856 t1 = IVEC2IS(dt_ij);
2857 t2 = IVEC2IS(dt_kj);
2858 t3 = IVEC2IS(dt_lj);
2862 t3 = pbc_rvec_sub(pbc, x[al], x[aj], dx_jl);
2869 rvec_inc(fshift[t1], f_i);
2870 rvec_inc(fshift[CENTRAL], f_j);
2871 rvec_inc(fshift[t2], f_k);
2872 rvec_inc(fshift[t3], f_l);
2880 real cbtdihs(int nbonds,
2881 const t_iatom forceatoms[], const t_iparams forceparams[],
2882 const rvec x[], rvec f[], rvec fshift[],
2883 const t_pbc *pbc, const t_graph *g,
2884 real gmx_unused lambda, real gmx_unused *dvdlambda,
2885 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2886 int gmx_unused *global_atom_index)
2888 int type, ai, aj, ak, al, i, d;
2892 rvec f_i, f_j, f_k, f_l;
2893 ivec jt, dt_ij, dt_kj, dt_lj;
2895 rvec delta_ante, delta_crnt, delta_post;
2896 rvec f_phi_ai, f_phi_aj, f_phi_ak, f_phi_al;
2897 rvec f_theta_ante_ai, f_theta_ante_aj, f_theta_ante_ak;
2898 rvec f_theta_post_aj, f_theta_post_ak, f_theta_post_al;
2904 for (i = 0; (i < nbonds); )
2906 type = forceatoms[i++];
2907 ai = forceatoms[i++];
2908 aj = forceatoms[i++];
2909 ak = forceatoms[i++];
2910 al = forceatoms[i++];
2913 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], vec_temp);
2914 pbc_rvec_sub(pbc, x[aj], x[ai], delta_ante);
2915 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], vec_temp);
2916 pbc_rvec_sub(pbc, x[ak], x[aj], delta_crnt);
2917 t3 = pbc_rvec_sub(pbc, x[ak], x[al], vec_temp);
2918 pbc_rvec_sub(pbc, x[al], x[ak], delta_post);
2920 /* \brief Compute factors for CBT potential
2921 * The combined bending-torsion potential goes to zero in a very smooth manner, eliminating the numerical
2922 * instabilities, when three coarse-grained particles align and the dihedral angle and standard
2923 * dihedral potentials cannot be calculated. The CBT potential is calculated using the formula:
2924 * \f[V_{\rm CBT}(\theta_{i-1}, \theta_i, \phi_i) = k_{\phi} \sin^3\theta_{i-1} \sin^3\theta_{i}
2925 * \sum_{n=0}^4 { a_n \cos^n\phi_i}\f] ({eq:CBT} and ref \cite{MonicaGoga2013} from the manual).
2926 * It contains in its expression not only the dihedral angle \f$\phi\f$
2927 * but also \f[\theta_{i-1}\f] (theta_ante bellow) and \f[\theta_{i}\f] (theta_post bellow)
2928 * --- the adjacent bending angles.
2929 * For more explanations see comments file "restcbt.h". */
2931 compute_factors_cbtdihs(type, forceparams, delta_ante, delta_crnt, delta_post,
2932 f_phi_ai, f_phi_aj, f_phi_ak, f_phi_al,
2933 f_theta_ante_ai, f_theta_ante_aj, f_theta_ante_ak,
2934 f_theta_post_aj, f_theta_post_ak, f_theta_post_al,
2938 /* Acumulate the resuts per beads */
2939 for (d = 0; d < DIM; d++)
2941 f_i[d] = f_phi_ai[d] + f_theta_ante_ai[d];
2942 f_j[d] = f_phi_aj[d] + f_theta_ante_aj[d] + f_theta_post_aj[d];
2943 f_k[d] = f_phi_ak[d] + f_theta_ante_ak[d] + f_theta_post_ak[d];
2944 f_l[d] = f_phi_al[d] + f_theta_post_al[d];
2947 /* Compute the potential energy */
2952 /* Updating the forces */
2953 rvec_inc(f[ai], f_i);
2954 rvec_inc(f[aj], f_j);
2955 rvec_inc(f[ak], f_k);
2956 rvec_inc(f[al], f_l);
2959 /* Updating the fshift forces for the pressure coupling */
2962 copy_ivec(SHIFT_IVEC(g, aj), jt);
2963 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
2964 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
2965 ivec_sub(SHIFT_IVEC(g, al), jt, dt_lj);
2966 t1 = IVEC2IS(dt_ij);
2967 t2 = IVEC2IS(dt_kj);
2968 t3 = IVEC2IS(dt_lj);
2972 t3 = pbc_rvec_sub(pbc, x[al], x[aj], dx_jl);
2979 rvec_inc(fshift[t1], f_i);
2980 rvec_inc(fshift[CENTRAL], f_j);
2981 rvec_inc(fshift[t2], f_k);
2982 rvec_inc(fshift[t3], f_l);
2988 real rbdihs(int nbonds,
2989 const t_iatom forceatoms[], const t_iparams forceparams[],
2990 const rvec x[], rvec f[], rvec fshift[],
2991 const t_pbc *pbc, const t_graph *g,
2992 real lambda, real *dvdlambda,
2993 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2994 int gmx_unused *global_atom_index)
2996 const real c0 = 0.0, c1 = 1.0, c2 = 2.0, c3 = 3.0, c4 = 4.0, c5 = 5.0;
2997 int type, ai, aj, ak, al, i, j;
2999 rvec r_ij, r_kj, r_kl, m, n;
3000 real parmA[NR_RBDIHS];
3001 real parmB[NR_RBDIHS];
3002 real parm[NR_RBDIHS];
3003 real cos_phi, phi, rbp, rbpBA;
3004 real v, sign, ddphi, sin_phi;
3006 real L1 = 1.0-lambda;
3010 for (i = 0; (i < nbonds); )
3012 type = forceatoms[i++];
3013 ai = forceatoms[i++];
3014 aj = forceatoms[i++];
3015 ak = forceatoms[i++];
3016 al = forceatoms[i++];
3018 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
3019 &sign, &t1, &t2, &t3); /* 84 */
3021 /* Change to polymer convention */
3028 phi -= M_PI; /* 1 */
3032 /* Beware of accuracy loss, cannot use 1-sqrt(cos^2) ! */
3035 for (j = 0; (j < NR_RBDIHS); j++)
3037 parmA[j] = forceparams[type].rbdihs.rbcA[j];
3038 parmB[j] = forceparams[type].rbdihs.rbcB[j];
3039 parm[j] = L1*parmA[j]+lambda*parmB[j];
3041 /* Calculate cosine powers */
3042 /* Calculate the energy */
3043 /* Calculate the derivative */
3046 dvdl_term += (parmB[0]-parmA[0]);
3051 rbpBA = parmB[1]-parmA[1];
3052 ddphi += rbp*cosfac;
3055 dvdl_term += cosfac*rbpBA;
3057 rbpBA = parmB[2]-parmA[2];
3058 ddphi += c2*rbp*cosfac;
3061 dvdl_term += cosfac*rbpBA;
3063 rbpBA = parmB[3]-parmA[3];
3064 ddphi += c3*rbp*cosfac;
3067 dvdl_term += cosfac*rbpBA;
3069 rbpBA = parmB[4]-parmA[4];
3070 ddphi += c4*rbp*cosfac;
3073 dvdl_term += cosfac*rbpBA;
3075 rbpBA = parmB[5]-parmA[5];
3076 ddphi += c5*rbp*cosfac;
3079 dvdl_term += cosfac*rbpBA;
3081 ddphi = -ddphi*sin_phi; /* 11 */
3083 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
3084 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
3087 *dvdlambda += dvdl_term;
3092 int cmap_setup_grid_index(int ip, int grid_spacing, int *ipm1, int *ipp1, int *ipp2)
3098 ip = ip + grid_spacing - 1;
3100 else if (ip > grid_spacing)
3102 ip = ip - grid_spacing - 1;
3111 im1 = grid_spacing - 1;
3113 else if (ip == grid_spacing-2)
3117 else if (ip == grid_spacing-1)
3131 real cmap_dihs(int nbonds,
3132 const t_iatom forceatoms[], const t_iparams forceparams[],
3133 const gmx_cmap_t *cmap_grid,
3134 const rvec x[], rvec f[], rvec fshift[],
3135 const t_pbc *pbc, const t_graph *g,
3136 real gmx_unused lambda, real gmx_unused *dvdlambda,
3137 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3138 int gmx_unused *global_atom_index)
3140 int i, j, k, n, idx;
3141 int ai, aj, ak, al, am;
3142 int a1i, a1j, a1k, a1l, a2i, a2j, a2k, a2l;
3144 int t11, t21, t31, t12, t22, t32;
3145 int iphi1, ip1m1, ip1p1, ip1p2;
3146 int iphi2, ip2m1, ip2p1, ip2p2;
3148 int pos1, pos2, pos3, pos4, tmp;
3150 real ty[4], ty1[4], ty2[4], ty12[4], tc[16], tx[16];
3151 real phi1, psi1, cos_phi1, sin_phi1, sign1, xphi1;
3152 real phi2, psi2, cos_phi2, sin_phi2, sign2, xphi2;
3153 real dx, xx, tt, tu, e, df1, df2, ddf1, ddf2, ddf12, vtot;
3154 real ra21, rb21, rg21, rg1, rgr1, ra2r1, rb2r1, rabr1;
3155 real ra22, rb22, rg22, rg2, rgr2, ra2r2, rb2r2, rabr2;
3156 real fg1, hg1, fga1, hgb1, gaa1, gbb1;
3157 real fg2, hg2, fga2, hgb2, gaa2, gbb2;
3160 rvec r1_ij, r1_kj, r1_kl, m1, n1;
3161 rvec r2_ij, r2_kj, r2_kl, m2, n2;
3162 rvec f1_i, f1_j, f1_k, f1_l;
3163 rvec f2_i, f2_j, f2_k, f2_l;
3164 rvec a1, b1, a2, b2;
3165 rvec f1, g1, h1, f2, g2, h2;
3166 rvec dtf1, dtg1, dth1, dtf2, dtg2, dth2;
3167 ivec jt1, dt1_ij, dt1_kj, dt1_lj;
3168 ivec jt2, dt2_ij, dt2_kj, dt2_lj;
3172 int loop_index[4][4] = {
3179 /* Total CMAP energy */
3182 for (n = 0; n < nbonds; )
3184 /* Five atoms are involved in the two torsions */
3185 type = forceatoms[n++];
3186 ai = forceatoms[n++];
3187 aj = forceatoms[n++];
3188 ak = forceatoms[n++];
3189 al = forceatoms[n++];
3190 am = forceatoms[n++];
3192 /* Which CMAP type is this */
3193 cmapA = forceparams[type].cmap.cmapA;
3194 cmapd = cmap_grid->cmapdata[cmapA].cmap;
3202 phi1 = dih_angle(x[a1i], x[a1j], x[a1k], x[a1l], pbc, r1_ij, r1_kj, r1_kl, m1, n1,
3203 &sign1, &t11, &t21, &t31); /* 84 */
3205 cos_phi1 = cos(phi1);
3207 a1[0] = r1_ij[1]*r1_kj[2]-r1_ij[2]*r1_kj[1];
3208 a1[1] = r1_ij[2]*r1_kj[0]-r1_ij[0]*r1_kj[2];
3209 a1[2] = r1_ij[0]*r1_kj[1]-r1_ij[1]*r1_kj[0]; /* 9 */
3211 b1[0] = r1_kl[1]*r1_kj[2]-r1_kl[2]*r1_kj[1];
3212 b1[1] = r1_kl[2]*r1_kj[0]-r1_kl[0]*r1_kj[2];
3213 b1[2] = r1_kl[0]*r1_kj[1]-r1_kl[1]*r1_kj[0]; /* 9 */
3215 tmp = pbc_rvec_sub(pbc, x[a1l], x[a1k], h1);
3217 ra21 = iprod(a1, a1); /* 5 */
3218 rb21 = iprod(b1, b1); /* 5 */
3219 rg21 = iprod(r1_kj, r1_kj); /* 5 */
3225 rabr1 = sqrt(ra2r1*rb2r1);
3227 sin_phi1 = rg1 * rabr1 * iprod(a1, h1) * (-1);
3229 if (cos_phi1 < -0.5 || cos_phi1 > 0.5)
3231 phi1 = asin(sin_phi1);
3241 phi1 = -M_PI - phi1;
3247 phi1 = acos(cos_phi1);
3255 xphi1 = phi1 + M_PI; /* 1 */
3257 /* Second torsion */
3263 phi2 = dih_angle(x[a2i], x[a2j], x[a2k], x[a2l], pbc, r2_ij, r2_kj, r2_kl, m2, n2,
3264 &sign2, &t12, &t22, &t32); /* 84 */
3266 cos_phi2 = cos(phi2);
3268 a2[0] = r2_ij[1]*r2_kj[2]-r2_ij[2]*r2_kj[1];
3269 a2[1] = r2_ij[2]*r2_kj[0]-r2_ij[0]*r2_kj[2];
3270 a2[2] = r2_ij[0]*r2_kj[1]-r2_ij[1]*r2_kj[0]; /* 9 */
3272 b2[0] = r2_kl[1]*r2_kj[2]-r2_kl[2]*r2_kj[1];
3273 b2[1] = r2_kl[2]*r2_kj[0]-r2_kl[0]*r2_kj[2];
3274 b2[2] = r2_kl[0]*r2_kj[1]-r2_kl[1]*r2_kj[0]; /* 9 */
3276 tmp = pbc_rvec_sub(pbc, x[a2l], x[a2k], h2);
3278 ra22 = iprod(a2, a2); /* 5 */
3279 rb22 = iprod(b2, b2); /* 5 */
3280 rg22 = iprod(r2_kj, r2_kj); /* 5 */
3286 rabr2 = sqrt(ra2r2*rb2r2);
3288 sin_phi2 = rg2 * rabr2 * iprod(a2, h2) * (-1);
3290 if (cos_phi2 < -0.5 || cos_phi2 > 0.5)
3292 phi2 = asin(sin_phi2);
3302 phi2 = -M_PI - phi2;
3308 phi2 = acos(cos_phi2);
3316 xphi2 = phi2 + M_PI; /* 1 */
3318 /* Range mangling */
3321 xphi1 = xphi1 + 2*M_PI;
3323 else if (xphi1 >= 2*M_PI)
3325 xphi1 = xphi1 - 2*M_PI;
3330 xphi2 = xphi2 + 2*M_PI;
3332 else if (xphi2 >= 2*M_PI)
3334 xphi2 = xphi2 - 2*M_PI;
3337 /* Number of grid points */
3338 dx = 2*M_PI / cmap_grid->grid_spacing;
3340 /* Where on the grid are we */
3341 iphi1 = (int)(xphi1/dx);
3342 iphi2 = (int)(xphi2/dx);
3344 iphi1 = cmap_setup_grid_index(iphi1, cmap_grid->grid_spacing, &ip1m1, &ip1p1, &ip1p2);
3345 iphi2 = cmap_setup_grid_index(iphi2, cmap_grid->grid_spacing, &ip2m1, &ip2p1, &ip2p2);
3347 pos1 = iphi1*cmap_grid->grid_spacing+iphi2;
3348 pos2 = ip1p1*cmap_grid->grid_spacing+iphi2;
3349 pos3 = ip1p1*cmap_grid->grid_spacing+ip2p1;
3350 pos4 = iphi1*cmap_grid->grid_spacing+ip2p1;
3352 ty[0] = cmapd[pos1*4];
3353 ty[1] = cmapd[pos2*4];
3354 ty[2] = cmapd[pos3*4];
3355 ty[3] = cmapd[pos4*4];
3357 ty1[0] = cmapd[pos1*4+1];
3358 ty1[1] = cmapd[pos2*4+1];
3359 ty1[2] = cmapd[pos3*4+1];
3360 ty1[3] = cmapd[pos4*4+1];
3362 ty2[0] = cmapd[pos1*4+2];
3363 ty2[1] = cmapd[pos2*4+2];
3364 ty2[2] = cmapd[pos3*4+2];
3365 ty2[3] = cmapd[pos4*4+2];
3367 ty12[0] = cmapd[pos1*4+3];
3368 ty12[1] = cmapd[pos2*4+3];
3369 ty12[2] = cmapd[pos3*4+3];
3370 ty12[3] = cmapd[pos4*4+3];
3372 /* Switch to degrees */
3373 dx = 360.0 / cmap_grid->grid_spacing;
3374 xphi1 = xphi1 * RAD2DEG;
3375 xphi2 = xphi2 * RAD2DEG;
3377 for (i = 0; i < 4; i++) /* 16 */
3380 tx[i+4] = ty1[i]*dx;
3381 tx[i+8] = ty2[i]*dx;
3382 tx[i+12] = ty12[i]*dx*dx;
3386 for (i = 0; i < 4; i++) /* 1056 */
3388 for (j = 0; j < 4; j++)
3391 for (k = 0; k < 16; k++)
3393 xx = xx + cmap_coeff_matrix[k*16+idx]*tx[k];
3401 tt = (xphi1-iphi1*dx)/dx;
3402 tu = (xphi2-iphi2*dx)/dx;
3411 for (i = 3; i >= 0; i--)
3413 l1 = loop_index[i][3];
3414 l2 = loop_index[i][2];
3415 l3 = loop_index[i][1];
3417 e = tt * e + ((tc[i*4+3]*tu+tc[i*4+2])*tu + tc[i*4+1])*tu+tc[i*4];
3418 df1 = tu * df1 + (3.0*tc[l1]*tt+2.0*tc[l2])*tt+tc[l3];
3419 df2 = tt * df2 + (3.0*tc[i*4+3]*tu+2.0*tc[i*4+2])*tu+tc[i*4+1];
3420 ddf1 = tu * ddf1 + 2.0*3.0*tc[l1]*tt+2.0*tc[l2];
3421 ddf2 = tt * ddf2 + 2.0*3.0*tc[4*i+3]*tu+2.0*tc[4*i+2];
3424 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) +
3425 3.0*tu*tu*(tc[7]+2.0*tc[11]*tt+3.0*tc[15]*tt*tt);
3430 ddf1 = ddf1 * fac * fac;
3431 ddf2 = ddf2 * fac * fac;
3432 ddf12 = ddf12 * fac * fac;
3437 /* Do forces - first torsion */
3438 fg1 = iprod(r1_ij, r1_kj);
3439 hg1 = iprod(r1_kl, r1_kj);
3440 fga1 = fg1*ra2r1*rgr1;
3441 hgb1 = hg1*rb2r1*rgr1;
3445 for (i = 0; i < DIM; i++)
3447 dtf1[i] = gaa1 * a1[i];
3448 dtg1[i] = fga1 * a1[i] - hgb1 * b1[i];
3449 dth1[i] = gbb1 * b1[i];
3451 f1[i] = df1 * dtf1[i];
3452 g1[i] = df1 * dtg1[i];
3453 h1[i] = df1 * dth1[i];
3456 f1_j[i] = -f1[i] - g1[i];
3457 f1_k[i] = h1[i] + g1[i];
3460 f[a1i][i] = f[a1i][i] + f1_i[i];
3461 f[a1j][i] = f[a1j][i] + f1_j[i]; /* - f1[i] - g1[i] */
3462 f[a1k][i] = f[a1k][i] + f1_k[i]; /* h1[i] + g1[i] */
3463 f[a1l][i] = f[a1l][i] + f1_l[i]; /* h1[i] */
3466 /* Do forces - second torsion */
3467 fg2 = iprod(r2_ij, r2_kj);
3468 hg2 = iprod(r2_kl, r2_kj);
3469 fga2 = fg2*ra2r2*rgr2;
3470 hgb2 = hg2*rb2r2*rgr2;
3474 for (i = 0; i < DIM; i++)
3476 dtf2[i] = gaa2 * a2[i];
3477 dtg2[i] = fga2 * a2[i] - hgb2 * b2[i];
3478 dth2[i] = gbb2 * b2[i];
3480 f2[i] = df2 * dtf2[i];
3481 g2[i] = df2 * dtg2[i];
3482 h2[i] = df2 * dth2[i];
3485 f2_j[i] = -f2[i] - g2[i];
3486 f2_k[i] = h2[i] + g2[i];
3489 f[a2i][i] = f[a2i][i] + f2_i[i]; /* f2[i] */
3490 f[a2j][i] = f[a2j][i] + f2_j[i]; /* - f2[i] - g2[i] */
3491 f[a2k][i] = f[a2k][i] + f2_k[i]; /* h2[i] + g2[i] */
3492 f[a2l][i] = f[a2l][i] + f2_l[i]; /* - h2[i] */
3498 copy_ivec(SHIFT_IVEC(g, a1j), jt1);
3499 ivec_sub(SHIFT_IVEC(g, a1i), jt1, dt1_ij);
3500 ivec_sub(SHIFT_IVEC(g, a1k), jt1, dt1_kj);
3501 ivec_sub(SHIFT_IVEC(g, a1l), jt1, dt1_lj);
3502 t11 = IVEC2IS(dt1_ij);
3503 t21 = IVEC2IS(dt1_kj);
3504 t31 = IVEC2IS(dt1_lj);
3506 copy_ivec(SHIFT_IVEC(g, a2j), jt2);
3507 ivec_sub(SHIFT_IVEC(g, a2i), jt2, dt2_ij);
3508 ivec_sub(SHIFT_IVEC(g, a2k), jt2, dt2_kj);
3509 ivec_sub(SHIFT_IVEC(g, a2l), jt2, dt2_lj);
3510 t12 = IVEC2IS(dt2_ij);
3511 t22 = IVEC2IS(dt2_kj);
3512 t32 = IVEC2IS(dt2_lj);
3516 t31 = pbc_rvec_sub(pbc, x[a1l], x[a1j], h1);
3517 t32 = pbc_rvec_sub(pbc, x[a2l], x[a2j], h2);
3525 rvec_inc(fshift[t11], f1_i);
3526 rvec_inc(fshift[CENTRAL], f1_j);
3527 rvec_inc(fshift[t21], f1_k);
3528 rvec_inc(fshift[t31], f1_l);
3530 rvec_inc(fshift[t21], f2_i);
3531 rvec_inc(fshift[CENTRAL], f2_j);
3532 rvec_inc(fshift[t22], f2_k);
3533 rvec_inc(fshift[t32], f2_l);
3540 /***********************************************************
3542 * G R O M O S 9 6 F U N C T I O N S
3544 ***********************************************************/
3545 real g96harmonic(real kA, real kB, real xA, real xB, real x, real lambda,
3548 const real half = 0.5;
3549 real L1, kk, x0, dx, dx2;
3550 real v, f, dvdlambda;
3553 kk = L1*kA+lambda*kB;
3554 x0 = L1*xA+lambda*xB;
3561 dvdlambda = half*(kB-kA)*dx2 + (xA-xB)*kk*dx;
3568 /* That was 21 flops */
3571 real g96bonds(int nbonds,
3572 const t_iatom forceatoms[], const t_iparams forceparams[],
3573 const rvec x[], rvec f[], rvec fshift[],
3574 const t_pbc *pbc, const t_graph *g,
3575 real lambda, real *dvdlambda,
3576 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3577 int gmx_unused *global_atom_index)
3579 int i, m, ki, ai, aj, type;
3580 real dr2, fbond, vbond, fij, vtot;
3585 for (i = 0; (i < nbonds); )
3587 type = forceatoms[i++];
3588 ai = forceatoms[i++];
3589 aj = forceatoms[i++];
3591 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
3592 dr2 = iprod(dx, dx); /* 5 */
3594 *dvdlambda += g96harmonic(forceparams[type].harmonic.krA,
3595 forceparams[type].harmonic.krB,
3596 forceparams[type].harmonic.rA,
3597 forceparams[type].harmonic.rB,
3598 dr2, lambda, &vbond, &fbond);
3600 vtot += 0.5*vbond; /* 1*/
3604 fprintf(debug, "G96-BONDS: dr = %10g vbond = %10g fbond = %10g\n",
3605 sqrt(dr2), vbond, fbond);
3611 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
3614 for (m = 0; (m < DIM); m++) /* 15 */
3619 fshift[ki][m] += fij;
3620 fshift[CENTRAL][m] -= fij;
3626 real g96bond_angle(const rvec xi, const rvec xj, const rvec xk, const t_pbc *pbc,
3627 rvec r_ij, rvec r_kj,
3629 /* Return value is the angle between the bonds i-j and j-k */
3633 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
3634 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
3636 costh = cos_angle(r_ij, r_kj); /* 25 */
3641 real g96angles(int nbonds,
3642 const t_iatom forceatoms[], const t_iparams forceparams[],
3643 const rvec x[], rvec f[], rvec fshift[],
3644 const t_pbc *pbc, const t_graph *g,
3645 real lambda, real *dvdlambda,
3646 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3647 int gmx_unused *global_atom_index)
3649 int i, ai, aj, ak, type, m, t1, t2;
3651 real cos_theta, dVdt, va, vtot;
3652 real rij_1, rij_2, rkj_1, rkj_2, rijrkj_1;
3654 ivec jt, dt_ij, dt_kj;
3657 for (i = 0; (i < nbonds); )
3659 type = forceatoms[i++];
3660 ai = forceatoms[i++];
3661 aj = forceatoms[i++];
3662 ak = forceatoms[i++];
3664 cos_theta = g96bond_angle(x[ai], x[aj], x[ak], pbc, r_ij, r_kj, &t1, &t2);
3666 *dvdlambda += g96harmonic(forceparams[type].harmonic.krA,
3667 forceparams[type].harmonic.krB,
3668 forceparams[type].harmonic.rA,
3669 forceparams[type].harmonic.rB,
3670 cos_theta, lambda, &va, &dVdt);
3673 rij_1 = gmx_invsqrt(iprod(r_ij, r_ij));
3674 rkj_1 = gmx_invsqrt(iprod(r_kj, r_kj));
3675 rij_2 = rij_1*rij_1;
3676 rkj_2 = rkj_1*rkj_1;
3677 rijrkj_1 = rij_1*rkj_1; /* 23 */
3682 fprintf(debug, "G96ANGLES: costheta = %10g vth = %10g dV/dct = %10g\n",
3683 cos_theta, va, dVdt);
3686 for (m = 0; (m < DIM); m++) /* 42 */
3688 f_i[m] = dVdt*(r_kj[m]*rijrkj_1 - r_ij[m]*rij_2*cos_theta);
3689 f_k[m] = dVdt*(r_ij[m]*rijrkj_1 - r_kj[m]*rkj_2*cos_theta);
3690 f_j[m] = -f_i[m]-f_k[m];
3698 copy_ivec(SHIFT_IVEC(g, aj), jt);
3700 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3701 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3702 t1 = IVEC2IS(dt_ij);
3703 t2 = IVEC2IS(dt_kj);
3705 rvec_inc(fshift[t1], f_i);
3706 rvec_inc(fshift[CENTRAL], f_j);
3707 rvec_inc(fshift[t2], f_k); /* 9 */
3713 real cross_bond_bond(int nbonds,
3714 const t_iatom forceatoms[], const t_iparams forceparams[],
3715 const rvec x[], rvec f[], rvec fshift[],
3716 const t_pbc *pbc, const t_graph *g,
3717 real gmx_unused lambda, real gmx_unused *dvdlambda,
3718 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3719 int gmx_unused *global_atom_index)
3721 /* Potential from Lawrence and Skimmer, Chem. Phys. Lett. 372 (2003)
3724 int i, ai, aj, ak, type, m, t1, t2;
3726 real vtot, vrr, s1, s2, r1, r2, r1e, r2e, krr;
3728 ivec jt, dt_ij, dt_kj;
3731 for (i = 0; (i < nbonds); )
3733 type = forceatoms[i++];
3734 ai = forceatoms[i++];
3735 aj = forceatoms[i++];
3736 ak = forceatoms[i++];
3737 r1e = forceparams[type].cross_bb.r1e;
3738 r2e = forceparams[type].cross_bb.r2e;
3739 krr = forceparams[type].cross_bb.krr;
3741 /* Compute distance vectors ... */
3742 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
3743 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
3745 /* ... and their lengths */
3749 /* Deviations from ideality */
3753 /* Energy (can be negative!) */
3758 svmul(-krr*s2/r1, r_ij, f_i);
3759 svmul(-krr*s1/r2, r_kj, f_k);
3761 for (m = 0; (m < DIM); m++) /* 12 */
3763 f_j[m] = -f_i[m] - f_k[m];
3772 copy_ivec(SHIFT_IVEC(g, aj), jt);
3774 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3775 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3776 t1 = IVEC2IS(dt_ij);
3777 t2 = IVEC2IS(dt_kj);
3779 rvec_inc(fshift[t1], f_i);
3780 rvec_inc(fshift[CENTRAL], f_j);
3781 rvec_inc(fshift[t2], f_k); /* 9 */
3787 real cross_bond_angle(int nbonds,
3788 const t_iatom forceatoms[], const t_iparams forceparams[],
3789 const rvec x[], rvec f[], rvec fshift[],
3790 const t_pbc *pbc, const t_graph *g,
3791 real gmx_unused lambda, real gmx_unused *dvdlambda,
3792 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3793 int gmx_unused *global_atom_index)
3795 /* Potential from Lawrence and Skimmer, Chem. Phys. Lett. 372 (2003)
3798 int i, ai, aj, ak, type, m, t1, t2, t3;
3799 rvec r_ij, r_kj, r_ik;
3800 real vtot, vrt, s1, s2, s3, r1, r2, r3, r1e, r2e, r3e, krt, k1, k2, k3;
3802 ivec jt, dt_ij, dt_kj;
3805 for (i = 0; (i < nbonds); )
3807 type = forceatoms[i++];
3808 ai = forceatoms[i++];
3809 aj = forceatoms[i++];
3810 ak = forceatoms[i++];
3811 r1e = forceparams[type].cross_ba.r1e;
3812 r2e = forceparams[type].cross_ba.r2e;
3813 r3e = forceparams[type].cross_ba.r3e;
3814 krt = forceparams[type].cross_ba.krt;
3816 /* Compute distance vectors ... */
3817 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
3818 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
3819 t3 = pbc_rvec_sub(pbc, x[ai], x[ak], r_ik);
3821 /* ... and their lengths */
3826 /* Deviations from ideality */
3831 /* Energy (can be negative!) */
3832 vrt = krt*s3*(s1+s2);
3838 k3 = -krt*(s1+s2)/r3;
3839 for (m = 0; (m < DIM); m++)
3841 f_i[m] = k1*r_ij[m] + k3*r_ik[m];
3842 f_k[m] = k2*r_kj[m] - k3*r_ik[m];
3843 f_j[m] = -f_i[m] - f_k[m];
3846 for (m = 0; (m < DIM); m++) /* 12 */
3856 copy_ivec(SHIFT_IVEC(g, aj), jt);
3858 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3859 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3860 t1 = IVEC2IS(dt_ij);
3861 t2 = IVEC2IS(dt_kj);
3863 rvec_inc(fshift[t1], f_i);
3864 rvec_inc(fshift[CENTRAL], f_j);
3865 rvec_inc(fshift[t2], f_k); /* 9 */
3871 static real bonded_tab(const char *type, int table_nr,
3872 const bondedtable_t *table, real kA, real kB, real r,
3873 real lambda, real *V, real *F)
3875 real k, tabscale, *VFtab, rt, eps, eps2, Yt, Ft, Geps, Heps2, Fp, VV, FF;
3877 real v, f, dvdlambda;
3879 k = (1.0 - lambda)*kA + lambda*kB;
3881 tabscale = table->scale;
3882 VFtab = table->data;
3888 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",
3889 type, table_nr, r, n0, n0+1, table->n);
3896 Geps = VFtab[nnn+2]*eps;
3897 Heps2 = VFtab[nnn+3]*eps2;
3898 Fp = Ft + Geps + Heps2;
3900 FF = Fp + Geps + 2.0*Heps2;
3902 *F = -k*FF*tabscale;
3904 dvdlambda = (kB - kA)*VV;
3908 /* That was 22 flops */
3911 real tab_bonds(int nbonds,
3912 const t_iatom forceatoms[], const t_iparams forceparams[],
3913 const rvec x[], rvec f[], rvec fshift[],
3914 const t_pbc *pbc, const t_graph *g,
3915 real lambda, real *dvdlambda,
3916 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
3917 int gmx_unused *global_atom_index)
3919 int i, m, ki, ai, aj, type, table;
3920 real dr, dr2, fbond, vbond, fij, vtot;
3925 for (i = 0; (i < nbonds); )
3927 type = forceatoms[i++];
3928 ai = forceatoms[i++];
3929 aj = forceatoms[i++];
3931 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
3932 dr2 = iprod(dx, dx); /* 5 */
3933 dr = dr2*gmx_invsqrt(dr2); /* 10 */
3935 table = forceparams[type].tab.table;
3937 *dvdlambda += bonded_tab("bond", table,
3938 &fcd->bondtab[table],
3939 forceparams[type].tab.kA,
3940 forceparams[type].tab.kB,
3941 dr, lambda, &vbond, &fbond); /* 22 */
3949 vtot += vbond; /* 1*/
3950 fbond *= gmx_invsqrt(dr2); /* 6 */
3954 fprintf(debug, "TABBONDS: dr = %10g vbond = %10g fbond = %10g\n",
3960 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
3963 for (m = 0; (m < DIM); m++) /* 15 */
3968 fshift[ki][m] += fij;
3969 fshift[CENTRAL][m] -= fij;
3975 real tab_angles(int nbonds,
3976 const t_iatom forceatoms[], const t_iparams forceparams[],
3977 const rvec x[], rvec f[], rvec fshift[],
3978 const t_pbc *pbc, const t_graph *g,
3979 real lambda, real *dvdlambda,
3980 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
3981 int gmx_unused *global_atom_index)
3983 int i, ai, aj, ak, t1, t2, type, table;
3985 real cos_theta, cos_theta2, theta, dVdt, va, vtot;
3986 ivec jt, dt_ij, dt_kj;
3989 for (i = 0; (i < nbonds); )
3991 type = forceatoms[i++];
3992 ai = forceatoms[i++];
3993 aj = forceatoms[i++];
3994 ak = forceatoms[i++];
3996 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
3997 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
3999 table = forceparams[type].tab.table;
4001 *dvdlambda += bonded_tab("angle", table,
4002 &fcd->angletab[table],
4003 forceparams[type].tab.kA,
4004 forceparams[type].tab.kB,
4005 theta, lambda, &va, &dVdt); /* 22 */
4008 cos_theta2 = sqr(cos_theta); /* 1 */
4017 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
4018 sth = st*cos_theta; /* 1 */
4022 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
4023 theta*RAD2DEG, va, dVdt);
4026 nrkj2 = iprod(r_kj, r_kj); /* 5 */
4027 nrij2 = iprod(r_ij, r_ij);
4029 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
4030 cii = sth/nrij2; /* 10 */
4031 ckk = sth/nrkj2; /* 10 */
4033 for (m = 0; (m < DIM); m++) /* 39 */
4035 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
4036 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
4037 f_j[m] = -f_i[m]-f_k[m];
4044 copy_ivec(SHIFT_IVEC(g, aj), jt);
4046 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
4047 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
4048 t1 = IVEC2IS(dt_ij);
4049 t2 = IVEC2IS(dt_kj);
4051 rvec_inc(fshift[t1], f_i);
4052 rvec_inc(fshift[CENTRAL], f_j);
4053 rvec_inc(fshift[t2], f_k);
4059 real tab_dihs(int nbonds,
4060 const t_iatom forceatoms[], const t_iparams forceparams[],
4061 const rvec x[], rvec f[], rvec fshift[],
4062 const t_pbc *pbc, const t_graph *g,
4063 real lambda, real *dvdlambda,
4064 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
4065 int gmx_unused *global_atom_index)
4067 int i, type, ai, aj, ak, al, table;
4069 rvec r_ij, r_kj, r_kl, m, n;
4070 real phi, sign, ddphi, vpd, vtot;
4073 for (i = 0; (i < nbonds); )
4075 type = forceatoms[i++];
4076 ai = forceatoms[i++];
4077 aj = forceatoms[i++];
4078 ak = forceatoms[i++];
4079 al = forceatoms[i++];
4081 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
4082 &sign, &t1, &t2, &t3); /* 84 */
4084 table = forceparams[type].tab.table;
4086 /* Hopefully phi+M_PI never results in values < 0 */
4087 *dvdlambda += bonded_tab("dihedral", table,
4088 &fcd->dihtab[table],
4089 forceparams[type].tab.kA,
4090 forceparams[type].tab.kB,
4091 phi+M_PI, lambda, &vpd, &ddphi);
4094 do_dih_fup(ai, aj, ak, al, -ddphi, r_ij, r_kj, r_kl, m, n,
4095 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
4098 fprintf(debug, "pdih: (%d,%d,%d,%d) phi=%g\n",
4099 ai, aj, ak, al, phi);
4106 /* Return if this is a potential calculated in bondfree.c,
4107 * i.e. an interaction that actually calculates a potential and
4108 * works on multiple atoms (not e.g. a connection or a position restraint).
4110 static gmx_inline gmx_bool ftype_is_bonded_potential(int ftype)
4113 (interaction_function[ftype].flags & IF_BOND) &&
4114 !(ftype == F_CONNBONDS || ftype == F_POSRES || ftype == F_FBPOSRES) &&
4115 (ftype < F_GB12 || ftype > F_GB14);
4118 static void divide_bondeds_over_threads(t_idef *idef, int nthreads)
4125 idef->nthreads = nthreads;
4127 if (F_NRE*(nthreads+1) > idef->il_thread_division_nalloc)
4129 idef->il_thread_division_nalloc = F_NRE*(nthreads+1);
4130 snew(idef->il_thread_division, idef->il_thread_division_nalloc);
4133 for (ftype = 0; ftype < F_NRE; ftype++)
4135 if (ftype_is_bonded_potential(ftype))
4137 nat1 = interaction_function[ftype].nratoms + 1;
4139 for (t = 0; t <= nthreads; t++)
4141 /* Divide the interactions equally over the threads.
4142 * When the different types of bonded interactions
4143 * are distributed roughly equally over the threads,
4144 * this should lead to well localized output into
4145 * the force buffer on each thread.
4146 * If this is not the case, a more advanced scheme
4147 * (not implemented yet) will do better.
4149 il_nr_thread = (((idef->il[ftype].nr/nat1)*t)/nthreads)*nat1;
4151 /* Ensure that distance restraint pairs with the same label
4152 * end up on the same thread.
4153 * This is slighlty tricky code, since the next for iteration
4154 * may have an initial il_nr_thread lower than the final value
4155 * in the previous iteration, but this will anyhow be increased
4156 * to the approriate value again by this while loop.
4158 while (ftype == F_DISRES &&
4160 il_nr_thread < idef->il[ftype].nr &&
4161 idef->iparams[idef->il[ftype].iatoms[il_nr_thread]].disres.label ==
4162 idef->iparams[idef->il[ftype].iatoms[il_nr_thread-nat1]].disres.label)
4164 il_nr_thread += nat1;
4167 idef->il_thread_division[ftype*(nthreads+1)+t] = il_nr_thread;
4174 calc_bonded_reduction_mask(const t_idef *idef,
4179 int ftype, nb, nat1, nb0, nb1, i, a;
4183 for (ftype = 0; ftype < F_NRE; ftype++)
4185 if (ftype_is_bonded_potential(ftype))
4187 nb = idef->il[ftype].nr;
4190 nat1 = interaction_function[ftype].nratoms + 1;
4192 /* Divide this interaction equally over the threads.
4193 * This is not stored: should match division in calc_bonds.
4195 nb0 = idef->il_thread_division[ftype*(nt+1)+t];
4196 nb1 = idef->il_thread_division[ftype*(nt+1)+t+1];
4198 for (i = nb0; i < nb1; i += nat1)
4200 for (a = 1; a < nat1; a++)
4202 mask |= (1U << (idef->il[ftype].iatoms[i+a]>>shift));
4212 void setup_bonded_threading(t_forcerec *fr, t_idef *idef)
4214 #define MAX_BLOCK_BITS 32
4218 assert(fr->nthreads >= 1);
4220 /* Divide the bonded interaction over the threads */
4221 divide_bondeds_over_threads(idef, fr->nthreads);
4223 if (fr->nthreads == 1)
4230 /* We divide the force array in a maximum of 32 blocks.
4231 * Minimum force block reduction size is 2^6=64.
4234 while (fr->natoms_force > (int)(MAX_BLOCK_BITS*(1U<<fr->red_ashift)))
4240 fprintf(debug, "bonded force buffer block atom shift %d bits\n",
4244 /* Determine to which blocks each thread's bonded force calculation
4245 * contributes. Store this is a mask for each thread.
4247 #pragma omp parallel for num_threads(fr->nthreads) schedule(static)
4248 for (t = 1; t < fr->nthreads; t++)
4250 fr->f_t[t].red_mask =
4251 calc_bonded_reduction_mask(idef, fr->red_ashift, t, fr->nthreads);
4254 /* Determine the maximum number of blocks we need to reduce over */
4257 for (t = 0; t < fr->nthreads; t++)
4260 for (b = 0; b < MAX_BLOCK_BITS; b++)
4262 if (fr->f_t[t].red_mask & (1U<<b))
4264 fr->red_nblock = max(fr->red_nblock, b+1);
4270 fprintf(debug, "thread %d flags %x count %d\n",
4271 t, fr->f_t[t].red_mask, c);
4277 fprintf(debug, "Number of blocks to reduce: %d of size %d\n",
4278 fr->red_nblock, 1<<fr->red_ashift);
4279 fprintf(debug, "Reduction density %.2f density/#thread %.2f\n",
4280 ctot*(1<<fr->red_ashift)/(double)fr->natoms_force,
4281 ctot*(1<<fr->red_ashift)/(double)(fr->natoms_force*fr->nthreads));
4285 static void zero_thread_forces(f_thread_t *f_t, int n,
4286 int nblock, int blocksize)
4288 int b, a0, a1, a, i, j;
4290 if (n > f_t->f_nalloc)
4292 f_t->f_nalloc = over_alloc_large(n);
4293 srenew(f_t->f, f_t->f_nalloc);
4296 if (f_t->red_mask != 0)
4298 for (b = 0; b < nblock; b++)
4300 if (f_t->red_mask && (1U<<b))
4303 a1 = min((b+1)*blocksize, n);
4304 for (a = a0; a < a1; a++)
4306 clear_rvec(f_t->f[a]);
4311 for (i = 0; i < SHIFTS; i++)
4313 clear_rvec(f_t->fshift[i]);
4315 for (i = 0; i < F_NRE; i++)
4319 for (i = 0; i < egNR; i++)
4321 for (j = 0; j < f_t->grpp.nener; j++)
4323 f_t->grpp.ener[i][j] = 0;
4326 for (i = 0; i < efptNR; i++)
4332 static void reduce_thread_force_buffer(int n, rvec *f,
4333 int nthreads, f_thread_t *f_t,
4334 int nblock, int block_size)
4336 /* The max thread number is arbitrary,
4337 * we used a fixed number to avoid memory management.
4338 * Using more than 16 threads is probably never useful performance wise.
4340 #define MAX_BONDED_THREADS 256
4343 if (nthreads > MAX_BONDED_THREADS)
4345 gmx_fatal(FARGS, "Can not reduce bonded forces on more than %d threads",
4346 MAX_BONDED_THREADS);
4349 /* This reduction can run on any number of threads,
4350 * independently of nthreads.
4352 #pragma omp parallel for num_threads(nthreads) schedule(static)
4353 for (b = 0; b < nblock; b++)
4355 rvec *fp[MAX_BONDED_THREADS];
4359 /* Determine which threads contribute to this block */
4361 for (ft = 1; ft < nthreads; ft++)
4363 if (f_t[ft].red_mask & (1U<<b))
4365 fp[nfb++] = f_t[ft].f;
4370 /* Reduce force buffers for threads that contribute */
4372 a1 = (b+1)*block_size;
4374 for (a = a0; a < a1; a++)
4376 for (fb = 0; fb < nfb; fb++)
4378 rvec_inc(f[a], fp[fb][a]);
4385 static void reduce_thread_forces(int n, rvec *f, rvec *fshift,
4386 real *ener, gmx_grppairener_t *grpp, real *dvdl,
4387 int nthreads, f_thread_t *f_t,
4388 int nblock, int block_size,
4389 gmx_bool bCalcEnerVir,
4394 /* Reduce the bonded force buffer */
4395 reduce_thread_force_buffer(n, f, nthreads, f_t, nblock, block_size);
4398 /* When necessary, reduce energy and virial using one thread only */
4403 for (i = 0; i < SHIFTS; i++)
4405 for (t = 1; t < nthreads; t++)
4407 rvec_inc(fshift[i], f_t[t].fshift[i]);
4410 for (i = 0; i < F_NRE; i++)
4412 for (t = 1; t < nthreads; t++)
4414 ener[i] += f_t[t].ener[i];
4417 for (i = 0; i < egNR; i++)
4419 for (j = 0; j < f_t[1].grpp.nener; j++)
4421 for (t = 1; t < nthreads; t++)
4424 grpp->ener[i][j] += f_t[t].grpp.ener[i][j];
4430 for (i = 0; i < efptNR; i++)
4433 for (t = 1; t < nthreads; t++)
4435 dvdl[i] += f_t[t].dvdl[i];
4442 static real calc_one_bond(FILE *fplog, int thread,
4443 int ftype, const t_idef *idef,
4444 rvec x[], rvec f[], rvec fshift[],
4446 const t_pbc *pbc, const t_graph *g,
4447 gmx_grppairener_t *grpp,
4449 real *lambda, real *dvdl,
4450 const t_mdatoms *md, t_fcdata *fcd,
4451 gmx_bool bCalcEnerVir,
4452 int *global_atom_index, gmx_bool bPrintSepPot)
4454 int nat1, nbonds, efptFTYPE;
4459 if (IS_RESTRAINT_TYPE(ftype))
4461 efptFTYPE = efptRESTRAINT;
4465 efptFTYPE = efptBONDED;
4468 nat1 = interaction_function[ftype].nratoms + 1;
4469 nbonds = idef->il[ftype].nr/nat1;
4470 iatoms = idef->il[ftype].iatoms;
4472 nb0 = idef->il_thread_division[ftype*(idef->nthreads+1)+thread];
4473 nbn = idef->il_thread_division[ftype*(idef->nthreads+1)+thread+1] - nb0;
4475 if (!IS_LISTED_LJ_C(ftype))
4477 if (ftype == F_CMAP)
4479 v = cmap_dihs(nbn, iatoms+nb0,
4480 idef->iparams, &idef->cmap_grid,
4481 (const rvec*)x, f, fshift,
4482 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
4483 md, fcd, global_atom_index);
4485 #ifdef GMX_SIMD_HAVE_REAL
4486 else if (ftype == F_ANGLES &&
4487 !bCalcEnerVir && fr->efep == efepNO)
4489 /* No energies, shift forces, dvdl */
4490 angles_noener_simd(nbn, idef->il[ftype].iatoms+nb0,
4493 pbc, g, lambda[efptFTYPE], md, fcd,
4498 else if (ftype == F_PDIHS &&
4499 !bCalcEnerVir && fr->efep == efepNO)
4501 /* No energies, shift forces, dvdl */
4502 #ifdef GMX_SIMD_HAVE_REAL
4507 (nbn, idef->il[ftype].iatoms+nb0,
4510 pbc, g, lambda[efptFTYPE], md, fcd,
4516 v = interaction_function[ftype].ifunc(nbn, iatoms+nb0,
4518 (const rvec*)x, f, fshift,
4519 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
4520 md, fcd, global_atom_index);
4524 fprintf(fplog, " %-23s #%4d V %12.5e dVdl %12.5e\n",
4525 interaction_function[ftype].longname,
4526 nbonds, v, lambda[efptFTYPE]);
4531 v = do_nonbonded_listed(ftype, nbn, iatoms+nb0, idef->iparams, (const rvec*)x, f, fshift,
4532 pbc, g, lambda, dvdl, md, fr, grpp, global_atom_index);
4536 fprintf(fplog, " %-5s + %-15s #%4d dVdl %12.5e\n",
4537 interaction_function[ftype].longname,
4538 interaction_function[F_LJ14].longname, nbonds, dvdl[efptVDW]);
4539 fprintf(fplog, " %-5s + %-15s #%4d dVdl %12.5e\n",
4540 interaction_function[ftype].longname,
4541 interaction_function[F_COUL14].longname, nbonds, dvdl[efptCOUL]);
4547 inc_nrnb(nrnb, interaction_function[ftype].nrnb_ind, nbonds);
4553 void calc_bonds(FILE *fplog, const gmx_multisim_t *ms,
4555 rvec x[], history_t *hist,
4556 rvec f[], t_forcerec *fr,
4557 const t_pbc *pbc, const t_graph *g,
4558 gmx_enerdata_t *enerd, t_nrnb *nrnb,
4560 const t_mdatoms *md,
4561 t_fcdata *fcd, int *global_atom_index,
4562 t_atomtypes gmx_unused *atype, gmx_genborn_t gmx_unused *born,
4564 gmx_bool bPrintSepPot, gmx_int64_t step)
4566 gmx_bool bCalcEnerVir;
4568 real v, dvdl[efptNR], dvdl_dum[efptNR]; /* The dummy array is to have a place to store the dhdl at other values
4569 of lambda, which will be thrown away in the end*/
4570 const t_pbc *pbc_null;
4574 assert(fr->nthreads == idef->nthreads);
4576 bCalcEnerVir = (force_flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY));
4578 for (i = 0; i < efptNR; i++)
4592 fprintf(fplog, "Step %s: bonded V and dVdl for this node\n",
4593 gmx_step_str(step, buf));
4599 p_graph(debug, "Bondage is fun", g);
4603 /* Do pre force calculation stuff which might require communication */
4604 if (idef->il[F_ORIRES].nr)
4606 enerd->term[F_ORIRESDEV] =
4607 calc_orires_dev(ms, idef->il[F_ORIRES].nr,
4608 idef->il[F_ORIRES].iatoms,
4609 idef->iparams, md, (const rvec*)x,
4610 pbc_null, fcd, hist);
4612 if (idef->il[F_DISRES].nr)
4614 calc_disres_R_6(idef->il[F_DISRES].nr,
4615 idef->il[F_DISRES].iatoms,
4616 idef->iparams, (const rvec*)x, pbc_null,
4619 if (fcd->disres.nsystems > 1)
4621 gmx_sum_sim(2*fcd->disres.nres, fcd->disres.Rt_6, ms);
4626 #pragma omp parallel for num_threads(fr->nthreads) schedule(static)
4627 for (thread = 0; thread < fr->nthreads; thread++)
4634 gmx_grppairener_t *grpp;
4639 fshift = fr->fshift;
4641 grpp = &enerd->grpp;
4646 zero_thread_forces(&fr->f_t[thread], fr->natoms_force,
4647 fr->red_nblock, 1<<fr->red_ashift);
4649 ft = fr->f_t[thread].f;
4650 fshift = fr->f_t[thread].fshift;
4651 epot = fr->f_t[thread].ener;
4652 grpp = &fr->f_t[thread].grpp;
4653 dvdlt = fr->f_t[thread].dvdl;
4655 /* Loop over all bonded force types to calculate the bonded forces */
4656 for (ftype = 0; (ftype < F_NRE); ftype++)
4658 if (idef->il[ftype].nr > 0 && ftype_is_bonded_potential(ftype))
4660 v = calc_one_bond(fplog, thread, ftype, idef, x,
4661 ft, fshift, fr, pbc_null, g, grpp,
4662 nrnb, lambda, dvdlt,
4663 md, fcd, bCalcEnerVir,
4664 global_atom_index, bPrintSepPot);
4669 if (fr->nthreads > 1)
4671 reduce_thread_forces(fr->natoms_force, f, fr->fshift,
4672 enerd->term, &enerd->grpp, dvdl,
4673 fr->nthreads, fr->f_t,
4674 fr->red_nblock, 1<<fr->red_ashift,
4676 force_flags & GMX_FORCE_DHDL);
4678 if (force_flags & GMX_FORCE_DHDL)
4680 for (i = 0; i < efptNR; i++)
4682 enerd->dvdl_nonlin[i] += dvdl[i];
4686 /* Copy the sum of violations for the distance restraints from fcd */
4689 enerd->term[F_DISRESVIOL] = fcd->disres.sumviol;
4694 void calc_bonds_lambda(FILE *fplog,
4698 const t_pbc *pbc, const t_graph *g,
4699 gmx_grppairener_t *grpp, real *epot, t_nrnb *nrnb,
4701 const t_mdatoms *md,
4703 int *global_atom_index)
4705 int i, ftype, nr_nonperturbed, nr;
4707 real dvdl_dum[efptNR];
4709 const t_pbc *pbc_null;
4721 /* Copy the whole idef, so we can modify the contents locally */
4723 idef_fe.nthreads = 1;
4724 snew(idef_fe.il_thread_division, F_NRE*(idef_fe.nthreads+1));
4726 /* We already have the forces, so we use temp buffers here */
4727 snew(f, fr->natoms_force);
4728 snew(fshift, SHIFTS);
4730 /* Loop over all bonded force types to calculate the bonded energies */
4731 for (ftype = 0; (ftype < F_NRE); ftype++)
4733 if (ftype_is_bonded_potential(ftype))
4735 /* Set the work range of thread 0 to the perturbed bondeds only */
4736 nr_nonperturbed = idef->il[ftype].nr_nonperturbed;
4737 nr = idef->il[ftype].nr;
4738 idef_fe.il_thread_division[ftype*2+0] = nr_nonperturbed;
4739 idef_fe.il_thread_division[ftype*2+1] = nr;
4741 /* This is only to get the flop count correct */
4742 idef_fe.il[ftype].nr = nr - nr_nonperturbed;
4744 if (nr - nr_nonperturbed > 0)
4746 v = calc_one_bond(fplog, 0, ftype, &idef_fe,
4747 x, f, fshift, fr, pbc_null, g,
4748 grpp, nrnb, lambda, dvdl_dum,
4750 global_atom_index, FALSE);
4759 sfree(idef_fe.il_thread_division);