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.
48 #include "gromacs/math/units.h"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/math/utilities.h"
51 #include "gromacs/legacyheaders/txtdump.h"
52 #include "gromacs/legacyheaders/ns.h"
53 #include "gromacs/legacyheaders/macros.h"
54 #include "gromacs/legacyheaders/names.h"
55 #include "gromacs/legacyheaders/disre.h"
56 #include "gromacs/legacyheaders/orires.h"
57 #include "gromacs/legacyheaders/force.h"
58 #include "gromacs/legacyheaders/nonbonded.h"
60 #include "gromacs/pbcutil/ishift.h"
61 #include "gromacs/pbcutil/mshift.h"
62 #include "gromacs/pbcutil/pbc.h"
63 #include "gromacs/simd/simd.h"
64 #include "gromacs/simd/simd_math.h"
65 #include "gromacs/simd/vector_operations.h"
66 #include "gromacs/utility/fatalerror.h"
67 #include "gromacs/utility/smalloc.h"
71 /* Find a better place for this? */
72 const int cmap_coeff_matrix[] = {
73 1, 0, -3, 2, 0, 0, 0, 0, -3, 0, 9, -6, 2, 0, -6, 4,
74 0, 0, 0, 0, 0, 0, 0, 0, 3, 0, -9, 6, -2, 0, 6, -4,
75 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 9, -6, 0, 0, -6, 4,
76 0, 0, 3, -2, 0, 0, 0, 0, 0, 0, -9, 6, 0, 0, 6, -4,
77 0, 0, 0, 0, 1, 0, -3, 2, -2, 0, 6, -4, 1, 0, -3, 2,
78 0, 0, 0, 0, 0, 0, 0, 0, -1, 0, 3, -2, 1, 0, -3, 2,
79 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -3, 2, 0, 0, 3, -2,
80 0, 0, 0, 0, 0, 0, 3, -2, 0, 0, -6, 4, 0, 0, 3, -2,
81 0, 1, -2, 1, 0, 0, 0, 0, 0, -3, 6, -3, 0, 2, -4, 2,
82 0, 0, 0, 0, 0, 0, 0, 0, 0, 3, -6, 3, 0, -2, 4, -2,
83 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -3, 3, 0, 0, 2, -2,
84 0, 0, -1, 1, 0, 0, 0, 0, 0, 0, 3, -3, 0, 0, -2, 2,
85 0, 0, 0, 0, 0, 1, -2, 1, 0, -2, 4, -2, 0, 1, -2, 1,
86 0, 0, 0, 0, 0, 0, 0, 0, 0, -1, 2, -1, 0, 1, -2, 1,
87 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, -1, 0, 0, -1, 1,
88 0, 0, 0, 0, 0, 0, -1, 1, 0, 0, 2, -2, 0, 0, -1, 1
92 /* TODO This function should go and live in nonbonded.c where it is
93 really needed. Here, it only supports giving a fatal error message
95 int glatnr(int *global_atom_index, int i)
99 if (global_atom_index == NULL)
105 atnr = global_atom_index[i] + 1;
111 /* TODO This kind of code appears in many places. Consolidate it */
112 static int pbc_rvec_sub(const t_pbc *pbc, const rvec xi, const rvec xj, rvec dx)
116 return pbc_dx_aiuc(pbc, xi, xj, dx);
120 rvec_sub(xi, xj, dx);
125 #ifdef GMX_SIMD_HAVE_REAL
127 /* SIMD PBC data structure, containing 1/boxdiag and the box vectors */
129 gmx_simd_real_t inv_bzz;
130 gmx_simd_real_t inv_byy;
131 gmx_simd_real_t inv_bxx;
140 /* Set the SIMD pbc data from a normal t_pbc struct */
141 static void set_pbc_simd(const t_pbc *pbc, pbc_simd_t *pbc_simd)
146 /* Setting inv_bdiag to 0 effectively turns off PBC */
147 clear_rvec(inv_bdiag);
150 for (d = 0; d < pbc->ndim_ePBC; d++)
152 inv_bdiag[d] = 1.0/pbc->box[d][d];
156 pbc_simd->inv_bzz = gmx_simd_set1_r(inv_bdiag[ZZ]);
157 pbc_simd->inv_byy = gmx_simd_set1_r(inv_bdiag[YY]);
158 pbc_simd->inv_bxx = gmx_simd_set1_r(inv_bdiag[XX]);
162 pbc_simd->bzx = gmx_simd_set1_r(pbc->box[ZZ][XX]);
163 pbc_simd->bzy = gmx_simd_set1_r(pbc->box[ZZ][YY]);
164 pbc_simd->bzz = gmx_simd_set1_r(pbc->box[ZZ][ZZ]);
165 pbc_simd->byx = gmx_simd_set1_r(pbc->box[YY][XX]);
166 pbc_simd->byy = gmx_simd_set1_r(pbc->box[YY][YY]);
167 pbc_simd->bxx = gmx_simd_set1_r(pbc->box[XX][XX]);
171 pbc_simd->bzx = gmx_simd_setzero_r();
172 pbc_simd->bzy = gmx_simd_setzero_r();
173 pbc_simd->bzz = gmx_simd_setzero_r();
174 pbc_simd->byx = gmx_simd_setzero_r();
175 pbc_simd->byy = gmx_simd_setzero_r();
176 pbc_simd->bxx = gmx_simd_setzero_r();
180 /* Correct distance vector *dx,*dy,*dz for PBC using SIMD */
181 static gmx_inline void
182 pbc_dx_simd(gmx_simd_real_t *dx, gmx_simd_real_t *dy, gmx_simd_real_t *dz,
183 const pbc_simd_t *pbc)
187 sh = gmx_simd_round_r(gmx_simd_mul_r(*dz, pbc->inv_bzz));
188 *dx = gmx_simd_fnmadd_r(sh, pbc->bzx, *dx);
189 *dy = gmx_simd_fnmadd_r(sh, pbc->bzy, *dy);
190 *dz = gmx_simd_fnmadd_r(sh, pbc->bzz, *dz);
192 sh = gmx_simd_round_r(gmx_simd_mul_r(*dy, pbc->inv_byy));
193 *dx = gmx_simd_fnmadd_r(sh, pbc->byx, *dx);
194 *dy = gmx_simd_fnmadd_r(sh, pbc->byy, *dy);
196 sh = gmx_simd_round_r(gmx_simd_mul_r(*dx, pbc->inv_bxx));
197 *dx = gmx_simd_fnmadd_r(sh, pbc->bxx, *dx);
200 #endif /* GMX_SIMD_HAVE_REAL */
203 * Morse potential bond by Frank Everdij
205 * Three parameters needed:
207 * b0 = equilibrium distance in nm
208 * be = beta in nm^-1 (actually, it's nu_e*Sqrt(2*pi*pi*mu/D_e))
209 * cb = well depth in kJ/mol
211 * Note: the potential is referenced to be +cb at infinite separation
212 * and zero at the equilibrium distance!
215 real morse_bonds(int nbonds,
216 const t_iatom forceatoms[], const t_iparams forceparams[],
217 const rvec x[], rvec f[], rvec fshift[],
218 const t_pbc *pbc, const t_graph *g,
219 real lambda, real *dvdlambda,
220 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
221 int gmx_unused *global_atom_index)
223 const real one = 1.0;
224 const real two = 2.0;
225 real dr, dr2, temp, omtemp, cbomtemp, fbond, vbond, fij, vtot;
226 real b0, be, cb, b0A, beA, cbA, b0B, beB, cbB, L1;
228 int i, m, ki, type, ai, aj;
232 for (i = 0; (i < nbonds); )
234 type = forceatoms[i++];
235 ai = forceatoms[i++];
236 aj = forceatoms[i++];
238 b0A = forceparams[type].morse.b0A;
239 beA = forceparams[type].morse.betaA;
240 cbA = forceparams[type].morse.cbA;
242 b0B = forceparams[type].morse.b0B;
243 beB = forceparams[type].morse.betaB;
244 cbB = forceparams[type].morse.cbB;
246 L1 = one-lambda; /* 1 */
247 b0 = L1*b0A + lambda*b0B; /* 3 */
248 be = L1*beA + lambda*beB; /* 3 */
249 cb = L1*cbA + lambda*cbB; /* 3 */
251 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
252 dr2 = iprod(dx, dx); /* 5 */
253 dr = dr2*gmx_invsqrt(dr2); /* 10 */
254 temp = exp(-be*(dr-b0)); /* 12 */
258 /* bonds are constrainted. This may _not_ include bond constraints if they are lambda dependent */
259 *dvdlambda += cbB-cbA;
263 omtemp = one-temp; /* 1 */
264 cbomtemp = cb*omtemp; /* 1 */
265 vbond = cbomtemp*omtemp; /* 1 */
266 fbond = -two*be*temp*cbomtemp*gmx_invsqrt(dr2); /* 9 */
267 vtot += vbond; /* 1 */
269 *dvdlambda += (cbB - cbA) * omtemp * omtemp - (2-2*omtemp)*omtemp * cb * ((b0B-b0A)*be - (beB-beA)*(dr-b0)); /* 15 */
273 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
277 for (m = 0; (m < DIM); m++) /* 15 */
282 fshift[ki][m] += fij;
283 fshift[CENTRAL][m] -= fij;
289 real cubic_bonds(int nbonds,
290 const t_iatom forceatoms[], const t_iparams forceparams[],
291 const rvec x[], rvec f[], rvec fshift[],
292 const t_pbc *pbc, const t_graph *g,
293 real gmx_unused lambda, real gmx_unused *dvdlambda,
294 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
295 int gmx_unused *global_atom_index)
297 const real three = 3.0;
298 const real two = 2.0;
300 real dr, dr2, dist, kdist, kdist2, fbond, vbond, fij, vtot;
302 int i, m, ki, type, ai, aj;
306 for (i = 0; (i < nbonds); )
308 type = forceatoms[i++];
309 ai = forceatoms[i++];
310 aj = forceatoms[i++];
312 b0 = forceparams[type].cubic.b0;
313 kb = forceparams[type].cubic.kb;
314 kcub = forceparams[type].cubic.kcub;
316 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
317 dr2 = iprod(dx, dx); /* 5 */
324 dr = dr2*gmx_invsqrt(dr2); /* 10 */
329 vbond = kdist2 + kcub*kdist2*dist;
330 fbond = -(two*kdist + three*kdist2*kcub)/dr;
332 vtot += vbond; /* 21 */
336 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
339 for (m = 0; (m < DIM); m++) /* 15 */
344 fshift[ki][m] += fij;
345 fshift[CENTRAL][m] -= fij;
351 real FENE_bonds(int nbonds,
352 const t_iatom forceatoms[], const t_iparams forceparams[],
353 const rvec x[], rvec f[], rvec fshift[],
354 const t_pbc *pbc, const t_graph *g,
355 real gmx_unused lambda, real gmx_unused *dvdlambda,
356 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
357 int *global_atom_index)
359 const real half = 0.5;
360 const real one = 1.0;
362 real dr2, bm2, omdr2obm2, fbond, vbond, fij, vtot;
364 int i, m, ki, type, ai, aj;
368 for (i = 0; (i < nbonds); )
370 type = forceatoms[i++];
371 ai = forceatoms[i++];
372 aj = forceatoms[i++];
374 bm = forceparams[type].fene.bm;
375 kb = forceparams[type].fene.kb;
377 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
378 dr2 = iprod(dx, dx); /* 5 */
390 "r^2 (%f) >= bm^2 (%f) in FENE bond between atoms %d and %d",
392 glatnr(global_atom_index, ai),
393 glatnr(global_atom_index, aj));
396 omdr2obm2 = one - dr2/bm2;
398 vbond = -half*kb*bm2*log(omdr2obm2);
399 fbond = -kb/omdr2obm2;
401 vtot += vbond; /* 35 */
405 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
408 for (m = 0; (m < DIM); m++) /* 15 */
413 fshift[ki][m] += fij;
414 fshift[CENTRAL][m] -= fij;
420 real harmonic(real kA, real kB, real xA, real xB, real x, real lambda,
423 const real half = 0.5;
424 real L1, kk, x0, dx, dx2;
425 real v, f, dvdlambda;
428 kk = L1*kA+lambda*kB;
429 x0 = L1*xA+lambda*xB;
436 dvdlambda = half*(kB-kA)*dx2 + (xA-xB)*kk*dx;
443 /* That was 19 flops */
447 real bonds(int nbonds,
448 const t_iatom forceatoms[], const t_iparams forceparams[],
449 const rvec x[], rvec f[], rvec fshift[],
450 const t_pbc *pbc, const t_graph *g,
451 real lambda, real *dvdlambda,
452 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
453 int gmx_unused *global_atom_index)
455 int i, m, ki, ai, aj, type;
456 real dr, dr2, fbond, vbond, fij, vtot;
461 for (i = 0; (i < nbonds); )
463 type = forceatoms[i++];
464 ai = forceatoms[i++];
465 aj = forceatoms[i++];
467 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
468 dr2 = iprod(dx, dx); /* 5 */
469 dr = dr2*gmx_invsqrt(dr2); /* 10 */
471 *dvdlambda += harmonic(forceparams[type].harmonic.krA,
472 forceparams[type].harmonic.krB,
473 forceparams[type].harmonic.rA,
474 forceparams[type].harmonic.rB,
475 dr, lambda, &vbond, &fbond); /* 19 */
483 vtot += vbond; /* 1*/
484 fbond *= gmx_invsqrt(dr2); /* 6 */
488 fprintf(debug, "BONDS: dr = %10g vbond = %10g fbond = %10g\n",
494 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
497 for (m = 0; (m < DIM); m++) /* 15 */
502 fshift[ki][m] += fij;
503 fshift[CENTRAL][m] -= fij;
509 real restraint_bonds(int nbonds,
510 const t_iatom forceatoms[], const t_iparams forceparams[],
511 const rvec x[], rvec f[], rvec fshift[],
512 const t_pbc *pbc, const t_graph *g,
513 real lambda, real *dvdlambda,
514 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
515 int gmx_unused *global_atom_index)
517 int i, m, ki, ai, aj, type;
518 real dr, dr2, fbond, vbond, fij, vtot;
520 real low, dlow, up1, dup1, up2, dup2, k, dk;
528 for (i = 0; (i < nbonds); )
530 type = forceatoms[i++];
531 ai = forceatoms[i++];
532 aj = forceatoms[i++];
534 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
535 dr2 = iprod(dx, dx); /* 5 */
536 dr = dr2*gmx_invsqrt(dr2); /* 10 */
538 low = L1*forceparams[type].restraint.lowA + lambda*forceparams[type].restraint.lowB;
539 dlow = -forceparams[type].restraint.lowA + forceparams[type].restraint.lowB;
540 up1 = L1*forceparams[type].restraint.up1A + lambda*forceparams[type].restraint.up1B;
541 dup1 = -forceparams[type].restraint.up1A + forceparams[type].restraint.up1B;
542 up2 = L1*forceparams[type].restraint.up2A + lambda*forceparams[type].restraint.up2B;
543 dup2 = -forceparams[type].restraint.up2A + forceparams[type].restraint.up2B;
544 k = L1*forceparams[type].restraint.kA + lambda*forceparams[type].restraint.kB;
545 dk = -forceparams[type].restraint.kA + forceparams[type].restraint.kB;
554 *dvdlambda += 0.5*dk*drh2 - k*dlow*drh;
567 *dvdlambda += 0.5*dk*drh2 - k*dup1*drh;
572 vbond = k*(up2 - up1)*(0.5*(up2 - up1) + drh);
573 fbond = -k*(up2 - up1);
574 *dvdlambda += dk*(up2 - up1)*(0.5*(up2 - up1) + drh)
575 + k*(dup2 - dup1)*(up2 - up1 + drh)
576 - k*(up2 - up1)*dup2;
584 vtot += vbond; /* 1*/
585 fbond *= gmx_invsqrt(dr2); /* 6 */
589 fprintf(debug, "BONDS: dr = %10g vbond = %10g fbond = %10g\n",
595 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
598 for (m = 0; (m < DIM); m++) /* 15 */
603 fshift[ki][m] += fij;
604 fshift[CENTRAL][m] -= fij;
611 real polarize(int nbonds,
612 const t_iatom forceatoms[], const t_iparams forceparams[],
613 const rvec x[], rvec f[], rvec fshift[],
614 const t_pbc *pbc, const t_graph *g,
615 real lambda, real *dvdlambda,
616 const t_mdatoms *md, t_fcdata gmx_unused *fcd,
617 int gmx_unused *global_atom_index)
619 int i, m, ki, ai, aj, type;
620 real dr, dr2, fbond, vbond, fij, vtot, ksh;
625 for (i = 0; (i < nbonds); )
627 type = forceatoms[i++];
628 ai = forceatoms[i++];
629 aj = forceatoms[i++];
630 ksh = sqr(md->chargeA[aj])*ONE_4PI_EPS0/forceparams[type].polarize.alpha;
633 fprintf(debug, "POL: local ai = %d aj = %d ksh = %.3f\n", ai, aj, ksh);
636 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
637 dr2 = iprod(dx, dx); /* 5 */
638 dr = dr2*gmx_invsqrt(dr2); /* 10 */
640 *dvdlambda += harmonic(ksh, ksh, 0, 0, dr, lambda, &vbond, &fbond); /* 19 */
647 vtot += vbond; /* 1*/
648 fbond *= gmx_invsqrt(dr2); /* 6 */
652 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
655 for (m = 0; (m < DIM); m++) /* 15 */
660 fshift[ki][m] += fij;
661 fshift[CENTRAL][m] -= fij;
667 real anharm_polarize(int nbonds,
668 const t_iatom forceatoms[], const t_iparams forceparams[],
669 const rvec x[], rvec f[], rvec fshift[],
670 const t_pbc *pbc, const t_graph *g,
671 real lambda, real *dvdlambda,
672 const t_mdatoms *md, t_fcdata gmx_unused *fcd,
673 int gmx_unused *global_atom_index)
675 int i, m, ki, ai, aj, type;
676 real dr, dr2, fbond, vbond, fij, vtot, ksh, khyp, drcut, ddr, ddr3;
681 for (i = 0; (i < nbonds); )
683 type = forceatoms[i++];
684 ai = forceatoms[i++];
685 aj = forceatoms[i++];
686 ksh = sqr(md->chargeA[aj])*ONE_4PI_EPS0/forceparams[type].anharm_polarize.alpha; /* 7*/
687 khyp = forceparams[type].anharm_polarize.khyp;
688 drcut = forceparams[type].anharm_polarize.drcut;
691 fprintf(debug, "POL: local ai = %d aj = %d ksh = %.3f\n", ai, aj, ksh);
694 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
695 dr2 = iprod(dx, dx); /* 5 */
696 dr = dr2*gmx_invsqrt(dr2); /* 10 */
698 *dvdlambda += harmonic(ksh, ksh, 0, 0, dr, lambda, &vbond, &fbond); /* 19 */
709 vbond += khyp*ddr*ddr3;
710 fbond -= 4*khyp*ddr3;
712 fbond *= gmx_invsqrt(dr2); /* 6 */
713 vtot += vbond; /* 1*/
717 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
720 for (m = 0; (m < DIM); m++) /* 15 */
725 fshift[ki][m] += fij;
726 fshift[CENTRAL][m] -= fij;
732 real water_pol(int nbonds,
733 const t_iatom forceatoms[], const t_iparams forceparams[],
734 const rvec x[], rvec f[], rvec gmx_unused fshift[],
735 const t_pbc gmx_unused *pbc, const t_graph gmx_unused *g,
736 real gmx_unused lambda, real gmx_unused *dvdlambda,
737 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
738 int gmx_unused *global_atom_index)
740 /* This routine implements anisotropic polarizibility for water, through
741 * a shell connected to a dummy with spring constant that differ in the
742 * three spatial dimensions in the molecular frame.
744 int i, m, aO, aH1, aH2, aD, aS, type, type0, ki;
746 rvec dOH1, dOH2, dHH, dOD, dDS, nW, kk, dx, kdx, proj;
750 real vtot, fij, r_HH, r_OD, r_nW, tx, ty, tz, qS;
755 type0 = forceatoms[0];
757 qS = md->chargeA[aS];
758 kk[XX] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_x;
759 kk[YY] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_y;
760 kk[ZZ] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_z;
761 r_HH = 1.0/forceparams[type0].wpol.rHH;
764 fprintf(debug, "WPOL: qS = %10.5f aS = %5d\n", qS, aS);
765 fprintf(debug, "WPOL: kk = %10.3f %10.3f %10.3f\n",
766 kk[XX], kk[YY], kk[ZZ]);
767 fprintf(debug, "WPOL: rOH = %10.3f rHH = %10.3f rOD = %10.3f\n",
768 forceparams[type0].wpol.rOH,
769 forceparams[type0].wpol.rHH,
770 forceparams[type0].wpol.rOD);
772 for (i = 0; (i < nbonds); i += 6)
774 type = forceatoms[i];
777 gmx_fatal(FARGS, "Sorry, type = %d, type0 = %d, file = %s, line = %d",
778 type, type0, __FILE__, __LINE__);
780 aO = forceatoms[i+1];
781 aH1 = forceatoms[i+2];
782 aH2 = forceatoms[i+3];
783 aD = forceatoms[i+4];
784 aS = forceatoms[i+5];
786 /* Compute vectors describing the water frame */
787 pbc_rvec_sub(pbc, x[aH1], x[aO], dOH1);
788 pbc_rvec_sub(pbc, x[aH2], x[aO], dOH2);
789 pbc_rvec_sub(pbc, x[aH2], x[aH1], dHH);
790 pbc_rvec_sub(pbc, x[aD], x[aO], dOD);
791 ki = pbc_rvec_sub(pbc, x[aS], x[aD], dDS);
792 cprod(dOH1, dOH2, nW);
794 /* Compute inverse length of normal vector
795 * (this one could be precomputed, but I'm too lazy now)
797 r_nW = gmx_invsqrt(iprod(nW, nW));
798 /* This is for precision, but does not make a big difference,
801 r_OD = gmx_invsqrt(iprod(dOD, dOD));
803 /* Normalize the vectors in the water frame */
805 svmul(r_HH, dHH, dHH);
806 svmul(r_OD, dOD, dOD);
808 /* Compute displacement of shell along components of the vector */
809 dx[ZZ] = iprod(dDS, dOD);
810 /* Compute projection on the XY plane: dDS - dx[ZZ]*dOD */
811 for (m = 0; (m < DIM); m++)
813 proj[m] = dDS[m]-dx[ZZ]*dOD[m];
816 /*dx[XX] = iprod(dDS,nW);
817 dx[YY] = iprod(dDS,dHH);*/
818 dx[XX] = iprod(proj, nW);
819 for (m = 0; (m < DIM); m++)
821 proj[m] -= dx[XX]*nW[m];
823 dx[YY] = iprod(proj, dHH);
828 fprintf(debug, "WPOL: dx2=%10g dy2=%10g dz2=%10g sum=%10g dDS^2=%10g\n",
829 sqr(dx[XX]), sqr(dx[YY]), sqr(dx[ZZ]), iprod(dx, dx), iprod(dDS, dDS));
830 fprintf(debug, "WPOL: dHH=(%10g,%10g,%10g)\n", dHH[XX], dHH[YY], dHH[ZZ]);
831 fprintf(debug, "WPOL: dOD=(%10g,%10g,%10g), 1/r_OD = %10g\n",
832 dOD[XX], dOD[YY], dOD[ZZ], 1/r_OD);
833 fprintf(debug, "WPOL: nW =(%10g,%10g,%10g), 1/r_nW = %10g\n",
834 nW[XX], nW[YY], nW[ZZ], 1/r_nW);
835 fprintf(debug, "WPOL: dx =%10g, dy =%10g, dz =%10g\n",
836 dx[XX], dx[YY], dx[ZZ]);
837 fprintf(debug, "WPOL: dDSx=%10g, dDSy=%10g, dDSz=%10g\n",
838 dDS[XX], dDS[YY], dDS[ZZ]);
841 /* Now compute the forces and energy */
842 kdx[XX] = kk[XX]*dx[XX];
843 kdx[YY] = kk[YY]*dx[YY];
844 kdx[ZZ] = kk[ZZ]*dx[ZZ];
845 vtot += iprod(dx, kdx);
849 ivec_sub(SHIFT_IVEC(g, aS), SHIFT_IVEC(g, aD), dt);
853 for (m = 0; (m < DIM); m++)
855 /* This is a tensor operation but written out for speed */
865 fshift[ki][m] += fij;
866 fshift[CENTRAL][m] -= fij;
871 fprintf(debug, "WPOL: vwpol=%g\n", 0.5*iprod(dx, kdx));
872 fprintf(debug, "WPOL: df = (%10g, %10g, %10g)\n", df[XX], df[YY], df[ZZ]);
880 static real do_1_thole(const rvec xi, const rvec xj, rvec fi, rvec fj,
881 const t_pbc *pbc, real qq,
882 rvec fshift[], real afac)
885 real r12sq, r12_1, r12bar, v0, v1, fscal, ebar, fff;
888 t = pbc_rvec_sub(pbc, xi, xj, r12); /* 3 */
890 r12sq = iprod(r12, r12); /* 5 */
891 r12_1 = gmx_invsqrt(r12sq); /* 5 */
892 r12bar = afac/r12_1; /* 5 */
893 v0 = qq*ONE_4PI_EPS0*r12_1; /* 2 */
894 ebar = exp(-r12bar); /* 5 */
895 v1 = (1-(1+0.5*r12bar)*ebar); /* 4 */
896 fscal = ((v0*r12_1)*v1 - v0*0.5*afac*ebar*(r12bar+1))*r12_1; /* 9 */
899 fprintf(debug, "THOLE: v0 = %.3f v1 = %.3f r12= % .3f r12bar = %.3f fscal = %.3f ebar = %.3f\n", v0, v1, 1/r12_1, r12bar, fscal, ebar);
902 for (m = 0; (m < DIM); m++)
908 fshift[CENTRAL][m] -= fff;
911 return v0*v1; /* 1 */
915 real thole_pol(int nbonds,
916 const t_iatom forceatoms[], const t_iparams forceparams[],
917 const rvec x[], rvec f[], rvec fshift[],
918 const t_pbc *pbc, const t_graph gmx_unused *g,
919 real gmx_unused lambda, real gmx_unused *dvdlambda,
920 const t_mdatoms *md, t_fcdata gmx_unused *fcd,
921 int gmx_unused *global_atom_index)
923 /* Interaction between two pairs of particles with opposite charge */
924 int i, type, a1, da1, a2, da2;
925 real q1, q2, qq, a, al1, al2, afac;
927 const real minusOneOnSix = -1.0/6.0;
929 for (i = 0; (i < nbonds); )
931 type = forceatoms[i++];
932 a1 = forceatoms[i++];
933 da1 = forceatoms[i++];
934 a2 = forceatoms[i++];
935 da2 = forceatoms[i++];
936 q1 = md->chargeA[da1];
937 q2 = md->chargeA[da2];
938 a = forceparams[type].thole.a;
939 al1 = forceparams[type].thole.alpha1;
940 al2 = forceparams[type].thole.alpha2;
942 afac = a*pow(al1*al2, minusOneOnSix);
943 V += do_1_thole(x[a1], x[a2], f[a1], f[a2], pbc, qq, fshift, afac);
944 V += do_1_thole(x[da1], x[a2], f[da1], f[a2], pbc, -qq, fshift, afac);
945 V += do_1_thole(x[a1], x[da2], f[a1], f[da2], pbc, -qq, fshift, afac);
946 V += do_1_thole(x[da1], x[da2], f[da1], f[da2], pbc, qq, fshift, afac);
952 real bond_angle(const rvec xi, const rvec xj, const rvec xk, const t_pbc *pbc,
953 rvec r_ij, rvec r_kj, real *costh,
955 /* Return value is the angle between the bonds i-j and j-k */
960 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
961 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
963 *costh = cos_angle(r_ij, r_kj); /* 25 */
964 th = acos(*costh); /* 10 */
969 real angles(int nbonds,
970 const t_iatom forceatoms[], const t_iparams forceparams[],
971 const rvec x[], rvec f[], rvec fshift[],
972 const t_pbc *pbc, const t_graph *g,
973 real lambda, real *dvdlambda,
974 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
975 int gmx_unused *global_atom_index)
977 int i, ai, aj, ak, t1, t2, type;
979 real cos_theta, cos_theta2, theta, dVdt, va, vtot;
980 ivec jt, dt_ij, dt_kj;
983 for (i = 0; i < nbonds; )
985 type = forceatoms[i++];
986 ai = forceatoms[i++];
987 aj = forceatoms[i++];
988 ak = forceatoms[i++];
990 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
991 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
993 *dvdlambda += harmonic(forceparams[type].harmonic.krA,
994 forceparams[type].harmonic.krB,
995 forceparams[type].harmonic.rA*DEG2RAD,
996 forceparams[type].harmonic.rB*DEG2RAD,
997 theta, lambda, &va, &dVdt); /* 21 */
1000 cos_theta2 = sqr(cos_theta);
1007 real nrkj_1, nrij_1;
1010 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1011 sth = st*cos_theta; /* 1 */
1015 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1016 theta*RAD2DEG, va, dVdt);
1019 nrij2 = iprod(r_ij, r_ij); /* 5 */
1020 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1022 nrij_1 = gmx_invsqrt(nrij2); /* 10 */
1023 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1025 cik = st*nrij_1*nrkj_1; /* 2 */
1026 cii = sth*nrij_1*nrij_1; /* 2 */
1027 ckk = sth*nrkj_1*nrkj_1; /* 2 */
1029 for (m = 0; m < DIM; m++)
1031 f_i[m] = -(cik*r_kj[m] - cii*r_ij[m]);
1032 f_k[m] = -(cik*r_ij[m] - ckk*r_kj[m]);
1033 f_j[m] = -f_i[m] - f_k[m];
1040 copy_ivec(SHIFT_IVEC(g, aj), jt);
1042 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1043 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1044 t1 = IVEC2IS(dt_ij);
1045 t2 = IVEC2IS(dt_kj);
1047 rvec_inc(fshift[t1], f_i);
1048 rvec_inc(fshift[CENTRAL], f_j);
1049 rvec_inc(fshift[t2], f_k);
1056 #ifdef GMX_SIMD_HAVE_REAL
1058 /* As angles, but using SIMD to calculate many dihedrals at once.
1059 * This routines does not calculate energies and shift forces.
1061 static gmx_inline void
1062 angles_noener_simd(int nbonds,
1063 const t_iatom forceatoms[], const t_iparams forceparams[],
1064 const rvec x[], rvec f[],
1065 const t_pbc *pbc, const t_graph gmx_unused *g,
1066 real gmx_unused lambda,
1067 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1068 int gmx_unused *global_atom_index)
1072 int type, ai[GMX_SIMD_REAL_WIDTH], aj[GMX_SIMD_REAL_WIDTH];
1073 int ak[GMX_SIMD_REAL_WIDTH];
1074 real coeff_array[2*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *coeff;
1075 real dr_array[2*DIM*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *dr;
1076 real f_buf_array[6*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *f_buf;
1077 gmx_simd_real_t k_S, theta0_S;
1078 gmx_simd_real_t rijx_S, rijy_S, rijz_S;
1079 gmx_simd_real_t rkjx_S, rkjy_S, rkjz_S;
1080 gmx_simd_real_t one_S;
1081 gmx_simd_real_t min_one_plus_eps_S;
1082 gmx_simd_real_t rij_rkj_S;
1083 gmx_simd_real_t nrij2_S, nrij_1_S;
1084 gmx_simd_real_t nrkj2_S, nrkj_1_S;
1085 gmx_simd_real_t cos_S, invsin_S;
1086 gmx_simd_real_t theta_S;
1087 gmx_simd_real_t st_S, sth_S;
1088 gmx_simd_real_t cik_S, cii_S, ckk_S;
1089 gmx_simd_real_t f_ix_S, f_iy_S, f_iz_S;
1090 gmx_simd_real_t f_kx_S, f_ky_S, f_kz_S;
1091 pbc_simd_t pbc_simd;
1093 /* Ensure register memory alignment */
1094 coeff = gmx_simd_align_r(coeff_array);
1095 dr = gmx_simd_align_r(dr_array);
1096 f_buf = gmx_simd_align_r(f_buf_array);
1098 set_pbc_simd(pbc, &pbc_simd);
1100 one_S = gmx_simd_set1_r(1.0);
1102 /* The smallest number > -1 */
1103 min_one_plus_eps_S = gmx_simd_set1_r(-1.0 + 2*GMX_REAL_EPS);
1105 /* nbonds is the number of angles times nfa1, here we step GMX_SIMD_REAL_WIDTH angles */
1106 for (i = 0; (i < nbonds); i += GMX_SIMD_REAL_WIDTH*nfa1)
1108 /* Collect atoms for GMX_SIMD_REAL_WIDTH angles.
1109 * iu indexes into forceatoms, we should not let iu go beyond nbonds.
1112 for (s = 0; s < GMX_SIMD_REAL_WIDTH; s++)
1114 type = forceatoms[iu];
1115 ai[s] = forceatoms[iu+1];
1116 aj[s] = forceatoms[iu+2];
1117 ak[s] = forceatoms[iu+3];
1119 coeff[s] = forceparams[type].harmonic.krA;
1120 coeff[GMX_SIMD_REAL_WIDTH+s] = forceparams[type].harmonic.rA*DEG2RAD;
1122 /* If you can't use pbc_dx_simd below for PBC, e.g. because
1123 * you can't round in SIMD, use pbc_rvec_sub here.
1125 /* Store the non PBC corrected distances packed and aligned */
1126 for (m = 0; m < DIM; m++)
1128 dr[s + m *GMX_SIMD_REAL_WIDTH] = x[ai[s]][m] - x[aj[s]][m];
1129 dr[s + (DIM+m)*GMX_SIMD_REAL_WIDTH] = x[ak[s]][m] - x[aj[s]][m];
1132 /* At the end fill the arrays with identical entries */
1133 if (iu + nfa1 < nbonds)
1139 k_S = gmx_simd_load_r(coeff);
1140 theta0_S = gmx_simd_load_r(coeff+GMX_SIMD_REAL_WIDTH);
1142 rijx_S = gmx_simd_load_r(dr + 0*GMX_SIMD_REAL_WIDTH);
1143 rijy_S = gmx_simd_load_r(dr + 1*GMX_SIMD_REAL_WIDTH);
1144 rijz_S = gmx_simd_load_r(dr + 2*GMX_SIMD_REAL_WIDTH);
1145 rkjx_S = gmx_simd_load_r(dr + 3*GMX_SIMD_REAL_WIDTH);
1146 rkjy_S = gmx_simd_load_r(dr + 4*GMX_SIMD_REAL_WIDTH);
1147 rkjz_S = gmx_simd_load_r(dr + 5*GMX_SIMD_REAL_WIDTH);
1149 pbc_dx_simd(&rijx_S, &rijy_S, &rijz_S, &pbc_simd);
1150 pbc_dx_simd(&rkjx_S, &rkjy_S, &rkjz_S, &pbc_simd);
1152 rij_rkj_S = gmx_simd_iprod_r(rijx_S, rijy_S, rijz_S,
1153 rkjx_S, rkjy_S, rkjz_S);
1155 nrij2_S = gmx_simd_norm2_r(rijx_S, rijy_S, rijz_S);
1156 nrkj2_S = gmx_simd_norm2_r(rkjx_S, rkjy_S, rkjz_S);
1158 nrij_1_S = gmx_simd_invsqrt_r(nrij2_S);
1159 nrkj_1_S = gmx_simd_invsqrt_r(nrkj2_S);
1161 cos_S = gmx_simd_mul_r(rij_rkj_S, gmx_simd_mul_r(nrij_1_S, nrkj_1_S));
1163 /* To allow for 180 degrees, we take the max of cos and -1 + 1bit,
1164 * so we can safely get the 1/sin from 1/sqrt(1 - cos^2).
1165 * This also ensures that rounding errors would cause the argument
1166 * of gmx_simd_acos_r to be < -1.
1167 * Note that we do not take precautions for cos(0)=1, so the outer
1168 * atoms in an angle should not be on top of each other.
1170 cos_S = gmx_simd_max_r(cos_S, min_one_plus_eps_S);
1172 theta_S = gmx_simd_acos_r(cos_S);
1174 invsin_S = gmx_simd_invsqrt_r(gmx_simd_sub_r(one_S, gmx_simd_mul_r(cos_S, cos_S)));
1176 st_S = gmx_simd_mul_r(gmx_simd_mul_r(k_S, gmx_simd_sub_r(theta0_S, theta_S)),
1178 sth_S = gmx_simd_mul_r(st_S, cos_S);
1180 cik_S = gmx_simd_mul_r(st_S, gmx_simd_mul_r(nrij_1_S, nrkj_1_S));
1181 cii_S = gmx_simd_mul_r(sth_S, gmx_simd_mul_r(nrij_1_S, nrij_1_S));
1182 ckk_S = gmx_simd_mul_r(sth_S, gmx_simd_mul_r(nrkj_1_S, nrkj_1_S));
1184 f_ix_S = gmx_simd_mul_r(cii_S, rijx_S);
1185 f_ix_S = gmx_simd_fnmadd_r(cik_S, rkjx_S, f_ix_S);
1186 f_iy_S = gmx_simd_mul_r(cii_S, rijy_S);
1187 f_iy_S = gmx_simd_fnmadd_r(cik_S, rkjy_S, f_iy_S);
1188 f_iz_S = gmx_simd_mul_r(cii_S, rijz_S);
1189 f_iz_S = gmx_simd_fnmadd_r(cik_S, rkjz_S, f_iz_S);
1190 f_kx_S = gmx_simd_mul_r(ckk_S, rkjx_S);
1191 f_kx_S = gmx_simd_fnmadd_r(cik_S, rijx_S, f_kx_S);
1192 f_ky_S = gmx_simd_mul_r(ckk_S, rkjy_S);
1193 f_ky_S = gmx_simd_fnmadd_r(cik_S, rijy_S, f_ky_S);
1194 f_kz_S = gmx_simd_mul_r(ckk_S, rkjz_S);
1195 f_kz_S = gmx_simd_fnmadd_r(cik_S, rijz_S, f_kz_S);
1197 gmx_simd_store_r(f_buf + 0*GMX_SIMD_REAL_WIDTH, f_ix_S);
1198 gmx_simd_store_r(f_buf + 1*GMX_SIMD_REAL_WIDTH, f_iy_S);
1199 gmx_simd_store_r(f_buf + 2*GMX_SIMD_REAL_WIDTH, f_iz_S);
1200 gmx_simd_store_r(f_buf + 3*GMX_SIMD_REAL_WIDTH, f_kx_S);
1201 gmx_simd_store_r(f_buf + 4*GMX_SIMD_REAL_WIDTH, f_ky_S);
1202 gmx_simd_store_r(f_buf + 5*GMX_SIMD_REAL_WIDTH, f_kz_S);
1208 for (m = 0; m < DIM; m++)
1210 f[ai[s]][m] += f_buf[s + m*GMX_SIMD_REAL_WIDTH];
1211 f[aj[s]][m] -= f_buf[s + m*GMX_SIMD_REAL_WIDTH] + f_buf[s + (DIM+m)*GMX_SIMD_REAL_WIDTH];
1212 f[ak[s]][m] += f_buf[s + (DIM+m)*GMX_SIMD_REAL_WIDTH];
1217 while (s < GMX_SIMD_REAL_WIDTH && iu < nbonds);
1221 #endif /* GMX_SIMD_HAVE_REAL */
1223 real linear_angles(int nbonds,
1224 const t_iatom forceatoms[], const t_iparams forceparams[],
1225 const rvec x[], rvec f[], rvec fshift[],
1226 const t_pbc *pbc, const t_graph *g,
1227 real lambda, real *dvdlambda,
1228 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1229 int gmx_unused *global_atom_index)
1231 int i, m, ai, aj, ak, t1, t2, type;
1233 real L1, kA, kB, aA, aB, dr, dr2, va, vtot, a, b, klin;
1234 ivec jt, dt_ij, dt_kj;
1235 rvec r_ij, r_kj, r_ik, dx;
1239 for (i = 0; (i < nbonds); )
1241 type = forceatoms[i++];
1242 ai = forceatoms[i++];
1243 aj = forceatoms[i++];
1244 ak = forceatoms[i++];
1246 kA = forceparams[type].linangle.klinA;
1247 kB = forceparams[type].linangle.klinB;
1248 klin = L1*kA + lambda*kB;
1250 aA = forceparams[type].linangle.aA;
1251 aB = forceparams[type].linangle.aB;
1252 a = L1*aA+lambda*aB;
1255 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
1256 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
1257 rvec_sub(r_ij, r_kj, r_ik);
1260 for (m = 0; (m < DIM); m++)
1262 dr = -a * r_ij[m] - b * r_kj[m];
1267 f_j[m] = -(f_i[m]+f_k[m]);
1273 *dvdlambda += 0.5*(kB-kA)*dr2 + klin*(aB-aA)*iprod(dx, r_ik);
1279 copy_ivec(SHIFT_IVEC(g, aj), jt);
1281 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1282 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1283 t1 = IVEC2IS(dt_ij);
1284 t2 = IVEC2IS(dt_kj);
1286 rvec_inc(fshift[t1], f_i);
1287 rvec_inc(fshift[CENTRAL], f_j);
1288 rvec_inc(fshift[t2], f_k);
1293 real urey_bradley(int nbonds,
1294 const t_iatom forceatoms[], const t_iparams forceparams[],
1295 const rvec x[], rvec f[], rvec fshift[],
1296 const t_pbc *pbc, const t_graph *g,
1297 real lambda, real *dvdlambda,
1298 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1299 int gmx_unused *global_atom_index)
1301 int i, m, ai, aj, ak, t1, t2, type, ki;
1302 rvec r_ij, r_kj, r_ik;
1303 real cos_theta, cos_theta2, theta;
1304 real dVdt, va, vtot, dr, dr2, vbond, fbond, fik;
1305 real kthA, th0A, kUBA, r13A, kthB, th0B, kUBB, r13B;
1306 ivec jt, dt_ij, dt_kj, dt_ik;
1309 for (i = 0; (i < nbonds); )
1311 type = forceatoms[i++];
1312 ai = forceatoms[i++];
1313 aj = forceatoms[i++];
1314 ak = forceatoms[i++];
1315 th0A = forceparams[type].u_b.thetaA*DEG2RAD;
1316 kthA = forceparams[type].u_b.kthetaA;
1317 r13A = forceparams[type].u_b.r13A;
1318 kUBA = forceparams[type].u_b.kUBA;
1319 th0B = forceparams[type].u_b.thetaB*DEG2RAD;
1320 kthB = forceparams[type].u_b.kthetaB;
1321 r13B = forceparams[type].u_b.r13B;
1322 kUBB = forceparams[type].u_b.kUBB;
1324 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
1325 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
1327 *dvdlambda += harmonic(kthA, kthB, th0A, th0B, theta, lambda, &va, &dVdt); /* 21 */
1330 ki = pbc_rvec_sub(pbc, x[ai], x[ak], r_ik); /* 3 */
1331 dr2 = iprod(r_ik, r_ik); /* 5 */
1332 dr = dr2*gmx_invsqrt(dr2); /* 10 */
1334 *dvdlambda += harmonic(kUBA, kUBB, r13A, r13B, dr, lambda, &vbond, &fbond); /* 19 */
1336 cos_theta2 = sqr(cos_theta); /* 1 */
1344 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1345 sth = st*cos_theta; /* 1 */
1349 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1350 theta*RAD2DEG, va, dVdt);
1353 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1354 nrij2 = iprod(r_ij, r_ij);
1356 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
1357 cii = sth/nrij2; /* 10 */
1358 ckk = sth/nrkj2; /* 10 */
1360 for (m = 0; (m < DIM); m++) /* 39 */
1362 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
1363 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
1364 f_j[m] = -f_i[m]-f_k[m];
1371 copy_ivec(SHIFT_IVEC(g, aj), jt);
1373 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1374 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1375 t1 = IVEC2IS(dt_ij);
1376 t2 = IVEC2IS(dt_kj);
1378 rvec_inc(fshift[t1], f_i);
1379 rvec_inc(fshift[CENTRAL], f_j);
1380 rvec_inc(fshift[t2], f_k);
1382 /* Time for the bond calculations */
1388 vtot += vbond; /* 1*/
1389 fbond *= gmx_invsqrt(dr2); /* 6 */
1393 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ak), dt_ik);
1394 ki = IVEC2IS(dt_ik);
1396 for (m = 0; (m < DIM); m++) /* 15 */
1398 fik = fbond*r_ik[m];
1401 fshift[ki][m] += fik;
1402 fshift[CENTRAL][m] -= fik;
1408 real quartic_angles(int nbonds,
1409 const t_iatom forceatoms[], const t_iparams forceparams[],
1410 const rvec x[], rvec f[], rvec fshift[],
1411 const t_pbc *pbc, const t_graph *g,
1412 real gmx_unused lambda, real gmx_unused *dvdlambda,
1413 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1414 int gmx_unused *global_atom_index)
1416 int i, j, ai, aj, ak, t1, t2, type;
1418 real cos_theta, cos_theta2, theta, dt, dVdt, va, dtp, c, vtot;
1419 ivec jt, dt_ij, dt_kj;
1422 for (i = 0; (i < nbonds); )
1424 type = forceatoms[i++];
1425 ai = forceatoms[i++];
1426 aj = forceatoms[i++];
1427 ak = forceatoms[i++];
1429 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
1430 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
1432 dt = theta - forceparams[type].qangle.theta*DEG2RAD; /* 2 */
1435 va = forceparams[type].qangle.c[0];
1437 for (j = 1; j <= 4; j++)
1439 c = forceparams[type].qangle.c[j];
1448 cos_theta2 = sqr(cos_theta); /* 1 */
1457 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1458 sth = st*cos_theta; /* 1 */
1462 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1463 theta*RAD2DEG, va, dVdt);
1466 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1467 nrij2 = iprod(r_ij, r_ij);
1469 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
1470 cii = sth/nrij2; /* 10 */
1471 ckk = sth/nrkj2; /* 10 */
1473 for (m = 0; (m < DIM); m++) /* 39 */
1475 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
1476 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
1477 f_j[m] = -f_i[m]-f_k[m];
1484 copy_ivec(SHIFT_IVEC(g, aj), jt);
1486 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1487 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1488 t1 = IVEC2IS(dt_ij);
1489 t2 = IVEC2IS(dt_kj);
1491 rvec_inc(fshift[t1], f_i);
1492 rvec_inc(fshift[CENTRAL], f_j);
1493 rvec_inc(fshift[t2], f_k);
1499 real dih_angle(const rvec xi, const rvec xj, const rvec xk, const rvec xl,
1501 rvec r_ij, rvec r_kj, rvec r_kl, rvec m, rvec n,
1502 real *sign, int *t1, int *t2, int *t3)
1506 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
1507 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
1508 *t3 = pbc_rvec_sub(pbc, xk, xl, r_kl); /* 3 */
1510 cprod(r_ij, r_kj, m); /* 9 */
1511 cprod(r_kj, r_kl, n); /* 9 */
1512 phi = gmx_angle(m, n); /* 49 (assuming 25 for atan2) */
1513 ipr = iprod(r_ij, n); /* 5 */
1514 (*sign) = (ipr < 0.0) ? -1.0 : 1.0;
1515 phi = (*sign)*phi; /* 1 */
1521 #ifdef GMX_SIMD_HAVE_REAL
1523 /* As dih_angle above, but calculates 4 dihedral angles at once using SIMD,
1524 * also calculates the pre-factor required for the dihedral force update.
1525 * Note that bv and buf should be register aligned.
1527 static gmx_inline void
1528 dih_angle_simd(const rvec *x,
1529 const int *ai, const int *aj, const int *ak, const int *al,
1530 const pbc_simd_t *pbc,
1532 gmx_simd_real_t *phi_S,
1533 gmx_simd_real_t *mx_S, gmx_simd_real_t *my_S, gmx_simd_real_t *mz_S,
1534 gmx_simd_real_t *nx_S, gmx_simd_real_t *ny_S, gmx_simd_real_t *nz_S,
1535 gmx_simd_real_t *nrkj_m2_S,
1536 gmx_simd_real_t *nrkj_n2_S,
1541 gmx_simd_real_t rijx_S, rijy_S, rijz_S;
1542 gmx_simd_real_t rkjx_S, rkjy_S, rkjz_S;
1543 gmx_simd_real_t rklx_S, rkly_S, rklz_S;
1544 gmx_simd_real_t cx_S, cy_S, cz_S;
1545 gmx_simd_real_t cn_S;
1546 gmx_simd_real_t s_S;
1547 gmx_simd_real_t ipr_S;
1548 gmx_simd_real_t iprm_S, iprn_S;
1549 gmx_simd_real_t nrkj2_S, nrkj_1_S, nrkj_2_S, nrkj_S;
1550 gmx_simd_real_t toler_S;
1551 gmx_simd_real_t p_S, q_S;
1552 gmx_simd_real_t nrkj2_min_S;
1553 gmx_simd_real_t real_eps_S;
1555 /* Used to avoid division by zero.
1556 * We take into acount that we multiply the result by real_eps_S.
1558 nrkj2_min_S = gmx_simd_set1_r(GMX_REAL_MIN/(2*GMX_REAL_EPS));
1560 /* The value of the last significant bit (GMX_REAL_EPS is half of that) */
1561 real_eps_S = gmx_simd_set1_r(2*GMX_REAL_EPS);
1563 for (s = 0; s < GMX_SIMD_REAL_WIDTH; s++)
1565 /* If you can't use pbc_dx_simd below for PBC, e.g. because
1566 * you can't round in SIMD, use pbc_rvec_sub here.
1568 for (m = 0; m < DIM; m++)
1570 dr[s + (0*DIM + m)*GMX_SIMD_REAL_WIDTH] = x[ai[s]][m] - x[aj[s]][m];
1571 dr[s + (1*DIM + m)*GMX_SIMD_REAL_WIDTH] = x[ak[s]][m] - x[aj[s]][m];
1572 dr[s + (2*DIM + m)*GMX_SIMD_REAL_WIDTH] = x[ak[s]][m] - x[al[s]][m];
1576 rijx_S = gmx_simd_load_r(dr + 0*GMX_SIMD_REAL_WIDTH);
1577 rijy_S = gmx_simd_load_r(dr + 1*GMX_SIMD_REAL_WIDTH);
1578 rijz_S = gmx_simd_load_r(dr + 2*GMX_SIMD_REAL_WIDTH);
1579 rkjx_S = gmx_simd_load_r(dr + 3*GMX_SIMD_REAL_WIDTH);
1580 rkjy_S = gmx_simd_load_r(dr + 4*GMX_SIMD_REAL_WIDTH);
1581 rkjz_S = gmx_simd_load_r(dr + 5*GMX_SIMD_REAL_WIDTH);
1582 rklx_S = gmx_simd_load_r(dr + 6*GMX_SIMD_REAL_WIDTH);
1583 rkly_S = gmx_simd_load_r(dr + 7*GMX_SIMD_REAL_WIDTH);
1584 rklz_S = gmx_simd_load_r(dr + 8*GMX_SIMD_REAL_WIDTH);
1586 pbc_dx_simd(&rijx_S, &rijy_S, &rijz_S, pbc);
1587 pbc_dx_simd(&rkjx_S, &rkjy_S, &rkjz_S, pbc);
1588 pbc_dx_simd(&rklx_S, &rkly_S, &rklz_S, pbc);
1590 gmx_simd_cprod_r(rijx_S, rijy_S, rijz_S,
1591 rkjx_S, rkjy_S, rkjz_S,
1594 gmx_simd_cprod_r(rkjx_S, rkjy_S, rkjz_S,
1595 rklx_S, rkly_S, rklz_S,
1598 gmx_simd_cprod_r(*mx_S, *my_S, *mz_S,
1599 *nx_S, *ny_S, *nz_S,
1600 &cx_S, &cy_S, &cz_S);
1602 cn_S = gmx_simd_sqrt_r(gmx_simd_norm2_r(cx_S, cy_S, cz_S));
1604 s_S = gmx_simd_iprod_r(*mx_S, *my_S, *mz_S, *nx_S, *ny_S, *nz_S);
1606 /* Determine the dihedral angle, the sign might need correction */
1607 *phi_S = gmx_simd_atan2_r(cn_S, s_S);
1609 ipr_S = gmx_simd_iprod_r(rijx_S, rijy_S, rijz_S,
1610 *nx_S, *ny_S, *nz_S);
1612 iprm_S = gmx_simd_norm2_r(*mx_S, *my_S, *mz_S);
1613 iprn_S = gmx_simd_norm2_r(*nx_S, *ny_S, *nz_S);
1615 nrkj2_S = gmx_simd_norm2_r(rkjx_S, rkjy_S, rkjz_S);
1617 /* Avoid division by zero. When zero, the result is multiplied by 0
1618 * anyhow, so the 3 max below do not affect the final result.
1620 nrkj2_S = gmx_simd_max_r(nrkj2_S, nrkj2_min_S);
1621 nrkj_1_S = gmx_simd_invsqrt_r(nrkj2_S);
1622 nrkj_2_S = gmx_simd_mul_r(nrkj_1_S, nrkj_1_S);
1623 nrkj_S = gmx_simd_mul_r(nrkj2_S, nrkj_1_S);
1625 toler_S = gmx_simd_mul_r(nrkj2_S, real_eps_S);
1627 /* Here the plain-C code uses a conditional, but we can't do that in SIMD.
1628 * So we take a max with the tolerance instead. Since we multiply with
1629 * m or n later, the max does not affect the results.
1631 iprm_S = gmx_simd_max_r(iprm_S, toler_S);
1632 iprn_S = gmx_simd_max_r(iprn_S, toler_S);
1633 *nrkj_m2_S = gmx_simd_mul_r(nrkj_S, gmx_simd_inv_r(iprm_S));
1634 *nrkj_n2_S = gmx_simd_mul_r(nrkj_S, gmx_simd_inv_r(iprn_S));
1636 /* Set sign of phi_S with the sign of ipr_S; phi_S is currently positive */
1637 *phi_S = gmx_simd_xor_sign_r(*phi_S, ipr_S);
1638 p_S = gmx_simd_iprod_r(rijx_S, rijy_S, rijz_S,
1639 rkjx_S, rkjy_S, rkjz_S);
1640 p_S = gmx_simd_mul_r(p_S, nrkj_2_S);
1642 q_S = gmx_simd_iprod_r(rklx_S, rkly_S, rklz_S,
1643 rkjx_S, rkjy_S, rkjz_S);
1644 q_S = gmx_simd_mul_r(q_S, nrkj_2_S);
1646 gmx_simd_store_r(p, p_S);
1647 gmx_simd_store_r(q, q_S);
1650 #endif /* GMX_SIMD_HAVE_REAL */
1653 void do_dih_fup(int i, int j, int k, int l, real ddphi,
1654 rvec r_ij, rvec r_kj, rvec r_kl,
1655 rvec m, rvec n, rvec f[], rvec fshift[],
1656 const t_pbc *pbc, const t_graph *g,
1657 const rvec x[], int t1, int t2, int t3)
1660 rvec f_i, f_j, f_k, f_l;
1661 rvec uvec, vvec, svec, dx_jl;
1662 real iprm, iprn, nrkj, nrkj2, nrkj_1, nrkj_2;
1663 real a, b, p, q, toler;
1664 ivec jt, dt_ij, dt_kj, dt_lj;
1666 iprm = iprod(m, m); /* 5 */
1667 iprn = iprod(n, n); /* 5 */
1668 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1669 toler = nrkj2*GMX_REAL_EPS;
1670 if ((iprm > toler) && (iprn > toler))
1672 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1673 nrkj_2 = nrkj_1*nrkj_1; /* 1 */
1674 nrkj = nrkj2*nrkj_1; /* 1 */
1675 a = -ddphi*nrkj/iprm; /* 11 */
1676 svmul(a, m, f_i); /* 3 */
1677 b = ddphi*nrkj/iprn; /* 11 */
1678 svmul(b, n, f_l); /* 3 */
1679 p = iprod(r_ij, r_kj); /* 5 */
1680 p *= nrkj_2; /* 1 */
1681 q = iprod(r_kl, r_kj); /* 5 */
1682 q *= nrkj_2; /* 1 */
1683 svmul(p, f_i, uvec); /* 3 */
1684 svmul(q, f_l, vvec); /* 3 */
1685 rvec_sub(uvec, vvec, svec); /* 3 */
1686 rvec_sub(f_i, svec, f_j); /* 3 */
1687 rvec_add(f_l, svec, f_k); /* 3 */
1688 rvec_inc(f[i], f_i); /* 3 */
1689 rvec_dec(f[j], f_j); /* 3 */
1690 rvec_dec(f[k], f_k); /* 3 */
1691 rvec_inc(f[l], f_l); /* 3 */
1695 copy_ivec(SHIFT_IVEC(g, j), jt);
1696 ivec_sub(SHIFT_IVEC(g, i), jt, dt_ij);
1697 ivec_sub(SHIFT_IVEC(g, k), jt, dt_kj);
1698 ivec_sub(SHIFT_IVEC(g, l), jt, dt_lj);
1699 t1 = IVEC2IS(dt_ij);
1700 t2 = IVEC2IS(dt_kj);
1701 t3 = IVEC2IS(dt_lj);
1705 t3 = pbc_rvec_sub(pbc, x[l], x[j], dx_jl);
1712 rvec_inc(fshift[t1], f_i);
1713 rvec_dec(fshift[CENTRAL], f_j);
1714 rvec_dec(fshift[t2], f_k);
1715 rvec_inc(fshift[t3], f_l);
1720 /* As do_dih_fup above, but without shift forces */
1722 do_dih_fup_noshiftf(int i, int j, int k, int l, real ddphi,
1723 rvec r_ij, rvec r_kj, rvec r_kl,
1724 rvec m, rvec n, rvec f[])
1726 rvec f_i, f_j, f_k, f_l;
1727 rvec uvec, vvec, svec;
1728 real iprm, iprn, nrkj, nrkj2, nrkj_1, nrkj_2;
1729 real a, b, p, q, toler;
1731 iprm = iprod(m, m); /* 5 */
1732 iprn = iprod(n, n); /* 5 */
1733 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1734 toler = nrkj2*GMX_REAL_EPS;
1735 if ((iprm > toler) && (iprn > toler))
1737 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1738 nrkj_2 = nrkj_1*nrkj_1; /* 1 */
1739 nrkj = nrkj2*nrkj_1; /* 1 */
1740 a = -ddphi*nrkj/iprm; /* 11 */
1741 svmul(a, m, f_i); /* 3 */
1742 b = ddphi*nrkj/iprn; /* 11 */
1743 svmul(b, n, f_l); /* 3 */
1744 p = iprod(r_ij, r_kj); /* 5 */
1745 p *= nrkj_2; /* 1 */
1746 q = iprod(r_kl, r_kj); /* 5 */
1747 q *= nrkj_2; /* 1 */
1748 svmul(p, f_i, uvec); /* 3 */
1749 svmul(q, f_l, vvec); /* 3 */
1750 rvec_sub(uvec, vvec, svec); /* 3 */
1751 rvec_sub(f_i, svec, f_j); /* 3 */
1752 rvec_add(f_l, svec, f_k); /* 3 */
1753 rvec_inc(f[i], f_i); /* 3 */
1754 rvec_dec(f[j], f_j); /* 3 */
1755 rvec_dec(f[k], f_k); /* 3 */
1756 rvec_inc(f[l], f_l); /* 3 */
1760 /* As do_dih_fup_noshiftf above, but with pre-calculated pre-factors */
1761 static gmx_inline void
1762 do_dih_fup_noshiftf_precalc(int i, int j, int k, int l,
1764 real f_i_x, real f_i_y, real f_i_z,
1765 real mf_l_x, real mf_l_y, real mf_l_z,
1768 rvec f_i, f_j, f_k, f_l;
1769 rvec uvec, vvec, svec;
1777 svmul(p, f_i, uvec);
1778 svmul(q, f_l, vvec);
1779 rvec_sub(uvec, vvec, svec);
1780 rvec_sub(f_i, svec, f_j);
1781 rvec_add(f_l, svec, f_k);
1782 rvec_inc(f[i], f_i);
1783 rvec_dec(f[j], f_j);
1784 rvec_dec(f[k], f_k);
1785 rvec_inc(f[l], f_l);
1789 real dopdihs(real cpA, real cpB, real phiA, real phiB, int mult,
1790 real phi, real lambda, real *V, real *F)
1792 real v, dvdlambda, mdphi, v1, sdphi, ddphi;
1793 real L1 = 1.0 - lambda;
1794 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1795 real dph0 = (phiB - phiA)*DEG2RAD;
1796 real cp = L1*cpA + lambda*cpB;
1798 mdphi = mult*phi - ph0;
1800 ddphi = -cp*mult*sdphi;
1801 v1 = 1.0 + cos(mdphi);
1804 dvdlambda = (cpB - cpA)*v1 + cp*dph0*sdphi;
1811 /* That was 40 flops */
1815 dopdihs_noener(real cpA, real cpB, real phiA, real phiB, int mult,
1816 real phi, real lambda, real *F)
1818 real mdphi, sdphi, ddphi;
1819 real L1 = 1.0 - lambda;
1820 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1821 real cp = L1*cpA + lambda*cpB;
1823 mdphi = mult*phi - ph0;
1825 ddphi = -cp*mult*sdphi;
1829 /* That was 20 flops */
1833 dopdihs_mdphi(real cpA, real cpB, real phiA, real phiB, int mult,
1834 real phi, real lambda, real *cp, real *mdphi)
1836 real L1 = 1.0 - lambda;
1837 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1839 *cp = L1*cpA + lambda*cpB;
1841 *mdphi = mult*phi - ph0;
1844 static real dopdihs_min(real cpA, real cpB, real phiA, real phiB, int mult,
1845 real phi, real lambda, real *V, real *F)
1846 /* similar to dopdihs, except for a minus sign *
1847 * and a different treatment of mult/phi0 */
1849 real v, dvdlambda, mdphi, v1, sdphi, ddphi;
1850 real L1 = 1.0 - lambda;
1851 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1852 real dph0 = (phiB - phiA)*DEG2RAD;
1853 real cp = L1*cpA + lambda*cpB;
1855 mdphi = mult*(phi-ph0);
1857 ddphi = cp*mult*sdphi;
1858 v1 = 1.0-cos(mdphi);
1861 dvdlambda = (cpB-cpA)*v1 + cp*dph0*sdphi;
1868 /* That was 40 flops */
1871 real pdihs(int nbonds,
1872 const t_iatom forceatoms[], const t_iparams forceparams[],
1873 const rvec x[], rvec f[], rvec fshift[],
1874 const t_pbc *pbc, const t_graph *g,
1875 real lambda, real *dvdlambda,
1876 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1877 int gmx_unused *global_atom_index)
1879 int i, type, ai, aj, ak, al;
1881 rvec r_ij, r_kj, r_kl, m, n;
1882 real phi, sign, ddphi, vpd, vtot;
1886 for (i = 0; (i < nbonds); )
1888 type = forceatoms[i++];
1889 ai = forceatoms[i++];
1890 aj = forceatoms[i++];
1891 ak = forceatoms[i++];
1892 al = forceatoms[i++];
1894 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
1895 &sign, &t1, &t2, &t3); /* 84 */
1896 *dvdlambda += dopdihs(forceparams[type].pdihs.cpA,
1897 forceparams[type].pdihs.cpB,
1898 forceparams[type].pdihs.phiA,
1899 forceparams[type].pdihs.phiB,
1900 forceparams[type].pdihs.mult,
1901 phi, lambda, &vpd, &ddphi);
1904 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
1905 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
1908 fprintf(debug, "pdih: (%d,%d,%d,%d) phi=%g\n",
1909 ai, aj, ak, al, phi);
1916 void make_dp_periodic(real *dp) /* 1 flop? */
1918 /* dp cannot be outside (-pi,pi) */
1923 else if (*dp < -M_PI)
1930 /* As pdihs above, but without calculating energies and shift forces */
1932 pdihs_noener(int nbonds,
1933 const t_iatom forceatoms[], const t_iparams forceparams[],
1934 const rvec x[], rvec f[],
1935 const t_pbc gmx_unused *pbc, const t_graph gmx_unused *g,
1937 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1938 int gmx_unused *global_atom_index)
1940 int i, type, ai, aj, ak, al;
1942 rvec r_ij, r_kj, r_kl, m, n;
1943 real phi, sign, ddphi_tot, ddphi;
1945 for (i = 0; (i < nbonds); )
1947 ai = forceatoms[i+1];
1948 aj = forceatoms[i+2];
1949 ak = forceatoms[i+3];
1950 al = forceatoms[i+4];
1952 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
1953 &sign, &t1, &t2, &t3);
1957 /* Loop over dihedrals working on the same atoms,
1958 * so we avoid recalculating angles and force distributions.
1962 type = forceatoms[i];
1963 dopdihs_noener(forceparams[type].pdihs.cpA,
1964 forceparams[type].pdihs.cpB,
1965 forceparams[type].pdihs.phiA,
1966 forceparams[type].pdihs.phiB,
1967 forceparams[type].pdihs.mult,
1968 phi, lambda, &ddphi);
1973 while (i < nbonds &&
1974 forceatoms[i+1] == ai &&
1975 forceatoms[i+2] == aj &&
1976 forceatoms[i+3] == ak &&
1977 forceatoms[i+4] == al);
1979 do_dih_fup_noshiftf(ai, aj, ak, al, ddphi_tot, r_ij, r_kj, r_kl, m, n, f);
1984 #ifdef GMX_SIMD_HAVE_REAL
1986 /* As pdihs_noner above, but using SIMD to calculate many dihedrals at once */
1988 pdihs_noener_simd(int nbonds,
1989 const t_iatom forceatoms[], const t_iparams forceparams[],
1990 const rvec x[], rvec f[],
1991 const t_pbc *pbc, const t_graph gmx_unused *g,
1992 real gmx_unused lambda,
1993 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1994 int gmx_unused *global_atom_index)
1998 int type, ai[GMX_SIMD_REAL_WIDTH], aj[GMX_SIMD_REAL_WIDTH], ak[GMX_SIMD_REAL_WIDTH], al[GMX_SIMD_REAL_WIDTH];
1999 real dr_array[3*DIM*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *dr;
2000 real buf_array[7*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *buf;
2001 real *cp, *phi0, *mult, *p, *q;
2002 gmx_simd_real_t phi0_S, phi_S;
2003 gmx_simd_real_t mx_S, my_S, mz_S;
2004 gmx_simd_real_t nx_S, ny_S, nz_S;
2005 gmx_simd_real_t nrkj_m2_S, nrkj_n2_S;
2006 gmx_simd_real_t cp_S, mdphi_S, mult_S;
2007 gmx_simd_real_t sin_S, cos_S;
2008 gmx_simd_real_t mddphi_S;
2009 gmx_simd_real_t sf_i_S, msf_l_S;
2010 pbc_simd_t pbc_simd;
2012 /* Ensure SIMD register alignment */
2013 dr = gmx_simd_align_r(dr_array);
2014 buf = gmx_simd_align_r(buf_array);
2016 /* Extract aligned pointer for parameters and variables */
2017 cp = buf + 0*GMX_SIMD_REAL_WIDTH;
2018 phi0 = buf + 1*GMX_SIMD_REAL_WIDTH;
2019 mult = buf + 2*GMX_SIMD_REAL_WIDTH;
2020 p = buf + 3*GMX_SIMD_REAL_WIDTH;
2021 q = buf + 4*GMX_SIMD_REAL_WIDTH;
2023 set_pbc_simd(pbc, &pbc_simd);
2025 /* nbonds is the number of dihedrals times nfa1, here we step GMX_SIMD_REAL_WIDTH dihs */
2026 for (i = 0; (i < nbonds); i += GMX_SIMD_REAL_WIDTH*nfa1)
2028 /* Collect atoms quadruplets for GMX_SIMD_REAL_WIDTH dihedrals.
2029 * iu indexes into forceatoms, we should not let iu go beyond nbonds.
2032 for (s = 0; s < GMX_SIMD_REAL_WIDTH; s++)
2034 type = forceatoms[iu];
2035 ai[s] = forceatoms[iu+1];
2036 aj[s] = forceatoms[iu+2];
2037 ak[s] = forceatoms[iu+3];
2038 al[s] = forceatoms[iu+4];
2040 cp[s] = forceparams[type].pdihs.cpA;
2041 phi0[s] = forceparams[type].pdihs.phiA*DEG2RAD;
2042 mult[s] = forceparams[type].pdihs.mult;
2044 /* At the end fill the arrays with identical entries */
2045 if (iu + nfa1 < nbonds)
2051 /* Caclulate GMX_SIMD_REAL_WIDTH dihedral angles at once */
2052 dih_angle_simd(x, ai, aj, ak, al, &pbc_simd,
2055 &mx_S, &my_S, &mz_S,
2056 &nx_S, &ny_S, &nz_S,
2061 cp_S = gmx_simd_load_r(cp);
2062 phi0_S = gmx_simd_load_r(phi0);
2063 mult_S = gmx_simd_load_r(mult);
2065 mdphi_S = gmx_simd_sub_r(gmx_simd_mul_r(mult_S, phi_S), phi0_S);
2067 /* Calculate GMX_SIMD_REAL_WIDTH sines at once */
2068 gmx_simd_sincos_r(mdphi_S, &sin_S, &cos_S);
2069 mddphi_S = gmx_simd_mul_r(gmx_simd_mul_r(cp_S, mult_S), sin_S);
2070 sf_i_S = gmx_simd_mul_r(mddphi_S, nrkj_m2_S);
2071 msf_l_S = gmx_simd_mul_r(mddphi_S, nrkj_n2_S);
2073 /* After this m?_S will contain f[i] */
2074 mx_S = gmx_simd_mul_r(sf_i_S, mx_S);
2075 my_S = gmx_simd_mul_r(sf_i_S, my_S);
2076 mz_S = gmx_simd_mul_r(sf_i_S, mz_S);
2078 /* After this m?_S will contain -f[l] */
2079 nx_S = gmx_simd_mul_r(msf_l_S, nx_S);
2080 ny_S = gmx_simd_mul_r(msf_l_S, ny_S);
2081 nz_S = gmx_simd_mul_r(msf_l_S, nz_S);
2083 gmx_simd_store_r(dr + 0*GMX_SIMD_REAL_WIDTH, mx_S);
2084 gmx_simd_store_r(dr + 1*GMX_SIMD_REAL_WIDTH, my_S);
2085 gmx_simd_store_r(dr + 2*GMX_SIMD_REAL_WIDTH, mz_S);
2086 gmx_simd_store_r(dr + 3*GMX_SIMD_REAL_WIDTH, nx_S);
2087 gmx_simd_store_r(dr + 4*GMX_SIMD_REAL_WIDTH, ny_S);
2088 gmx_simd_store_r(dr + 5*GMX_SIMD_REAL_WIDTH, nz_S);
2094 do_dih_fup_noshiftf_precalc(ai[s], aj[s], ak[s], al[s],
2096 dr[ XX *GMX_SIMD_REAL_WIDTH+s],
2097 dr[ YY *GMX_SIMD_REAL_WIDTH+s],
2098 dr[ ZZ *GMX_SIMD_REAL_WIDTH+s],
2099 dr[(DIM+XX)*GMX_SIMD_REAL_WIDTH+s],
2100 dr[(DIM+YY)*GMX_SIMD_REAL_WIDTH+s],
2101 dr[(DIM+ZZ)*GMX_SIMD_REAL_WIDTH+s],
2106 while (s < GMX_SIMD_REAL_WIDTH && iu < nbonds);
2110 #endif /* GMX_SIMD_HAVE_REAL */
2113 real idihs(int nbonds,
2114 const t_iatom forceatoms[], const t_iparams forceparams[],
2115 const rvec x[], rvec f[], rvec fshift[],
2116 const t_pbc *pbc, const t_graph *g,
2117 real lambda, real *dvdlambda,
2118 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2119 int gmx_unused *global_atom_index)
2121 int i, type, ai, aj, ak, al;
2123 real phi, phi0, dphi0, ddphi, sign, vtot;
2124 rvec r_ij, r_kj, r_kl, m, n;
2125 real L1, kk, dp, dp2, kA, kB, pA, pB, dvdl_term;
2130 for (i = 0; (i < nbonds); )
2132 type = forceatoms[i++];
2133 ai = forceatoms[i++];
2134 aj = forceatoms[i++];
2135 ak = forceatoms[i++];
2136 al = forceatoms[i++];
2138 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2139 &sign, &t1, &t2, &t3); /* 84 */
2141 /* phi can jump if phi0 is close to Pi/-Pi, which will cause huge
2142 * force changes if we just apply a normal harmonic.
2143 * Instead, we first calculate phi-phi0 and take it modulo (-Pi,Pi).
2144 * This means we will never have the periodicity problem, unless
2145 * the dihedral is Pi away from phiO, which is very unlikely due to
2148 kA = forceparams[type].harmonic.krA;
2149 kB = forceparams[type].harmonic.krB;
2150 pA = forceparams[type].harmonic.rA;
2151 pB = forceparams[type].harmonic.rB;
2153 kk = L1*kA + lambda*kB;
2154 phi0 = (L1*pA + lambda*pB)*DEG2RAD;
2155 dphi0 = (pB - pA)*DEG2RAD;
2159 make_dp_periodic(&dp);
2166 dvdl_term += 0.5*(kB - kA)*dp2 - kk*dphi0*dp;
2168 do_dih_fup(ai, aj, ak, al, -ddphi, r_ij, r_kj, r_kl, m, n,
2169 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2174 fprintf(debug, "idih: (%d,%d,%d,%d) phi=%g\n",
2175 ai, aj, ak, al, phi);
2180 *dvdlambda += dvdl_term;
2185 /*! \brief returns dx, rdist, and dpdl for functions posres() and fbposres()
2187 static void posres_dx(const rvec x, const rvec pos0A, const rvec pos0B,
2188 const rvec comA_sc, const rvec comB_sc,
2190 t_pbc *pbc, int refcoord_scaling, int npbcdim,
2191 rvec dx, rvec rdist, rvec dpdl)
2194 real posA, posB, L1, ref = 0.;
2199 for (m = 0; m < DIM; m++)
2205 switch (refcoord_scaling)
2209 rdist[m] = L1*posA + lambda*posB;
2210 dpdl[m] = posB - posA;
2213 /* Box relative coordinates are stored for dimensions with pbc */
2214 posA *= pbc->box[m][m];
2215 posB *= pbc->box[m][m];
2216 assert(npbcdim <= DIM);
2217 for (d = m+1; d < npbcdim; d++)
2219 posA += pos0A[d]*pbc->box[d][m];
2220 posB += pos0B[d]*pbc->box[d][m];
2222 ref = L1*posA + lambda*posB;
2224 dpdl[m] = posB - posA;
2227 ref = L1*comA_sc[m] + lambda*comB_sc[m];
2228 rdist[m] = L1*posA + lambda*posB;
2229 dpdl[m] = comB_sc[m] - comA_sc[m] + posB - posA;
2232 gmx_fatal(FARGS, "No such scaling method implemented");
2237 ref = L1*posA + lambda*posB;
2239 dpdl[m] = posB - posA;
2242 /* We do pbc_dx with ref+rdist,
2243 * since with only ref we can be up to half a box vector wrong.
2245 pos[m] = ref + rdist[m];
2250 pbc_dx(pbc, x, pos, dx);
2254 rvec_sub(x, pos, dx);
2258 /*! \brief Computes forces and potential for flat-bottom cylindrical restraints.
2259 * Returns the flat-bottom potential. */
2260 static real do_fbposres_cylinder(int fbdim, rvec fm, rvec dx, real rfb, real kk, gmx_bool bInvert)
2263 real dr, dr2, invdr, v, rfb2;
2269 for (d = 0; d < DIM; d++)
2278 ( (dr2 > rfb2 && bInvert == FALSE ) || (dr2 < rfb2 && bInvert == TRUE ) )
2283 v = 0.5*kk*sqr(dr - rfb);
2284 for (d = 0; d < DIM; d++)
2288 fm[d] = -kk*(dr-rfb)*dx[d]*invdr; /* Force pointing to the center */
2296 /*! \brief Adds forces of flat-bottomed positions restraints to f[]
2297 * and fixes vir_diag. Returns the flat-bottomed potential. */
2298 real fbposres(int nbonds,
2299 const t_iatom forceatoms[], const t_iparams forceparams[],
2300 const rvec x[], rvec f[], rvec vir_diag,
2302 int refcoord_scaling, int ePBC, rvec com)
2303 /* compute flat-bottomed positions restraints */
2305 int i, ai, m, d, type, npbcdim = 0, fbdim;
2306 const t_iparams *pr;
2308 real dr, dr2, rfb, rfb2, fact;
2309 rvec com_sc, rdist, dx, dpdl, fm;
2312 npbcdim = ePBC2npbcdim(ePBC);
2314 if (refcoord_scaling == erscCOM)
2317 for (m = 0; m < npbcdim; m++)
2319 assert(npbcdim <= DIM);
2320 for (d = m; d < npbcdim; d++)
2322 com_sc[m] += com[d]*pbc->box[d][m];
2328 for (i = 0; (i < nbonds); )
2330 type = forceatoms[i++];
2331 ai = forceatoms[i++];
2332 pr = &forceparams[type];
2334 /* same calculation as for normal posres, but with identical A and B states, and lambda==0 */
2335 posres_dx(x[ai], forceparams[type].fbposres.pos0, forceparams[type].fbposres.pos0,
2336 com_sc, com_sc, 0.0,
2337 pbc, refcoord_scaling, npbcdim,
2343 kk = pr->fbposres.k;
2344 rfb = pr->fbposres.r;
2347 /* with rfb<0, push particle out of the sphere/cylinder/layer */
2355 switch (pr->fbposres.geom)
2357 case efbposresSPHERE:
2358 /* spherical flat-bottom posres */
2361 ( (dr2 > rfb2 && bInvert == FALSE ) || (dr2 < rfb2 && bInvert == TRUE ) )
2365 v = 0.5*kk*sqr(dr - rfb);
2366 fact = -kk*(dr-rfb)/dr; /* Force pointing to the center pos0 */
2367 svmul(fact, dx, fm);
2370 case efbposresCYLINDERX:
2371 /* cylindrical flat-bottom posres in y-z plane. fm[XX] = 0. */
2373 v = do_fbposres_cylinder(fbdim, fm, dx, rfb, kk, bInvert);
2375 case efbposresCYLINDERY:
2376 /* cylindrical flat-bottom posres in x-z plane. fm[YY] = 0. */
2378 v = do_fbposres_cylinder(fbdim, fm, dx, rfb, kk, bInvert);
2380 case efbposresCYLINDER:
2381 /* equivalent to efbposresCYLINDERZ for backwards compatibility */
2382 case efbposresCYLINDERZ:
2383 /* cylindrical flat-bottom posres in x-y plane. fm[ZZ] = 0. */
2385 v = do_fbposres_cylinder(fbdim, fm, dx, rfb, kk, bInvert);
2387 case efbposresX: /* fbdim=XX */
2388 case efbposresY: /* fbdim=YY */
2389 case efbposresZ: /* fbdim=ZZ */
2390 /* 1D flat-bottom potential */
2391 fbdim = pr->fbposres.geom - efbposresX;
2393 if ( ( dr > rfb && bInvert == FALSE ) || ( 0 < dr && dr < rfb && bInvert == TRUE ) )
2395 v = 0.5*kk*sqr(dr - rfb);
2396 fm[fbdim] = -kk*(dr - rfb);
2398 else if ( (dr < (-rfb) && bInvert == FALSE ) || ( (-rfb) < dr && dr < 0 && bInvert == TRUE ))
2400 v = 0.5*kk*sqr(dr + rfb);
2401 fm[fbdim] = -kk*(dr + rfb);
2408 for (m = 0; (m < DIM); m++)
2411 /* Here we correct for the pbc_dx which included rdist */
2412 vir_diag[m] -= 0.5*(dx[m] + rdist[m])*fm[m];
2420 real posres(int nbonds,
2421 const t_iatom forceatoms[], const t_iparams forceparams[],
2422 const rvec x[], rvec f[], rvec vir_diag,
2424 real lambda, real *dvdlambda,
2425 int refcoord_scaling, int ePBC, rvec comA, rvec comB)
2427 int i, ai, m, d, type, npbcdim = 0;
2428 const t_iparams *pr;
2431 rvec comA_sc, comB_sc, rdist, dpdl, dx;
2432 gmx_bool bForceValid = TRUE;
2434 if ((f == NULL) || (vir_diag == NULL)) /* should both be null together! */
2436 bForceValid = FALSE;
2439 npbcdim = ePBC2npbcdim(ePBC);
2441 if (refcoord_scaling == erscCOM)
2443 clear_rvec(comA_sc);
2444 clear_rvec(comB_sc);
2445 for (m = 0; m < npbcdim; m++)
2447 assert(npbcdim <= DIM);
2448 for (d = m; d < npbcdim; d++)
2450 comA_sc[m] += comA[d]*pbc->box[d][m];
2451 comB_sc[m] += comB[d]*pbc->box[d][m];
2459 for (i = 0; (i < nbonds); )
2461 type = forceatoms[i++];
2462 ai = forceatoms[i++];
2463 pr = &forceparams[type];
2465 /* return dx, rdist, and dpdl */
2466 posres_dx(x[ai], forceparams[type].posres.pos0A, forceparams[type].posres.pos0B,
2467 comA_sc, comB_sc, lambda,
2468 pbc, refcoord_scaling, npbcdim,
2471 for (m = 0; (m < DIM); m++)
2473 kk = L1*pr->posres.fcA[m] + lambda*pr->posres.fcB[m];
2475 vtot += 0.5*kk*dx[m]*dx[m];
2477 0.5*(pr->posres.fcB[m] - pr->posres.fcA[m])*dx[m]*dx[m]
2480 /* Here we correct for the pbc_dx which included rdist */
2484 vir_diag[m] -= 0.5*(dx[m] + rdist[m])*fm;
2492 static real low_angres(int nbonds,
2493 const t_iatom forceatoms[], const t_iparams forceparams[],
2494 const rvec x[], rvec f[], rvec fshift[],
2495 const t_pbc *pbc, const t_graph *g,
2496 real lambda, real *dvdlambda,
2499 int i, m, type, ai, aj, ak, al;
2501 real phi, cos_phi, cos_phi2, vid, vtot, dVdphi;
2502 rvec r_ij, r_kl, f_i, f_k = {0, 0, 0};
2503 real st, sth, nrij2, nrkl2, c, cij, ckl;
2506 t2 = 0; /* avoid warning with gcc-3.3. It is never used uninitialized */
2509 ak = al = 0; /* to avoid warnings */
2510 for (i = 0; i < nbonds; )
2512 type = forceatoms[i++];
2513 ai = forceatoms[i++];
2514 aj = forceatoms[i++];
2515 t1 = pbc_rvec_sub(pbc, x[aj], x[ai], r_ij); /* 3 */
2518 ak = forceatoms[i++];
2519 al = forceatoms[i++];
2520 t2 = pbc_rvec_sub(pbc, x[al], x[ak], r_kl); /* 3 */
2529 cos_phi = cos_angle(r_ij, r_kl); /* 25 */
2530 phi = acos(cos_phi); /* 10 */
2532 *dvdlambda += dopdihs_min(forceparams[type].pdihs.cpA,
2533 forceparams[type].pdihs.cpB,
2534 forceparams[type].pdihs.phiA,
2535 forceparams[type].pdihs.phiB,
2536 forceparams[type].pdihs.mult,
2537 phi, lambda, &vid, &dVdphi); /* 40 */
2541 cos_phi2 = sqr(cos_phi); /* 1 */
2544 st = -dVdphi*gmx_invsqrt(1 - cos_phi2); /* 12 */
2545 sth = st*cos_phi; /* 1 */
2546 nrij2 = iprod(r_ij, r_ij); /* 5 */
2547 nrkl2 = iprod(r_kl, r_kl); /* 5 */
2549 c = st*gmx_invsqrt(nrij2*nrkl2); /* 11 */
2550 cij = sth/nrij2; /* 10 */
2551 ckl = sth/nrkl2; /* 10 */
2553 for (m = 0; m < DIM; m++) /* 18+18 */
2555 f_i[m] = (c*r_kl[m]-cij*r_ij[m]);
2560 f_k[m] = (c*r_ij[m]-ckl*r_kl[m]);
2568 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
2571 rvec_inc(fshift[t1], f_i);
2572 rvec_dec(fshift[CENTRAL], f_i);
2577 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, al), dt);
2580 rvec_inc(fshift[t2], f_k);
2581 rvec_dec(fshift[CENTRAL], f_k);
2586 return vtot; /* 184 / 157 (bZAxis) total */
2589 real angres(int nbonds,
2590 const t_iatom forceatoms[], const t_iparams forceparams[],
2591 const rvec x[], rvec f[], rvec fshift[],
2592 const t_pbc *pbc, const t_graph *g,
2593 real lambda, real *dvdlambda,
2594 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2595 int gmx_unused *global_atom_index)
2597 return low_angres(nbonds, forceatoms, forceparams, x, f, fshift, pbc, g,
2598 lambda, dvdlambda, FALSE);
2601 real angresz(int nbonds,
2602 const t_iatom forceatoms[], const t_iparams forceparams[],
2603 const rvec x[], rvec f[], rvec fshift[],
2604 const t_pbc *pbc, const t_graph *g,
2605 real lambda, real *dvdlambda,
2606 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2607 int gmx_unused *global_atom_index)
2609 return low_angres(nbonds, forceatoms, forceparams, x, f, fshift, pbc, g,
2610 lambda, dvdlambda, TRUE);
2613 real dihres(int nbonds,
2614 const t_iatom forceatoms[], const t_iparams forceparams[],
2615 const rvec x[], rvec f[], rvec fshift[],
2616 const t_pbc *pbc, const t_graph *g,
2617 real lambda, real *dvdlambda,
2618 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2619 int gmx_unused *global_atom_index)
2622 int ai, aj, ak, al, i, k, type, t1, t2, t3;
2623 real phi0A, phi0B, dphiA, dphiB, kfacA, kfacB, phi0, dphi, kfac;
2624 real phi, ddphi, ddp, ddp2, dp, sign, d2r, L1;
2625 rvec r_ij, r_kj, r_kl, m, n;
2632 for (i = 0; (i < nbonds); )
2634 type = forceatoms[i++];
2635 ai = forceatoms[i++];
2636 aj = forceatoms[i++];
2637 ak = forceatoms[i++];
2638 al = forceatoms[i++];
2640 phi0A = forceparams[type].dihres.phiA*d2r;
2641 dphiA = forceparams[type].dihres.dphiA*d2r;
2642 kfacA = forceparams[type].dihres.kfacA;
2644 phi0B = forceparams[type].dihres.phiB*d2r;
2645 dphiB = forceparams[type].dihres.dphiB*d2r;
2646 kfacB = forceparams[type].dihres.kfacB;
2648 phi0 = L1*phi0A + lambda*phi0B;
2649 dphi = L1*dphiA + lambda*dphiB;
2650 kfac = L1*kfacA + lambda*kfacB;
2652 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2653 &sign, &t1, &t2, &t3);
2658 fprintf(debug, "dihres[%d]: %d %d %d %d : phi=%f, dphi=%f, kfac=%f\n",
2659 k++, ai, aj, ak, al, phi0, dphi, kfac);
2661 /* phi can jump if phi0 is close to Pi/-Pi, which will cause huge
2662 * force changes if we just apply a normal harmonic.
2663 * Instead, we first calculate phi-phi0 and take it modulo (-Pi,Pi).
2664 * This means we will never have the periodicity problem, unless
2665 * the dihedral is Pi away from phiO, which is very unlikely due to
2669 make_dp_periodic(&dp);
2675 else if (dp < -dphi)
2687 vtot += 0.5*kfac*ddp2;
2690 *dvdlambda += 0.5*(kfacB - kfacA)*ddp2;
2691 /* lambda dependence from changing restraint distances */
2694 *dvdlambda -= kfac*ddp*((dphiB - dphiA)+(phi0B - phi0A));
2698 *dvdlambda += kfac*ddp*((dphiB - dphiA)-(phi0B - phi0A));
2700 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
2701 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2708 real unimplemented(int gmx_unused nbonds,
2709 const t_iatom gmx_unused forceatoms[], const t_iparams gmx_unused forceparams[],
2710 const rvec gmx_unused x[], rvec gmx_unused f[], rvec gmx_unused fshift[],
2711 const t_pbc gmx_unused *pbc, const t_graph gmx_unused *g,
2712 real gmx_unused lambda, real gmx_unused *dvdlambda,
2713 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2714 int gmx_unused *global_atom_index)
2716 gmx_impl("*** you are using a not implemented function");
2718 return 0.0; /* To make the compiler happy */
2721 real restrangles(int nbonds,
2722 const t_iatom forceatoms[], const t_iparams forceparams[],
2723 const rvec x[], rvec f[], rvec fshift[],
2724 const t_pbc *pbc, const t_graph *g,
2725 real gmx_unused lambda, real gmx_unused *dvdlambda,
2726 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2727 int gmx_unused *global_atom_index)
2729 int i, d, ai, aj, ak, type, m;
2732 ivec jt, dt_ij, dt_kj;
2734 real prefactor, ratio_ante, ratio_post;
2735 rvec delta_ante, delta_post, vec_temp;
2738 for (i = 0; (i < nbonds); )
2740 type = forceatoms[i++];
2741 ai = forceatoms[i++];
2742 aj = forceatoms[i++];
2743 ak = forceatoms[i++];
2745 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], vec_temp);
2746 pbc_rvec_sub(pbc, x[aj], x[ai], delta_ante);
2747 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], delta_post);
2750 /* This function computes factors needed for restricted angle potential.
2751 * The restricted angle potential is used in coarse-grained simulations to avoid singularities
2752 * when three particles align and the dihedral angle and dihedral potential
2753 * cannot be calculated. This potential is calculated using the formula:
2754 real restrangles(int nbonds,
2755 const t_iatom forceatoms[],const t_iparams forceparams[],
2756 const rvec x[],rvec f[],rvec fshift[],
2757 const t_pbc *pbc,const t_graph *g,
2758 real gmx_unused lambda,real gmx_unused *dvdlambda,
2759 const t_mdatoms gmx_unused *md,t_fcdata gmx_unused *fcd,
2760 int gmx_unused *global_atom_index)
2762 int i, d, ai, aj, ak, type, m;
2766 ivec jt,dt_ij,dt_kj;
2768 real prefactor, ratio_ante, ratio_post;
2769 rvec delta_ante, delta_post, vec_temp;
2772 for(i=0; (i<nbonds); )
2774 type = forceatoms[i++];
2775 ai = forceatoms[i++];
2776 aj = forceatoms[i++];
2777 ak = forceatoms[i++];
2779 * \f[V_{\rm ReB}(\theta_i) = \frac{1}{2} k_{\theta} \frac{(\cos\theta_i - \cos\theta_0)^2}
2780 * {\sin^2\theta_i}\f] ({eq:ReB} and ref \cite{MonicaGoga2013} from the manual).
2781 * For more explanations see comments file "restcbt.h". */
2783 compute_factors_restangles(type, forceparams, delta_ante, delta_post,
2784 &prefactor, &ratio_ante, &ratio_post, &v);
2786 /* Forces are computed per component */
2787 for (d = 0; d < DIM; d++)
2789 f_i[d] = prefactor * (ratio_ante * delta_ante[d] - delta_post[d]);
2790 f_j[d] = prefactor * ((ratio_post + 1.0) * delta_post[d] - (ratio_ante + 1.0) * delta_ante[d]);
2791 f_k[d] = prefactor * (delta_ante[d] - ratio_post * delta_post[d]);
2794 /* Computation of potential energy */
2800 for (m = 0; (m < DIM); m++)
2809 copy_ivec(SHIFT_IVEC(g, aj), jt);
2810 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
2811 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
2812 t1 = IVEC2IS(dt_ij);
2813 t2 = IVEC2IS(dt_kj);
2816 rvec_inc(fshift[t1], f_i);
2817 rvec_inc(fshift[CENTRAL], f_j);
2818 rvec_inc(fshift[t2], f_k);
2824 real restrdihs(int nbonds,
2825 const t_iatom forceatoms[], const t_iparams forceparams[],
2826 const rvec x[], rvec f[], rvec fshift[],
2827 const t_pbc *pbc, const t_graph *g,
2828 real gmx_unused lambda, real gmx_unused *dvlambda,
2829 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2830 int gmx_unused *global_atom_index)
2832 int i, d, type, ai, aj, ak, al;
2833 rvec f_i, f_j, f_k, f_l;
2835 ivec jt, dt_ij, dt_kj, dt_lj;
2838 rvec delta_ante, delta_crnt, delta_post, vec_temp;
2839 real factor_phi_ai_ante, factor_phi_ai_crnt, factor_phi_ai_post;
2840 real factor_phi_aj_ante, factor_phi_aj_crnt, factor_phi_aj_post;
2841 real factor_phi_ak_ante, factor_phi_ak_crnt, factor_phi_ak_post;
2842 real factor_phi_al_ante, factor_phi_al_crnt, factor_phi_al_post;
2847 for (i = 0; (i < nbonds); )
2849 type = forceatoms[i++];
2850 ai = forceatoms[i++];
2851 aj = forceatoms[i++];
2852 ak = forceatoms[i++];
2853 al = forceatoms[i++];
2855 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], vec_temp);
2856 pbc_rvec_sub(pbc, x[aj], x[ai], delta_ante);
2857 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], delta_crnt);
2858 pbc_rvec_sub(pbc, x[ak], x[al], vec_temp);
2859 pbc_rvec_sub(pbc, x[al], x[ak], delta_post);
2861 /* This function computes factors needed for restricted angle potential.
2862 * The restricted angle potential is used in coarse-grained simulations to avoid singularities
2863 * when three particles align and the dihedral angle and dihedral potential cannot be calculated.
2864 * This potential is calculated using the formula:
2865 * \f[V_{\rm ReB}(\theta_i) = \frac{1}{2} k_{\theta}
2866 * \frac{(\cos\theta_i - \cos\theta_0)^2}{\sin^2\theta_i}\f]
2867 * ({eq:ReB} and ref \cite{MonicaGoga2013} from the manual).
2868 * For more explanations see comments file "restcbt.h" */
2870 compute_factors_restrdihs(type, forceparams,
2871 delta_ante, delta_crnt, delta_post,
2872 &factor_phi_ai_ante, &factor_phi_ai_crnt, &factor_phi_ai_post,
2873 &factor_phi_aj_ante, &factor_phi_aj_crnt, &factor_phi_aj_post,
2874 &factor_phi_ak_ante, &factor_phi_ak_crnt, &factor_phi_ak_post,
2875 &factor_phi_al_ante, &factor_phi_al_crnt, &factor_phi_al_post,
2876 &prefactor_phi, &v);
2879 /* Computation of forces per component */
2880 for (d = 0; d < DIM; d++)
2882 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]);
2883 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]);
2884 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]);
2885 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]);
2887 /* Computation of the energy */
2893 /* Updating the forces */
2895 rvec_inc(f[ai], f_i);
2896 rvec_inc(f[aj], f_j);
2897 rvec_inc(f[ak], f_k);
2898 rvec_inc(f[al], f_l);
2901 /* Updating the fshift forces for the pressure coupling */
2904 copy_ivec(SHIFT_IVEC(g, aj), jt);
2905 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
2906 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
2907 ivec_sub(SHIFT_IVEC(g, al), jt, dt_lj);
2908 t1 = IVEC2IS(dt_ij);
2909 t2 = IVEC2IS(dt_kj);
2910 t3 = IVEC2IS(dt_lj);
2914 t3 = pbc_rvec_sub(pbc, x[al], x[aj], dx_jl);
2921 rvec_inc(fshift[t1], f_i);
2922 rvec_inc(fshift[CENTRAL], f_j);
2923 rvec_inc(fshift[t2], f_k);
2924 rvec_inc(fshift[t3], f_l);
2932 real cbtdihs(int nbonds,
2933 const t_iatom forceatoms[], const t_iparams forceparams[],
2934 const rvec x[], rvec f[], rvec fshift[],
2935 const t_pbc *pbc, const t_graph *g,
2936 real gmx_unused lambda, real gmx_unused *dvdlambda,
2937 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2938 int gmx_unused *global_atom_index)
2940 int type, ai, aj, ak, al, i, d;
2944 rvec f_i, f_j, f_k, f_l;
2945 ivec jt, dt_ij, dt_kj, dt_lj;
2947 rvec delta_ante, delta_crnt, delta_post;
2948 rvec f_phi_ai, f_phi_aj, f_phi_ak, f_phi_al;
2949 rvec f_theta_ante_ai, f_theta_ante_aj, f_theta_ante_ak;
2950 rvec f_theta_post_aj, f_theta_post_ak, f_theta_post_al;
2956 for (i = 0; (i < nbonds); )
2958 type = forceatoms[i++];
2959 ai = forceatoms[i++];
2960 aj = forceatoms[i++];
2961 ak = forceatoms[i++];
2962 al = forceatoms[i++];
2965 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], vec_temp);
2966 pbc_rvec_sub(pbc, x[aj], x[ai], delta_ante);
2967 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], vec_temp);
2968 pbc_rvec_sub(pbc, x[ak], x[aj], delta_crnt);
2969 pbc_rvec_sub(pbc, x[ak], x[al], vec_temp);
2970 pbc_rvec_sub(pbc, x[al], x[ak], delta_post);
2972 /* \brief Compute factors for CBT potential
2973 * The combined bending-torsion potential goes to zero in a very smooth manner, eliminating the numerical
2974 * instabilities, when three coarse-grained particles align and the dihedral angle and standard
2975 * dihedral potentials cannot be calculated. The CBT potential is calculated using the formula:
2976 * \f[V_{\rm CBT}(\theta_{i-1}, \theta_i, \phi_i) = k_{\phi} \sin^3\theta_{i-1} \sin^3\theta_{i}
2977 * \sum_{n=0}^4 { a_n \cos^n\phi_i}\f] ({eq:CBT} and ref \cite{MonicaGoga2013} from the manual).
2978 * It contains in its expression not only the dihedral angle \f$\phi\f$
2979 * but also \f[\theta_{i-1}\f] (theta_ante bellow) and \f[\theta_{i}\f] (theta_post bellow)
2980 * --- the adjacent bending angles.
2981 * For more explanations see comments file "restcbt.h". */
2983 compute_factors_cbtdihs(type, forceparams, delta_ante, delta_crnt, delta_post,
2984 f_phi_ai, f_phi_aj, f_phi_ak, f_phi_al,
2985 f_theta_ante_ai, f_theta_ante_aj, f_theta_ante_ak,
2986 f_theta_post_aj, f_theta_post_ak, f_theta_post_al,
2990 /* Acumulate the resuts per beads */
2991 for (d = 0; d < DIM; d++)
2993 f_i[d] = f_phi_ai[d] + f_theta_ante_ai[d];
2994 f_j[d] = f_phi_aj[d] + f_theta_ante_aj[d] + f_theta_post_aj[d];
2995 f_k[d] = f_phi_ak[d] + f_theta_ante_ak[d] + f_theta_post_ak[d];
2996 f_l[d] = f_phi_al[d] + f_theta_post_al[d];
2999 /* Compute the potential energy */
3004 /* Updating the forces */
3005 rvec_inc(f[ai], f_i);
3006 rvec_inc(f[aj], f_j);
3007 rvec_inc(f[ak], f_k);
3008 rvec_inc(f[al], f_l);
3011 /* Updating the fshift forces for the pressure coupling */
3014 copy_ivec(SHIFT_IVEC(g, aj), jt);
3015 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3016 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3017 ivec_sub(SHIFT_IVEC(g, al), jt, dt_lj);
3018 t1 = IVEC2IS(dt_ij);
3019 t2 = IVEC2IS(dt_kj);
3020 t3 = IVEC2IS(dt_lj);
3024 t3 = pbc_rvec_sub(pbc, x[al], x[aj], dx_jl);
3031 rvec_inc(fshift[t1], f_i);
3032 rvec_inc(fshift[CENTRAL], f_j);
3033 rvec_inc(fshift[t2], f_k);
3034 rvec_inc(fshift[t3], f_l);
3040 real rbdihs(int nbonds,
3041 const t_iatom forceatoms[], const t_iparams forceparams[],
3042 const rvec x[], rvec f[], rvec fshift[],
3043 const t_pbc *pbc, const t_graph *g,
3044 real lambda, real *dvdlambda,
3045 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3046 int gmx_unused *global_atom_index)
3048 const real c0 = 0.0, c1 = 1.0, c2 = 2.0, c3 = 3.0, c4 = 4.0, c5 = 5.0;
3049 int type, ai, aj, ak, al, i, j;
3051 rvec r_ij, r_kj, r_kl, m, n;
3052 real parmA[NR_RBDIHS];
3053 real parmB[NR_RBDIHS];
3054 real parm[NR_RBDIHS];
3055 real cos_phi, phi, rbp, rbpBA;
3056 real v, sign, ddphi, sin_phi;
3058 real L1 = 1.0-lambda;
3062 for (i = 0; (i < nbonds); )
3064 type = forceatoms[i++];
3065 ai = forceatoms[i++];
3066 aj = forceatoms[i++];
3067 ak = forceatoms[i++];
3068 al = forceatoms[i++];
3070 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
3071 &sign, &t1, &t2, &t3); /* 84 */
3073 /* Change to polymer convention */
3080 phi -= M_PI; /* 1 */
3084 /* Beware of accuracy loss, cannot use 1-sqrt(cos^2) ! */
3087 for (j = 0; (j < NR_RBDIHS); j++)
3089 parmA[j] = forceparams[type].rbdihs.rbcA[j];
3090 parmB[j] = forceparams[type].rbdihs.rbcB[j];
3091 parm[j] = L1*parmA[j]+lambda*parmB[j];
3093 /* Calculate cosine powers */
3094 /* Calculate the energy */
3095 /* Calculate the derivative */
3098 dvdl_term += (parmB[0]-parmA[0]);
3103 rbpBA = parmB[1]-parmA[1];
3104 ddphi += rbp*cosfac;
3107 dvdl_term += cosfac*rbpBA;
3109 rbpBA = parmB[2]-parmA[2];
3110 ddphi += c2*rbp*cosfac;
3113 dvdl_term += cosfac*rbpBA;
3115 rbpBA = parmB[3]-parmA[3];
3116 ddphi += c3*rbp*cosfac;
3119 dvdl_term += cosfac*rbpBA;
3121 rbpBA = parmB[4]-parmA[4];
3122 ddphi += c4*rbp*cosfac;
3125 dvdl_term += cosfac*rbpBA;
3127 rbpBA = parmB[5]-parmA[5];
3128 ddphi += c5*rbp*cosfac;
3131 dvdl_term += cosfac*rbpBA;
3133 ddphi = -ddphi*sin_phi; /* 11 */
3135 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
3136 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
3139 *dvdlambda += dvdl_term;
3144 int cmap_setup_grid_index(int ip, int grid_spacing, int *ipm1, int *ipp1, int *ipp2)
3150 ip = ip + grid_spacing - 1;
3152 else if (ip > grid_spacing)
3154 ip = ip - grid_spacing - 1;
3163 im1 = grid_spacing - 1;
3165 else if (ip == grid_spacing-2)
3169 else if (ip == grid_spacing-1)
3183 real cmap_dihs(int nbonds,
3184 const t_iatom forceatoms[], const t_iparams forceparams[],
3185 const gmx_cmap_t *cmap_grid,
3186 const rvec x[], rvec f[], rvec fshift[],
3187 const t_pbc *pbc, const t_graph *g,
3188 real gmx_unused lambda, real gmx_unused *dvdlambda,
3189 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3190 int gmx_unused *global_atom_index)
3192 int i, j, k, n, idx;
3193 int ai, aj, ak, al, am;
3194 int a1i, a1j, a1k, a1l, a2i, a2j, a2k, a2l;
3196 int t11, t21, t31, t12, t22, t32;
3197 int iphi1, ip1m1, ip1p1, ip1p2;
3198 int iphi2, ip2m1, ip2p1, ip2p2;
3200 int pos1, pos2, pos3, pos4;
3202 real ty[4], ty1[4], ty2[4], ty12[4], tc[16], tx[16];
3203 real phi1, cos_phi1, sin_phi1, sign1, xphi1;
3204 real phi2, cos_phi2, sin_phi2, sign2, xphi2;
3205 real dx, xx, tt, tu, e, df1, df2, vtot;
3206 real ra21, rb21, rg21, rg1, rgr1, ra2r1, rb2r1, rabr1;
3207 real ra22, rb22, rg22, rg2, rgr2, ra2r2, rb2r2, rabr2;
3208 real fg1, hg1, fga1, hgb1, gaa1, gbb1;
3209 real fg2, hg2, fga2, hgb2, gaa2, gbb2;
3212 rvec r1_ij, r1_kj, r1_kl, m1, n1;
3213 rvec r2_ij, r2_kj, r2_kl, m2, n2;
3214 rvec f1_i, f1_j, f1_k, f1_l;
3215 rvec f2_i, f2_j, f2_k, f2_l;
3216 rvec a1, b1, a2, b2;
3217 rvec f1, g1, h1, f2, g2, h2;
3218 rvec dtf1, dtg1, dth1, dtf2, dtg2, dth2;
3219 ivec jt1, dt1_ij, dt1_kj, dt1_lj;
3220 ivec jt2, dt2_ij, dt2_kj, dt2_lj;
3224 int loop_index[4][4] = {
3231 /* Total CMAP energy */
3234 for (n = 0; n < nbonds; )
3236 /* Five atoms are involved in the two torsions */
3237 type = forceatoms[n++];
3238 ai = forceatoms[n++];
3239 aj = forceatoms[n++];
3240 ak = forceatoms[n++];
3241 al = forceatoms[n++];
3242 am = forceatoms[n++];
3244 /* Which CMAP type is this */
3245 cmapA = forceparams[type].cmap.cmapA;
3246 cmapd = cmap_grid->cmapdata[cmapA].cmap;
3254 phi1 = dih_angle(x[a1i], x[a1j], x[a1k], x[a1l], pbc, r1_ij, r1_kj, r1_kl, m1, n1,
3255 &sign1, &t11, &t21, &t31); /* 84 */
3257 cos_phi1 = cos(phi1);
3259 a1[0] = r1_ij[1]*r1_kj[2]-r1_ij[2]*r1_kj[1];
3260 a1[1] = r1_ij[2]*r1_kj[0]-r1_ij[0]*r1_kj[2];
3261 a1[2] = r1_ij[0]*r1_kj[1]-r1_ij[1]*r1_kj[0]; /* 9 */
3263 b1[0] = r1_kl[1]*r1_kj[2]-r1_kl[2]*r1_kj[1];
3264 b1[1] = r1_kl[2]*r1_kj[0]-r1_kl[0]*r1_kj[2];
3265 b1[2] = r1_kl[0]*r1_kj[1]-r1_kl[1]*r1_kj[0]; /* 9 */
3267 pbc_rvec_sub(pbc, x[a1l], x[a1k], h1);
3269 ra21 = iprod(a1, a1); /* 5 */
3270 rb21 = iprod(b1, b1); /* 5 */
3271 rg21 = iprod(r1_kj, r1_kj); /* 5 */
3277 rabr1 = sqrt(ra2r1*rb2r1);
3279 sin_phi1 = rg1 * rabr1 * iprod(a1, h1) * (-1);
3281 if (cos_phi1 < -0.5 || cos_phi1 > 0.5)
3283 phi1 = asin(sin_phi1);
3293 phi1 = -M_PI - phi1;
3299 phi1 = acos(cos_phi1);
3307 xphi1 = phi1 + M_PI; /* 1 */
3309 /* Second torsion */
3315 phi2 = dih_angle(x[a2i], x[a2j], x[a2k], x[a2l], pbc, r2_ij, r2_kj, r2_kl, m2, n2,
3316 &sign2, &t12, &t22, &t32); /* 84 */
3318 cos_phi2 = cos(phi2);
3320 a2[0] = r2_ij[1]*r2_kj[2]-r2_ij[2]*r2_kj[1];
3321 a2[1] = r2_ij[2]*r2_kj[0]-r2_ij[0]*r2_kj[2];
3322 a2[2] = r2_ij[0]*r2_kj[1]-r2_ij[1]*r2_kj[0]; /* 9 */
3324 b2[0] = r2_kl[1]*r2_kj[2]-r2_kl[2]*r2_kj[1];
3325 b2[1] = r2_kl[2]*r2_kj[0]-r2_kl[0]*r2_kj[2];
3326 b2[2] = r2_kl[0]*r2_kj[1]-r2_kl[1]*r2_kj[0]; /* 9 */
3328 pbc_rvec_sub(pbc, x[a2l], x[a2k], h2);
3330 ra22 = iprod(a2, a2); /* 5 */
3331 rb22 = iprod(b2, b2); /* 5 */
3332 rg22 = iprod(r2_kj, r2_kj); /* 5 */
3338 rabr2 = sqrt(ra2r2*rb2r2);
3340 sin_phi2 = rg2 * rabr2 * iprod(a2, h2) * (-1);
3342 if (cos_phi2 < -0.5 || cos_phi2 > 0.5)
3344 phi2 = asin(sin_phi2);
3354 phi2 = -M_PI - phi2;
3360 phi2 = acos(cos_phi2);
3368 xphi2 = phi2 + M_PI; /* 1 */
3370 /* Range mangling */
3373 xphi1 = xphi1 + 2*M_PI;
3375 else if (xphi1 >= 2*M_PI)
3377 xphi1 = xphi1 - 2*M_PI;
3382 xphi2 = xphi2 + 2*M_PI;
3384 else if (xphi2 >= 2*M_PI)
3386 xphi2 = xphi2 - 2*M_PI;
3389 /* Number of grid points */
3390 dx = 2*M_PI / cmap_grid->grid_spacing;
3392 /* Where on the grid are we */
3393 iphi1 = static_cast<int>(xphi1/dx);
3394 iphi2 = static_cast<int>(xphi2/dx);
3396 iphi1 = cmap_setup_grid_index(iphi1, cmap_grid->grid_spacing, &ip1m1, &ip1p1, &ip1p2);
3397 iphi2 = cmap_setup_grid_index(iphi2, cmap_grid->grid_spacing, &ip2m1, &ip2p1, &ip2p2);
3399 pos1 = iphi1*cmap_grid->grid_spacing+iphi2;
3400 pos2 = ip1p1*cmap_grid->grid_spacing+iphi2;
3401 pos3 = ip1p1*cmap_grid->grid_spacing+ip2p1;
3402 pos4 = iphi1*cmap_grid->grid_spacing+ip2p1;
3404 ty[0] = cmapd[pos1*4];
3405 ty[1] = cmapd[pos2*4];
3406 ty[2] = cmapd[pos3*4];
3407 ty[3] = cmapd[pos4*4];
3409 ty1[0] = cmapd[pos1*4+1];
3410 ty1[1] = cmapd[pos2*4+1];
3411 ty1[2] = cmapd[pos3*4+1];
3412 ty1[3] = cmapd[pos4*4+1];
3414 ty2[0] = cmapd[pos1*4+2];
3415 ty2[1] = cmapd[pos2*4+2];
3416 ty2[2] = cmapd[pos3*4+2];
3417 ty2[3] = cmapd[pos4*4+2];
3419 ty12[0] = cmapd[pos1*4+3];
3420 ty12[1] = cmapd[pos2*4+3];
3421 ty12[2] = cmapd[pos3*4+3];
3422 ty12[3] = cmapd[pos4*4+3];
3424 /* Switch to degrees */
3425 dx = 360.0 / cmap_grid->grid_spacing;
3426 xphi1 = xphi1 * RAD2DEG;
3427 xphi2 = xphi2 * RAD2DEG;
3429 for (i = 0; i < 4; i++) /* 16 */
3432 tx[i+4] = ty1[i]*dx;
3433 tx[i+8] = ty2[i]*dx;
3434 tx[i+12] = ty12[i]*dx*dx;
3438 for (i = 0; i < 4; i++) /* 1056 */
3440 for (j = 0; j < 4; j++)
3443 for (k = 0; k < 16; k++)
3445 xx = xx + cmap_coeff_matrix[k*16+idx]*tx[k];
3453 tt = (xphi1-iphi1*dx)/dx;
3454 tu = (xphi2-iphi2*dx)/dx;
3460 for (i = 3; i >= 0; i--)
3462 l1 = loop_index[i][3];
3463 l2 = loop_index[i][2];
3464 l3 = loop_index[i][1];
3466 e = tt * e + ((tc[i*4+3]*tu+tc[i*4+2])*tu + tc[i*4+1])*tu+tc[i*4];
3467 df1 = tu * df1 + (3.0*tc[l1]*tt+2.0*tc[l2])*tt+tc[l3];
3468 df2 = tt * df2 + (3.0*tc[i*4+3]*tu+2.0*tc[i*4+2])*tu+tc[i*4+1];
3478 /* Do forces - first torsion */
3479 fg1 = iprod(r1_ij, r1_kj);
3480 hg1 = iprod(r1_kl, r1_kj);
3481 fga1 = fg1*ra2r1*rgr1;
3482 hgb1 = hg1*rb2r1*rgr1;
3486 for (i = 0; i < DIM; i++)
3488 dtf1[i] = gaa1 * a1[i];
3489 dtg1[i] = fga1 * a1[i] - hgb1 * b1[i];
3490 dth1[i] = gbb1 * b1[i];
3492 f1[i] = df1 * dtf1[i];
3493 g1[i] = df1 * dtg1[i];
3494 h1[i] = df1 * dth1[i];
3497 f1_j[i] = -f1[i] - g1[i];
3498 f1_k[i] = h1[i] + g1[i];
3501 f[a1i][i] = f[a1i][i] + f1_i[i];
3502 f[a1j][i] = f[a1j][i] + f1_j[i]; /* - f1[i] - g1[i] */
3503 f[a1k][i] = f[a1k][i] + f1_k[i]; /* h1[i] + g1[i] */
3504 f[a1l][i] = f[a1l][i] + f1_l[i]; /* h1[i] */
3507 /* Do forces - second torsion */
3508 fg2 = iprod(r2_ij, r2_kj);
3509 hg2 = iprod(r2_kl, r2_kj);
3510 fga2 = fg2*ra2r2*rgr2;
3511 hgb2 = hg2*rb2r2*rgr2;
3515 for (i = 0; i < DIM; i++)
3517 dtf2[i] = gaa2 * a2[i];
3518 dtg2[i] = fga2 * a2[i] - hgb2 * b2[i];
3519 dth2[i] = gbb2 * b2[i];
3521 f2[i] = df2 * dtf2[i];
3522 g2[i] = df2 * dtg2[i];
3523 h2[i] = df2 * dth2[i];
3526 f2_j[i] = -f2[i] - g2[i];
3527 f2_k[i] = h2[i] + g2[i];
3530 f[a2i][i] = f[a2i][i] + f2_i[i]; /* f2[i] */
3531 f[a2j][i] = f[a2j][i] + f2_j[i]; /* - f2[i] - g2[i] */
3532 f[a2k][i] = f[a2k][i] + f2_k[i]; /* h2[i] + g2[i] */
3533 f[a2l][i] = f[a2l][i] + f2_l[i]; /* - h2[i] */
3539 copy_ivec(SHIFT_IVEC(g, a1j), jt1);
3540 ivec_sub(SHIFT_IVEC(g, a1i), jt1, dt1_ij);
3541 ivec_sub(SHIFT_IVEC(g, a1k), jt1, dt1_kj);
3542 ivec_sub(SHIFT_IVEC(g, a1l), jt1, dt1_lj);
3543 t11 = IVEC2IS(dt1_ij);
3544 t21 = IVEC2IS(dt1_kj);
3545 t31 = IVEC2IS(dt1_lj);
3547 copy_ivec(SHIFT_IVEC(g, a2j), jt2);
3548 ivec_sub(SHIFT_IVEC(g, a2i), jt2, dt2_ij);
3549 ivec_sub(SHIFT_IVEC(g, a2k), jt2, dt2_kj);
3550 ivec_sub(SHIFT_IVEC(g, a2l), jt2, dt2_lj);
3551 t12 = IVEC2IS(dt2_ij);
3552 t22 = IVEC2IS(dt2_kj);
3553 t32 = IVEC2IS(dt2_lj);
3557 t31 = pbc_rvec_sub(pbc, x[a1l], x[a1j], h1);
3558 t32 = pbc_rvec_sub(pbc, x[a2l], x[a2j], h2);
3566 rvec_inc(fshift[t11], f1_i);
3567 rvec_inc(fshift[CENTRAL], f1_j);
3568 rvec_inc(fshift[t21], f1_k);
3569 rvec_inc(fshift[t31], f1_l);
3571 rvec_inc(fshift[t21], f2_i);
3572 rvec_inc(fshift[CENTRAL], f2_j);
3573 rvec_inc(fshift[t22], f2_k);
3574 rvec_inc(fshift[t32], f2_l);
3581 /***********************************************************
3583 * G R O M O S 9 6 F U N C T I O N S
3585 ***********************************************************/
3586 real g96harmonic(real kA, real kB, real xA, real xB, real x, real lambda,
3589 const real half = 0.5;
3590 real L1, kk, x0, dx, dx2;
3591 real v, f, dvdlambda;
3594 kk = L1*kA+lambda*kB;
3595 x0 = L1*xA+lambda*xB;
3602 dvdlambda = half*(kB-kA)*dx2 + (xA-xB)*kk*dx;
3609 /* That was 21 flops */
3612 real g96bonds(int nbonds,
3613 const t_iatom forceatoms[], const t_iparams forceparams[],
3614 const rvec x[], rvec f[], rvec fshift[],
3615 const t_pbc *pbc, const t_graph *g,
3616 real lambda, real *dvdlambda,
3617 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3618 int gmx_unused *global_atom_index)
3620 int i, m, ki, ai, aj, type;
3621 real dr2, fbond, vbond, fij, vtot;
3626 for (i = 0; (i < nbonds); )
3628 type = forceatoms[i++];
3629 ai = forceatoms[i++];
3630 aj = forceatoms[i++];
3632 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
3633 dr2 = iprod(dx, dx); /* 5 */
3635 *dvdlambda += g96harmonic(forceparams[type].harmonic.krA,
3636 forceparams[type].harmonic.krB,
3637 forceparams[type].harmonic.rA,
3638 forceparams[type].harmonic.rB,
3639 dr2, lambda, &vbond, &fbond);
3641 vtot += 0.5*vbond; /* 1*/
3645 fprintf(debug, "G96-BONDS: dr = %10g vbond = %10g fbond = %10g\n",
3646 sqrt(dr2), vbond, fbond);
3652 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
3655 for (m = 0; (m < DIM); m++) /* 15 */
3660 fshift[ki][m] += fij;
3661 fshift[CENTRAL][m] -= fij;
3667 real g96bond_angle(const rvec xi, const rvec xj, const rvec xk, const t_pbc *pbc,
3668 rvec r_ij, rvec r_kj,
3670 /* Return value is the angle between the bonds i-j and j-k */
3674 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
3675 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
3677 costh = cos_angle(r_ij, r_kj); /* 25 */
3682 real g96angles(int nbonds,
3683 const t_iatom forceatoms[], const t_iparams forceparams[],
3684 const rvec x[], rvec f[], rvec fshift[],
3685 const t_pbc *pbc, const t_graph *g,
3686 real lambda, real *dvdlambda,
3687 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3688 int gmx_unused *global_atom_index)
3690 int i, ai, aj, ak, type, m, t1, t2;
3692 real cos_theta, dVdt, va, vtot;
3693 real rij_1, rij_2, rkj_1, rkj_2, rijrkj_1;
3695 ivec jt, dt_ij, dt_kj;
3698 for (i = 0; (i < nbonds); )
3700 type = forceatoms[i++];
3701 ai = forceatoms[i++];
3702 aj = forceatoms[i++];
3703 ak = forceatoms[i++];
3705 cos_theta = g96bond_angle(x[ai], x[aj], x[ak], pbc, r_ij, r_kj, &t1, &t2);
3707 *dvdlambda += g96harmonic(forceparams[type].harmonic.krA,
3708 forceparams[type].harmonic.krB,
3709 forceparams[type].harmonic.rA,
3710 forceparams[type].harmonic.rB,
3711 cos_theta, lambda, &va, &dVdt);
3714 rij_1 = gmx_invsqrt(iprod(r_ij, r_ij));
3715 rkj_1 = gmx_invsqrt(iprod(r_kj, r_kj));
3716 rij_2 = rij_1*rij_1;
3717 rkj_2 = rkj_1*rkj_1;
3718 rijrkj_1 = rij_1*rkj_1; /* 23 */
3723 fprintf(debug, "G96ANGLES: costheta = %10g vth = %10g dV/dct = %10g\n",
3724 cos_theta, va, dVdt);
3727 for (m = 0; (m < DIM); m++) /* 42 */
3729 f_i[m] = dVdt*(r_kj[m]*rijrkj_1 - r_ij[m]*rij_2*cos_theta);
3730 f_k[m] = dVdt*(r_ij[m]*rijrkj_1 - r_kj[m]*rkj_2*cos_theta);
3731 f_j[m] = -f_i[m]-f_k[m];
3739 copy_ivec(SHIFT_IVEC(g, aj), jt);
3741 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3742 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3743 t1 = IVEC2IS(dt_ij);
3744 t2 = IVEC2IS(dt_kj);
3746 rvec_inc(fshift[t1], f_i);
3747 rvec_inc(fshift[CENTRAL], f_j);
3748 rvec_inc(fshift[t2], f_k); /* 9 */
3754 real cross_bond_bond(int nbonds,
3755 const t_iatom forceatoms[], const t_iparams forceparams[],
3756 const rvec x[], rvec f[], rvec fshift[],
3757 const t_pbc *pbc, const t_graph *g,
3758 real gmx_unused lambda, real gmx_unused *dvdlambda,
3759 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3760 int gmx_unused *global_atom_index)
3762 /* Potential from Lawrence and Skimmer, Chem. Phys. Lett. 372 (2003)
3765 int i, ai, aj, ak, type, m, t1, t2;
3767 real vtot, vrr, s1, s2, r1, r2, r1e, r2e, krr;
3769 ivec jt, dt_ij, dt_kj;
3772 for (i = 0; (i < nbonds); )
3774 type = forceatoms[i++];
3775 ai = forceatoms[i++];
3776 aj = forceatoms[i++];
3777 ak = forceatoms[i++];
3778 r1e = forceparams[type].cross_bb.r1e;
3779 r2e = forceparams[type].cross_bb.r2e;
3780 krr = forceparams[type].cross_bb.krr;
3782 /* Compute distance vectors ... */
3783 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
3784 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
3786 /* ... and their lengths */
3790 /* Deviations from ideality */
3794 /* Energy (can be negative!) */
3799 svmul(-krr*s2/r1, r_ij, f_i);
3800 svmul(-krr*s1/r2, r_kj, f_k);
3802 for (m = 0; (m < DIM); m++) /* 12 */
3804 f_j[m] = -f_i[m] - f_k[m];
3813 copy_ivec(SHIFT_IVEC(g, aj), jt);
3815 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3816 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3817 t1 = IVEC2IS(dt_ij);
3818 t2 = IVEC2IS(dt_kj);
3820 rvec_inc(fshift[t1], f_i);
3821 rvec_inc(fshift[CENTRAL], f_j);
3822 rvec_inc(fshift[t2], f_k); /* 9 */
3828 real cross_bond_angle(int nbonds,
3829 const t_iatom forceatoms[], const t_iparams forceparams[],
3830 const rvec x[], rvec f[], rvec fshift[],
3831 const t_pbc *pbc, const t_graph *g,
3832 real gmx_unused lambda, real gmx_unused *dvdlambda,
3833 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3834 int gmx_unused *global_atom_index)
3836 /* Potential from Lawrence and Skimmer, Chem. Phys. Lett. 372 (2003)
3839 int i, ai, aj, ak, type, m, t1, t2;
3840 rvec r_ij, r_kj, r_ik;
3841 real vtot, vrt, s1, s2, s3, r1, r2, r3, r1e, r2e, r3e, krt, k1, k2, k3;
3843 ivec jt, dt_ij, dt_kj;
3846 for (i = 0; (i < nbonds); )
3848 type = forceatoms[i++];
3849 ai = forceatoms[i++];
3850 aj = forceatoms[i++];
3851 ak = forceatoms[i++];
3852 r1e = forceparams[type].cross_ba.r1e;
3853 r2e = forceparams[type].cross_ba.r2e;
3854 r3e = forceparams[type].cross_ba.r3e;
3855 krt = forceparams[type].cross_ba.krt;
3857 /* Compute distance vectors ... */
3858 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
3859 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
3860 pbc_rvec_sub(pbc, x[ai], x[ak], r_ik);
3862 /* ... and their lengths */
3867 /* Deviations from ideality */
3872 /* Energy (can be negative!) */
3873 vrt = krt*s3*(s1+s2);
3879 k3 = -krt*(s1+s2)/r3;
3880 for (m = 0; (m < DIM); m++)
3882 f_i[m] = k1*r_ij[m] + k3*r_ik[m];
3883 f_k[m] = k2*r_kj[m] - k3*r_ik[m];
3884 f_j[m] = -f_i[m] - f_k[m];
3887 for (m = 0; (m < DIM); m++) /* 12 */
3897 copy_ivec(SHIFT_IVEC(g, aj), jt);
3899 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3900 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3901 t1 = IVEC2IS(dt_ij);
3902 t2 = IVEC2IS(dt_kj);
3904 rvec_inc(fshift[t1], f_i);
3905 rvec_inc(fshift[CENTRAL], f_j);
3906 rvec_inc(fshift[t2], f_k); /* 9 */
3912 static real bonded_tab(const char *type, int table_nr,
3913 const bondedtable_t *table, real kA, real kB, real r,
3914 real lambda, real *V, real *F)
3916 real k, tabscale, *VFtab, rt, eps, eps2, Yt, Ft, Geps, Heps2, Fp, VV, FF;
3920 k = (1.0 - lambda)*kA + lambda*kB;
3922 tabscale = table->scale;
3923 VFtab = table->data;
3926 n0 = static_cast<int>(rt);
3929 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",
3930 type, table_nr, r, n0, n0+1, table->n);
3937 Geps = VFtab[nnn+2]*eps;
3938 Heps2 = VFtab[nnn+3]*eps2;
3939 Fp = Ft + Geps + Heps2;
3941 FF = Fp + Geps + 2.0*Heps2;
3943 *F = -k*FF*tabscale;
3945 dvdlambda = (kB - kA)*VV;
3949 /* That was 22 flops */
3952 real tab_bonds(int nbonds,
3953 const t_iatom forceatoms[], const t_iparams forceparams[],
3954 const rvec x[], rvec f[], rvec fshift[],
3955 const t_pbc *pbc, const t_graph *g,
3956 real lambda, real *dvdlambda,
3957 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
3958 int gmx_unused *global_atom_index)
3960 int i, m, ki, ai, aj, type, table;
3961 real dr, dr2, fbond, vbond, fij, vtot;
3966 for (i = 0; (i < nbonds); )
3968 type = forceatoms[i++];
3969 ai = forceatoms[i++];
3970 aj = forceatoms[i++];
3972 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
3973 dr2 = iprod(dx, dx); /* 5 */
3974 dr = dr2*gmx_invsqrt(dr2); /* 10 */
3976 table = forceparams[type].tab.table;
3978 *dvdlambda += bonded_tab("bond", table,
3979 &fcd->bondtab[table],
3980 forceparams[type].tab.kA,
3981 forceparams[type].tab.kB,
3982 dr, lambda, &vbond, &fbond); /* 22 */
3990 vtot += vbond; /* 1*/
3991 fbond *= gmx_invsqrt(dr2); /* 6 */
3995 fprintf(debug, "TABBONDS: dr = %10g vbond = %10g fbond = %10g\n",
4001 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
4004 for (m = 0; (m < DIM); m++) /* 15 */
4009 fshift[ki][m] += fij;
4010 fshift[CENTRAL][m] -= fij;
4016 real tab_angles(int nbonds,
4017 const t_iatom forceatoms[], const t_iparams forceparams[],
4018 const rvec x[], rvec f[], rvec fshift[],
4019 const t_pbc *pbc, const t_graph *g,
4020 real lambda, real *dvdlambda,
4021 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
4022 int gmx_unused *global_atom_index)
4024 int i, ai, aj, ak, t1, t2, type, table;
4026 real cos_theta, cos_theta2, theta, dVdt, va, vtot;
4027 ivec jt, dt_ij, dt_kj;
4030 for (i = 0; (i < nbonds); )
4032 type = forceatoms[i++];
4033 ai = forceatoms[i++];
4034 aj = forceatoms[i++];
4035 ak = forceatoms[i++];
4037 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
4038 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
4040 table = forceparams[type].tab.table;
4042 *dvdlambda += bonded_tab("angle", table,
4043 &fcd->angletab[table],
4044 forceparams[type].tab.kA,
4045 forceparams[type].tab.kB,
4046 theta, lambda, &va, &dVdt); /* 22 */
4049 cos_theta2 = sqr(cos_theta); /* 1 */
4058 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
4059 sth = st*cos_theta; /* 1 */
4063 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
4064 theta*RAD2DEG, va, dVdt);
4067 nrkj2 = iprod(r_kj, r_kj); /* 5 */
4068 nrij2 = iprod(r_ij, r_ij);
4070 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
4071 cii = sth/nrij2; /* 10 */
4072 ckk = sth/nrkj2; /* 10 */
4074 for (m = 0; (m < DIM); m++) /* 39 */
4076 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
4077 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
4078 f_j[m] = -f_i[m]-f_k[m];
4085 copy_ivec(SHIFT_IVEC(g, aj), jt);
4087 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
4088 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
4089 t1 = IVEC2IS(dt_ij);
4090 t2 = IVEC2IS(dt_kj);
4092 rvec_inc(fshift[t1], f_i);
4093 rvec_inc(fshift[CENTRAL], f_j);
4094 rvec_inc(fshift[t2], f_k);
4100 real tab_dihs(int nbonds,
4101 const t_iatom forceatoms[], const t_iparams forceparams[],
4102 const rvec x[], rvec f[], rvec fshift[],
4103 const t_pbc *pbc, const t_graph *g,
4104 real lambda, real *dvdlambda,
4105 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
4106 int gmx_unused *global_atom_index)
4108 int i, type, ai, aj, ak, al, table;
4110 rvec r_ij, r_kj, r_kl, m, n;
4111 real phi, sign, ddphi, vpd, vtot;
4114 for (i = 0; (i < nbonds); )
4116 type = forceatoms[i++];
4117 ai = forceatoms[i++];
4118 aj = forceatoms[i++];
4119 ak = forceatoms[i++];
4120 al = forceatoms[i++];
4122 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
4123 &sign, &t1, &t2, &t3); /* 84 */
4125 table = forceparams[type].tab.table;
4127 /* Hopefully phi+M_PI never results in values < 0 */
4128 *dvdlambda += bonded_tab("dihedral", table,
4129 &fcd->dihtab[table],
4130 forceparams[type].tab.kA,
4131 forceparams[type].tab.kB,
4132 phi+M_PI, lambda, &vpd, &ddphi);
4135 do_dih_fup(ai, aj, ak, al, -ddphi, r_ij, r_kj, r_kl, m, n,
4136 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
4139 fprintf(debug, "pdih: (%d,%d,%d,%d) phi=%g\n",
4140 ai, aj, ak, al, phi);
4147 /* TODO This function could go away when idef is not a big bucket of
4149 gmx_bool ftype_is_bonded_potential(int ftype)
4152 (interaction_function[ftype].flags & IF_BOND) &&
4153 !(ftype == F_CONNBONDS || ftype == F_POSRES || ftype == F_FBPOSRES) &&
4154 (ftype < F_GB12 || ftype > F_GB14);
4157 static void zero_thread_forces(f_thread_t *f_t, int n,
4158 int nblock, int blocksize)
4160 int b, a0, a1, a, i, j;
4162 if (n > f_t->f_nalloc)
4164 f_t->f_nalloc = over_alloc_large(n);
4165 srenew(f_t->f, f_t->f_nalloc);
4168 if (f_t->red_mask != 0)
4170 for (b = 0; b < nblock; b++)
4172 if (f_t->red_mask && (1U<<b))
4175 a1 = std::min((b+1)*blocksize, n);
4176 for (a = a0; a < a1; a++)
4178 clear_rvec(f_t->f[a]);
4183 for (i = 0; i < SHIFTS; i++)
4185 clear_rvec(f_t->fshift[i]);
4187 for (i = 0; i < F_NRE; i++)
4191 for (i = 0; i < egNR; i++)
4193 for (j = 0; j < f_t->grpp.nener; j++)
4195 f_t->grpp.ener[i][j] = 0;
4198 for (i = 0; i < efptNR; i++)
4204 static void reduce_thread_force_buffer(int n, rvec *f,
4205 int nthreads, f_thread_t *f_t,
4206 int nblock, int block_size)
4208 /* The max thread number is arbitrary,
4209 * we used a fixed number to avoid memory management.
4210 * Using more than 16 threads is probably never useful performance wise.
4212 #define MAX_BONDED_THREADS 256
4215 if (nthreads > MAX_BONDED_THREADS)
4217 gmx_fatal(FARGS, "Can not reduce bonded forces on more than %d threads",
4218 MAX_BONDED_THREADS);
4221 /* This reduction can run on any number of threads,
4222 * independently of nthreads.
4224 #pragma omp parallel for num_threads(nthreads) schedule(static)
4225 for (b = 0; b < nblock; b++)
4227 rvec *fp[MAX_BONDED_THREADS];
4231 /* Determine which threads contribute to this block */
4233 for (ft = 1; ft < nthreads; ft++)
4235 if (f_t[ft].red_mask & (1U<<b))
4237 fp[nfb++] = f_t[ft].f;
4242 /* Reduce force buffers for threads that contribute */
4244 a1 = (b+1)*block_size;
4245 a1 = std::min(a1, n);
4246 for (a = a0; a < a1; a++)
4248 for (fb = 0; fb < nfb; fb++)
4250 rvec_inc(f[a], fp[fb][a]);
4257 static void reduce_thread_forces(int n, rvec *f, rvec *fshift,
4258 real *ener, gmx_grppairener_t *grpp, real *dvdl,
4259 int nthreads, f_thread_t *f_t,
4260 int nblock, int block_size,
4261 gmx_bool bCalcEnerVir,
4266 /* Reduce the bonded force buffer */
4267 reduce_thread_force_buffer(n, f, nthreads, f_t, nblock, block_size);
4270 /* When necessary, reduce energy and virial using one thread only */
4275 for (i = 0; i < SHIFTS; i++)
4277 for (t = 1; t < nthreads; t++)
4279 rvec_inc(fshift[i], f_t[t].fshift[i]);
4282 for (i = 0; i < F_NRE; i++)
4284 for (t = 1; t < nthreads; t++)
4286 ener[i] += f_t[t].ener[i];
4289 for (i = 0; i < egNR; i++)
4291 for (j = 0; j < f_t[1].grpp.nener; j++)
4293 for (t = 1; t < nthreads; t++)
4296 grpp->ener[i][j] += f_t[t].grpp.ener[i][j];
4302 for (i = 0; i < efptNR; i++)
4305 for (t = 1; t < nthreads; t++)
4307 dvdl[i] += f_t[t].dvdl[i];
4314 static real calc_one_bond(int thread,
4315 int ftype, const t_idef *idef,
4316 const rvec x[], rvec f[], rvec fshift[],
4318 const t_pbc *pbc, const t_graph *g,
4319 gmx_grppairener_t *grpp,
4321 real *lambda, real *dvdl,
4322 const t_mdatoms *md, t_fcdata *fcd,
4323 gmx_bool bCalcEnerVir,
4324 int *global_atom_index)
4326 int nat1, nbonds, efptFTYPE;
4331 if (IS_RESTRAINT_TYPE(ftype))
4333 efptFTYPE = efptRESTRAINT;
4337 efptFTYPE = efptBONDED;
4340 nat1 = interaction_function[ftype].nratoms + 1;
4341 nbonds = idef->il[ftype].nr/nat1;
4342 iatoms = idef->il[ftype].iatoms;
4344 nb0 = idef->il_thread_division[ftype*(idef->nthreads+1)+thread];
4345 nbn = idef->il_thread_division[ftype*(idef->nthreads+1)+thread+1] - nb0;
4347 if (!IS_LISTED_LJ_C(ftype))
4349 if (ftype == F_CMAP)
4351 v = cmap_dihs(nbn, iatoms+nb0,
4352 idef->iparams, &idef->cmap_grid,
4354 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
4355 md, fcd, global_atom_index);
4357 #ifdef GMX_SIMD_HAVE_REAL
4358 else if (ftype == F_ANGLES &&
4359 !bCalcEnerVir && fr->efep == efepNO)
4361 /* No energies, shift forces, dvdl */
4362 angles_noener_simd(nbn, idef->il[ftype].iatoms+nb0,
4365 pbc, g, lambda[efptFTYPE], md, fcd,
4370 else if (ftype == F_PDIHS &&
4371 !bCalcEnerVir && fr->efep == efepNO)
4373 /* No energies, shift forces, dvdl */
4374 #ifdef GMX_SIMD_HAVE_REAL
4379 (nbn, idef->il[ftype].iatoms+nb0,
4382 pbc, g, lambda[efptFTYPE], md, fcd,
4388 v = interaction_function[ftype].ifunc(nbn, iatoms+nb0,
4391 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
4392 md, fcd, global_atom_index);
4397 v = do_nonbonded_listed(ftype, nbn, iatoms+nb0, idef->iparams, x, f, fshift,
4398 pbc, g, lambda, dvdl, md, fr, grpp, global_atom_index);
4403 inc_nrnb(nrnb, interaction_function[ftype].nrnb_ind, nbonds);
4409 void calc_bonds(const gmx_multisim_t *ms,
4411 const rvec x[], history_t *hist,
4412 rvec f[], t_forcerec *fr,
4413 const t_pbc *pbc, const t_graph *g,
4414 gmx_enerdata_t *enerd, t_nrnb *nrnb,
4416 const t_mdatoms *md,
4417 t_fcdata *fcd, int *global_atom_index,
4420 gmx_bool bCalcEnerVir;
4422 real dvdl[efptNR]; /* The dummy array is to have a place to store the dhdl at other values
4423 of lambda, which will be thrown away in the end*/
4424 const t_pbc *pbc_null;
4427 assert(fr->nthreads == idef->nthreads);
4429 bCalcEnerVir = (force_flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY));
4431 for (i = 0; i < efptNR; i++)
4447 p_graph(debug, "Bondage is fun", g);
4451 /* Do pre force calculation stuff which might require communication */
4452 if (idef->il[F_ORIRES].nr)
4454 enerd->term[F_ORIRESDEV] =
4455 calc_orires_dev(ms, idef->il[F_ORIRES].nr,
4456 idef->il[F_ORIRES].iatoms,
4457 idef->iparams, md, x,
4458 pbc_null, fcd, hist);
4460 if (idef->il[F_DISRES].nr)
4462 calc_disres_R_6(idef->il[F_DISRES].nr,
4463 idef->il[F_DISRES].iatoms,
4464 idef->iparams, x, pbc_null,
4467 if (fcd->disres.nsystems > 1)
4469 gmx_sum_sim(2*fcd->disres.nres, fcd->disres.Rt_6, ms);
4474 #pragma omp parallel for num_threads(fr->nthreads) schedule(static)
4475 for (thread = 0; thread < fr->nthreads; thread++)
4482 gmx_grppairener_t *grpp;
4487 fshift = fr->fshift;
4489 grpp = &enerd->grpp;
4494 zero_thread_forces(&fr->f_t[thread], fr->natoms_force,
4495 fr->red_nblock, 1<<fr->red_ashift);
4497 ft = fr->f_t[thread].f;
4498 fshift = fr->f_t[thread].fshift;
4499 epot = fr->f_t[thread].ener;
4500 grpp = &fr->f_t[thread].grpp;
4501 dvdlt = fr->f_t[thread].dvdl;
4503 /* Loop over all bonded force types to calculate the bonded forces */
4504 for (ftype = 0; (ftype < F_NRE); ftype++)
4506 if (idef->il[ftype].nr > 0 && ftype_is_bonded_potential(ftype))
4508 v = calc_one_bond(thread, ftype, idef, x,
4509 ft, fshift, fr, pbc_null, g, grpp,
4510 nrnb, lambda, dvdlt,
4511 md, fcd, bCalcEnerVir,
4517 if (fr->nthreads > 1)
4519 reduce_thread_forces(fr->natoms_force, f, fr->fshift,
4520 enerd->term, &enerd->grpp, dvdl,
4521 fr->nthreads, fr->f_t,
4522 fr->red_nblock, 1<<fr->red_ashift,
4524 force_flags & GMX_FORCE_DHDL);
4526 if (force_flags & GMX_FORCE_DHDL)
4528 for (i = 0; i < efptNR; i++)
4530 enerd->dvdl_nonlin[i] += dvdl[i];
4534 /* Copy the sum of violations for the distance restraints from fcd */
4537 enerd->term[F_DISRESVIOL] = fcd->disres.sumviol;
4542 void calc_bonds_lambda(const t_idef *idef,
4545 const t_pbc *pbc, const t_graph *g,
4546 gmx_grppairener_t *grpp, real *epot, t_nrnb *nrnb,
4548 const t_mdatoms *md,
4550 int *global_atom_index)
4552 int ftype, nr_nonperturbed, nr;
4554 real dvdl_dum[efptNR];
4556 const t_pbc *pbc_null;
4568 /* Copy the whole idef, so we can modify the contents locally */
4570 idef_fe.nthreads = 1;
4571 snew(idef_fe.il_thread_division, F_NRE*(idef_fe.nthreads+1));
4573 /* We already have the forces, so we use temp buffers here */
4574 snew(f, fr->natoms_force);
4575 snew(fshift, SHIFTS);
4577 /* Loop over all bonded force types to calculate the bonded energies */
4578 for (ftype = 0; (ftype < F_NRE); ftype++)
4580 if (ftype_is_bonded_potential(ftype))
4582 /* Set the work range of thread 0 to the perturbed bondeds only */
4583 nr_nonperturbed = idef->il[ftype].nr_nonperturbed;
4584 nr = idef->il[ftype].nr;
4585 idef_fe.il_thread_division[ftype*2+0] = nr_nonperturbed;
4586 idef_fe.il_thread_division[ftype*2+1] = nr;
4588 /* This is only to get the flop count correct */
4589 idef_fe.il[ftype].nr = nr - nr_nonperturbed;
4591 if (nr - nr_nonperturbed > 0)
4593 v = calc_one_bond(0, ftype, &idef_fe,
4594 x, f, fshift, fr, pbc_null, g,
4595 grpp, nrnb, lambda, dvdl_dum,
4606 sfree(idef_fe.il_thread_division);