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.
39 * \brief This file defines functions necessary for mdrun and tools to
40 * compute energies and forces for bonded interactions.
42 * \author Mark Abraham <mark.j.abraham@gmail.com>
44 * \ingroup module_listed-forces
58 #include "gromacs/legacyheaders/disre.h"
59 #include "gromacs/legacyheaders/force.h"
60 #include "gromacs/legacyheaders/macros.h"
61 #include "gromacs/legacyheaders/names.h"
62 #include "gromacs/legacyheaders/nonbonded.h"
63 #include "gromacs/legacyheaders/ns.h"
64 #include "gromacs/legacyheaders/orires.h"
65 #include "gromacs/legacyheaders/txtdump.h"
66 #include "gromacs/math/units.h"
67 #include "gromacs/math/utilities.h"
68 #include "gromacs/math/vec.h"
69 #include "gromacs/pbcutil/ishift.h"
70 #include "gromacs/pbcutil/mshift.h"
71 #include "gromacs/pbcutil/pbc.h"
72 #include "gromacs/simd/simd.h"
73 #include "gromacs/simd/simd_math.h"
74 #include "gromacs/simd/vector_operations.h"
75 #include "gromacs/utility/fatalerror.h"
76 #include "gromacs/utility/smalloc.h"
80 /*! \brief Mysterious CMAP coefficient matrix */
81 const int cmap_coeff_matrix[] = {
82 1, 0, -3, 2, 0, 0, 0, 0, -3, 0, 9, -6, 2, 0, -6, 4,
83 0, 0, 0, 0, 0, 0, 0, 0, 3, 0, -9, 6, -2, 0, 6, -4,
84 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 9, -6, 0, 0, -6, 4,
85 0, 0, 3, -2, 0, 0, 0, 0, 0, 0, -9, 6, 0, 0, 6, -4,
86 0, 0, 0, 0, 1, 0, -3, 2, -2, 0, 6, -4, 1, 0, -3, 2,
87 0, 0, 0, 0, 0, 0, 0, 0, -1, 0, 3, -2, 1, 0, -3, 2,
88 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -3, 2, 0, 0, 3, -2,
89 0, 0, 0, 0, 0, 0, 3, -2, 0, 0, -6, 4, 0, 0, 3, -2,
90 0, 1, -2, 1, 0, 0, 0, 0, 0, -3, 6, -3, 0, 2, -4, 2,
91 0, 0, 0, 0, 0, 0, 0, 0, 0, 3, -6, 3, 0, -2, 4, -2,
92 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -3, 3, 0, 0, 2, -2,
93 0, 0, -1, 1, 0, 0, 0, 0, 0, 0, 3, -3, 0, 0, -2, 2,
94 0, 0, 0, 0, 0, 1, -2, 1, 0, -2, 4, -2, 0, 1, -2, 1,
95 0, 0, 0, 0, 0, 0, 0, 0, 0, -1, 2, -1, 0, 1, -2, 1,
96 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, -1, 0, 0, -1, 1,
97 0, 0, 0, 0, 0, 0, -1, 1, 0, 0, 2, -2, 0, 0, -1, 1
101 /* TODO This function should go and live in nonbonded.c where it is
102 really needed. Here, it only supports giving a fatal error message
104 int glatnr(int *global_atom_index, int i)
108 if (global_atom_index == NULL)
114 atnr = global_atom_index[i] + 1;
120 /*! \brief Compute dx = xi - xj, modulo PBC if non-NULL
122 * \todo This kind of code appears in many places. Consolidate it */
123 static int pbc_rvec_sub(const t_pbc *pbc, const rvec xi, const rvec xj, rvec dx)
127 return pbc_dx_aiuc(pbc, xi, xj, dx);
131 rvec_sub(xi, xj, dx);
136 #ifdef GMX_SIMD_HAVE_REAL
138 /* SIMD PBC data structure, containing 1/boxdiag and the box vectors */
140 gmx_simd_real_t inv_bzz;
141 gmx_simd_real_t inv_byy;
142 gmx_simd_real_t inv_bxx;
151 /*! \brief Set the SIMD pbc data from a normal t_pbc struct */
152 static void set_pbc_simd(const t_pbc *pbc, pbc_simd_t *pbc_simd)
157 /* Setting inv_bdiag to 0 effectively turns off PBC */
158 clear_rvec(inv_bdiag);
161 for (d = 0; d < pbc->ndim_ePBC; d++)
163 inv_bdiag[d] = 1.0/pbc->box[d][d];
167 pbc_simd->inv_bzz = gmx_simd_set1_r(inv_bdiag[ZZ]);
168 pbc_simd->inv_byy = gmx_simd_set1_r(inv_bdiag[YY]);
169 pbc_simd->inv_bxx = gmx_simd_set1_r(inv_bdiag[XX]);
173 pbc_simd->bzx = gmx_simd_set1_r(pbc->box[ZZ][XX]);
174 pbc_simd->bzy = gmx_simd_set1_r(pbc->box[ZZ][YY]);
175 pbc_simd->bzz = gmx_simd_set1_r(pbc->box[ZZ][ZZ]);
176 pbc_simd->byx = gmx_simd_set1_r(pbc->box[YY][XX]);
177 pbc_simd->byy = gmx_simd_set1_r(pbc->box[YY][YY]);
178 pbc_simd->bxx = gmx_simd_set1_r(pbc->box[XX][XX]);
182 pbc_simd->bzx = gmx_simd_setzero_r();
183 pbc_simd->bzy = gmx_simd_setzero_r();
184 pbc_simd->bzz = gmx_simd_setzero_r();
185 pbc_simd->byx = gmx_simd_setzero_r();
186 pbc_simd->byy = gmx_simd_setzero_r();
187 pbc_simd->bxx = gmx_simd_setzero_r();
191 /*! \brief Correct distance vector *dx,*dy,*dz for PBC using SIMD */
192 static gmx_inline void
193 pbc_dx_simd(gmx_simd_real_t *dx, gmx_simd_real_t *dy, gmx_simd_real_t *dz,
194 const pbc_simd_t *pbc)
198 sh = gmx_simd_round_r(gmx_simd_mul_r(*dz, pbc->inv_bzz));
199 *dx = gmx_simd_fnmadd_r(sh, pbc->bzx, *dx);
200 *dy = gmx_simd_fnmadd_r(sh, pbc->bzy, *dy);
201 *dz = gmx_simd_fnmadd_r(sh, pbc->bzz, *dz);
203 sh = gmx_simd_round_r(gmx_simd_mul_r(*dy, pbc->inv_byy));
204 *dx = gmx_simd_fnmadd_r(sh, pbc->byx, *dx);
205 *dy = gmx_simd_fnmadd_r(sh, pbc->byy, *dy);
207 sh = gmx_simd_round_r(gmx_simd_mul_r(*dx, pbc->inv_bxx));
208 *dx = gmx_simd_fnmadd_r(sh, pbc->bxx, *dx);
211 #endif /* GMX_SIMD_HAVE_REAL */
213 /*! \brief Morse potential bond
215 * By Frank Everdij. Three parameters needed:
217 * b0 = equilibrium distance in nm
218 * be = beta in nm^-1 (actually, it's nu_e*Sqrt(2*pi*pi*mu/D_e))
219 * cb = well depth in kJ/mol
221 * Note: the potential is referenced to be +cb at infinite separation
222 * and zero at the equilibrium distance!
224 real morse_bonds(int nbonds,
225 const t_iatom forceatoms[], const t_iparams forceparams[],
226 const rvec x[], rvec f[], rvec fshift[],
227 const t_pbc *pbc, const t_graph *g,
228 real lambda, real *dvdlambda,
229 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
230 int gmx_unused *global_atom_index)
232 const real one = 1.0;
233 const real two = 2.0;
234 real dr, dr2, temp, omtemp, cbomtemp, fbond, vbond, fij, vtot;
235 real b0, be, cb, b0A, beA, cbA, b0B, beB, cbB, L1;
237 int i, m, ki, type, ai, aj;
241 for (i = 0; (i < nbonds); )
243 type = forceatoms[i++];
244 ai = forceatoms[i++];
245 aj = forceatoms[i++];
247 b0A = forceparams[type].morse.b0A;
248 beA = forceparams[type].morse.betaA;
249 cbA = forceparams[type].morse.cbA;
251 b0B = forceparams[type].morse.b0B;
252 beB = forceparams[type].morse.betaB;
253 cbB = forceparams[type].morse.cbB;
255 L1 = one-lambda; /* 1 */
256 b0 = L1*b0A + lambda*b0B; /* 3 */
257 be = L1*beA + lambda*beB; /* 3 */
258 cb = L1*cbA + lambda*cbB; /* 3 */
260 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
261 dr2 = iprod(dx, dx); /* 5 */
262 dr = dr2*gmx_invsqrt(dr2); /* 10 */
263 temp = exp(-be*(dr-b0)); /* 12 */
267 /* bonds are constrainted. This may _not_ include bond constraints if they are lambda dependent */
268 *dvdlambda += cbB-cbA;
272 omtemp = one-temp; /* 1 */
273 cbomtemp = cb*omtemp; /* 1 */
274 vbond = cbomtemp*omtemp; /* 1 */
275 fbond = -two*be*temp*cbomtemp*gmx_invsqrt(dr2); /* 9 */
276 vtot += vbond; /* 1 */
278 *dvdlambda += (cbB - cbA) * omtemp * omtemp - (2-2*omtemp)*omtemp * cb * ((b0B-b0A)*be - (beB-beA)*(dr-b0)); /* 15 */
282 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
286 for (m = 0; (m < DIM); m++) /* 15 */
291 fshift[ki][m] += fij;
292 fshift[CENTRAL][m] -= fij;
299 real cubic_bonds(int nbonds,
300 const t_iatom forceatoms[], const t_iparams forceparams[],
301 const rvec x[], rvec f[], rvec fshift[],
302 const t_pbc *pbc, const t_graph *g,
303 real gmx_unused lambda, real gmx_unused *dvdlambda,
304 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
305 int gmx_unused *global_atom_index)
307 const real three = 3.0;
308 const real two = 2.0;
310 real dr, dr2, dist, kdist, kdist2, fbond, vbond, fij, vtot;
312 int i, m, ki, type, ai, aj;
316 for (i = 0; (i < nbonds); )
318 type = forceatoms[i++];
319 ai = forceatoms[i++];
320 aj = forceatoms[i++];
322 b0 = forceparams[type].cubic.b0;
323 kb = forceparams[type].cubic.kb;
324 kcub = forceparams[type].cubic.kcub;
326 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
327 dr2 = iprod(dx, dx); /* 5 */
334 dr = dr2*gmx_invsqrt(dr2); /* 10 */
339 vbond = kdist2 + kcub*kdist2*dist;
340 fbond = -(two*kdist + three*kdist2*kcub)/dr;
342 vtot += vbond; /* 21 */
346 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
349 for (m = 0; (m < DIM); m++) /* 15 */
354 fshift[ki][m] += fij;
355 fshift[CENTRAL][m] -= fij;
361 real FENE_bonds(int nbonds,
362 const t_iatom forceatoms[], const t_iparams forceparams[],
363 const rvec x[], rvec f[], rvec fshift[],
364 const t_pbc *pbc, const t_graph *g,
365 real gmx_unused lambda, real gmx_unused *dvdlambda,
366 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
367 int *global_atom_index)
369 const real half = 0.5;
370 const real one = 1.0;
372 real dr2, bm2, omdr2obm2, fbond, vbond, fij, vtot;
374 int i, m, ki, type, ai, aj;
378 for (i = 0; (i < nbonds); )
380 type = forceatoms[i++];
381 ai = forceatoms[i++];
382 aj = forceatoms[i++];
384 bm = forceparams[type].fene.bm;
385 kb = forceparams[type].fene.kb;
387 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
388 dr2 = iprod(dx, dx); /* 5 */
400 "r^2 (%f) >= bm^2 (%f) in FENE bond between atoms %d and %d",
402 glatnr(global_atom_index, ai),
403 glatnr(global_atom_index, aj));
406 omdr2obm2 = one - dr2/bm2;
408 vbond = -half*kb*bm2*log(omdr2obm2);
409 fbond = -kb/omdr2obm2;
411 vtot += vbond; /* 35 */
415 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
418 for (m = 0; (m < DIM); m++) /* 15 */
423 fshift[ki][m] += fij;
424 fshift[CENTRAL][m] -= fij;
430 real harmonic(real kA, real kB, real xA, real xB, real x, real lambda,
433 const real half = 0.5;
434 real L1, kk, x0, dx, dx2;
435 real v, f, dvdlambda;
438 kk = L1*kA+lambda*kB;
439 x0 = L1*xA+lambda*xB;
446 dvdlambda = half*(kB-kA)*dx2 + (xA-xB)*kk*dx;
453 /* That was 19 flops */
457 real bonds(int nbonds,
458 const t_iatom forceatoms[], const t_iparams forceparams[],
459 const rvec x[], rvec f[], rvec fshift[],
460 const t_pbc *pbc, const t_graph *g,
461 real lambda, real *dvdlambda,
462 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
463 int gmx_unused *global_atom_index)
465 int i, m, ki, ai, aj, type;
466 real dr, dr2, fbond, vbond, fij, vtot;
471 for (i = 0; (i < nbonds); )
473 type = forceatoms[i++];
474 ai = forceatoms[i++];
475 aj = forceatoms[i++];
477 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
478 dr2 = iprod(dx, dx); /* 5 */
479 dr = dr2*gmx_invsqrt(dr2); /* 10 */
481 *dvdlambda += harmonic(forceparams[type].harmonic.krA,
482 forceparams[type].harmonic.krB,
483 forceparams[type].harmonic.rA,
484 forceparams[type].harmonic.rB,
485 dr, lambda, &vbond, &fbond); /* 19 */
493 vtot += vbond; /* 1*/
494 fbond *= gmx_invsqrt(dr2); /* 6 */
498 fprintf(debug, "BONDS: dr = %10g vbond = %10g fbond = %10g\n",
504 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
507 for (m = 0; (m < DIM); m++) /* 15 */
512 fshift[ki][m] += fij;
513 fshift[CENTRAL][m] -= fij;
519 real restraint_bonds(int nbonds,
520 const t_iatom forceatoms[], const t_iparams forceparams[],
521 const rvec x[], rvec f[], rvec fshift[],
522 const t_pbc *pbc, const t_graph *g,
523 real lambda, real *dvdlambda,
524 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
525 int gmx_unused *global_atom_index)
527 int i, m, ki, ai, aj, type;
528 real dr, dr2, fbond, vbond, fij, vtot;
530 real low, dlow, up1, dup1, up2, dup2, k, dk;
538 for (i = 0; (i < nbonds); )
540 type = forceatoms[i++];
541 ai = forceatoms[i++];
542 aj = forceatoms[i++];
544 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
545 dr2 = iprod(dx, dx); /* 5 */
546 dr = dr2*gmx_invsqrt(dr2); /* 10 */
548 low = L1*forceparams[type].restraint.lowA + lambda*forceparams[type].restraint.lowB;
549 dlow = -forceparams[type].restraint.lowA + forceparams[type].restraint.lowB;
550 up1 = L1*forceparams[type].restraint.up1A + lambda*forceparams[type].restraint.up1B;
551 dup1 = -forceparams[type].restraint.up1A + forceparams[type].restraint.up1B;
552 up2 = L1*forceparams[type].restraint.up2A + lambda*forceparams[type].restraint.up2B;
553 dup2 = -forceparams[type].restraint.up2A + forceparams[type].restraint.up2B;
554 k = L1*forceparams[type].restraint.kA + lambda*forceparams[type].restraint.kB;
555 dk = -forceparams[type].restraint.kA + forceparams[type].restraint.kB;
564 *dvdlambda += 0.5*dk*drh2 - k*dlow*drh;
577 *dvdlambda += 0.5*dk*drh2 - k*dup1*drh;
582 vbond = k*(up2 - up1)*(0.5*(up2 - up1) + drh);
583 fbond = -k*(up2 - up1);
584 *dvdlambda += dk*(up2 - up1)*(0.5*(up2 - up1) + drh)
585 + k*(dup2 - dup1)*(up2 - up1 + drh)
586 - k*(up2 - up1)*dup2;
594 vtot += vbond; /* 1*/
595 fbond *= gmx_invsqrt(dr2); /* 6 */
599 fprintf(debug, "BONDS: dr = %10g vbond = %10g fbond = %10g\n",
605 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
608 for (m = 0; (m < DIM); m++) /* 15 */
613 fshift[ki][m] += fij;
614 fshift[CENTRAL][m] -= fij;
621 real polarize(int nbonds,
622 const t_iatom forceatoms[], const t_iparams forceparams[],
623 const rvec x[], rvec f[], rvec fshift[],
624 const t_pbc *pbc, const t_graph *g,
625 real lambda, real *dvdlambda,
626 const t_mdatoms *md, t_fcdata gmx_unused *fcd,
627 int gmx_unused *global_atom_index)
629 int i, m, ki, ai, aj, type;
630 real dr, dr2, fbond, vbond, fij, vtot, ksh;
635 for (i = 0; (i < nbonds); )
637 type = forceatoms[i++];
638 ai = forceatoms[i++];
639 aj = forceatoms[i++];
640 ksh = sqr(md->chargeA[aj])*ONE_4PI_EPS0/forceparams[type].polarize.alpha;
643 fprintf(debug, "POL: local ai = %d aj = %d ksh = %.3f\n", ai, aj, ksh);
646 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
647 dr2 = iprod(dx, dx); /* 5 */
648 dr = dr2*gmx_invsqrt(dr2); /* 10 */
650 *dvdlambda += harmonic(ksh, ksh, 0, 0, dr, lambda, &vbond, &fbond); /* 19 */
657 vtot += vbond; /* 1*/
658 fbond *= gmx_invsqrt(dr2); /* 6 */
662 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
665 for (m = 0; (m < DIM); m++) /* 15 */
670 fshift[ki][m] += fij;
671 fshift[CENTRAL][m] -= fij;
677 real anharm_polarize(int nbonds,
678 const t_iatom forceatoms[], const t_iparams forceparams[],
679 const rvec x[], rvec f[], rvec fshift[],
680 const t_pbc *pbc, const t_graph *g,
681 real lambda, real *dvdlambda,
682 const t_mdatoms *md, t_fcdata gmx_unused *fcd,
683 int gmx_unused *global_atom_index)
685 int i, m, ki, ai, aj, type;
686 real dr, dr2, fbond, vbond, fij, vtot, ksh, khyp, drcut, ddr, ddr3;
691 for (i = 0; (i < nbonds); )
693 type = forceatoms[i++];
694 ai = forceatoms[i++];
695 aj = forceatoms[i++];
696 ksh = sqr(md->chargeA[aj])*ONE_4PI_EPS0/forceparams[type].anharm_polarize.alpha; /* 7*/
697 khyp = forceparams[type].anharm_polarize.khyp;
698 drcut = forceparams[type].anharm_polarize.drcut;
701 fprintf(debug, "POL: local ai = %d aj = %d ksh = %.3f\n", ai, aj, ksh);
704 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
705 dr2 = iprod(dx, dx); /* 5 */
706 dr = dr2*gmx_invsqrt(dr2); /* 10 */
708 *dvdlambda += harmonic(ksh, ksh, 0, 0, dr, lambda, &vbond, &fbond); /* 19 */
719 vbond += khyp*ddr*ddr3;
720 fbond -= 4*khyp*ddr3;
722 fbond *= gmx_invsqrt(dr2); /* 6 */
723 vtot += vbond; /* 1*/
727 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
730 for (m = 0; (m < DIM); m++) /* 15 */
735 fshift[ki][m] += fij;
736 fshift[CENTRAL][m] -= fij;
742 real water_pol(int nbonds,
743 const t_iatom forceatoms[], const t_iparams forceparams[],
744 const rvec x[], rvec f[], rvec gmx_unused fshift[],
745 const t_pbc gmx_unused *pbc, const t_graph gmx_unused *g,
746 real gmx_unused lambda, real gmx_unused *dvdlambda,
747 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
748 int gmx_unused *global_atom_index)
750 /* This routine implements anisotropic polarizibility for water, through
751 * a shell connected to a dummy with spring constant that differ in the
752 * three spatial dimensions in the molecular frame.
754 int i, m, aO, aH1, aH2, aD, aS, type, type0, ki;
756 rvec dOH1, dOH2, dHH, dOD, dDS, nW, kk, dx, kdx, proj;
760 real vtot, fij, r_HH, r_OD, r_nW, tx, ty, tz, qS;
765 type0 = forceatoms[0];
767 qS = md->chargeA[aS];
768 kk[XX] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_x;
769 kk[YY] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_y;
770 kk[ZZ] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_z;
771 r_HH = 1.0/forceparams[type0].wpol.rHH;
774 fprintf(debug, "WPOL: qS = %10.5f aS = %5d\n", qS, aS);
775 fprintf(debug, "WPOL: kk = %10.3f %10.3f %10.3f\n",
776 kk[XX], kk[YY], kk[ZZ]);
777 fprintf(debug, "WPOL: rOH = %10.3f rHH = %10.3f rOD = %10.3f\n",
778 forceparams[type0].wpol.rOH,
779 forceparams[type0].wpol.rHH,
780 forceparams[type0].wpol.rOD);
782 for (i = 0; (i < nbonds); i += 6)
784 type = forceatoms[i];
787 gmx_fatal(FARGS, "Sorry, type = %d, type0 = %d, file = %s, line = %d",
788 type, type0, __FILE__, __LINE__);
790 aO = forceatoms[i+1];
791 aH1 = forceatoms[i+2];
792 aH2 = forceatoms[i+3];
793 aD = forceatoms[i+4];
794 aS = forceatoms[i+5];
796 /* Compute vectors describing the water frame */
797 pbc_rvec_sub(pbc, x[aH1], x[aO], dOH1);
798 pbc_rvec_sub(pbc, x[aH2], x[aO], dOH2);
799 pbc_rvec_sub(pbc, x[aH2], x[aH1], dHH);
800 pbc_rvec_sub(pbc, x[aD], x[aO], dOD);
801 ki = pbc_rvec_sub(pbc, x[aS], x[aD], dDS);
802 cprod(dOH1, dOH2, nW);
804 /* Compute inverse length of normal vector
805 * (this one could be precomputed, but I'm too lazy now)
807 r_nW = gmx_invsqrt(iprod(nW, nW));
808 /* This is for precision, but does not make a big difference,
811 r_OD = gmx_invsqrt(iprod(dOD, dOD));
813 /* Normalize the vectors in the water frame */
815 svmul(r_HH, dHH, dHH);
816 svmul(r_OD, dOD, dOD);
818 /* Compute displacement of shell along components of the vector */
819 dx[ZZ] = iprod(dDS, dOD);
820 /* Compute projection on the XY plane: dDS - dx[ZZ]*dOD */
821 for (m = 0; (m < DIM); m++)
823 proj[m] = dDS[m]-dx[ZZ]*dOD[m];
826 /*dx[XX] = iprod(dDS,nW);
827 dx[YY] = iprod(dDS,dHH);*/
828 dx[XX] = iprod(proj, nW);
829 for (m = 0; (m < DIM); m++)
831 proj[m] -= dx[XX]*nW[m];
833 dx[YY] = iprod(proj, dHH);
838 fprintf(debug, "WPOL: dx2=%10g dy2=%10g dz2=%10g sum=%10g dDS^2=%10g\n",
839 sqr(dx[XX]), sqr(dx[YY]), sqr(dx[ZZ]), iprod(dx, dx), iprod(dDS, dDS));
840 fprintf(debug, "WPOL: dHH=(%10g,%10g,%10g)\n", dHH[XX], dHH[YY], dHH[ZZ]);
841 fprintf(debug, "WPOL: dOD=(%10g,%10g,%10g), 1/r_OD = %10g\n",
842 dOD[XX], dOD[YY], dOD[ZZ], 1/r_OD);
843 fprintf(debug, "WPOL: nW =(%10g,%10g,%10g), 1/r_nW = %10g\n",
844 nW[XX], nW[YY], nW[ZZ], 1/r_nW);
845 fprintf(debug, "WPOL: dx =%10g, dy =%10g, dz =%10g\n",
846 dx[XX], dx[YY], dx[ZZ]);
847 fprintf(debug, "WPOL: dDSx=%10g, dDSy=%10g, dDSz=%10g\n",
848 dDS[XX], dDS[YY], dDS[ZZ]);
851 /* Now compute the forces and energy */
852 kdx[XX] = kk[XX]*dx[XX];
853 kdx[YY] = kk[YY]*dx[YY];
854 kdx[ZZ] = kk[ZZ]*dx[ZZ];
855 vtot += iprod(dx, kdx);
859 ivec_sub(SHIFT_IVEC(g, aS), SHIFT_IVEC(g, aD), dt);
863 for (m = 0; (m < DIM); m++)
865 /* This is a tensor operation but written out for speed */
875 fshift[ki][m] += fij;
876 fshift[CENTRAL][m] -= fij;
881 fprintf(debug, "WPOL: vwpol=%g\n", 0.5*iprod(dx, kdx));
882 fprintf(debug, "WPOL: df = (%10g, %10g, %10g)\n", df[XX], df[YY], df[ZZ]);
890 static real do_1_thole(const rvec xi, const rvec xj, rvec fi, rvec fj,
891 const t_pbc *pbc, real qq,
892 rvec fshift[], real afac)
895 real r12sq, r12_1, r12bar, v0, v1, fscal, ebar, fff;
898 t = pbc_rvec_sub(pbc, xi, xj, r12); /* 3 */
900 r12sq = iprod(r12, r12); /* 5 */
901 r12_1 = gmx_invsqrt(r12sq); /* 5 */
902 r12bar = afac/r12_1; /* 5 */
903 v0 = qq*ONE_4PI_EPS0*r12_1; /* 2 */
904 ebar = exp(-r12bar); /* 5 */
905 v1 = (1-(1+0.5*r12bar)*ebar); /* 4 */
906 fscal = ((v0*r12_1)*v1 - v0*0.5*afac*ebar*(r12bar+1))*r12_1; /* 9 */
909 fprintf(debug, "THOLE: v0 = %.3f v1 = %.3f r12= % .3f r12bar = %.3f fscal = %.3f ebar = %.3f\n", v0, v1, 1/r12_1, r12bar, fscal, ebar);
912 for (m = 0; (m < DIM); m++)
918 fshift[CENTRAL][m] -= fff;
921 return v0*v1; /* 1 */
925 real thole_pol(int nbonds,
926 const t_iatom forceatoms[], const t_iparams forceparams[],
927 const rvec x[], rvec f[], rvec fshift[],
928 const t_pbc *pbc, const t_graph gmx_unused *g,
929 real gmx_unused lambda, real gmx_unused *dvdlambda,
930 const t_mdatoms *md, t_fcdata gmx_unused *fcd,
931 int gmx_unused *global_atom_index)
933 /* Interaction between two pairs of particles with opposite charge */
934 int i, type, a1, da1, a2, da2;
935 real q1, q2, qq, a, al1, al2, afac;
937 const real minusOneOnSix = -1.0/6.0;
939 for (i = 0; (i < nbonds); )
941 type = forceatoms[i++];
942 a1 = forceatoms[i++];
943 da1 = forceatoms[i++];
944 a2 = forceatoms[i++];
945 da2 = forceatoms[i++];
946 q1 = md->chargeA[da1];
947 q2 = md->chargeA[da2];
948 a = forceparams[type].thole.a;
949 al1 = forceparams[type].thole.alpha1;
950 al2 = forceparams[type].thole.alpha2;
952 afac = a*pow(al1*al2, minusOneOnSix);
953 V += do_1_thole(x[a1], x[a2], f[a1], f[a2], pbc, qq, fshift, afac);
954 V += do_1_thole(x[da1], x[a2], f[da1], f[a2], pbc, -qq, fshift, afac);
955 V += do_1_thole(x[a1], x[da2], f[a1], f[da2], pbc, -qq, fshift, afac);
956 V += do_1_thole(x[da1], x[da2], f[da1], f[da2], pbc, qq, fshift, afac);
962 real bond_angle(const rvec xi, const rvec xj, const rvec xk, const t_pbc *pbc,
963 rvec r_ij, rvec r_kj, real *costh,
965 /* Return value is the angle between the bonds i-j and j-k */
970 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
971 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
973 *costh = cos_angle(r_ij, r_kj); /* 25 */
974 th = acos(*costh); /* 10 */
979 real angles(int nbonds,
980 const t_iatom forceatoms[], const t_iparams forceparams[],
981 const rvec x[], rvec f[], rvec fshift[],
982 const t_pbc *pbc, const t_graph *g,
983 real lambda, real *dvdlambda,
984 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
985 int gmx_unused *global_atom_index)
987 int i, ai, aj, ak, t1, t2, type;
989 real cos_theta, cos_theta2, theta, dVdt, va, vtot;
990 ivec jt, dt_ij, dt_kj;
993 for (i = 0; i < nbonds; )
995 type = forceatoms[i++];
996 ai = forceatoms[i++];
997 aj = forceatoms[i++];
998 ak = forceatoms[i++];
1000 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
1001 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
1003 *dvdlambda += harmonic(forceparams[type].harmonic.krA,
1004 forceparams[type].harmonic.krB,
1005 forceparams[type].harmonic.rA*DEG2RAD,
1006 forceparams[type].harmonic.rB*DEG2RAD,
1007 theta, lambda, &va, &dVdt); /* 21 */
1010 cos_theta2 = sqr(cos_theta);
1017 real nrkj_1, nrij_1;
1020 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1021 sth = st*cos_theta; /* 1 */
1025 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1026 theta*RAD2DEG, va, dVdt);
1029 nrij2 = iprod(r_ij, r_ij); /* 5 */
1030 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1032 nrij_1 = gmx_invsqrt(nrij2); /* 10 */
1033 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1035 cik = st*nrij_1*nrkj_1; /* 2 */
1036 cii = sth*nrij_1*nrij_1; /* 2 */
1037 ckk = sth*nrkj_1*nrkj_1; /* 2 */
1039 for (m = 0; m < DIM; m++)
1041 f_i[m] = -(cik*r_kj[m] - cii*r_ij[m]);
1042 f_k[m] = -(cik*r_ij[m] - ckk*r_kj[m]);
1043 f_j[m] = -f_i[m] - f_k[m];
1050 copy_ivec(SHIFT_IVEC(g, aj), jt);
1052 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1053 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1054 t1 = IVEC2IS(dt_ij);
1055 t2 = IVEC2IS(dt_kj);
1057 rvec_inc(fshift[t1], f_i);
1058 rvec_inc(fshift[CENTRAL], f_j);
1059 rvec_inc(fshift[t2], f_k);
1066 #ifdef GMX_SIMD_HAVE_REAL
1068 /* As angles, but using SIMD to calculate many dihedrals at once.
1069 * This routines does not calculate energies and shift forces.
1071 static gmx_inline void
1072 angles_noener_simd(int nbonds,
1073 const t_iatom forceatoms[], const t_iparams forceparams[],
1074 const rvec x[], rvec f[],
1075 const t_pbc *pbc, const t_graph gmx_unused *g,
1076 real gmx_unused lambda,
1077 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1078 int gmx_unused *global_atom_index)
1082 int type, ai[GMX_SIMD_REAL_WIDTH], aj[GMX_SIMD_REAL_WIDTH];
1083 int ak[GMX_SIMD_REAL_WIDTH];
1084 real coeff_array[2*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *coeff;
1085 real dr_array[2*DIM*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *dr;
1086 real f_buf_array[6*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *f_buf;
1087 gmx_simd_real_t k_S, theta0_S;
1088 gmx_simd_real_t rijx_S, rijy_S, rijz_S;
1089 gmx_simd_real_t rkjx_S, rkjy_S, rkjz_S;
1090 gmx_simd_real_t one_S;
1091 gmx_simd_real_t min_one_plus_eps_S;
1092 gmx_simd_real_t rij_rkj_S;
1093 gmx_simd_real_t nrij2_S, nrij_1_S;
1094 gmx_simd_real_t nrkj2_S, nrkj_1_S;
1095 gmx_simd_real_t cos_S, invsin_S;
1096 gmx_simd_real_t theta_S;
1097 gmx_simd_real_t st_S, sth_S;
1098 gmx_simd_real_t cik_S, cii_S, ckk_S;
1099 gmx_simd_real_t f_ix_S, f_iy_S, f_iz_S;
1100 gmx_simd_real_t f_kx_S, f_ky_S, f_kz_S;
1101 pbc_simd_t pbc_simd;
1103 /* Ensure register memory alignment */
1104 coeff = gmx_simd_align_r(coeff_array);
1105 dr = gmx_simd_align_r(dr_array);
1106 f_buf = gmx_simd_align_r(f_buf_array);
1108 set_pbc_simd(pbc, &pbc_simd);
1110 one_S = gmx_simd_set1_r(1.0);
1112 /* The smallest number > -1 */
1113 min_one_plus_eps_S = gmx_simd_set1_r(-1.0 + 2*GMX_REAL_EPS);
1115 /* nbonds is the number of angles times nfa1, here we step GMX_SIMD_REAL_WIDTH angles */
1116 for (i = 0; (i < nbonds); i += GMX_SIMD_REAL_WIDTH*nfa1)
1118 /* Collect atoms for GMX_SIMD_REAL_WIDTH angles.
1119 * iu indexes into forceatoms, we should not let iu go beyond nbonds.
1122 for (s = 0; s < GMX_SIMD_REAL_WIDTH; s++)
1124 type = forceatoms[iu];
1125 ai[s] = forceatoms[iu+1];
1126 aj[s] = forceatoms[iu+2];
1127 ak[s] = forceatoms[iu+3];
1129 coeff[s] = forceparams[type].harmonic.krA;
1130 coeff[GMX_SIMD_REAL_WIDTH+s] = forceparams[type].harmonic.rA*DEG2RAD;
1132 /* If you can't use pbc_dx_simd below for PBC, e.g. because
1133 * you can't round in SIMD, use pbc_rvec_sub here.
1135 /* Store the non PBC corrected distances packed and aligned */
1136 for (m = 0; m < DIM; m++)
1138 dr[s + m *GMX_SIMD_REAL_WIDTH] = x[ai[s]][m] - x[aj[s]][m];
1139 dr[s + (DIM+m)*GMX_SIMD_REAL_WIDTH] = x[ak[s]][m] - x[aj[s]][m];
1142 /* At the end fill the arrays with identical entries */
1143 if (iu + nfa1 < nbonds)
1149 k_S = gmx_simd_load_r(coeff);
1150 theta0_S = gmx_simd_load_r(coeff+GMX_SIMD_REAL_WIDTH);
1152 rijx_S = gmx_simd_load_r(dr + 0*GMX_SIMD_REAL_WIDTH);
1153 rijy_S = gmx_simd_load_r(dr + 1*GMX_SIMD_REAL_WIDTH);
1154 rijz_S = gmx_simd_load_r(dr + 2*GMX_SIMD_REAL_WIDTH);
1155 rkjx_S = gmx_simd_load_r(dr + 3*GMX_SIMD_REAL_WIDTH);
1156 rkjy_S = gmx_simd_load_r(dr + 4*GMX_SIMD_REAL_WIDTH);
1157 rkjz_S = gmx_simd_load_r(dr + 5*GMX_SIMD_REAL_WIDTH);
1159 pbc_dx_simd(&rijx_S, &rijy_S, &rijz_S, &pbc_simd);
1160 pbc_dx_simd(&rkjx_S, &rkjy_S, &rkjz_S, &pbc_simd);
1162 rij_rkj_S = gmx_simd_iprod_r(rijx_S, rijy_S, rijz_S,
1163 rkjx_S, rkjy_S, rkjz_S);
1165 nrij2_S = gmx_simd_norm2_r(rijx_S, rijy_S, rijz_S);
1166 nrkj2_S = gmx_simd_norm2_r(rkjx_S, rkjy_S, rkjz_S);
1168 nrij_1_S = gmx_simd_invsqrt_r(nrij2_S);
1169 nrkj_1_S = gmx_simd_invsqrt_r(nrkj2_S);
1171 cos_S = gmx_simd_mul_r(rij_rkj_S, gmx_simd_mul_r(nrij_1_S, nrkj_1_S));
1173 /* To allow for 180 degrees, we take the max of cos and -1 + 1bit,
1174 * so we can safely get the 1/sin from 1/sqrt(1 - cos^2).
1175 * This also ensures that rounding errors would cause the argument
1176 * of gmx_simd_acos_r to be < -1.
1177 * Note that we do not take precautions for cos(0)=1, so the outer
1178 * atoms in an angle should not be on top of each other.
1180 cos_S = gmx_simd_max_r(cos_S, min_one_plus_eps_S);
1182 theta_S = gmx_simd_acos_r(cos_S);
1184 invsin_S = gmx_simd_invsqrt_r(gmx_simd_sub_r(one_S, gmx_simd_mul_r(cos_S, cos_S)));
1186 st_S = gmx_simd_mul_r(gmx_simd_mul_r(k_S, gmx_simd_sub_r(theta0_S, theta_S)),
1188 sth_S = gmx_simd_mul_r(st_S, cos_S);
1190 cik_S = gmx_simd_mul_r(st_S, gmx_simd_mul_r(nrij_1_S, nrkj_1_S));
1191 cii_S = gmx_simd_mul_r(sth_S, gmx_simd_mul_r(nrij_1_S, nrij_1_S));
1192 ckk_S = gmx_simd_mul_r(sth_S, gmx_simd_mul_r(nrkj_1_S, nrkj_1_S));
1194 f_ix_S = gmx_simd_mul_r(cii_S, rijx_S);
1195 f_ix_S = gmx_simd_fnmadd_r(cik_S, rkjx_S, f_ix_S);
1196 f_iy_S = gmx_simd_mul_r(cii_S, rijy_S);
1197 f_iy_S = gmx_simd_fnmadd_r(cik_S, rkjy_S, f_iy_S);
1198 f_iz_S = gmx_simd_mul_r(cii_S, rijz_S);
1199 f_iz_S = gmx_simd_fnmadd_r(cik_S, rkjz_S, f_iz_S);
1200 f_kx_S = gmx_simd_mul_r(ckk_S, rkjx_S);
1201 f_kx_S = gmx_simd_fnmadd_r(cik_S, rijx_S, f_kx_S);
1202 f_ky_S = gmx_simd_mul_r(ckk_S, rkjy_S);
1203 f_ky_S = gmx_simd_fnmadd_r(cik_S, rijy_S, f_ky_S);
1204 f_kz_S = gmx_simd_mul_r(ckk_S, rkjz_S);
1205 f_kz_S = gmx_simd_fnmadd_r(cik_S, rijz_S, f_kz_S);
1207 gmx_simd_store_r(f_buf + 0*GMX_SIMD_REAL_WIDTH, f_ix_S);
1208 gmx_simd_store_r(f_buf + 1*GMX_SIMD_REAL_WIDTH, f_iy_S);
1209 gmx_simd_store_r(f_buf + 2*GMX_SIMD_REAL_WIDTH, f_iz_S);
1210 gmx_simd_store_r(f_buf + 3*GMX_SIMD_REAL_WIDTH, f_kx_S);
1211 gmx_simd_store_r(f_buf + 4*GMX_SIMD_REAL_WIDTH, f_ky_S);
1212 gmx_simd_store_r(f_buf + 5*GMX_SIMD_REAL_WIDTH, f_kz_S);
1218 for (m = 0; m < DIM; m++)
1220 f[ai[s]][m] += f_buf[s + m*GMX_SIMD_REAL_WIDTH];
1221 f[aj[s]][m] -= f_buf[s + m*GMX_SIMD_REAL_WIDTH] + f_buf[s + (DIM+m)*GMX_SIMD_REAL_WIDTH];
1222 f[ak[s]][m] += f_buf[s + (DIM+m)*GMX_SIMD_REAL_WIDTH];
1227 while (s < GMX_SIMD_REAL_WIDTH && iu < nbonds);
1231 #endif /* GMX_SIMD_HAVE_REAL */
1233 real linear_angles(int nbonds,
1234 const t_iatom forceatoms[], const t_iparams forceparams[],
1235 const rvec x[], rvec f[], rvec fshift[],
1236 const t_pbc *pbc, const t_graph *g,
1237 real lambda, real *dvdlambda,
1238 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1239 int gmx_unused *global_atom_index)
1241 int i, m, ai, aj, ak, t1, t2, type;
1243 real L1, kA, kB, aA, aB, dr, dr2, va, vtot, a, b, klin;
1244 ivec jt, dt_ij, dt_kj;
1245 rvec r_ij, r_kj, r_ik, dx;
1249 for (i = 0; (i < nbonds); )
1251 type = forceatoms[i++];
1252 ai = forceatoms[i++];
1253 aj = forceatoms[i++];
1254 ak = forceatoms[i++];
1256 kA = forceparams[type].linangle.klinA;
1257 kB = forceparams[type].linangle.klinB;
1258 klin = L1*kA + lambda*kB;
1260 aA = forceparams[type].linangle.aA;
1261 aB = forceparams[type].linangle.aB;
1262 a = L1*aA+lambda*aB;
1265 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
1266 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
1267 rvec_sub(r_ij, r_kj, r_ik);
1270 for (m = 0; (m < DIM); m++)
1272 dr = -a * r_ij[m] - b * r_kj[m];
1277 f_j[m] = -(f_i[m]+f_k[m]);
1283 *dvdlambda += 0.5*(kB-kA)*dr2 + klin*(aB-aA)*iprod(dx, r_ik);
1289 copy_ivec(SHIFT_IVEC(g, aj), jt);
1291 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1292 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1293 t1 = IVEC2IS(dt_ij);
1294 t2 = IVEC2IS(dt_kj);
1296 rvec_inc(fshift[t1], f_i);
1297 rvec_inc(fshift[CENTRAL], f_j);
1298 rvec_inc(fshift[t2], f_k);
1303 real urey_bradley(int nbonds,
1304 const t_iatom forceatoms[], const t_iparams forceparams[],
1305 const rvec x[], rvec f[], rvec fshift[],
1306 const t_pbc *pbc, const t_graph *g,
1307 real lambda, real *dvdlambda,
1308 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1309 int gmx_unused *global_atom_index)
1311 int i, m, ai, aj, ak, t1, t2, type, ki;
1312 rvec r_ij, r_kj, r_ik;
1313 real cos_theta, cos_theta2, theta;
1314 real dVdt, va, vtot, dr, dr2, vbond, fbond, fik;
1315 real kthA, th0A, kUBA, r13A, kthB, th0B, kUBB, r13B;
1316 ivec jt, dt_ij, dt_kj, dt_ik;
1319 for (i = 0; (i < nbonds); )
1321 type = forceatoms[i++];
1322 ai = forceatoms[i++];
1323 aj = forceatoms[i++];
1324 ak = forceatoms[i++];
1325 th0A = forceparams[type].u_b.thetaA*DEG2RAD;
1326 kthA = forceparams[type].u_b.kthetaA;
1327 r13A = forceparams[type].u_b.r13A;
1328 kUBA = forceparams[type].u_b.kUBA;
1329 th0B = forceparams[type].u_b.thetaB*DEG2RAD;
1330 kthB = forceparams[type].u_b.kthetaB;
1331 r13B = forceparams[type].u_b.r13B;
1332 kUBB = forceparams[type].u_b.kUBB;
1334 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
1335 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
1337 *dvdlambda += harmonic(kthA, kthB, th0A, th0B, theta, lambda, &va, &dVdt); /* 21 */
1340 ki = pbc_rvec_sub(pbc, x[ai], x[ak], r_ik); /* 3 */
1341 dr2 = iprod(r_ik, r_ik); /* 5 */
1342 dr = dr2*gmx_invsqrt(dr2); /* 10 */
1344 *dvdlambda += harmonic(kUBA, kUBB, r13A, r13B, dr, lambda, &vbond, &fbond); /* 19 */
1346 cos_theta2 = sqr(cos_theta); /* 1 */
1354 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1355 sth = st*cos_theta; /* 1 */
1359 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1360 theta*RAD2DEG, va, dVdt);
1363 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1364 nrij2 = iprod(r_ij, r_ij);
1366 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
1367 cii = sth/nrij2; /* 10 */
1368 ckk = sth/nrkj2; /* 10 */
1370 for (m = 0; (m < DIM); m++) /* 39 */
1372 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
1373 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
1374 f_j[m] = -f_i[m]-f_k[m];
1381 copy_ivec(SHIFT_IVEC(g, aj), jt);
1383 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1384 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1385 t1 = IVEC2IS(dt_ij);
1386 t2 = IVEC2IS(dt_kj);
1388 rvec_inc(fshift[t1], f_i);
1389 rvec_inc(fshift[CENTRAL], f_j);
1390 rvec_inc(fshift[t2], f_k);
1392 /* Time for the bond calculations */
1398 vtot += vbond; /* 1*/
1399 fbond *= gmx_invsqrt(dr2); /* 6 */
1403 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ak), dt_ik);
1404 ki = IVEC2IS(dt_ik);
1406 for (m = 0; (m < DIM); m++) /* 15 */
1408 fik = fbond*r_ik[m];
1411 fshift[ki][m] += fik;
1412 fshift[CENTRAL][m] -= fik;
1418 real quartic_angles(int nbonds,
1419 const t_iatom forceatoms[], const t_iparams forceparams[],
1420 const rvec x[], rvec f[], rvec fshift[],
1421 const t_pbc *pbc, const t_graph *g,
1422 real gmx_unused lambda, real gmx_unused *dvdlambda,
1423 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1424 int gmx_unused *global_atom_index)
1426 int i, j, ai, aj, ak, t1, t2, type;
1428 real cos_theta, cos_theta2, theta, dt, dVdt, va, dtp, c, vtot;
1429 ivec jt, dt_ij, dt_kj;
1432 for (i = 0; (i < nbonds); )
1434 type = forceatoms[i++];
1435 ai = forceatoms[i++];
1436 aj = forceatoms[i++];
1437 ak = forceatoms[i++];
1439 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
1440 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
1442 dt = theta - forceparams[type].qangle.theta*DEG2RAD; /* 2 */
1445 va = forceparams[type].qangle.c[0];
1447 for (j = 1; j <= 4; j++)
1449 c = forceparams[type].qangle.c[j];
1458 cos_theta2 = sqr(cos_theta); /* 1 */
1467 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1468 sth = st*cos_theta; /* 1 */
1472 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1473 theta*RAD2DEG, va, dVdt);
1476 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1477 nrij2 = iprod(r_ij, r_ij);
1479 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
1480 cii = sth/nrij2; /* 10 */
1481 ckk = sth/nrkj2; /* 10 */
1483 for (m = 0; (m < DIM); m++) /* 39 */
1485 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
1486 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
1487 f_j[m] = -f_i[m]-f_k[m];
1494 copy_ivec(SHIFT_IVEC(g, aj), jt);
1496 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1497 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1498 t1 = IVEC2IS(dt_ij);
1499 t2 = IVEC2IS(dt_kj);
1501 rvec_inc(fshift[t1], f_i);
1502 rvec_inc(fshift[CENTRAL], f_j);
1503 rvec_inc(fshift[t2], f_k);
1509 real dih_angle(const rvec xi, const rvec xj, const rvec xk, const rvec xl,
1511 rvec r_ij, rvec r_kj, rvec r_kl, rvec m, rvec n,
1512 real *sign, int *t1, int *t2, int *t3)
1516 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
1517 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
1518 *t3 = pbc_rvec_sub(pbc, xk, xl, r_kl); /* 3 */
1520 cprod(r_ij, r_kj, m); /* 9 */
1521 cprod(r_kj, r_kl, n); /* 9 */
1522 phi = gmx_angle(m, n); /* 49 (assuming 25 for atan2) */
1523 ipr = iprod(r_ij, n); /* 5 */
1524 (*sign) = (ipr < 0.0) ? -1.0 : 1.0;
1525 phi = (*sign)*phi; /* 1 */
1531 #ifdef GMX_SIMD_HAVE_REAL
1533 /* As dih_angle above, but calculates 4 dihedral angles at once using SIMD,
1534 * also calculates the pre-factor required for the dihedral force update.
1535 * Note that bv and buf should be register aligned.
1537 static gmx_inline void
1538 dih_angle_simd(const rvec *x,
1539 const int *ai, const int *aj, const int *ak, const int *al,
1540 const pbc_simd_t *pbc,
1542 gmx_simd_real_t *phi_S,
1543 gmx_simd_real_t *mx_S, gmx_simd_real_t *my_S, gmx_simd_real_t *mz_S,
1544 gmx_simd_real_t *nx_S, gmx_simd_real_t *ny_S, gmx_simd_real_t *nz_S,
1545 gmx_simd_real_t *nrkj_m2_S,
1546 gmx_simd_real_t *nrkj_n2_S,
1551 gmx_simd_real_t rijx_S, rijy_S, rijz_S;
1552 gmx_simd_real_t rkjx_S, rkjy_S, rkjz_S;
1553 gmx_simd_real_t rklx_S, rkly_S, rklz_S;
1554 gmx_simd_real_t cx_S, cy_S, cz_S;
1555 gmx_simd_real_t cn_S;
1556 gmx_simd_real_t s_S;
1557 gmx_simd_real_t ipr_S;
1558 gmx_simd_real_t iprm_S, iprn_S;
1559 gmx_simd_real_t nrkj2_S, nrkj_1_S, nrkj_2_S, nrkj_S;
1560 gmx_simd_real_t toler_S;
1561 gmx_simd_real_t p_S, q_S;
1562 gmx_simd_real_t nrkj2_min_S;
1563 gmx_simd_real_t real_eps_S;
1565 /* Used to avoid division by zero.
1566 * We take into acount that we multiply the result by real_eps_S.
1568 nrkj2_min_S = gmx_simd_set1_r(GMX_REAL_MIN/(2*GMX_REAL_EPS));
1570 /* The value of the last significant bit (GMX_REAL_EPS is half of that) */
1571 real_eps_S = gmx_simd_set1_r(2*GMX_REAL_EPS);
1573 for (s = 0; s < GMX_SIMD_REAL_WIDTH; s++)
1575 /* If you can't use pbc_dx_simd below for PBC, e.g. because
1576 * you can't round in SIMD, use pbc_rvec_sub here.
1578 for (m = 0; m < DIM; m++)
1580 dr[s + (0*DIM + m)*GMX_SIMD_REAL_WIDTH] = x[ai[s]][m] - x[aj[s]][m];
1581 dr[s + (1*DIM + m)*GMX_SIMD_REAL_WIDTH] = x[ak[s]][m] - x[aj[s]][m];
1582 dr[s + (2*DIM + m)*GMX_SIMD_REAL_WIDTH] = x[ak[s]][m] - x[al[s]][m];
1586 rijx_S = gmx_simd_load_r(dr + 0*GMX_SIMD_REAL_WIDTH);
1587 rijy_S = gmx_simd_load_r(dr + 1*GMX_SIMD_REAL_WIDTH);
1588 rijz_S = gmx_simd_load_r(dr + 2*GMX_SIMD_REAL_WIDTH);
1589 rkjx_S = gmx_simd_load_r(dr + 3*GMX_SIMD_REAL_WIDTH);
1590 rkjy_S = gmx_simd_load_r(dr + 4*GMX_SIMD_REAL_WIDTH);
1591 rkjz_S = gmx_simd_load_r(dr + 5*GMX_SIMD_REAL_WIDTH);
1592 rklx_S = gmx_simd_load_r(dr + 6*GMX_SIMD_REAL_WIDTH);
1593 rkly_S = gmx_simd_load_r(dr + 7*GMX_SIMD_REAL_WIDTH);
1594 rklz_S = gmx_simd_load_r(dr + 8*GMX_SIMD_REAL_WIDTH);
1596 pbc_dx_simd(&rijx_S, &rijy_S, &rijz_S, pbc);
1597 pbc_dx_simd(&rkjx_S, &rkjy_S, &rkjz_S, pbc);
1598 pbc_dx_simd(&rklx_S, &rkly_S, &rklz_S, pbc);
1600 gmx_simd_cprod_r(rijx_S, rijy_S, rijz_S,
1601 rkjx_S, rkjy_S, rkjz_S,
1604 gmx_simd_cprod_r(rkjx_S, rkjy_S, rkjz_S,
1605 rklx_S, rkly_S, rklz_S,
1608 gmx_simd_cprod_r(*mx_S, *my_S, *mz_S,
1609 *nx_S, *ny_S, *nz_S,
1610 &cx_S, &cy_S, &cz_S);
1612 cn_S = gmx_simd_sqrt_r(gmx_simd_norm2_r(cx_S, cy_S, cz_S));
1614 s_S = gmx_simd_iprod_r(*mx_S, *my_S, *mz_S, *nx_S, *ny_S, *nz_S);
1616 /* Determine the dihedral angle, the sign might need correction */
1617 *phi_S = gmx_simd_atan2_r(cn_S, s_S);
1619 ipr_S = gmx_simd_iprod_r(rijx_S, rijy_S, rijz_S,
1620 *nx_S, *ny_S, *nz_S);
1622 iprm_S = gmx_simd_norm2_r(*mx_S, *my_S, *mz_S);
1623 iprn_S = gmx_simd_norm2_r(*nx_S, *ny_S, *nz_S);
1625 nrkj2_S = gmx_simd_norm2_r(rkjx_S, rkjy_S, rkjz_S);
1627 /* Avoid division by zero. When zero, the result is multiplied by 0
1628 * anyhow, so the 3 max below do not affect the final result.
1630 nrkj2_S = gmx_simd_max_r(nrkj2_S, nrkj2_min_S);
1631 nrkj_1_S = gmx_simd_invsqrt_r(nrkj2_S);
1632 nrkj_2_S = gmx_simd_mul_r(nrkj_1_S, nrkj_1_S);
1633 nrkj_S = gmx_simd_mul_r(nrkj2_S, nrkj_1_S);
1635 toler_S = gmx_simd_mul_r(nrkj2_S, real_eps_S);
1637 /* Here the plain-C code uses a conditional, but we can't do that in SIMD.
1638 * So we take a max with the tolerance instead. Since we multiply with
1639 * m or n later, the max does not affect the results.
1641 iprm_S = gmx_simd_max_r(iprm_S, toler_S);
1642 iprn_S = gmx_simd_max_r(iprn_S, toler_S);
1643 *nrkj_m2_S = gmx_simd_mul_r(nrkj_S, gmx_simd_inv_r(iprm_S));
1644 *nrkj_n2_S = gmx_simd_mul_r(nrkj_S, gmx_simd_inv_r(iprn_S));
1646 /* Set sign of phi_S with the sign of ipr_S; phi_S is currently positive */
1647 *phi_S = gmx_simd_xor_sign_r(*phi_S, ipr_S);
1648 p_S = gmx_simd_iprod_r(rijx_S, rijy_S, rijz_S,
1649 rkjx_S, rkjy_S, rkjz_S);
1650 p_S = gmx_simd_mul_r(p_S, nrkj_2_S);
1652 q_S = gmx_simd_iprod_r(rklx_S, rkly_S, rklz_S,
1653 rkjx_S, rkjy_S, rkjz_S);
1654 q_S = gmx_simd_mul_r(q_S, nrkj_2_S);
1656 gmx_simd_store_r(p, p_S);
1657 gmx_simd_store_r(q, q_S);
1660 #endif /* GMX_SIMD_HAVE_REAL */
1663 void do_dih_fup(int i, int j, int k, int l, real ddphi,
1664 rvec r_ij, rvec r_kj, rvec r_kl,
1665 rvec m, rvec n, rvec f[], rvec fshift[],
1666 const t_pbc *pbc, const t_graph *g,
1667 const rvec x[], int t1, int t2, int t3)
1670 rvec f_i, f_j, f_k, f_l;
1671 rvec uvec, vvec, svec, dx_jl;
1672 real iprm, iprn, nrkj, nrkj2, nrkj_1, nrkj_2;
1673 real a, b, p, q, toler;
1674 ivec jt, dt_ij, dt_kj, dt_lj;
1676 iprm = iprod(m, m); /* 5 */
1677 iprn = iprod(n, n); /* 5 */
1678 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1679 toler = nrkj2*GMX_REAL_EPS;
1680 if ((iprm > toler) && (iprn > toler))
1682 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1683 nrkj_2 = nrkj_1*nrkj_1; /* 1 */
1684 nrkj = nrkj2*nrkj_1; /* 1 */
1685 a = -ddphi*nrkj/iprm; /* 11 */
1686 svmul(a, m, f_i); /* 3 */
1687 b = ddphi*nrkj/iprn; /* 11 */
1688 svmul(b, n, f_l); /* 3 */
1689 p = iprod(r_ij, r_kj); /* 5 */
1690 p *= nrkj_2; /* 1 */
1691 q = iprod(r_kl, r_kj); /* 5 */
1692 q *= nrkj_2; /* 1 */
1693 svmul(p, f_i, uvec); /* 3 */
1694 svmul(q, f_l, vvec); /* 3 */
1695 rvec_sub(uvec, vvec, svec); /* 3 */
1696 rvec_sub(f_i, svec, f_j); /* 3 */
1697 rvec_add(f_l, svec, f_k); /* 3 */
1698 rvec_inc(f[i], f_i); /* 3 */
1699 rvec_dec(f[j], f_j); /* 3 */
1700 rvec_dec(f[k], f_k); /* 3 */
1701 rvec_inc(f[l], f_l); /* 3 */
1705 copy_ivec(SHIFT_IVEC(g, j), jt);
1706 ivec_sub(SHIFT_IVEC(g, i), jt, dt_ij);
1707 ivec_sub(SHIFT_IVEC(g, k), jt, dt_kj);
1708 ivec_sub(SHIFT_IVEC(g, l), jt, dt_lj);
1709 t1 = IVEC2IS(dt_ij);
1710 t2 = IVEC2IS(dt_kj);
1711 t3 = IVEC2IS(dt_lj);
1715 t3 = pbc_rvec_sub(pbc, x[l], x[j], dx_jl);
1722 rvec_inc(fshift[t1], f_i);
1723 rvec_dec(fshift[CENTRAL], f_j);
1724 rvec_dec(fshift[t2], f_k);
1725 rvec_inc(fshift[t3], f_l);
1730 /* As do_dih_fup above, but without shift forces */
1732 do_dih_fup_noshiftf(int i, int j, int k, int l, real ddphi,
1733 rvec r_ij, rvec r_kj, rvec r_kl,
1734 rvec m, rvec n, rvec f[])
1736 rvec f_i, f_j, f_k, f_l;
1737 rvec uvec, vvec, svec;
1738 real iprm, iprn, nrkj, nrkj2, nrkj_1, nrkj_2;
1739 real a, b, p, q, toler;
1741 iprm = iprod(m, m); /* 5 */
1742 iprn = iprod(n, n); /* 5 */
1743 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1744 toler = nrkj2*GMX_REAL_EPS;
1745 if ((iprm > toler) && (iprn > toler))
1747 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1748 nrkj_2 = nrkj_1*nrkj_1; /* 1 */
1749 nrkj = nrkj2*nrkj_1; /* 1 */
1750 a = -ddphi*nrkj/iprm; /* 11 */
1751 svmul(a, m, f_i); /* 3 */
1752 b = ddphi*nrkj/iprn; /* 11 */
1753 svmul(b, n, f_l); /* 3 */
1754 p = iprod(r_ij, r_kj); /* 5 */
1755 p *= nrkj_2; /* 1 */
1756 q = iprod(r_kl, r_kj); /* 5 */
1757 q *= nrkj_2; /* 1 */
1758 svmul(p, f_i, uvec); /* 3 */
1759 svmul(q, f_l, vvec); /* 3 */
1760 rvec_sub(uvec, vvec, svec); /* 3 */
1761 rvec_sub(f_i, svec, f_j); /* 3 */
1762 rvec_add(f_l, svec, f_k); /* 3 */
1763 rvec_inc(f[i], f_i); /* 3 */
1764 rvec_dec(f[j], f_j); /* 3 */
1765 rvec_dec(f[k], f_k); /* 3 */
1766 rvec_inc(f[l], f_l); /* 3 */
1770 /* As do_dih_fup_noshiftf above, but with pre-calculated pre-factors */
1771 static gmx_inline void
1772 do_dih_fup_noshiftf_precalc(int i, int j, int k, int l,
1774 real f_i_x, real f_i_y, real f_i_z,
1775 real mf_l_x, real mf_l_y, real mf_l_z,
1778 rvec f_i, f_j, f_k, f_l;
1779 rvec uvec, vvec, svec;
1787 svmul(p, f_i, uvec);
1788 svmul(q, f_l, vvec);
1789 rvec_sub(uvec, vvec, svec);
1790 rvec_sub(f_i, svec, f_j);
1791 rvec_add(f_l, svec, f_k);
1792 rvec_inc(f[i], f_i);
1793 rvec_dec(f[j], f_j);
1794 rvec_dec(f[k], f_k);
1795 rvec_inc(f[l], f_l);
1799 real dopdihs(real cpA, real cpB, real phiA, real phiB, int mult,
1800 real phi, real lambda, real *V, real *F)
1802 real v, dvdlambda, mdphi, v1, sdphi, ddphi;
1803 real L1 = 1.0 - lambda;
1804 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1805 real dph0 = (phiB - phiA)*DEG2RAD;
1806 real cp = L1*cpA + lambda*cpB;
1808 mdphi = mult*phi - ph0;
1810 ddphi = -cp*mult*sdphi;
1811 v1 = 1.0 + cos(mdphi);
1814 dvdlambda = (cpB - cpA)*v1 + cp*dph0*sdphi;
1821 /* That was 40 flops */
1825 dopdihs_noener(real cpA, real cpB, real phiA, real phiB, int mult,
1826 real phi, real lambda, real *F)
1828 real mdphi, sdphi, ddphi;
1829 real L1 = 1.0 - lambda;
1830 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1831 real cp = L1*cpA + lambda*cpB;
1833 mdphi = mult*phi - ph0;
1835 ddphi = -cp*mult*sdphi;
1839 /* That was 20 flops */
1843 dopdihs_mdphi(real cpA, real cpB, real phiA, real phiB, int mult,
1844 real phi, real lambda, real *cp, real *mdphi)
1846 real L1 = 1.0 - lambda;
1847 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1849 *cp = L1*cpA + lambda*cpB;
1851 *mdphi = mult*phi - ph0;
1854 static real dopdihs_min(real cpA, real cpB, real phiA, real phiB, int mult,
1855 real phi, real lambda, real *V, real *F)
1856 /* similar to dopdihs, except for a minus sign *
1857 * and a different treatment of mult/phi0 */
1859 real v, dvdlambda, mdphi, v1, sdphi, ddphi;
1860 real L1 = 1.0 - lambda;
1861 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1862 real dph0 = (phiB - phiA)*DEG2RAD;
1863 real cp = L1*cpA + lambda*cpB;
1865 mdphi = mult*(phi-ph0);
1867 ddphi = cp*mult*sdphi;
1868 v1 = 1.0-cos(mdphi);
1871 dvdlambda = (cpB-cpA)*v1 + cp*dph0*sdphi;
1878 /* That was 40 flops */
1881 real pdihs(int nbonds,
1882 const t_iatom forceatoms[], const t_iparams forceparams[],
1883 const rvec x[], rvec f[], rvec fshift[],
1884 const t_pbc *pbc, const t_graph *g,
1885 real lambda, real *dvdlambda,
1886 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1887 int gmx_unused *global_atom_index)
1889 int i, type, ai, aj, ak, al;
1891 rvec r_ij, r_kj, r_kl, m, n;
1892 real phi, sign, ddphi, vpd, vtot;
1896 for (i = 0; (i < nbonds); )
1898 type = forceatoms[i++];
1899 ai = forceatoms[i++];
1900 aj = forceatoms[i++];
1901 ak = forceatoms[i++];
1902 al = forceatoms[i++];
1904 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
1905 &sign, &t1, &t2, &t3); /* 84 */
1906 *dvdlambda += dopdihs(forceparams[type].pdihs.cpA,
1907 forceparams[type].pdihs.cpB,
1908 forceparams[type].pdihs.phiA,
1909 forceparams[type].pdihs.phiB,
1910 forceparams[type].pdihs.mult,
1911 phi, lambda, &vpd, &ddphi);
1914 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
1915 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
1918 fprintf(debug, "pdih: (%d,%d,%d,%d) phi=%g\n",
1919 ai, aj, ak, al, phi);
1926 void make_dp_periodic(real *dp) /* 1 flop? */
1928 /* dp cannot be outside (-pi,pi) */
1933 else if (*dp < -M_PI)
1940 /* As pdihs above, but without calculating energies and shift forces */
1942 pdihs_noener(int nbonds,
1943 const t_iatom forceatoms[], const t_iparams forceparams[],
1944 const rvec x[], rvec f[],
1945 const t_pbc gmx_unused *pbc, const t_graph gmx_unused *g,
1947 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1948 int gmx_unused *global_atom_index)
1950 int i, type, ai, aj, ak, al;
1952 rvec r_ij, r_kj, r_kl, m, n;
1953 real phi, sign, ddphi_tot, ddphi;
1955 for (i = 0; (i < nbonds); )
1957 ai = forceatoms[i+1];
1958 aj = forceatoms[i+2];
1959 ak = forceatoms[i+3];
1960 al = forceatoms[i+4];
1962 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
1963 &sign, &t1, &t2, &t3);
1967 /* Loop over dihedrals working on the same atoms,
1968 * so we avoid recalculating angles and force distributions.
1972 type = forceatoms[i];
1973 dopdihs_noener(forceparams[type].pdihs.cpA,
1974 forceparams[type].pdihs.cpB,
1975 forceparams[type].pdihs.phiA,
1976 forceparams[type].pdihs.phiB,
1977 forceparams[type].pdihs.mult,
1978 phi, lambda, &ddphi);
1983 while (i < nbonds &&
1984 forceatoms[i+1] == ai &&
1985 forceatoms[i+2] == aj &&
1986 forceatoms[i+3] == ak &&
1987 forceatoms[i+4] == al);
1989 do_dih_fup_noshiftf(ai, aj, ak, al, ddphi_tot, r_ij, r_kj, r_kl, m, n, f);
1994 #ifdef GMX_SIMD_HAVE_REAL
1996 /* As pdihs_noner above, but using SIMD to calculate many dihedrals at once */
1998 pdihs_noener_simd(int nbonds,
1999 const t_iatom forceatoms[], const t_iparams forceparams[],
2000 const rvec x[], rvec f[],
2001 const t_pbc *pbc, const t_graph gmx_unused *g,
2002 real gmx_unused lambda,
2003 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2004 int gmx_unused *global_atom_index)
2008 int type, ai[GMX_SIMD_REAL_WIDTH], aj[GMX_SIMD_REAL_WIDTH], ak[GMX_SIMD_REAL_WIDTH], al[GMX_SIMD_REAL_WIDTH];
2009 real dr_array[3*DIM*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *dr;
2010 real buf_array[7*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *buf;
2011 real *cp, *phi0, *mult, *p, *q;
2012 gmx_simd_real_t phi0_S, phi_S;
2013 gmx_simd_real_t mx_S, my_S, mz_S;
2014 gmx_simd_real_t nx_S, ny_S, nz_S;
2015 gmx_simd_real_t nrkj_m2_S, nrkj_n2_S;
2016 gmx_simd_real_t cp_S, mdphi_S, mult_S;
2017 gmx_simd_real_t sin_S, cos_S;
2018 gmx_simd_real_t mddphi_S;
2019 gmx_simd_real_t sf_i_S, msf_l_S;
2020 pbc_simd_t pbc_simd;
2022 /* Ensure SIMD register alignment */
2023 dr = gmx_simd_align_r(dr_array);
2024 buf = gmx_simd_align_r(buf_array);
2026 /* Extract aligned pointer for parameters and variables */
2027 cp = buf + 0*GMX_SIMD_REAL_WIDTH;
2028 phi0 = buf + 1*GMX_SIMD_REAL_WIDTH;
2029 mult = buf + 2*GMX_SIMD_REAL_WIDTH;
2030 p = buf + 3*GMX_SIMD_REAL_WIDTH;
2031 q = buf + 4*GMX_SIMD_REAL_WIDTH;
2033 set_pbc_simd(pbc, &pbc_simd);
2035 /* nbonds is the number of dihedrals times nfa1, here we step GMX_SIMD_REAL_WIDTH dihs */
2036 for (i = 0; (i < nbonds); i += GMX_SIMD_REAL_WIDTH*nfa1)
2038 /* Collect atoms quadruplets for GMX_SIMD_REAL_WIDTH dihedrals.
2039 * iu indexes into forceatoms, we should not let iu go beyond nbonds.
2042 for (s = 0; s < GMX_SIMD_REAL_WIDTH; s++)
2044 type = forceatoms[iu];
2045 ai[s] = forceatoms[iu+1];
2046 aj[s] = forceatoms[iu+2];
2047 ak[s] = forceatoms[iu+3];
2048 al[s] = forceatoms[iu+4];
2050 cp[s] = forceparams[type].pdihs.cpA;
2051 phi0[s] = forceparams[type].pdihs.phiA*DEG2RAD;
2052 mult[s] = forceparams[type].pdihs.mult;
2054 /* At the end fill the arrays with identical entries */
2055 if (iu + nfa1 < nbonds)
2061 /* Caclulate GMX_SIMD_REAL_WIDTH dihedral angles at once */
2062 dih_angle_simd(x, ai, aj, ak, al, &pbc_simd,
2065 &mx_S, &my_S, &mz_S,
2066 &nx_S, &ny_S, &nz_S,
2071 cp_S = gmx_simd_load_r(cp);
2072 phi0_S = gmx_simd_load_r(phi0);
2073 mult_S = gmx_simd_load_r(mult);
2075 mdphi_S = gmx_simd_sub_r(gmx_simd_mul_r(mult_S, phi_S), phi0_S);
2077 /* Calculate GMX_SIMD_REAL_WIDTH sines at once */
2078 gmx_simd_sincos_r(mdphi_S, &sin_S, &cos_S);
2079 mddphi_S = gmx_simd_mul_r(gmx_simd_mul_r(cp_S, mult_S), sin_S);
2080 sf_i_S = gmx_simd_mul_r(mddphi_S, nrkj_m2_S);
2081 msf_l_S = gmx_simd_mul_r(mddphi_S, nrkj_n2_S);
2083 /* After this m?_S will contain f[i] */
2084 mx_S = gmx_simd_mul_r(sf_i_S, mx_S);
2085 my_S = gmx_simd_mul_r(sf_i_S, my_S);
2086 mz_S = gmx_simd_mul_r(sf_i_S, mz_S);
2088 /* After this m?_S will contain -f[l] */
2089 nx_S = gmx_simd_mul_r(msf_l_S, nx_S);
2090 ny_S = gmx_simd_mul_r(msf_l_S, ny_S);
2091 nz_S = gmx_simd_mul_r(msf_l_S, nz_S);
2093 gmx_simd_store_r(dr + 0*GMX_SIMD_REAL_WIDTH, mx_S);
2094 gmx_simd_store_r(dr + 1*GMX_SIMD_REAL_WIDTH, my_S);
2095 gmx_simd_store_r(dr + 2*GMX_SIMD_REAL_WIDTH, mz_S);
2096 gmx_simd_store_r(dr + 3*GMX_SIMD_REAL_WIDTH, nx_S);
2097 gmx_simd_store_r(dr + 4*GMX_SIMD_REAL_WIDTH, ny_S);
2098 gmx_simd_store_r(dr + 5*GMX_SIMD_REAL_WIDTH, nz_S);
2104 do_dih_fup_noshiftf_precalc(ai[s], aj[s], ak[s], al[s],
2106 dr[ XX *GMX_SIMD_REAL_WIDTH+s],
2107 dr[ YY *GMX_SIMD_REAL_WIDTH+s],
2108 dr[ ZZ *GMX_SIMD_REAL_WIDTH+s],
2109 dr[(DIM+XX)*GMX_SIMD_REAL_WIDTH+s],
2110 dr[(DIM+YY)*GMX_SIMD_REAL_WIDTH+s],
2111 dr[(DIM+ZZ)*GMX_SIMD_REAL_WIDTH+s],
2116 while (s < GMX_SIMD_REAL_WIDTH && iu < nbonds);
2120 #endif /* GMX_SIMD_HAVE_REAL */
2123 real idihs(int nbonds,
2124 const t_iatom forceatoms[], const t_iparams forceparams[],
2125 const rvec x[], rvec f[], rvec fshift[],
2126 const t_pbc *pbc, const t_graph *g,
2127 real lambda, real *dvdlambda,
2128 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2129 int gmx_unused *global_atom_index)
2131 int i, type, ai, aj, ak, al;
2133 real phi, phi0, dphi0, ddphi, sign, vtot;
2134 rvec r_ij, r_kj, r_kl, m, n;
2135 real L1, kk, dp, dp2, kA, kB, pA, pB, dvdl_term;
2140 for (i = 0; (i < nbonds); )
2142 type = forceatoms[i++];
2143 ai = forceatoms[i++];
2144 aj = forceatoms[i++];
2145 ak = forceatoms[i++];
2146 al = forceatoms[i++];
2148 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2149 &sign, &t1, &t2, &t3); /* 84 */
2151 /* phi can jump if phi0 is close to Pi/-Pi, which will cause huge
2152 * force changes if we just apply a normal harmonic.
2153 * Instead, we first calculate phi-phi0 and take it modulo (-Pi,Pi).
2154 * This means we will never have the periodicity problem, unless
2155 * the dihedral is Pi away from phiO, which is very unlikely due to
2158 kA = forceparams[type].harmonic.krA;
2159 kB = forceparams[type].harmonic.krB;
2160 pA = forceparams[type].harmonic.rA;
2161 pB = forceparams[type].harmonic.rB;
2163 kk = L1*kA + lambda*kB;
2164 phi0 = (L1*pA + lambda*pB)*DEG2RAD;
2165 dphi0 = (pB - pA)*DEG2RAD;
2169 make_dp_periodic(&dp);
2176 dvdl_term += 0.5*(kB - kA)*dp2 - kk*dphi0*dp;
2178 do_dih_fup(ai, aj, ak, al, -ddphi, r_ij, r_kj, r_kl, m, n,
2179 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2184 fprintf(debug, "idih: (%d,%d,%d,%d) phi=%g\n",
2185 ai, aj, ak, al, phi);
2190 *dvdlambda += dvdl_term;
2195 /*! \brief returns dx, rdist, and dpdl for functions posres() and fbposres()
2197 static void posres_dx(const rvec x, const rvec pos0A, const rvec pos0B,
2198 const rvec comA_sc, const rvec comB_sc,
2200 t_pbc *pbc, int refcoord_scaling, int npbcdim,
2201 rvec dx, rvec rdist, rvec dpdl)
2204 real posA, posB, L1, ref = 0.;
2209 for (m = 0; m < DIM; m++)
2215 switch (refcoord_scaling)
2219 rdist[m] = L1*posA + lambda*posB;
2220 dpdl[m] = posB - posA;
2223 /* Box relative coordinates are stored for dimensions with pbc */
2224 posA *= pbc->box[m][m];
2225 posB *= pbc->box[m][m];
2226 assert(npbcdim <= DIM);
2227 for (d = m+1; d < npbcdim; d++)
2229 posA += pos0A[d]*pbc->box[d][m];
2230 posB += pos0B[d]*pbc->box[d][m];
2232 ref = L1*posA + lambda*posB;
2234 dpdl[m] = posB - posA;
2237 ref = L1*comA_sc[m] + lambda*comB_sc[m];
2238 rdist[m] = L1*posA + lambda*posB;
2239 dpdl[m] = comB_sc[m] - comA_sc[m] + posB - posA;
2242 gmx_fatal(FARGS, "No such scaling method implemented");
2247 ref = L1*posA + lambda*posB;
2249 dpdl[m] = posB - posA;
2252 /* We do pbc_dx with ref+rdist,
2253 * since with only ref we can be up to half a box vector wrong.
2255 pos[m] = ref + rdist[m];
2260 pbc_dx(pbc, x, pos, dx);
2264 rvec_sub(x, pos, dx);
2268 /*! \brief Computes forces and potential for flat-bottom cylindrical restraints.
2269 * Returns the flat-bottom potential. */
2270 static real do_fbposres_cylinder(int fbdim, rvec fm, rvec dx, real rfb, real kk, gmx_bool bInvert)
2273 real dr, dr2, invdr, v, rfb2;
2279 for (d = 0; d < DIM; d++)
2288 ( (dr2 > rfb2 && bInvert == FALSE ) || (dr2 < rfb2 && bInvert == TRUE ) )
2293 v = 0.5*kk*sqr(dr - rfb);
2294 for (d = 0; d < DIM; d++)
2298 fm[d] = -kk*(dr-rfb)*dx[d]*invdr; /* Force pointing to the center */
2306 /*! \brief Adds forces of flat-bottomed positions restraints to f[]
2307 * and fixes vir_diag. Returns the flat-bottomed potential. */
2308 real fbposres(int nbonds,
2309 const t_iatom forceatoms[], const t_iparams forceparams[],
2310 const rvec x[], rvec f[], rvec vir_diag,
2312 int refcoord_scaling, int ePBC, rvec com)
2313 /* compute flat-bottomed positions restraints */
2315 int i, ai, m, d, type, npbcdim = 0, fbdim;
2316 const t_iparams *pr;
2318 real dr, dr2, rfb, rfb2, fact;
2319 rvec com_sc, rdist, dx, dpdl, fm;
2322 npbcdim = ePBC2npbcdim(ePBC);
2324 if (refcoord_scaling == erscCOM)
2327 for (m = 0; m < npbcdim; m++)
2329 assert(npbcdim <= DIM);
2330 for (d = m; d < npbcdim; d++)
2332 com_sc[m] += com[d]*pbc->box[d][m];
2338 for (i = 0; (i < nbonds); )
2340 type = forceatoms[i++];
2341 ai = forceatoms[i++];
2342 pr = &forceparams[type];
2344 /* same calculation as for normal posres, but with identical A and B states, and lambda==0 */
2345 posres_dx(x[ai], forceparams[type].fbposres.pos0, forceparams[type].fbposres.pos0,
2346 com_sc, com_sc, 0.0,
2347 pbc, refcoord_scaling, npbcdim,
2353 kk = pr->fbposres.k;
2354 rfb = pr->fbposres.r;
2357 /* with rfb<0, push particle out of the sphere/cylinder/layer */
2365 switch (pr->fbposres.geom)
2367 case efbposresSPHERE:
2368 /* spherical flat-bottom posres */
2371 ( (dr2 > rfb2 && bInvert == FALSE ) || (dr2 < rfb2 && bInvert == TRUE ) )
2375 v = 0.5*kk*sqr(dr - rfb);
2376 fact = -kk*(dr-rfb)/dr; /* Force pointing to the center pos0 */
2377 svmul(fact, dx, fm);
2380 case efbposresCYLINDERX:
2381 /* cylindrical flat-bottom posres in y-z plane. fm[XX] = 0. */
2383 v = do_fbposres_cylinder(fbdim, fm, dx, rfb, kk, bInvert);
2385 case efbposresCYLINDERY:
2386 /* cylindrical flat-bottom posres in x-z plane. fm[YY] = 0. */
2388 v = do_fbposres_cylinder(fbdim, fm, dx, rfb, kk, bInvert);
2390 case efbposresCYLINDER:
2391 /* equivalent to efbposresCYLINDERZ for backwards compatibility */
2392 case efbposresCYLINDERZ:
2393 /* cylindrical flat-bottom posres in x-y plane. fm[ZZ] = 0. */
2395 v = do_fbposres_cylinder(fbdim, fm, dx, rfb, kk, bInvert);
2397 case efbposresX: /* fbdim=XX */
2398 case efbposresY: /* fbdim=YY */
2399 case efbposresZ: /* fbdim=ZZ */
2400 /* 1D flat-bottom potential */
2401 fbdim = pr->fbposres.geom - efbposresX;
2403 if ( ( dr > rfb && bInvert == FALSE ) || ( 0 < dr && dr < rfb && bInvert == TRUE ) )
2405 v = 0.5*kk*sqr(dr - rfb);
2406 fm[fbdim] = -kk*(dr - rfb);
2408 else if ( (dr < (-rfb) && bInvert == FALSE ) || ( (-rfb) < dr && dr < 0 && bInvert == TRUE ))
2410 v = 0.5*kk*sqr(dr + rfb);
2411 fm[fbdim] = -kk*(dr + rfb);
2418 for (m = 0; (m < DIM); m++)
2421 /* Here we correct for the pbc_dx which included rdist */
2422 vir_diag[m] -= 0.5*(dx[m] + rdist[m])*fm[m];
2430 real posres(int nbonds,
2431 const t_iatom forceatoms[], const t_iparams forceparams[],
2432 const rvec x[], rvec f[], rvec vir_diag,
2434 real lambda, real *dvdlambda,
2435 int refcoord_scaling, int ePBC, rvec comA, rvec comB)
2437 int i, ai, m, d, type, npbcdim = 0;
2438 const t_iparams *pr;
2441 rvec comA_sc, comB_sc, rdist, dpdl, dx;
2442 gmx_bool bForceValid = TRUE;
2444 if ((f == NULL) || (vir_diag == NULL)) /* should both be null together! */
2446 bForceValid = FALSE;
2449 npbcdim = ePBC2npbcdim(ePBC);
2451 if (refcoord_scaling == erscCOM)
2453 clear_rvec(comA_sc);
2454 clear_rvec(comB_sc);
2455 for (m = 0; m < npbcdim; m++)
2457 assert(npbcdim <= DIM);
2458 for (d = m; d < npbcdim; d++)
2460 comA_sc[m] += comA[d]*pbc->box[d][m];
2461 comB_sc[m] += comB[d]*pbc->box[d][m];
2469 for (i = 0; (i < nbonds); )
2471 type = forceatoms[i++];
2472 ai = forceatoms[i++];
2473 pr = &forceparams[type];
2475 /* return dx, rdist, and dpdl */
2476 posres_dx(x[ai], forceparams[type].posres.pos0A, forceparams[type].posres.pos0B,
2477 comA_sc, comB_sc, lambda,
2478 pbc, refcoord_scaling, npbcdim,
2481 for (m = 0; (m < DIM); m++)
2483 kk = L1*pr->posres.fcA[m] + lambda*pr->posres.fcB[m];
2485 vtot += 0.5*kk*dx[m]*dx[m];
2487 0.5*(pr->posres.fcB[m] - pr->posres.fcA[m])*dx[m]*dx[m]
2490 /* Here we correct for the pbc_dx which included rdist */
2494 vir_diag[m] -= 0.5*(dx[m] + rdist[m])*fm;
2502 static real low_angres(int nbonds,
2503 const t_iatom forceatoms[], const t_iparams forceparams[],
2504 const rvec x[], rvec f[], rvec fshift[],
2505 const t_pbc *pbc, const t_graph *g,
2506 real lambda, real *dvdlambda,
2509 int i, m, type, ai, aj, ak, al;
2511 real phi, cos_phi, cos_phi2, vid, vtot, dVdphi;
2512 rvec r_ij, r_kl, f_i, f_k = {0, 0, 0};
2513 real st, sth, nrij2, nrkl2, c, cij, ckl;
2516 t2 = 0; /* avoid warning with gcc-3.3. It is never used uninitialized */
2519 ak = al = 0; /* to avoid warnings */
2520 for (i = 0; i < nbonds; )
2522 type = forceatoms[i++];
2523 ai = forceatoms[i++];
2524 aj = forceatoms[i++];
2525 t1 = pbc_rvec_sub(pbc, x[aj], x[ai], r_ij); /* 3 */
2528 ak = forceatoms[i++];
2529 al = forceatoms[i++];
2530 t2 = pbc_rvec_sub(pbc, x[al], x[ak], r_kl); /* 3 */
2539 cos_phi = cos_angle(r_ij, r_kl); /* 25 */
2540 phi = acos(cos_phi); /* 10 */
2542 *dvdlambda += dopdihs_min(forceparams[type].pdihs.cpA,
2543 forceparams[type].pdihs.cpB,
2544 forceparams[type].pdihs.phiA,
2545 forceparams[type].pdihs.phiB,
2546 forceparams[type].pdihs.mult,
2547 phi, lambda, &vid, &dVdphi); /* 40 */
2551 cos_phi2 = sqr(cos_phi); /* 1 */
2554 st = -dVdphi*gmx_invsqrt(1 - cos_phi2); /* 12 */
2555 sth = st*cos_phi; /* 1 */
2556 nrij2 = iprod(r_ij, r_ij); /* 5 */
2557 nrkl2 = iprod(r_kl, r_kl); /* 5 */
2559 c = st*gmx_invsqrt(nrij2*nrkl2); /* 11 */
2560 cij = sth/nrij2; /* 10 */
2561 ckl = sth/nrkl2; /* 10 */
2563 for (m = 0; m < DIM; m++) /* 18+18 */
2565 f_i[m] = (c*r_kl[m]-cij*r_ij[m]);
2570 f_k[m] = (c*r_ij[m]-ckl*r_kl[m]);
2578 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
2581 rvec_inc(fshift[t1], f_i);
2582 rvec_dec(fshift[CENTRAL], f_i);
2587 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, al), dt);
2590 rvec_inc(fshift[t2], f_k);
2591 rvec_dec(fshift[CENTRAL], f_k);
2596 return vtot; /* 184 / 157 (bZAxis) total */
2599 real angres(int nbonds,
2600 const t_iatom forceatoms[], const t_iparams forceparams[],
2601 const rvec x[], rvec f[], rvec fshift[],
2602 const t_pbc *pbc, const t_graph *g,
2603 real lambda, real *dvdlambda,
2604 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2605 int gmx_unused *global_atom_index)
2607 return low_angres(nbonds, forceatoms, forceparams, x, f, fshift, pbc, g,
2608 lambda, dvdlambda, FALSE);
2611 real angresz(int nbonds,
2612 const t_iatom forceatoms[], const t_iparams forceparams[],
2613 const rvec x[], rvec f[], rvec fshift[],
2614 const t_pbc *pbc, const t_graph *g,
2615 real lambda, real *dvdlambda,
2616 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2617 int gmx_unused *global_atom_index)
2619 return low_angres(nbonds, forceatoms, forceparams, x, f, fshift, pbc, g,
2620 lambda, dvdlambda, TRUE);
2623 real dihres(int nbonds,
2624 const t_iatom forceatoms[], const t_iparams forceparams[],
2625 const rvec x[], rvec f[], rvec fshift[],
2626 const t_pbc *pbc, const t_graph *g,
2627 real lambda, real *dvdlambda,
2628 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2629 int gmx_unused *global_atom_index)
2632 int ai, aj, ak, al, i, k, type, t1, t2, t3;
2633 real phi0A, phi0B, dphiA, dphiB, kfacA, kfacB, phi0, dphi, kfac;
2634 real phi, ddphi, ddp, ddp2, dp, sign, d2r, L1;
2635 rvec r_ij, r_kj, r_kl, m, n;
2642 for (i = 0; (i < nbonds); )
2644 type = forceatoms[i++];
2645 ai = forceatoms[i++];
2646 aj = forceatoms[i++];
2647 ak = forceatoms[i++];
2648 al = forceatoms[i++];
2650 phi0A = forceparams[type].dihres.phiA*d2r;
2651 dphiA = forceparams[type].dihres.dphiA*d2r;
2652 kfacA = forceparams[type].dihres.kfacA;
2654 phi0B = forceparams[type].dihres.phiB*d2r;
2655 dphiB = forceparams[type].dihres.dphiB*d2r;
2656 kfacB = forceparams[type].dihres.kfacB;
2658 phi0 = L1*phi0A + lambda*phi0B;
2659 dphi = L1*dphiA + lambda*dphiB;
2660 kfac = L1*kfacA + lambda*kfacB;
2662 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2663 &sign, &t1, &t2, &t3);
2668 fprintf(debug, "dihres[%d]: %d %d %d %d : phi=%f, dphi=%f, kfac=%f\n",
2669 k++, ai, aj, ak, al, phi0, dphi, kfac);
2671 /* phi can jump if phi0 is close to Pi/-Pi, which will cause huge
2672 * force changes if we just apply a normal harmonic.
2673 * Instead, we first calculate phi-phi0 and take it modulo (-Pi,Pi).
2674 * This means we will never have the periodicity problem, unless
2675 * the dihedral is Pi away from phiO, which is very unlikely due to
2679 make_dp_periodic(&dp);
2685 else if (dp < -dphi)
2697 vtot += 0.5*kfac*ddp2;
2700 *dvdlambda += 0.5*(kfacB - kfacA)*ddp2;
2701 /* lambda dependence from changing restraint distances */
2704 *dvdlambda -= kfac*ddp*((dphiB - dphiA)+(phi0B - phi0A));
2708 *dvdlambda += kfac*ddp*((dphiB - dphiA)-(phi0B - phi0A));
2710 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
2711 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2718 real unimplemented(int gmx_unused nbonds,
2719 const t_iatom gmx_unused forceatoms[], const t_iparams gmx_unused forceparams[],
2720 const rvec gmx_unused x[], rvec gmx_unused f[], rvec gmx_unused fshift[],
2721 const t_pbc gmx_unused *pbc, const t_graph gmx_unused *g,
2722 real gmx_unused lambda, real gmx_unused *dvdlambda,
2723 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2724 int gmx_unused *global_atom_index)
2726 gmx_impl("*** you are using a not implemented function");
2728 return 0.0; /* To make the compiler happy */
2731 real restrangles(int nbonds,
2732 const t_iatom forceatoms[], const t_iparams forceparams[],
2733 const rvec x[], rvec f[], rvec fshift[],
2734 const t_pbc *pbc, const t_graph *g,
2735 real gmx_unused lambda, real gmx_unused *dvdlambda,
2736 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2737 int gmx_unused *global_atom_index)
2739 int i, d, ai, aj, ak, type, m;
2742 ivec jt, dt_ij, dt_kj;
2744 real prefactor, ratio_ante, ratio_post;
2745 rvec delta_ante, delta_post, vec_temp;
2748 for (i = 0; (i < nbonds); )
2750 type = forceatoms[i++];
2751 ai = forceatoms[i++];
2752 aj = forceatoms[i++];
2753 ak = forceatoms[i++];
2755 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], vec_temp);
2756 pbc_rvec_sub(pbc, x[aj], x[ai], delta_ante);
2757 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], delta_post);
2760 /* This function computes factors needed for restricted angle potential.
2761 * The restricted angle potential is used in coarse-grained simulations to avoid singularities
2762 * when three particles align and the dihedral angle and dihedral potential
2763 * cannot be calculated. This potential is calculated using the formula:
2764 real restrangles(int nbonds,
2765 const t_iatom forceatoms[],const t_iparams forceparams[],
2766 const rvec x[],rvec f[],rvec fshift[],
2767 const t_pbc *pbc,const t_graph *g,
2768 real gmx_unused lambda,real gmx_unused *dvdlambda,
2769 const t_mdatoms gmx_unused *md,t_fcdata gmx_unused *fcd,
2770 int gmx_unused *global_atom_index)
2772 int i, d, ai, aj, ak, type, m;
2776 ivec jt,dt_ij,dt_kj;
2778 real prefactor, ratio_ante, ratio_post;
2779 rvec delta_ante, delta_post, vec_temp;
2782 for(i=0; (i<nbonds); )
2784 type = forceatoms[i++];
2785 ai = forceatoms[i++];
2786 aj = forceatoms[i++];
2787 ak = forceatoms[i++];
2789 * \f[V_{\rm ReB}(\theta_i) = \frac{1}{2} k_{\theta} \frac{(\cos\theta_i - \cos\theta_0)^2}
2790 * {\sin^2\theta_i}\f] ({eq:ReB} and ref \cite{MonicaGoga2013} from the manual).
2791 * For more explanations see comments file "restcbt.h". */
2793 compute_factors_restangles(type, forceparams, delta_ante, delta_post,
2794 &prefactor, &ratio_ante, &ratio_post, &v);
2796 /* Forces are computed per component */
2797 for (d = 0; d < DIM; d++)
2799 f_i[d] = prefactor * (ratio_ante * delta_ante[d] - delta_post[d]);
2800 f_j[d] = prefactor * ((ratio_post + 1.0) * delta_post[d] - (ratio_ante + 1.0) * delta_ante[d]);
2801 f_k[d] = prefactor * (delta_ante[d] - ratio_post * delta_post[d]);
2804 /* Computation of potential energy */
2810 for (m = 0; (m < DIM); m++)
2819 copy_ivec(SHIFT_IVEC(g, aj), jt);
2820 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
2821 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
2822 t1 = IVEC2IS(dt_ij);
2823 t2 = IVEC2IS(dt_kj);
2826 rvec_inc(fshift[t1], f_i);
2827 rvec_inc(fshift[CENTRAL], f_j);
2828 rvec_inc(fshift[t2], f_k);
2834 real restrdihs(int nbonds,
2835 const t_iatom forceatoms[], const t_iparams forceparams[],
2836 const rvec x[], rvec f[], rvec fshift[],
2837 const t_pbc *pbc, const t_graph *g,
2838 real gmx_unused lambda, real gmx_unused *dvlambda,
2839 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2840 int gmx_unused *global_atom_index)
2842 int i, d, type, ai, aj, ak, al;
2843 rvec f_i, f_j, f_k, f_l;
2845 ivec jt, dt_ij, dt_kj, dt_lj;
2848 rvec delta_ante, delta_crnt, delta_post, vec_temp;
2849 real factor_phi_ai_ante, factor_phi_ai_crnt, factor_phi_ai_post;
2850 real factor_phi_aj_ante, factor_phi_aj_crnt, factor_phi_aj_post;
2851 real factor_phi_ak_ante, factor_phi_ak_crnt, factor_phi_ak_post;
2852 real factor_phi_al_ante, factor_phi_al_crnt, factor_phi_al_post;
2857 for (i = 0; (i < nbonds); )
2859 type = forceatoms[i++];
2860 ai = forceatoms[i++];
2861 aj = forceatoms[i++];
2862 ak = forceatoms[i++];
2863 al = forceatoms[i++];
2865 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], vec_temp);
2866 pbc_rvec_sub(pbc, x[aj], x[ai], delta_ante);
2867 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], delta_crnt);
2868 pbc_rvec_sub(pbc, x[ak], x[al], vec_temp);
2869 pbc_rvec_sub(pbc, x[al], x[ak], delta_post);
2871 /* This function computes factors needed for restricted angle potential.
2872 * The restricted angle potential is used in coarse-grained simulations to avoid singularities
2873 * when three particles align and the dihedral angle and dihedral potential cannot be calculated.
2874 * This potential is calculated using the formula:
2875 * \f[V_{\rm ReB}(\theta_i) = \frac{1}{2} k_{\theta}
2876 * \frac{(\cos\theta_i - \cos\theta_0)^2}{\sin^2\theta_i}\f]
2877 * ({eq:ReB} and ref \cite{MonicaGoga2013} from the manual).
2878 * For more explanations see comments file "restcbt.h" */
2880 compute_factors_restrdihs(type, forceparams,
2881 delta_ante, delta_crnt, delta_post,
2882 &factor_phi_ai_ante, &factor_phi_ai_crnt, &factor_phi_ai_post,
2883 &factor_phi_aj_ante, &factor_phi_aj_crnt, &factor_phi_aj_post,
2884 &factor_phi_ak_ante, &factor_phi_ak_crnt, &factor_phi_ak_post,
2885 &factor_phi_al_ante, &factor_phi_al_crnt, &factor_phi_al_post,
2886 &prefactor_phi, &v);
2889 /* Computation of forces per component */
2890 for (d = 0; d < DIM; d++)
2892 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]);
2893 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]);
2894 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]);
2895 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]);
2897 /* Computation of the energy */
2903 /* Updating the forces */
2905 rvec_inc(f[ai], f_i);
2906 rvec_inc(f[aj], f_j);
2907 rvec_inc(f[ak], f_k);
2908 rvec_inc(f[al], f_l);
2911 /* Updating the fshift forces for the pressure coupling */
2914 copy_ivec(SHIFT_IVEC(g, aj), jt);
2915 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
2916 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
2917 ivec_sub(SHIFT_IVEC(g, al), jt, dt_lj);
2918 t1 = IVEC2IS(dt_ij);
2919 t2 = IVEC2IS(dt_kj);
2920 t3 = IVEC2IS(dt_lj);
2924 t3 = pbc_rvec_sub(pbc, x[al], x[aj], dx_jl);
2931 rvec_inc(fshift[t1], f_i);
2932 rvec_inc(fshift[CENTRAL], f_j);
2933 rvec_inc(fshift[t2], f_k);
2934 rvec_inc(fshift[t3], f_l);
2942 real cbtdihs(int nbonds,
2943 const t_iatom forceatoms[], const t_iparams forceparams[],
2944 const rvec x[], rvec f[], rvec fshift[],
2945 const t_pbc *pbc, const t_graph *g,
2946 real gmx_unused lambda, real gmx_unused *dvdlambda,
2947 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2948 int gmx_unused *global_atom_index)
2950 int type, ai, aj, ak, al, i, d;
2954 rvec f_i, f_j, f_k, f_l;
2955 ivec jt, dt_ij, dt_kj, dt_lj;
2957 rvec delta_ante, delta_crnt, delta_post;
2958 rvec f_phi_ai, f_phi_aj, f_phi_ak, f_phi_al;
2959 rvec f_theta_ante_ai, f_theta_ante_aj, f_theta_ante_ak;
2960 rvec f_theta_post_aj, f_theta_post_ak, f_theta_post_al;
2966 for (i = 0; (i < nbonds); )
2968 type = forceatoms[i++];
2969 ai = forceatoms[i++];
2970 aj = forceatoms[i++];
2971 ak = forceatoms[i++];
2972 al = forceatoms[i++];
2975 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], vec_temp);
2976 pbc_rvec_sub(pbc, x[aj], x[ai], delta_ante);
2977 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], vec_temp);
2978 pbc_rvec_sub(pbc, x[ak], x[aj], delta_crnt);
2979 pbc_rvec_sub(pbc, x[ak], x[al], vec_temp);
2980 pbc_rvec_sub(pbc, x[al], x[ak], delta_post);
2982 /* \brief Compute factors for CBT potential
2983 * The combined bending-torsion potential goes to zero in a very smooth manner, eliminating the numerical
2984 * instabilities, when three coarse-grained particles align and the dihedral angle and standard
2985 * dihedral potentials cannot be calculated. The CBT potential is calculated using the formula:
2986 * \f[V_{\rm CBT}(\theta_{i-1}, \theta_i, \phi_i) = k_{\phi} \sin^3\theta_{i-1} \sin^3\theta_{i}
2987 * \sum_{n=0}^4 { a_n \cos^n\phi_i}\f] ({eq:CBT} and ref \cite{MonicaGoga2013} from the manual).
2988 * It contains in its expression not only the dihedral angle \f$\phi\f$
2989 * but also \f[\theta_{i-1}\f] (theta_ante bellow) and \f[\theta_{i}\f] (theta_post bellow)
2990 * --- the adjacent bending angles.
2991 * For more explanations see comments file "restcbt.h". */
2993 compute_factors_cbtdihs(type, forceparams, delta_ante, delta_crnt, delta_post,
2994 f_phi_ai, f_phi_aj, f_phi_ak, f_phi_al,
2995 f_theta_ante_ai, f_theta_ante_aj, f_theta_ante_ak,
2996 f_theta_post_aj, f_theta_post_ak, f_theta_post_al,
3000 /* Acumulate the resuts per beads */
3001 for (d = 0; d < DIM; d++)
3003 f_i[d] = f_phi_ai[d] + f_theta_ante_ai[d];
3004 f_j[d] = f_phi_aj[d] + f_theta_ante_aj[d] + f_theta_post_aj[d];
3005 f_k[d] = f_phi_ak[d] + f_theta_ante_ak[d] + f_theta_post_ak[d];
3006 f_l[d] = f_phi_al[d] + f_theta_post_al[d];
3009 /* Compute the potential energy */
3014 /* Updating the forces */
3015 rvec_inc(f[ai], f_i);
3016 rvec_inc(f[aj], f_j);
3017 rvec_inc(f[ak], f_k);
3018 rvec_inc(f[al], f_l);
3021 /* Updating the fshift forces for the pressure coupling */
3024 copy_ivec(SHIFT_IVEC(g, aj), jt);
3025 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3026 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3027 ivec_sub(SHIFT_IVEC(g, al), jt, dt_lj);
3028 t1 = IVEC2IS(dt_ij);
3029 t2 = IVEC2IS(dt_kj);
3030 t3 = IVEC2IS(dt_lj);
3034 t3 = pbc_rvec_sub(pbc, x[al], x[aj], dx_jl);
3041 rvec_inc(fshift[t1], f_i);
3042 rvec_inc(fshift[CENTRAL], f_j);
3043 rvec_inc(fshift[t2], f_k);
3044 rvec_inc(fshift[t3], f_l);
3050 real rbdihs(int nbonds,
3051 const t_iatom forceatoms[], const t_iparams forceparams[],
3052 const rvec x[], rvec f[], rvec fshift[],
3053 const t_pbc *pbc, const t_graph *g,
3054 real lambda, real *dvdlambda,
3055 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3056 int gmx_unused *global_atom_index)
3058 const real c0 = 0.0, c1 = 1.0, c2 = 2.0, c3 = 3.0, c4 = 4.0, c5 = 5.0;
3059 int type, ai, aj, ak, al, i, j;
3061 rvec r_ij, r_kj, r_kl, m, n;
3062 real parmA[NR_RBDIHS];
3063 real parmB[NR_RBDIHS];
3064 real parm[NR_RBDIHS];
3065 real cos_phi, phi, rbp, rbpBA;
3066 real v, sign, ddphi, sin_phi;
3068 real L1 = 1.0-lambda;
3072 for (i = 0; (i < nbonds); )
3074 type = forceatoms[i++];
3075 ai = forceatoms[i++];
3076 aj = forceatoms[i++];
3077 ak = forceatoms[i++];
3078 al = forceatoms[i++];
3080 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
3081 &sign, &t1, &t2, &t3); /* 84 */
3083 /* Change to polymer convention */
3090 phi -= M_PI; /* 1 */
3094 /* Beware of accuracy loss, cannot use 1-sqrt(cos^2) ! */
3097 for (j = 0; (j < NR_RBDIHS); j++)
3099 parmA[j] = forceparams[type].rbdihs.rbcA[j];
3100 parmB[j] = forceparams[type].rbdihs.rbcB[j];
3101 parm[j] = L1*parmA[j]+lambda*parmB[j];
3103 /* Calculate cosine powers */
3104 /* Calculate the energy */
3105 /* Calculate the derivative */
3108 dvdl_term += (parmB[0]-parmA[0]);
3113 rbpBA = parmB[1]-parmA[1];
3114 ddphi += rbp*cosfac;
3117 dvdl_term += cosfac*rbpBA;
3119 rbpBA = parmB[2]-parmA[2];
3120 ddphi += c2*rbp*cosfac;
3123 dvdl_term += cosfac*rbpBA;
3125 rbpBA = parmB[3]-parmA[3];
3126 ddphi += c3*rbp*cosfac;
3129 dvdl_term += cosfac*rbpBA;
3131 rbpBA = parmB[4]-parmA[4];
3132 ddphi += c4*rbp*cosfac;
3135 dvdl_term += cosfac*rbpBA;
3137 rbpBA = parmB[5]-parmA[5];
3138 ddphi += c5*rbp*cosfac;
3141 dvdl_term += cosfac*rbpBA;
3143 ddphi = -ddphi*sin_phi; /* 11 */
3145 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
3146 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
3149 *dvdlambda += dvdl_term;
3156 /*! \brief Mysterious undocumented function */
3158 cmap_setup_grid_index(int ip, int grid_spacing, int *ipm1, int *ipp1, int *ipp2)
3164 ip = ip + grid_spacing - 1;
3166 else if (ip > grid_spacing)
3168 ip = ip - grid_spacing - 1;
3177 im1 = grid_spacing - 1;
3179 else if (ip == grid_spacing-2)
3183 else if (ip == grid_spacing-1)
3197 /*! \brief Compute CMAP dihedral energies and forces */
3199 cmap_dihs(int nbonds,
3200 const t_iatom forceatoms[], const t_iparams forceparams[],
3201 const gmx_cmap_t *cmap_grid,
3202 const rvec x[], rvec f[], rvec fshift[],
3203 const t_pbc *pbc, const t_graph *g,
3204 real gmx_unused lambda, real gmx_unused *dvdlambda,
3205 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3206 int gmx_unused *global_atom_index)
3208 int i, j, k, n, idx;
3209 int ai, aj, ak, al, am;
3210 int a1i, a1j, a1k, a1l, a2i, a2j, a2k, a2l;
3212 int t11, t21, t31, t12, t22, t32;
3213 int iphi1, ip1m1, ip1p1, ip1p2;
3214 int iphi2, ip2m1, ip2p1, ip2p2;
3216 int pos1, pos2, pos3, pos4;
3218 real ty[4], ty1[4], ty2[4], ty12[4], tc[16], tx[16];
3219 real phi1, cos_phi1, sin_phi1, sign1, xphi1;
3220 real phi2, cos_phi2, sin_phi2, sign2, xphi2;
3221 real dx, xx, tt, tu, e, df1, df2, vtot;
3222 real ra21, rb21, rg21, rg1, rgr1, ra2r1, rb2r1, rabr1;
3223 real ra22, rb22, rg22, rg2, rgr2, ra2r2, rb2r2, rabr2;
3224 real fg1, hg1, fga1, hgb1, gaa1, gbb1;
3225 real fg2, hg2, fga2, hgb2, gaa2, gbb2;
3228 rvec r1_ij, r1_kj, r1_kl, m1, n1;
3229 rvec r2_ij, r2_kj, r2_kl, m2, n2;
3230 rvec f1_i, f1_j, f1_k, f1_l;
3231 rvec f2_i, f2_j, f2_k, f2_l;
3232 rvec a1, b1, a2, b2;
3233 rvec f1, g1, h1, f2, g2, h2;
3234 rvec dtf1, dtg1, dth1, dtf2, dtg2, dth2;
3235 ivec jt1, dt1_ij, dt1_kj, dt1_lj;
3236 ivec jt2, dt2_ij, dt2_kj, dt2_lj;
3240 int loop_index[4][4] = {
3247 /* Total CMAP energy */
3250 for (n = 0; n < nbonds; )
3252 /* Five atoms are involved in the two torsions */
3253 type = forceatoms[n++];
3254 ai = forceatoms[n++];
3255 aj = forceatoms[n++];
3256 ak = forceatoms[n++];
3257 al = forceatoms[n++];
3258 am = forceatoms[n++];
3260 /* Which CMAP type is this */
3261 cmapA = forceparams[type].cmap.cmapA;
3262 cmapd = cmap_grid->cmapdata[cmapA].cmap;
3270 phi1 = dih_angle(x[a1i], x[a1j], x[a1k], x[a1l], pbc, r1_ij, r1_kj, r1_kl, m1, n1,
3271 &sign1, &t11, &t21, &t31); /* 84 */
3273 cos_phi1 = cos(phi1);
3275 a1[0] = r1_ij[1]*r1_kj[2]-r1_ij[2]*r1_kj[1];
3276 a1[1] = r1_ij[2]*r1_kj[0]-r1_ij[0]*r1_kj[2];
3277 a1[2] = r1_ij[0]*r1_kj[1]-r1_ij[1]*r1_kj[0]; /* 9 */
3279 b1[0] = r1_kl[1]*r1_kj[2]-r1_kl[2]*r1_kj[1];
3280 b1[1] = r1_kl[2]*r1_kj[0]-r1_kl[0]*r1_kj[2];
3281 b1[2] = r1_kl[0]*r1_kj[1]-r1_kl[1]*r1_kj[0]; /* 9 */
3283 pbc_rvec_sub(pbc, x[a1l], x[a1k], h1);
3285 ra21 = iprod(a1, a1); /* 5 */
3286 rb21 = iprod(b1, b1); /* 5 */
3287 rg21 = iprod(r1_kj, r1_kj); /* 5 */
3293 rabr1 = sqrt(ra2r1*rb2r1);
3295 sin_phi1 = rg1 * rabr1 * iprod(a1, h1) * (-1);
3297 if (cos_phi1 < -0.5 || cos_phi1 > 0.5)
3299 phi1 = asin(sin_phi1);
3309 phi1 = -M_PI - phi1;
3315 phi1 = acos(cos_phi1);
3323 xphi1 = phi1 + M_PI; /* 1 */
3325 /* Second torsion */
3331 phi2 = dih_angle(x[a2i], x[a2j], x[a2k], x[a2l], pbc, r2_ij, r2_kj, r2_kl, m2, n2,
3332 &sign2, &t12, &t22, &t32); /* 84 */
3334 cos_phi2 = cos(phi2);
3336 a2[0] = r2_ij[1]*r2_kj[2]-r2_ij[2]*r2_kj[1];
3337 a2[1] = r2_ij[2]*r2_kj[0]-r2_ij[0]*r2_kj[2];
3338 a2[2] = r2_ij[0]*r2_kj[1]-r2_ij[1]*r2_kj[0]; /* 9 */
3340 b2[0] = r2_kl[1]*r2_kj[2]-r2_kl[2]*r2_kj[1];
3341 b2[1] = r2_kl[2]*r2_kj[0]-r2_kl[0]*r2_kj[2];
3342 b2[2] = r2_kl[0]*r2_kj[1]-r2_kl[1]*r2_kj[0]; /* 9 */
3344 pbc_rvec_sub(pbc, x[a2l], x[a2k], h2);
3346 ra22 = iprod(a2, a2); /* 5 */
3347 rb22 = iprod(b2, b2); /* 5 */
3348 rg22 = iprod(r2_kj, r2_kj); /* 5 */
3354 rabr2 = sqrt(ra2r2*rb2r2);
3356 sin_phi2 = rg2 * rabr2 * iprod(a2, h2) * (-1);
3358 if (cos_phi2 < -0.5 || cos_phi2 > 0.5)
3360 phi2 = asin(sin_phi2);
3370 phi2 = -M_PI - phi2;
3376 phi2 = acos(cos_phi2);
3384 xphi2 = phi2 + M_PI; /* 1 */
3386 /* Range mangling */
3389 xphi1 = xphi1 + 2*M_PI;
3391 else if (xphi1 >= 2*M_PI)
3393 xphi1 = xphi1 - 2*M_PI;
3398 xphi2 = xphi2 + 2*M_PI;
3400 else if (xphi2 >= 2*M_PI)
3402 xphi2 = xphi2 - 2*M_PI;
3405 /* Number of grid points */
3406 dx = 2*M_PI / cmap_grid->grid_spacing;
3408 /* Where on the grid are we */
3409 iphi1 = static_cast<int>(xphi1/dx);
3410 iphi2 = static_cast<int>(xphi2/dx);
3412 iphi1 = cmap_setup_grid_index(iphi1, cmap_grid->grid_spacing, &ip1m1, &ip1p1, &ip1p2);
3413 iphi2 = cmap_setup_grid_index(iphi2, cmap_grid->grid_spacing, &ip2m1, &ip2p1, &ip2p2);
3415 pos1 = iphi1*cmap_grid->grid_spacing+iphi2;
3416 pos2 = ip1p1*cmap_grid->grid_spacing+iphi2;
3417 pos3 = ip1p1*cmap_grid->grid_spacing+ip2p1;
3418 pos4 = iphi1*cmap_grid->grid_spacing+ip2p1;
3420 ty[0] = cmapd[pos1*4];
3421 ty[1] = cmapd[pos2*4];
3422 ty[2] = cmapd[pos3*4];
3423 ty[3] = cmapd[pos4*4];
3425 ty1[0] = cmapd[pos1*4+1];
3426 ty1[1] = cmapd[pos2*4+1];
3427 ty1[2] = cmapd[pos3*4+1];
3428 ty1[3] = cmapd[pos4*4+1];
3430 ty2[0] = cmapd[pos1*4+2];
3431 ty2[1] = cmapd[pos2*4+2];
3432 ty2[2] = cmapd[pos3*4+2];
3433 ty2[3] = cmapd[pos4*4+2];
3435 ty12[0] = cmapd[pos1*4+3];
3436 ty12[1] = cmapd[pos2*4+3];
3437 ty12[2] = cmapd[pos3*4+3];
3438 ty12[3] = cmapd[pos4*4+3];
3440 /* Switch to degrees */
3441 dx = 360.0 / cmap_grid->grid_spacing;
3442 xphi1 = xphi1 * RAD2DEG;
3443 xphi2 = xphi2 * RAD2DEG;
3445 for (i = 0; i < 4; i++) /* 16 */
3448 tx[i+4] = ty1[i]*dx;
3449 tx[i+8] = ty2[i]*dx;
3450 tx[i+12] = ty12[i]*dx*dx;
3454 for (i = 0; i < 4; i++) /* 1056 */
3456 for (j = 0; j < 4; j++)
3459 for (k = 0; k < 16; k++)
3461 xx = xx + cmap_coeff_matrix[k*16+idx]*tx[k];
3469 tt = (xphi1-iphi1*dx)/dx;
3470 tu = (xphi2-iphi2*dx)/dx;
3476 for (i = 3; i >= 0; i--)
3478 l1 = loop_index[i][3];
3479 l2 = loop_index[i][2];
3480 l3 = loop_index[i][1];
3482 e = tt * e + ((tc[i*4+3]*tu+tc[i*4+2])*tu + tc[i*4+1])*tu+tc[i*4];
3483 df1 = tu * df1 + (3.0*tc[l1]*tt+2.0*tc[l2])*tt+tc[l3];
3484 df2 = tt * df2 + (3.0*tc[i*4+3]*tu+2.0*tc[i*4+2])*tu+tc[i*4+1];
3494 /* Do forces - first torsion */
3495 fg1 = iprod(r1_ij, r1_kj);
3496 hg1 = iprod(r1_kl, r1_kj);
3497 fga1 = fg1*ra2r1*rgr1;
3498 hgb1 = hg1*rb2r1*rgr1;
3502 for (i = 0; i < DIM; i++)
3504 dtf1[i] = gaa1 * a1[i];
3505 dtg1[i] = fga1 * a1[i] - hgb1 * b1[i];
3506 dth1[i] = gbb1 * b1[i];
3508 f1[i] = df1 * dtf1[i];
3509 g1[i] = df1 * dtg1[i];
3510 h1[i] = df1 * dth1[i];
3513 f1_j[i] = -f1[i] - g1[i];
3514 f1_k[i] = h1[i] + g1[i];
3517 f[a1i][i] = f[a1i][i] + f1_i[i];
3518 f[a1j][i] = f[a1j][i] + f1_j[i]; /* - f1[i] - g1[i] */
3519 f[a1k][i] = f[a1k][i] + f1_k[i]; /* h1[i] + g1[i] */
3520 f[a1l][i] = f[a1l][i] + f1_l[i]; /* h1[i] */
3523 /* Do forces - second torsion */
3524 fg2 = iprod(r2_ij, r2_kj);
3525 hg2 = iprod(r2_kl, r2_kj);
3526 fga2 = fg2*ra2r2*rgr2;
3527 hgb2 = hg2*rb2r2*rgr2;
3531 for (i = 0; i < DIM; i++)
3533 dtf2[i] = gaa2 * a2[i];
3534 dtg2[i] = fga2 * a2[i] - hgb2 * b2[i];
3535 dth2[i] = gbb2 * b2[i];
3537 f2[i] = df2 * dtf2[i];
3538 g2[i] = df2 * dtg2[i];
3539 h2[i] = df2 * dth2[i];
3542 f2_j[i] = -f2[i] - g2[i];
3543 f2_k[i] = h2[i] + g2[i];
3546 f[a2i][i] = f[a2i][i] + f2_i[i]; /* f2[i] */
3547 f[a2j][i] = f[a2j][i] + f2_j[i]; /* - f2[i] - g2[i] */
3548 f[a2k][i] = f[a2k][i] + f2_k[i]; /* h2[i] + g2[i] */
3549 f[a2l][i] = f[a2l][i] + f2_l[i]; /* - h2[i] */
3555 copy_ivec(SHIFT_IVEC(g, a1j), jt1);
3556 ivec_sub(SHIFT_IVEC(g, a1i), jt1, dt1_ij);
3557 ivec_sub(SHIFT_IVEC(g, a1k), jt1, dt1_kj);
3558 ivec_sub(SHIFT_IVEC(g, a1l), jt1, dt1_lj);
3559 t11 = IVEC2IS(dt1_ij);
3560 t21 = IVEC2IS(dt1_kj);
3561 t31 = IVEC2IS(dt1_lj);
3563 copy_ivec(SHIFT_IVEC(g, a2j), jt2);
3564 ivec_sub(SHIFT_IVEC(g, a2i), jt2, dt2_ij);
3565 ivec_sub(SHIFT_IVEC(g, a2k), jt2, dt2_kj);
3566 ivec_sub(SHIFT_IVEC(g, a2l), jt2, dt2_lj);
3567 t12 = IVEC2IS(dt2_ij);
3568 t22 = IVEC2IS(dt2_kj);
3569 t32 = IVEC2IS(dt2_lj);
3573 t31 = pbc_rvec_sub(pbc, x[a1l], x[a1j], h1);
3574 t32 = pbc_rvec_sub(pbc, x[a2l], x[a2j], h2);
3582 rvec_inc(fshift[t11], f1_i);
3583 rvec_inc(fshift[CENTRAL], f1_j);
3584 rvec_inc(fshift[t21], f1_k);
3585 rvec_inc(fshift[t31], f1_l);
3587 rvec_inc(fshift[t21], f2_i);
3588 rvec_inc(fshift[CENTRAL], f2_j);
3589 rvec_inc(fshift[t22], f2_k);
3590 rvec_inc(fshift[t32], f2_l);
3597 /***********************************************************
3599 * G R O M O S 9 6 F U N C T I O N S
3601 ***********************************************************/
3602 real g96harmonic(real kA, real kB, real xA, real xB, real x, real lambda,
3605 const real half = 0.5;
3606 real L1, kk, x0, dx, dx2;
3607 real v, f, dvdlambda;
3610 kk = L1*kA+lambda*kB;
3611 x0 = L1*xA+lambda*xB;
3618 dvdlambda = half*(kB-kA)*dx2 + (xA-xB)*kk*dx;
3625 /* That was 21 flops */
3628 real g96bonds(int nbonds,
3629 const t_iatom forceatoms[], const t_iparams forceparams[],
3630 const rvec x[], rvec f[], rvec fshift[],
3631 const t_pbc *pbc, const t_graph *g,
3632 real lambda, real *dvdlambda,
3633 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3634 int gmx_unused *global_atom_index)
3636 int i, m, ki, ai, aj, type;
3637 real dr2, fbond, vbond, fij, vtot;
3642 for (i = 0; (i < nbonds); )
3644 type = forceatoms[i++];
3645 ai = forceatoms[i++];
3646 aj = forceatoms[i++];
3648 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
3649 dr2 = iprod(dx, dx); /* 5 */
3651 *dvdlambda += g96harmonic(forceparams[type].harmonic.krA,
3652 forceparams[type].harmonic.krB,
3653 forceparams[type].harmonic.rA,
3654 forceparams[type].harmonic.rB,
3655 dr2, lambda, &vbond, &fbond);
3657 vtot += 0.5*vbond; /* 1*/
3661 fprintf(debug, "G96-BONDS: dr = %10g vbond = %10g fbond = %10g\n",
3662 sqrt(dr2), vbond, fbond);
3668 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
3671 for (m = 0; (m < DIM); m++) /* 15 */
3676 fshift[ki][m] += fij;
3677 fshift[CENTRAL][m] -= fij;
3683 real g96bond_angle(const rvec xi, const rvec xj, const rvec xk, const t_pbc *pbc,
3684 rvec r_ij, rvec r_kj,
3686 /* Return value is the angle between the bonds i-j and j-k */
3690 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
3691 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
3693 costh = cos_angle(r_ij, r_kj); /* 25 */
3698 real g96angles(int nbonds,
3699 const t_iatom forceatoms[], const t_iparams forceparams[],
3700 const rvec x[], rvec f[], rvec fshift[],
3701 const t_pbc *pbc, const t_graph *g,
3702 real lambda, real *dvdlambda,
3703 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3704 int gmx_unused *global_atom_index)
3706 int i, ai, aj, ak, type, m, t1, t2;
3708 real cos_theta, dVdt, va, vtot;
3709 real rij_1, rij_2, rkj_1, rkj_2, rijrkj_1;
3711 ivec jt, dt_ij, dt_kj;
3714 for (i = 0; (i < nbonds); )
3716 type = forceatoms[i++];
3717 ai = forceatoms[i++];
3718 aj = forceatoms[i++];
3719 ak = forceatoms[i++];
3721 cos_theta = g96bond_angle(x[ai], x[aj], x[ak], pbc, r_ij, r_kj, &t1, &t2);
3723 *dvdlambda += g96harmonic(forceparams[type].harmonic.krA,
3724 forceparams[type].harmonic.krB,
3725 forceparams[type].harmonic.rA,
3726 forceparams[type].harmonic.rB,
3727 cos_theta, lambda, &va, &dVdt);
3730 rij_1 = gmx_invsqrt(iprod(r_ij, r_ij));
3731 rkj_1 = gmx_invsqrt(iprod(r_kj, r_kj));
3732 rij_2 = rij_1*rij_1;
3733 rkj_2 = rkj_1*rkj_1;
3734 rijrkj_1 = rij_1*rkj_1; /* 23 */
3739 fprintf(debug, "G96ANGLES: costheta = %10g vth = %10g dV/dct = %10g\n",
3740 cos_theta, va, dVdt);
3743 for (m = 0; (m < DIM); m++) /* 42 */
3745 f_i[m] = dVdt*(r_kj[m]*rijrkj_1 - r_ij[m]*rij_2*cos_theta);
3746 f_k[m] = dVdt*(r_ij[m]*rijrkj_1 - r_kj[m]*rkj_2*cos_theta);
3747 f_j[m] = -f_i[m]-f_k[m];
3755 copy_ivec(SHIFT_IVEC(g, aj), jt);
3757 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3758 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3759 t1 = IVEC2IS(dt_ij);
3760 t2 = IVEC2IS(dt_kj);
3762 rvec_inc(fshift[t1], f_i);
3763 rvec_inc(fshift[CENTRAL], f_j);
3764 rvec_inc(fshift[t2], f_k); /* 9 */
3770 real cross_bond_bond(int nbonds,
3771 const t_iatom forceatoms[], const t_iparams forceparams[],
3772 const rvec x[], rvec f[], rvec fshift[],
3773 const t_pbc *pbc, const t_graph *g,
3774 real gmx_unused lambda, real gmx_unused *dvdlambda,
3775 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3776 int gmx_unused *global_atom_index)
3778 /* Potential from Lawrence and Skimmer, Chem. Phys. Lett. 372 (2003)
3781 int i, ai, aj, ak, type, m, t1, t2;
3783 real vtot, vrr, s1, s2, r1, r2, r1e, r2e, krr;
3785 ivec jt, dt_ij, dt_kj;
3788 for (i = 0; (i < nbonds); )
3790 type = forceatoms[i++];
3791 ai = forceatoms[i++];
3792 aj = forceatoms[i++];
3793 ak = forceatoms[i++];
3794 r1e = forceparams[type].cross_bb.r1e;
3795 r2e = forceparams[type].cross_bb.r2e;
3796 krr = forceparams[type].cross_bb.krr;
3798 /* Compute distance vectors ... */
3799 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
3800 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
3802 /* ... and their lengths */
3806 /* Deviations from ideality */
3810 /* Energy (can be negative!) */
3815 svmul(-krr*s2/r1, r_ij, f_i);
3816 svmul(-krr*s1/r2, r_kj, f_k);
3818 for (m = 0; (m < DIM); m++) /* 12 */
3820 f_j[m] = -f_i[m] - f_k[m];
3829 copy_ivec(SHIFT_IVEC(g, aj), jt);
3831 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3832 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3833 t1 = IVEC2IS(dt_ij);
3834 t2 = IVEC2IS(dt_kj);
3836 rvec_inc(fshift[t1], f_i);
3837 rvec_inc(fshift[CENTRAL], f_j);
3838 rvec_inc(fshift[t2], f_k); /* 9 */
3844 real cross_bond_angle(int nbonds,
3845 const t_iatom forceatoms[], const t_iparams forceparams[],
3846 const rvec x[], rvec f[], rvec fshift[],
3847 const t_pbc *pbc, const t_graph *g,
3848 real gmx_unused lambda, real gmx_unused *dvdlambda,
3849 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3850 int gmx_unused *global_atom_index)
3852 /* Potential from Lawrence and Skimmer, Chem. Phys. Lett. 372 (2003)
3855 int i, ai, aj, ak, type, m, t1, t2;
3856 rvec r_ij, r_kj, r_ik;
3857 real vtot, vrt, s1, s2, s3, r1, r2, r3, r1e, r2e, r3e, krt, k1, k2, k3;
3859 ivec jt, dt_ij, dt_kj;
3862 for (i = 0; (i < nbonds); )
3864 type = forceatoms[i++];
3865 ai = forceatoms[i++];
3866 aj = forceatoms[i++];
3867 ak = forceatoms[i++];
3868 r1e = forceparams[type].cross_ba.r1e;
3869 r2e = forceparams[type].cross_ba.r2e;
3870 r3e = forceparams[type].cross_ba.r3e;
3871 krt = forceparams[type].cross_ba.krt;
3873 /* Compute distance vectors ... */
3874 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
3875 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
3876 pbc_rvec_sub(pbc, x[ai], x[ak], r_ik);
3878 /* ... and their lengths */
3883 /* Deviations from ideality */
3888 /* Energy (can be negative!) */
3889 vrt = krt*s3*(s1+s2);
3895 k3 = -krt*(s1+s2)/r3;
3896 for (m = 0; (m < DIM); m++)
3898 f_i[m] = k1*r_ij[m] + k3*r_ik[m];
3899 f_k[m] = k2*r_kj[m] - k3*r_ik[m];
3900 f_j[m] = -f_i[m] - f_k[m];
3903 for (m = 0; (m < DIM); m++) /* 12 */
3913 copy_ivec(SHIFT_IVEC(g, aj), jt);
3915 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3916 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3917 t1 = IVEC2IS(dt_ij);
3918 t2 = IVEC2IS(dt_kj);
3920 rvec_inc(fshift[t1], f_i);
3921 rvec_inc(fshift[CENTRAL], f_j);
3922 rvec_inc(fshift[t2], f_k); /* 9 */
3928 static real bonded_tab(const char *type, int table_nr,
3929 const bondedtable_t *table, real kA, real kB, real r,
3930 real lambda, real *V, real *F)
3932 real k, tabscale, *VFtab, rt, eps, eps2, Yt, Ft, Geps, Heps2, Fp, VV, FF;
3936 k = (1.0 - lambda)*kA + lambda*kB;
3938 tabscale = table->scale;
3939 VFtab = table->data;
3942 n0 = static_cast<int>(rt);
3945 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",
3946 type, table_nr, r, n0, n0+1, table->n);
3953 Geps = VFtab[nnn+2]*eps;
3954 Heps2 = VFtab[nnn+3]*eps2;
3955 Fp = Ft + Geps + Heps2;
3957 FF = Fp + Geps + 2.0*Heps2;
3959 *F = -k*FF*tabscale;
3961 dvdlambda = (kB - kA)*VV;
3965 /* That was 22 flops */
3968 real tab_bonds(int nbonds,
3969 const t_iatom forceatoms[], const t_iparams forceparams[],
3970 const rvec x[], rvec f[], rvec fshift[],
3971 const t_pbc *pbc, const t_graph *g,
3972 real lambda, real *dvdlambda,
3973 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
3974 int gmx_unused *global_atom_index)
3976 int i, m, ki, ai, aj, type, table;
3977 real dr, dr2, fbond, vbond, fij, vtot;
3982 for (i = 0; (i < nbonds); )
3984 type = forceatoms[i++];
3985 ai = forceatoms[i++];
3986 aj = forceatoms[i++];
3988 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
3989 dr2 = iprod(dx, dx); /* 5 */
3990 dr = dr2*gmx_invsqrt(dr2); /* 10 */
3992 table = forceparams[type].tab.table;
3994 *dvdlambda += bonded_tab("bond", table,
3995 &fcd->bondtab[table],
3996 forceparams[type].tab.kA,
3997 forceparams[type].tab.kB,
3998 dr, lambda, &vbond, &fbond); /* 22 */
4006 vtot += vbond; /* 1*/
4007 fbond *= gmx_invsqrt(dr2); /* 6 */
4011 fprintf(debug, "TABBONDS: dr = %10g vbond = %10g fbond = %10g\n",
4017 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
4020 for (m = 0; (m < DIM); m++) /* 15 */
4025 fshift[ki][m] += fij;
4026 fshift[CENTRAL][m] -= fij;
4032 real tab_angles(int nbonds,
4033 const t_iatom forceatoms[], const t_iparams forceparams[],
4034 const rvec x[], rvec f[], rvec fshift[],
4035 const t_pbc *pbc, const t_graph *g,
4036 real lambda, real *dvdlambda,
4037 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
4038 int gmx_unused *global_atom_index)
4040 int i, ai, aj, ak, t1, t2, type, table;
4042 real cos_theta, cos_theta2, theta, dVdt, va, vtot;
4043 ivec jt, dt_ij, dt_kj;
4046 for (i = 0; (i < nbonds); )
4048 type = forceatoms[i++];
4049 ai = forceatoms[i++];
4050 aj = forceatoms[i++];
4051 ak = forceatoms[i++];
4053 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
4054 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
4056 table = forceparams[type].tab.table;
4058 *dvdlambda += bonded_tab("angle", table,
4059 &fcd->angletab[table],
4060 forceparams[type].tab.kA,
4061 forceparams[type].tab.kB,
4062 theta, lambda, &va, &dVdt); /* 22 */
4065 cos_theta2 = sqr(cos_theta); /* 1 */
4074 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
4075 sth = st*cos_theta; /* 1 */
4079 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
4080 theta*RAD2DEG, va, dVdt);
4083 nrkj2 = iprod(r_kj, r_kj); /* 5 */
4084 nrij2 = iprod(r_ij, r_ij);
4086 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
4087 cii = sth/nrij2; /* 10 */
4088 ckk = sth/nrkj2; /* 10 */
4090 for (m = 0; (m < DIM); m++) /* 39 */
4092 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
4093 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
4094 f_j[m] = -f_i[m]-f_k[m];
4101 copy_ivec(SHIFT_IVEC(g, aj), jt);
4103 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
4104 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
4105 t1 = IVEC2IS(dt_ij);
4106 t2 = IVEC2IS(dt_kj);
4108 rvec_inc(fshift[t1], f_i);
4109 rvec_inc(fshift[CENTRAL], f_j);
4110 rvec_inc(fshift[t2], f_k);
4116 real tab_dihs(int nbonds,
4117 const t_iatom forceatoms[], const t_iparams forceparams[],
4118 const rvec x[], rvec f[], rvec fshift[],
4119 const t_pbc *pbc, const t_graph *g,
4120 real lambda, real *dvdlambda,
4121 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
4122 int gmx_unused *global_atom_index)
4124 int i, type, ai, aj, ak, al, table;
4126 rvec r_ij, r_kj, r_kl, m, n;
4127 real phi, sign, ddphi, vpd, vtot;
4130 for (i = 0; (i < nbonds); )
4132 type = forceatoms[i++];
4133 ai = forceatoms[i++];
4134 aj = forceatoms[i++];
4135 ak = forceatoms[i++];
4136 al = forceatoms[i++];
4138 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
4139 &sign, &t1, &t2, &t3); /* 84 */
4141 table = forceparams[type].tab.table;
4143 /* Hopefully phi+M_PI never results in values < 0 */
4144 *dvdlambda += bonded_tab("dihedral", table,
4145 &fcd->dihtab[table],
4146 forceparams[type].tab.kA,
4147 forceparams[type].tab.kB,
4148 phi+M_PI, lambda, &vpd, &ddphi);
4151 do_dih_fup(ai, aj, ak, al, -ddphi, r_ij, r_kj, r_kl, m, n,
4152 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
4155 fprintf(debug, "pdih: (%d,%d,%d,%d) phi=%g\n",
4156 ai, aj, ak, al, phi);
4166 ftype_is_bonded_potential(int ftype)
4169 (interaction_function[ftype].flags & IF_BOND) &&
4170 !(ftype == F_CONNBONDS || ftype == F_POSRES || ftype == F_FBPOSRES) &&
4171 (ftype < F_GB12 || ftype > F_GB14);
4174 /*! \brief Zero thread-local force-output buffers */
4175 static void zero_thread_forces(f_thread_t *f_t, int n,
4176 int nblock, int blocksize)
4178 int b, a0, a1, a, i, j;
4180 if (n > f_t->f_nalloc)
4182 f_t->f_nalloc = over_alloc_large(n);
4183 srenew(f_t->f, f_t->f_nalloc);
4186 if (f_t->red_mask != 0)
4188 for (b = 0; b < nblock; b++)
4190 if (f_t->red_mask && (1U<<b))
4193 a1 = std::min((b+1)*blocksize, n);
4194 for (a = a0; a < a1; a++)
4196 clear_rvec(f_t->f[a]);
4201 for (i = 0; i < SHIFTS; i++)
4203 clear_rvec(f_t->fshift[i]);
4205 for (i = 0; i < F_NRE; i++)
4209 for (i = 0; i < egNR; i++)
4211 for (j = 0; j < f_t->grpp.nener; j++)
4213 f_t->grpp.ener[i][j] = 0;
4216 for (i = 0; i < efptNR; i++)
4222 /*! \brief The max thread number is arbitrary, we used a fixed number
4223 * to avoid memory management. Using more than 16 threads is probably
4224 * never useful performance wise. */
4225 #define MAX_BONDED_THREADS 256
4227 /*! \brief Reduce thread-local force buffers */
4228 static void reduce_thread_force_buffer(int n, rvec *f,
4229 int nthreads, f_thread_t *f_t,
4230 int nblock, int block_size)
4234 if (nthreads > MAX_BONDED_THREADS)
4236 gmx_fatal(FARGS, "Can not reduce bonded forces on more than %d threads",
4237 MAX_BONDED_THREADS);
4240 /* This reduction can run on any number of threads,
4241 * independently of nthreads.
4243 #pragma omp parallel for num_threads(nthreads) schedule(static)
4244 for (b = 0; b < nblock; b++)
4246 rvec *fp[MAX_BONDED_THREADS];
4250 /* Determine which threads contribute to this block */
4252 for (ft = 1; ft < nthreads; ft++)
4254 if (f_t[ft].red_mask & (1U<<b))
4256 fp[nfb++] = f_t[ft].f;
4261 /* Reduce force buffers for threads that contribute */
4263 a1 = (b+1)*block_size;
4264 a1 = std::min(a1, n);
4265 for (a = a0; a < a1; a++)
4267 for (fb = 0; fb < nfb; fb++)
4269 rvec_inc(f[a], fp[fb][a]);
4276 /*! \brief Reduce thread-local forces */
4277 static void reduce_thread_forces(int n, rvec *f, rvec *fshift,
4278 real *ener, gmx_grppairener_t *grpp, real *dvdl,
4279 int nthreads, f_thread_t *f_t,
4280 int nblock, int block_size,
4281 gmx_bool bCalcEnerVir,
4286 /* Reduce the bonded force buffer */
4287 reduce_thread_force_buffer(n, f, nthreads, f_t, nblock, block_size);
4290 /* When necessary, reduce energy and virial using one thread only */
4295 for (i = 0; i < SHIFTS; i++)
4297 for (t = 1; t < nthreads; t++)
4299 rvec_inc(fshift[i], f_t[t].fshift[i]);
4302 for (i = 0; i < F_NRE; i++)
4304 for (t = 1; t < nthreads; t++)
4306 ener[i] += f_t[t].ener[i];
4309 for (i = 0; i < egNR; i++)
4311 for (j = 0; j < f_t[1].grpp.nener; j++)
4313 for (t = 1; t < nthreads; t++)
4316 grpp->ener[i][j] += f_t[t].grpp.ener[i][j];
4322 for (i = 0; i < efptNR; i++)
4325 for (t = 1; t < nthreads; t++)
4327 dvdl[i] += f_t[t].dvdl[i];
4334 /*! \brief Calculate one element of the list of bonded interactions
4336 static real calc_one_bond(int thread,
4337 int ftype, const t_idef *idef,
4338 const rvec x[], rvec f[], rvec fshift[],
4340 const t_pbc *pbc, const t_graph *g,
4341 gmx_grppairener_t *grpp,
4343 real *lambda, real *dvdl,
4344 const t_mdatoms *md, t_fcdata *fcd,
4345 gmx_bool bCalcEnerVir,
4346 int *global_atom_index)
4348 int nat1, nbonds, efptFTYPE;
4353 if (IS_RESTRAINT_TYPE(ftype))
4355 efptFTYPE = efptRESTRAINT;
4359 efptFTYPE = efptBONDED;
4362 nat1 = interaction_function[ftype].nratoms + 1;
4363 nbonds = idef->il[ftype].nr/nat1;
4364 iatoms = idef->il[ftype].iatoms;
4366 nb0 = idef->il_thread_division[ftype*(idef->nthreads+1)+thread];
4367 nbn = idef->il_thread_division[ftype*(idef->nthreads+1)+thread+1] - nb0;
4369 if (!IS_LISTED_LJ_C(ftype))
4371 if (ftype == F_CMAP)
4373 v = cmap_dihs(nbn, iatoms+nb0,
4374 idef->iparams, &idef->cmap_grid,
4376 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
4377 md, fcd, global_atom_index);
4379 #ifdef GMX_SIMD_HAVE_REAL
4380 else if (ftype == F_ANGLES &&
4381 !bCalcEnerVir && fr->efep == efepNO)
4383 /* No energies, shift forces, dvdl */
4384 angles_noener_simd(nbn, idef->il[ftype].iatoms+nb0,
4387 pbc, g, lambda[efptFTYPE], md, fcd,
4392 else if (ftype == F_PDIHS &&
4393 !bCalcEnerVir && fr->efep == efepNO)
4395 /* No energies, shift forces, dvdl */
4396 #ifdef GMX_SIMD_HAVE_REAL
4401 (nbn, idef->il[ftype].iatoms+nb0,
4404 pbc, g, lambda[efptFTYPE], md, fcd,
4410 v = interaction_function[ftype].ifunc(nbn, iatoms+nb0,
4413 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
4414 md, fcd, global_atom_index);
4419 v = do_nonbonded_listed(ftype, nbn, iatoms+nb0, idef->iparams, x, f, fshift,
4420 pbc, g, lambda, dvdl, md, fr, grpp, global_atom_index);
4425 inc_nrnb(nrnb, interaction_function[ftype].nrnb_ind, nbonds);
4431 void calc_bonds(const gmx_multisim_t *ms,
4433 const rvec x[], history_t *hist,
4434 rvec f[], t_forcerec *fr,
4435 const struct t_pbc *pbc, const struct t_graph *g,
4436 gmx_enerdata_t *enerd, t_nrnb *nrnb,
4438 const t_mdatoms *md,
4439 t_fcdata *fcd, int *global_atom_index,
4442 gmx_bool bCalcEnerVir;
4444 real dvdl[efptNR]; /* The dummy array is to have a place to store the dhdl at other values
4445 of lambda, which will be thrown away in the end*/
4446 const t_pbc *pbc_null;
4449 assert(fr->nthreads == idef->nthreads);
4451 bCalcEnerVir = (force_flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY));
4453 for (i = 0; i < efptNR; i++)
4469 p_graph(debug, "Bondage is fun", g);
4473 /* Do pre force calculation stuff which might require communication */
4474 if (idef->il[F_ORIRES].nr)
4476 enerd->term[F_ORIRESDEV] =
4477 calc_orires_dev(ms, idef->il[F_ORIRES].nr,
4478 idef->il[F_ORIRES].iatoms,
4479 idef->iparams, md, x,
4480 pbc_null, fcd, hist);
4482 if (idef->il[F_DISRES].nr)
4484 calc_disres_R_6(idef->il[F_DISRES].nr,
4485 idef->il[F_DISRES].iatoms,
4486 idef->iparams, x, pbc_null,
4489 if (fcd->disres.nsystems > 1)
4491 gmx_sum_sim(2*fcd->disres.nres, fcd->disres.Rt_6, ms);
4496 #pragma omp parallel for num_threads(fr->nthreads) schedule(static)
4497 for (thread = 0; thread < fr->nthreads; thread++)
4504 gmx_grppairener_t *grpp;
4509 fshift = fr->fshift;
4511 grpp = &enerd->grpp;
4516 zero_thread_forces(&fr->f_t[thread], fr->natoms_force,
4517 fr->red_nblock, 1<<fr->red_ashift);
4519 ft = fr->f_t[thread].f;
4520 fshift = fr->f_t[thread].fshift;
4521 epot = fr->f_t[thread].ener;
4522 grpp = &fr->f_t[thread].grpp;
4523 dvdlt = fr->f_t[thread].dvdl;
4525 /* Loop over all bonded force types to calculate the bonded forces */
4526 for (ftype = 0; (ftype < F_NRE); ftype++)
4528 if (idef->il[ftype].nr > 0 && ftype_is_bonded_potential(ftype))
4530 v = calc_one_bond(thread, ftype, idef, x,
4531 ft, fshift, fr, pbc_null, g, grpp,
4532 nrnb, lambda, dvdlt,
4533 md, fcd, bCalcEnerVir,
4539 if (fr->nthreads > 1)
4541 reduce_thread_forces(fr->natoms_force, f, fr->fshift,
4542 enerd->term, &enerd->grpp, dvdl,
4543 fr->nthreads, fr->f_t,
4544 fr->red_nblock, 1<<fr->red_ashift,
4546 force_flags & GMX_FORCE_DHDL);
4548 if (force_flags & GMX_FORCE_DHDL)
4550 for (i = 0; i < efptNR; i++)
4552 enerd->dvdl_nonlin[i] += dvdl[i];
4556 /* Copy the sum of violations for the distance restraints from fcd */
4559 enerd->term[F_DISRESVIOL] = fcd->disres.sumviol;
4564 void calc_bonds_lambda(const t_idef *idef,
4567 const struct t_pbc *pbc, const struct t_graph *g,
4568 gmx_grppairener_t *grpp, real *epot, t_nrnb *nrnb,
4570 const t_mdatoms *md,
4572 int *global_atom_index)
4574 int ftype, nr_nonperturbed, nr;
4576 real dvdl_dum[efptNR];
4578 const t_pbc *pbc_null;
4590 /* Copy the whole idef, so we can modify the contents locally */
4592 idef_fe.nthreads = 1;
4593 snew(idef_fe.il_thread_division, F_NRE*(idef_fe.nthreads+1));
4595 /* We already have the forces, so we use temp buffers here */
4596 snew(f, fr->natoms_force);
4597 snew(fshift, SHIFTS);
4599 /* Loop over all bonded force types to calculate the bonded energies */
4600 for (ftype = 0; (ftype < F_NRE); ftype++)
4602 if (ftype_is_bonded_potential(ftype))
4604 /* Set the work range of thread 0 to the perturbed bondeds only */
4605 nr_nonperturbed = idef->il[ftype].nr_nonperturbed;
4606 nr = idef->il[ftype].nr;
4607 idef_fe.il_thread_division[ftype*2+0] = nr_nonperturbed;
4608 idef_fe.il_thread_division[ftype*2+1] = nr;
4610 /* This is only to get the flop count correct */
4611 idef_fe.il[ftype].nr = nr - nr_nonperturbed;
4613 if (nr - nr_nonperturbed > 0)
4615 v = calc_one_bond(0, ftype, &idef_fe,
4616 x, f, fshift, fr, pbc_null, g,
4617 grpp, nrnb, lambda, dvdl_dum,
4628 sfree(idef_fe.il_thread_division);