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.
49 #include "gromacs/bonded/restcbt.h"
50 #include "gromacs/legacyheaders/disre.h"
51 #include "gromacs/legacyheaders/force.h"
52 #include "gromacs/legacyheaders/macros.h"
53 #include "gromacs/legacyheaders/names.h"
54 #include "gromacs/legacyheaders/nonbonded.h"
55 #include "gromacs/legacyheaders/ns.h"
56 #include "gromacs/legacyheaders/orires.h"
57 #include "gromacs/legacyheaders/txtdump.h"
58 #include "gromacs/math/units.h"
59 #include "gromacs/math/utilities.h"
60 #include "gromacs/math/vec.h"
61 #include "gromacs/pbcutil/ishift.h"
62 #include "gromacs/pbcutil/mshift.h"
63 #include "gromacs/pbcutil/pbc.h"
64 #include "gromacs/simd/simd.h"
65 #include "gromacs/simd/simd_math.h"
66 #include "gromacs/simd/vector_operations.h"
67 #include "gromacs/utility/fatalerror.h"
68 #include "gromacs/utility/smalloc.h"
70 /* Find a better place for this? */
71 const int cmap_coeff_matrix[] = {
72 1, 0, -3, 2, 0, 0, 0, 0, -3, 0, 9, -6, 2, 0, -6, 4,
73 0, 0, 0, 0, 0, 0, 0, 0, 3, 0, -9, 6, -2, 0, 6, -4,
74 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 9, -6, 0, 0, -6, 4,
75 0, 0, 3, -2, 0, 0, 0, 0, 0, 0, -9, 6, 0, 0, 6, -4,
76 0, 0, 0, 0, 1, 0, -3, 2, -2, 0, 6, -4, 1, 0, -3, 2,
77 0, 0, 0, 0, 0, 0, 0, 0, -1, 0, 3, -2, 1, 0, -3, 2,
78 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -3, 2, 0, 0, 3, -2,
79 0, 0, 0, 0, 0, 0, 3, -2, 0, 0, -6, 4, 0, 0, 3, -2,
80 0, 1, -2, 1, 0, 0, 0, 0, 0, -3, 6, -3, 0, 2, -4, 2,
81 0, 0, 0, 0, 0, 0, 0, 0, 0, 3, -6, 3, 0, -2, 4, -2,
82 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -3, 3, 0, 0, 2, -2,
83 0, 0, -1, 1, 0, 0, 0, 0, 0, 0, 3, -3, 0, 0, -2, 2,
84 0, 0, 0, 0, 0, 1, -2, 1, 0, -2, 4, -2, 0, 1, -2, 1,
85 0, 0, 0, 0, 0, 0, 0, 0, 0, -1, 2, -1, 0, 1, -2, 1,
86 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, -1, 0, 0, -1, 1,
87 0, 0, 0, 0, 0, 0, -1, 1, 0, 0, 2, -2, 0, 0, -1, 1
91 /* TODO This function should go and live in nonbonded.c where it is
92 really needed. Here, it only supports giving a fatal error message
94 int glatnr(int *global_atom_index, int i)
98 if (global_atom_index == NULL)
104 atnr = global_atom_index[i] + 1;
110 /* TODO This kind of code appears in many places. Consolidate it */
111 static int pbc_rvec_sub(const t_pbc *pbc, const rvec xi, const rvec xj, rvec dx)
115 return pbc_dx_aiuc(pbc, xi, xj, dx);
119 rvec_sub(xi, xj, dx);
124 #ifdef GMX_SIMD_HAVE_REAL
126 /* SIMD PBC data structure, containing 1/boxdiag and the box vectors */
128 gmx_simd_real_t inv_bzz;
129 gmx_simd_real_t inv_byy;
130 gmx_simd_real_t inv_bxx;
139 /* Set the SIMD pbc data from a normal t_pbc struct */
140 static void set_pbc_simd(const t_pbc *pbc, pbc_simd_t *pbc_simd)
145 /* Setting inv_bdiag to 0 effectively turns off PBC */
146 clear_rvec(inv_bdiag);
149 for (d = 0; d < pbc->ndim_ePBC; d++)
151 inv_bdiag[d] = 1.0/pbc->box[d][d];
155 pbc_simd->inv_bzz = gmx_simd_set1_r(inv_bdiag[ZZ]);
156 pbc_simd->inv_byy = gmx_simd_set1_r(inv_bdiag[YY]);
157 pbc_simd->inv_bxx = gmx_simd_set1_r(inv_bdiag[XX]);
161 pbc_simd->bzx = gmx_simd_set1_r(pbc->box[ZZ][XX]);
162 pbc_simd->bzy = gmx_simd_set1_r(pbc->box[ZZ][YY]);
163 pbc_simd->bzz = gmx_simd_set1_r(pbc->box[ZZ][ZZ]);
164 pbc_simd->byx = gmx_simd_set1_r(pbc->box[YY][XX]);
165 pbc_simd->byy = gmx_simd_set1_r(pbc->box[YY][YY]);
166 pbc_simd->bxx = gmx_simd_set1_r(pbc->box[XX][XX]);
170 pbc_simd->bzx = gmx_simd_setzero_r();
171 pbc_simd->bzy = gmx_simd_setzero_r();
172 pbc_simd->bzz = gmx_simd_setzero_r();
173 pbc_simd->byx = gmx_simd_setzero_r();
174 pbc_simd->byy = gmx_simd_setzero_r();
175 pbc_simd->bxx = gmx_simd_setzero_r();
179 /* Correct distance vector *dx,*dy,*dz for PBC using SIMD */
180 static gmx_inline void
181 pbc_dx_simd(gmx_simd_real_t *dx, gmx_simd_real_t *dy, gmx_simd_real_t *dz,
182 const pbc_simd_t *pbc)
186 sh = gmx_simd_round_r(gmx_simd_mul_r(*dz, pbc->inv_bzz));
187 *dx = gmx_simd_fnmadd_r(sh, pbc->bzx, *dx);
188 *dy = gmx_simd_fnmadd_r(sh, pbc->bzy, *dy);
189 *dz = gmx_simd_fnmadd_r(sh, pbc->bzz, *dz);
191 sh = gmx_simd_round_r(gmx_simd_mul_r(*dy, pbc->inv_byy));
192 *dx = gmx_simd_fnmadd_r(sh, pbc->byx, *dx);
193 *dy = gmx_simd_fnmadd_r(sh, pbc->byy, *dy);
195 sh = gmx_simd_round_r(gmx_simd_mul_r(*dx, pbc->inv_bxx));
196 *dx = gmx_simd_fnmadd_r(sh, pbc->bxx, *dx);
199 #endif /* GMX_SIMD_HAVE_REAL */
202 * Morse potential bond by Frank Everdij
204 * Three parameters needed:
206 * b0 = equilibrium distance in nm
207 * be = beta in nm^-1 (actually, it's nu_e*Sqrt(2*pi*pi*mu/D_e))
208 * cb = well depth in kJ/mol
210 * Note: the potential is referenced to be +cb at infinite separation
211 * and zero at the equilibrium distance!
214 real morse_bonds(int nbonds,
215 const t_iatom forceatoms[], const t_iparams forceparams[],
216 const rvec x[], rvec f[], rvec fshift[],
217 const t_pbc *pbc, const t_graph *g,
218 real lambda, real *dvdlambda,
219 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
220 int gmx_unused *global_atom_index)
222 const real one = 1.0;
223 const real two = 2.0;
224 real dr, dr2, temp, omtemp, cbomtemp, fbond, vbond, fij, vtot;
225 real b0, be, cb, b0A, beA, cbA, b0B, beB, cbB, L1;
227 int i, m, ki, type, ai, aj;
231 for (i = 0; (i < nbonds); )
233 type = forceatoms[i++];
234 ai = forceatoms[i++];
235 aj = forceatoms[i++];
237 b0A = forceparams[type].morse.b0A;
238 beA = forceparams[type].morse.betaA;
239 cbA = forceparams[type].morse.cbA;
241 b0B = forceparams[type].morse.b0B;
242 beB = forceparams[type].morse.betaB;
243 cbB = forceparams[type].morse.cbB;
245 L1 = one-lambda; /* 1 */
246 b0 = L1*b0A + lambda*b0B; /* 3 */
247 be = L1*beA + lambda*beB; /* 3 */
248 cb = L1*cbA + lambda*cbB; /* 3 */
250 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
251 dr2 = iprod(dx, dx); /* 5 */
252 dr = dr2*gmx_invsqrt(dr2); /* 10 */
253 temp = exp(-be*(dr-b0)); /* 12 */
257 /* bonds are constrainted. This may _not_ include bond constraints if they are lambda dependent */
258 *dvdlambda += cbB-cbA;
262 omtemp = one-temp; /* 1 */
263 cbomtemp = cb*omtemp; /* 1 */
264 vbond = cbomtemp*omtemp; /* 1 */
265 fbond = -two*be*temp*cbomtemp*gmx_invsqrt(dr2); /* 9 */
266 vtot += vbond; /* 1 */
268 *dvdlambda += (cbB - cbA) * omtemp * omtemp - (2-2*omtemp)*omtemp * cb * ((b0B-b0A)*be - (beB-beA)*(dr-b0)); /* 15 */
272 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
276 for (m = 0; (m < DIM); m++) /* 15 */
281 fshift[ki][m] += fij;
282 fshift[CENTRAL][m] -= fij;
288 real cubic_bonds(int nbonds,
289 const t_iatom forceatoms[], const t_iparams forceparams[],
290 const rvec x[], rvec f[], rvec fshift[],
291 const t_pbc *pbc, const t_graph *g,
292 real gmx_unused lambda, real gmx_unused *dvdlambda,
293 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
294 int gmx_unused *global_atom_index)
296 const real three = 3.0;
297 const real two = 2.0;
299 real dr, dr2, dist, kdist, kdist2, fbond, vbond, fij, vtot;
301 int i, m, ki, type, ai, aj;
305 for (i = 0; (i < nbonds); )
307 type = forceatoms[i++];
308 ai = forceatoms[i++];
309 aj = forceatoms[i++];
311 b0 = forceparams[type].cubic.b0;
312 kb = forceparams[type].cubic.kb;
313 kcub = forceparams[type].cubic.kcub;
315 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
316 dr2 = iprod(dx, dx); /* 5 */
323 dr = dr2*gmx_invsqrt(dr2); /* 10 */
328 vbond = kdist2 + kcub*kdist2*dist;
329 fbond = -(two*kdist + three*kdist2*kcub)/dr;
331 vtot += vbond; /* 21 */
335 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
338 for (m = 0; (m < DIM); m++) /* 15 */
343 fshift[ki][m] += fij;
344 fshift[CENTRAL][m] -= fij;
350 real FENE_bonds(int nbonds,
351 const t_iatom forceatoms[], const t_iparams forceparams[],
352 const rvec x[], rvec f[], rvec fshift[],
353 const t_pbc *pbc, const t_graph *g,
354 real gmx_unused lambda, real gmx_unused *dvdlambda,
355 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
356 int *global_atom_index)
358 const real half = 0.5;
359 const real one = 1.0;
361 real dr2, bm2, omdr2obm2, fbond, vbond, fij, vtot;
363 int i, m, ki, type, ai, aj;
367 for (i = 0; (i < nbonds); )
369 type = forceatoms[i++];
370 ai = forceatoms[i++];
371 aj = forceatoms[i++];
373 bm = forceparams[type].fene.bm;
374 kb = forceparams[type].fene.kb;
376 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
377 dr2 = iprod(dx, dx); /* 5 */
389 "r^2 (%f) >= bm^2 (%f) in FENE bond between atoms %d and %d",
391 glatnr(global_atom_index, ai),
392 glatnr(global_atom_index, aj));
395 omdr2obm2 = one - dr2/bm2;
397 vbond = -half*kb*bm2*log(omdr2obm2);
398 fbond = -kb/omdr2obm2;
400 vtot += vbond; /* 35 */
404 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
407 for (m = 0; (m < DIM); m++) /* 15 */
412 fshift[ki][m] += fij;
413 fshift[CENTRAL][m] -= fij;
419 real harmonic(real kA, real kB, real xA, real xB, real x, real lambda,
422 const real half = 0.5;
423 real L1, kk, x0, dx, dx2;
424 real v, f, dvdlambda;
427 kk = L1*kA+lambda*kB;
428 x0 = L1*xA+lambda*xB;
435 dvdlambda = half*(kB-kA)*dx2 + (xA-xB)*kk*dx;
442 /* That was 19 flops */
446 real bonds(int nbonds,
447 const t_iatom forceatoms[], const t_iparams forceparams[],
448 const rvec x[], rvec f[], rvec fshift[],
449 const t_pbc *pbc, const t_graph *g,
450 real lambda, real *dvdlambda,
451 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
452 int gmx_unused *global_atom_index)
454 int i, m, ki, ai, aj, type;
455 real dr, dr2, fbond, vbond, fij, vtot;
460 for (i = 0; (i < nbonds); )
462 type = forceatoms[i++];
463 ai = forceatoms[i++];
464 aj = forceatoms[i++];
466 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
467 dr2 = iprod(dx, dx); /* 5 */
468 dr = dr2*gmx_invsqrt(dr2); /* 10 */
470 *dvdlambda += harmonic(forceparams[type].harmonic.krA,
471 forceparams[type].harmonic.krB,
472 forceparams[type].harmonic.rA,
473 forceparams[type].harmonic.rB,
474 dr, lambda, &vbond, &fbond); /* 19 */
482 vtot += vbond; /* 1*/
483 fbond *= gmx_invsqrt(dr2); /* 6 */
487 fprintf(debug, "BONDS: dr = %10g vbond = %10g fbond = %10g\n",
493 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
496 for (m = 0; (m < DIM); m++) /* 15 */
501 fshift[ki][m] += fij;
502 fshift[CENTRAL][m] -= fij;
508 real restraint_bonds(int nbonds,
509 const t_iatom forceatoms[], const t_iparams forceparams[],
510 const rvec x[], rvec f[], rvec fshift[],
511 const t_pbc *pbc, const t_graph *g,
512 real lambda, real *dvdlambda,
513 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
514 int gmx_unused *global_atom_index)
516 int i, m, ki, ai, aj, type;
517 real dr, dr2, fbond, vbond, fij, vtot;
519 real low, dlow, up1, dup1, up2, dup2, k, dk;
527 for (i = 0; (i < nbonds); )
529 type = forceatoms[i++];
530 ai = forceatoms[i++];
531 aj = forceatoms[i++];
533 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
534 dr2 = iprod(dx, dx); /* 5 */
535 dr = dr2*gmx_invsqrt(dr2); /* 10 */
537 low = L1*forceparams[type].restraint.lowA + lambda*forceparams[type].restraint.lowB;
538 dlow = -forceparams[type].restraint.lowA + forceparams[type].restraint.lowB;
539 up1 = L1*forceparams[type].restraint.up1A + lambda*forceparams[type].restraint.up1B;
540 dup1 = -forceparams[type].restraint.up1A + forceparams[type].restraint.up1B;
541 up2 = L1*forceparams[type].restraint.up2A + lambda*forceparams[type].restraint.up2B;
542 dup2 = -forceparams[type].restraint.up2A + forceparams[type].restraint.up2B;
543 k = L1*forceparams[type].restraint.kA + lambda*forceparams[type].restraint.kB;
544 dk = -forceparams[type].restraint.kA + forceparams[type].restraint.kB;
553 *dvdlambda += 0.5*dk*drh2 - k*dlow*drh;
566 *dvdlambda += 0.5*dk*drh2 - k*dup1*drh;
571 vbond = k*(up2 - up1)*(0.5*(up2 - up1) + drh);
572 fbond = -k*(up2 - up1);
573 *dvdlambda += dk*(up2 - up1)*(0.5*(up2 - up1) + drh)
574 + k*(dup2 - dup1)*(up2 - up1 + drh)
575 - k*(up2 - up1)*dup2;
583 vtot += vbond; /* 1*/
584 fbond *= gmx_invsqrt(dr2); /* 6 */
588 fprintf(debug, "BONDS: dr = %10g vbond = %10g fbond = %10g\n",
594 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
597 for (m = 0; (m < DIM); m++) /* 15 */
602 fshift[ki][m] += fij;
603 fshift[CENTRAL][m] -= fij;
610 real polarize(int nbonds,
611 const t_iatom forceatoms[], const t_iparams forceparams[],
612 const rvec x[], rvec f[], rvec fshift[],
613 const t_pbc *pbc, const t_graph *g,
614 real lambda, real *dvdlambda,
615 const t_mdatoms *md, t_fcdata gmx_unused *fcd,
616 int gmx_unused *global_atom_index)
618 int i, m, ki, ai, aj, type;
619 real dr, dr2, fbond, vbond, fij, vtot, ksh;
624 for (i = 0; (i < nbonds); )
626 type = forceatoms[i++];
627 ai = forceatoms[i++];
628 aj = forceatoms[i++];
629 ksh = sqr(md->chargeA[aj])*ONE_4PI_EPS0/forceparams[type].polarize.alpha;
632 fprintf(debug, "POL: local ai = %d aj = %d ksh = %.3f\n", ai, aj, ksh);
635 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
636 dr2 = iprod(dx, dx); /* 5 */
637 dr = dr2*gmx_invsqrt(dr2); /* 10 */
639 *dvdlambda += harmonic(ksh, ksh, 0, 0, dr, lambda, &vbond, &fbond); /* 19 */
646 vtot += vbond; /* 1*/
647 fbond *= gmx_invsqrt(dr2); /* 6 */
651 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
654 for (m = 0; (m < DIM); m++) /* 15 */
659 fshift[ki][m] += fij;
660 fshift[CENTRAL][m] -= fij;
666 real anharm_polarize(int nbonds,
667 const t_iatom forceatoms[], const t_iparams forceparams[],
668 const rvec x[], rvec f[], rvec fshift[],
669 const t_pbc *pbc, const t_graph *g,
670 real lambda, real *dvdlambda,
671 const t_mdatoms *md, t_fcdata gmx_unused *fcd,
672 int gmx_unused *global_atom_index)
674 int i, m, ki, ai, aj, type;
675 real dr, dr2, fbond, vbond, fij, vtot, ksh, khyp, drcut, ddr, ddr3;
680 for (i = 0; (i < nbonds); )
682 type = forceatoms[i++];
683 ai = forceatoms[i++];
684 aj = forceatoms[i++];
685 ksh = sqr(md->chargeA[aj])*ONE_4PI_EPS0/forceparams[type].anharm_polarize.alpha; /* 7*/
686 khyp = forceparams[type].anharm_polarize.khyp;
687 drcut = forceparams[type].anharm_polarize.drcut;
690 fprintf(debug, "POL: local ai = %d aj = %d ksh = %.3f\n", ai, aj, ksh);
693 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
694 dr2 = iprod(dx, dx); /* 5 */
695 dr = dr2*gmx_invsqrt(dr2); /* 10 */
697 *dvdlambda += harmonic(ksh, ksh, 0, 0, dr, lambda, &vbond, &fbond); /* 19 */
708 vbond += khyp*ddr*ddr3;
709 fbond -= 4*khyp*ddr3;
711 fbond *= gmx_invsqrt(dr2); /* 6 */
712 vtot += vbond; /* 1*/
716 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
719 for (m = 0; (m < DIM); m++) /* 15 */
724 fshift[ki][m] += fij;
725 fshift[CENTRAL][m] -= fij;
731 real water_pol(int nbonds,
732 const t_iatom forceatoms[], const t_iparams forceparams[],
733 const rvec x[], rvec f[], rvec gmx_unused fshift[],
734 const t_pbc gmx_unused *pbc, const t_graph gmx_unused *g,
735 real gmx_unused lambda, real gmx_unused *dvdlambda,
736 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
737 int gmx_unused *global_atom_index)
739 /* This routine implements anisotropic polarizibility for water, through
740 * a shell connected to a dummy with spring constant that differ in the
741 * three spatial dimensions in the molecular frame.
743 int i, m, aO, aH1, aH2, aD, aS, type, type0, ki;
745 rvec dOH1, dOH2, dHH, dOD, dDS, nW, kk, dx, kdx, proj;
749 real vtot, fij, r_HH, r_OD, r_nW, tx, ty, tz, qS;
754 type0 = forceatoms[0];
756 qS = md->chargeA[aS];
757 kk[XX] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_x;
758 kk[YY] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_y;
759 kk[ZZ] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_z;
760 r_HH = 1.0/forceparams[type0].wpol.rHH;
763 fprintf(debug, "WPOL: qS = %10.5f aS = %5d\n", qS, aS);
764 fprintf(debug, "WPOL: kk = %10.3f %10.3f %10.3f\n",
765 kk[XX], kk[YY], kk[ZZ]);
766 fprintf(debug, "WPOL: rOH = %10.3f rHH = %10.3f rOD = %10.3f\n",
767 forceparams[type0].wpol.rOH,
768 forceparams[type0].wpol.rHH,
769 forceparams[type0].wpol.rOD);
771 for (i = 0; (i < nbonds); i += 6)
773 type = forceatoms[i];
776 gmx_fatal(FARGS, "Sorry, type = %d, type0 = %d, file = %s, line = %d",
777 type, type0, __FILE__, __LINE__);
779 aO = forceatoms[i+1];
780 aH1 = forceatoms[i+2];
781 aH2 = forceatoms[i+3];
782 aD = forceatoms[i+4];
783 aS = forceatoms[i+5];
785 /* Compute vectors describing the water frame */
786 pbc_rvec_sub(pbc, x[aH1], x[aO], dOH1);
787 pbc_rvec_sub(pbc, x[aH2], x[aO], dOH2);
788 pbc_rvec_sub(pbc, x[aH2], x[aH1], dHH);
789 pbc_rvec_sub(pbc, x[aD], x[aO], dOD);
790 ki = pbc_rvec_sub(pbc, x[aS], x[aD], dDS);
791 cprod(dOH1, dOH2, nW);
793 /* Compute inverse length of normal vector
794 * (this one could be precomputed, but I'm too lazy now)
796 r_nW = gmx_invsqrt(iprod(nW, nW));
797 /* This is for precision, but does not make a big difference,
800 r_OD = gmx_invsqrt(iprod(dOD, dOD));
802 /* Normalize the vectors in the water frame */
804 svmul(r_HH, dHH, dHH);
805 svmul(r_OD, dOD, dOD);
807 /* Compute displacement of shell along components of the vector */
808 dx[ZZ] = iprod(dDS, dOD);
809 /* Compute projection on the XY plane: dDS - dx[ZZ]*dOD */
810 for (m = 0; (m < DIM); m++)
812 proj[m] = dDS[m]-dx[ZZ]*dOD[m];
815 /*dx[XX] = iprod(dDS,nW);
816 dx[YY] = iprod(dDS,dHH);*/
817 dx[XX] = iprod(proj, nW);
818 for (m = 0; (m < DIM); m++)
820 proj[m] -= dx[XX]*nW[m];
822 dx[YY] = iprod(proj, dHH);
827 fprintf(debug, "WPOL: dx2=%10g dy2=%10g dz2=%10g sum=%10g dDS^2=%10g\n",
828 sqr(dx[XX]), sqr(dx[YY]), sqr(dx[ZZ]), iprod(dx, dx), iprod(dDS, dDS));
829 fprintf(debug, "WPOL: dHH=(%10g,%10g,%10g)\n", dHH[XX], dHH[YY], dHH[ZZ]);
830 fprintf(debug, "WPOL: dOD=(%10g,%10g,%10g), 1/r_OD = %10g\n",
831 dOD[XX], dOD[YY], dOD[ZZ], 1/r_OD);
832 fprintf(debug, "WPOL: nW =(%10g,%10g,%10g), 1/r_nW = %10g\n",
833 nW[XX], nW[YY], nW[ZZ], 1/r_nW);
834 fprintf(debug, "WPOL: dx =%10g, dy =%10g, dz =%10g\n",
835 dx[XX], dx[YY], dx[ZZ]);
836 fprintf(debug, "WPOL: dDSx=%10g, dDSy=%10g, dDSz=%10g\n",
837 dDS[XX], dDS[YY], dDS[ZZ]);
840 /* Now compute the forces and energy */
841 kdx[XX] = kk[XX]*dx[XX];
842 kdx[YY] = kk[YY]*dx[YY];
843 kdx[ZZ] = kk[ZZ]*dx[ZZ];
844 vtot += iprod(dx, kdx);
848 ivec_sub(SHIFT_IVEC(g, aS), SHIFT_IVEC(g, aD), dt);
852 for (m = 0; (m < DIM); m++)
854 /* This is a tensor operation but written out for speed */
864 fshift[ki][m] += fij;
865 fshift[CENTRAL][m] -= fij;
870 fprintf(debug, "WPOL: vwpol=%g\n", 0.5*iprod(dx, kdx));
871 fprintf(debug, "WPOL: df = (%10g, %10g, %10g)\n", df[XX], df[YY], df[ZZ]);
879 static real do_1_thole(const rvec xi, const rvec xj, rvec fi, rvec fj,
880 const t_pbc *pbc, real qq,
881 rvec fshift[], real afac)
884 real r12sq, r12_1, r12bar, v0, v1, fscal, ebar, fff;
887 t = pbc_rvec_sub(pbc, xi, xj, r12); /* 3 */
889 r12sq = iprod(r12, r12); /* 5 */
890 r12_1 = gmx_invsqrt(r12sq); /* 5 */
891 r12bar = afac/r12_1; /* 5 */
892 v0 = qq*ONE_4PI_EPS0*r12_1; /* 2 */
893 ebar = exp(-r12bar); /* 5 */
894 v1 = (1-(1+0.5*r12bar)*ebar); /* 4 */
895 fscal = ((v0*r12_1)*v1 - v0*0.5*afac*ebar*(r12bar+1))*r12_1; /* 9 */
898 fprintf(debug, "THOLE: v0 = %.3f v1 = %.3f r12= % .3f r12bar = %.3f fscal = %.3f ebar = %.3f\n", v0, v1, 1/r12_1, r12bar, fscal, ebar);
901 for (m = 0; (m < DIM); m++)
907 fshift[CENTRAL][m] -= fff;
910 return v0*v1; /* 1 */
914 real thole_pol(int nbonds,
915 const t_iatom forceatoms[], const t_iparams forceparams[],
916 const rvec x[], rvec f[], rvec fshift[],
917 const t_pbc *pbc, const t_graph gmx_unused *g,
918 real gmx_unused lambda, real gmx_unused *dvdlambda,
919 const t_mdatoms *md, t_fcdata gmx_unused *fcd,
920 int gmx_unused *global_atom_index)
922 /* Interaction between two pairs of particles with opposite charge */
923 int i, type, a1, da1, a2, da2;
924 real q1, q2, qq, a, al1, al2, afac;
926 const real minusOneOnSix = -1.0/6.0;
928 for (i = 0; (i < nbonds); )
930 type = forceatoms[i++];
931 a1 = forceatoms[i++];
932 da1 = forceatoms[i++];
933 a2 = forceatoms[i++];
934 da2 = forceatoms[i++];
935 q1 = md->chargeA[da1];
936 q2 = md->chargeA[da2];
937 a = forceparams[type].thole.a;
938 al1 = forceparams[type].thole.alpha1;
939 al2 = forceparams[type].thole.alpha2;
941 afac = a*pow(al1*al2, minusOneOnSix);
942 V += do_1_thole(x[a1], x[a2], f[a1], f[a2], pbc, qq, fshift, afac);
943 V += do_1_thole(x[da1], x[a2], f[da1], f[a2], pbc, -qq, fshift, afac);
944 V += do_1_thole(x[a1], x[da2], f[a1], f[da2], pbc, -qq, fshift, afac);
945 V += do_1_thole(x[da1], x[da2], f[da1], f[da2], pbc, qq, fshift, afac);
951 real bond_angle(const rvec xi, const rvec xj, const rvec xk, const t_pbc *pbc,
952 rvec r_ij, rvec r_kj, real *costh,
954 /* Return value is the angle between the bonds i-j and j-k */
959 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
960 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
962 *costh = cos_angle(r_ij, r_kj); /* 25 */
963 th = acos(*costh); /* 10 */
968 real angles(int nbonds,
969 const t_iatom forceatoms[], const t_iparams forceparams[],
970 const rvec x[], rvec f[], rvec fshift[],
971 const t_pbc *pbc, const t_graph *g,
972 real lambda, real *dvdlambda,
973 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
974 int gmx_unused *global_atom_index)
976 int i, ai, aj, ak, t1, t2, type;
978 real cos_theta, cos_theta2, theta, dVdt, va, vtot;
979 ivec jt, dt_ij, dt_kj;
982 for (i = 0; i < nbonds; )
984 type = forceatoms[i++];
985 ai = forceatoms[i++];
986 aj = forceatoms[i++];
987 ak = forceatoms[i++];
989 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
990 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
992 *dvdlambda += harmonic(forceparams[type].harmonic.krA,
993 forceparams[type].harmonic.krB,
994 forceparams[type].harmonic.rA*DEG2RAD,
995 forceparams[type].harmonic.rB*DEG2RAD,
996 theta, lambda, &va, &dVdt); /* 21 */
999 cos_theta2 = sqr(cos_theta);
1006 real nrkj_1, nrij_1;
1009 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1010 sth = st*cos_theta; /* 1 */
1014 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1015 theta*RAD2DEG, va, dVdt);
1018 nrij2 = iprod(r_ij, r_ij); /* 5 */
1019 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1021 nrij_1 = gmx_invsqrt(nrij2); /* 10 */
1022 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1024 cik = st*nrij_1*nrkj_1; /* 2 */
1025 cii = sth*nrij_1*nrij_1; /* 2 */
1026 ckk = sth*nrkj_1*nrkj_1; /* 2 */
1028 for (m = 0; m < DIM; m++)
1030 f_i[m] = -(cik*r_kj[m] - cii*r_ij[m]);
1031 f_k[m] = -(cik*r_ij[m] - ckk*r_kj[m]);
1032 f_j[m] = -f_i[m] - f_k[m];
1039 copy_ivec(SHIFT_IVEC(g, aj), jt);
1041 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1042 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1043 t1 = IVEC2IS(dt_ij);
1044 t2 = IVEC2IS(dt_kj);
1046 rvec_inc(fshift[t1], f_i);
1047 rvec_inc(fshift[CENTRAL], f_j);
1048 rvec_inc(fshift[t2], f_k);
1055 #ifdef GMX_SIMD_HAVE_REAL
1057 /* As angles, but using SIMD to calculate many dihedrals at once.
1058 * This routines does not calculate energies and shift forces.
1060 static gmx_inline void
1061 angles_noener_simd(int nbonds,
1062 const t_iatom forceatoms[], const t_iparams forceparams[],
1063 const rvec x[], rvec f[],
1064 const t_pbc *pbc, const t_graph gmx_unused *g,
1065 real gmx_unused lambda,
1066 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1067 int gmx_unused *global_atom_index)
1071 int type, ai[GMX_SIMD_REAL_WIDTH], aj[GMX_SIMD_REAL_WIDTH];
1072 int ak[GMX_SIMD_REAL_WIDTH];
1073 real coeff_array[2*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *coeff;
1074 real dr_array[2*DIM*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *dr;
1075 real f_buf_array[6*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *f_buf;
1076 gmx_simd_real_t k_S, theta0_S;
1077 gmx_simd_real_t rijx_S, rijy_S, rijz_S;
1078 gmx_simd_real_t rkjx_S, rkjy_S, rkjz_S;
1079 gmx_simd_real_t one_S;
1080 gmx_simd_real_t min_one_plus_eps_S;
1081 gmx_simd_real_t rij_rkj_S;
1082 gmx_simd_real_t nrij2_S, nrij_1_S;
1083 gmx_simd_real_t nrkj2_S, nrkj_1_S;
1084 gmx_simd_real_t cos_S, invsin_S;
1085 gmx_simd_real_t theta_S;
1086 gmx_simd_real_t st_S, sth_S;
1087 gmx_simd_real_t cik_S, cii_S, ckk_S;
1088 gmx_simd_real_t f_ix_S, f_iy_S, f_iz_S;
1089 gmx_simd_real_t f_kx_S, f_ky_S, f_kz_S;
1090 pbc_simd_t pbc_simd;
1092 /* Ensure register memory alignment */
1093 coeff = gmx_simd_align_r(coeff_array);
1094 dr = gmx_simd_align_r(dr_array);
1095 f_buf = gmx_simd_align_r(f_buf_array);
1097 set_pbc_simd(pbc, &pbc_simd);
1099 one_S = gmx_simd_set1_r(1.0);
1101 /* The smallest number > -1 */
1102 min_one_plus_eps_S = gmx_simd_set1_r(-1.0 + 2*GMX_REAL_EPS);
1104 /* nbonds is the number of angles times nfa1, here we step GMX_SIMD_REAL_WIDTH angles */
1105 for (i = 0; (i < nbonds); i += GMX_SIMD_REAL_WIDTH*nfa1)
1107 /* Collect atoms for GMX_SIMD_REAL_WIDTH angles.
1108 * iu indexes into forceatoms, we should not let iu go beyond nbonds.
1111 for (s = 0; s < GMX_SIMD_REAL_WIDTH; s++)
1113 type = forceatoms[iu];
1114 ai[s] = forceatoms[iu+1];
1115 aj[s] = forceatoms[iu+2];
1116 ak[s] = forceatoms[iu+3];
1118 coeff[s] = forceparams[type].harmonic.krA;
1119 coeff[GMX_SIMD_REAL_WIDTH+s] = forceparams[type].harmonic.rA*DEG2RAD;
1121 /* If you can't use pbc_dx_simd below for PBC, e.g. because
1122 * you can't round in SIMD, use pbc_rvec_sub here.
1124 /* Store the non PBC corrected distances packed and aligned */
1125 for (m = 0; m < DIM; m++)
1127 dr[s + m *GMX_SIMD_REAL_WIDTH] = x[ai[s]][m] - x[aj[s]][m];
1128 dr[s + (DIM+m)*GMX_SIMD_REAL_WIDTH] = x[ak[s]][m] - x[aj[s]][m];
1131 /* At the end fill the arrays with identical entries */
1132 if (iu + nfa1 < nbonds)
1138 k_S = gmx_simd_load_r(coeff);
1139 theta0_S = gmx_simd_load_r(coeff+GMX_SIMD_REAL_WIDTH);
1141 rijx_S = gmx_simd_load_r(dr + 0*GMX_SIMD_REAL_WIDTH);
1142 rijy_S = gmx_simd_load_r(dr + 1*GMX_SIMD_REAL_WIDTH);
1143 rijz_S = gmx_simd_load_r(dr + 2*GMX_SIMD_REAL_WIDTH);
1144 rkjx_S = gmx_simd_load_r(dr + 3*GMX_SIMD_REAL_WIDTH);
1145 rkjy_S = gmx_simd_load_r(dr + 4*GMX_SIMD_REAL_WIDTH);
1146 rkjz_S = gmx_simd_load_r(dr + 5*GMX_SIMD_REAL_WIDTH);
1148 pbc_dx_simd(&rijx_S, &rijy_S, &rijz_S, &pbc_simd);
1149 pbc_dx_simd(&rkjx_S, &rkjy_S, &rkjz_S, &pbc_simd);
1151 rij_rkj_S = gmx_simd_iprod_r(rijx_S, rijy_S, rijz_S,
1152 rkjx_S, rkjy_S, rkjz_S);
1154 nrij2_S = gmx_simd_norm2_r(rijx_S, rijy_S, rijz_S);
1155 nrkj2_S = gmx_simd_norm2_r(rkjx_S, rkjy_S, rkjz_S);
1157 nrij_1_S = gmx_simd_invsqrt_r(nrij2_S);
1158 nrkj_1_S = gmx_simd_invsqrt_r(nrkj2_S);
1160 cos_S = gmx_simd_mul_r(rij_rkj_S, gmx_simd_mul_r(nrij_1_S, nrkj_1_S));
1162 /* To allow for 180 degrees, we take the max of cos and -1 + 1bit,
1163 * so we can safely get the 1/sin from 1/sqrt(1 - cos^2).
1164 * This also ensures that rounding errors would cause the argument
1165 * of gmx_simd_acos_r to be < -1.
1166 * Note that we do not take precautions for cos(0)=1, so the outer
1167 * atoms in an angle should not be on top of each other.
1169 cos_S = gmx_simd_max_r(cos_S, min_one_plus_eps_S);
1171 theta_S = gmx_simd_acos_r(cos_S);
1173 invsin_S = gmx_simd_invsqrt_r(gmx_simd_sub_r(one_S, gmx_simd_mul_r(cos_S, cos_S)));
1175 st_S = gmx_simd_mul_r(gmx_simd_mul_r(k_S, gmx_simd_sub_r(theta0_S, theta_S)),
1177 sth_S = gmx_simd_mul_r(st_S, cos_S);
1179 cik_S = gmx_simd_mul_r(st_S, gmx_simd_mul_r(nrij_1_S, nrkj_1_S));
1180 cii_S = gmx_simd_mul_r(sth_S, gmx_simd_mul_r(nrij_1_S, nrij_1_S));
1181 ckk_S = gmx_simd_mul_r(sth_S, gmx_simd_mul_r(nrkj_1_S, nrkj_1_S));
1183 f_ix_S = gmx_simd_mul_r(cii_S, rijx_S);
1184 f_ix_S = gmx_simd_fnmadd_r(cik_S, rkjx_S, f_ix_S);
1185 f_iy_S = gmx_simd_mul_r(cii_S, rijy_S);
1186 f_iy_S = gmx_simd_fnmadd_r(cik_S, rkjy_S, f_iy_S);
1187 f_iz_S = gmx_simd_mul_r(cii_S, rijz_S);
1188 f_iz_S = gmx_simd_fnmadd_r(cik_S, rkjz_S, f_iz_S);
1189 f_kx_S = gmx_simd_mul_r(ckk_S, rkjx_S);
1190 f_kx_S = gmx_simd_fnmadd_r(cik_S, rijx_S, f_kx_S);
1191 f_ky_S = gmx_simd_mul_r(ckk_S, rkjy_S);
1192 f_ky_S = gmx_simd_fnmadd_r(cik_S, rijy_S, f_ky_S);
1193 f_kz_S = gmx_simd_mul_r(ckk_S, rkjz_S);
1194 f_kz_S = gmx_simd_fnmadd_r(cik_S, rijz_S, f_kz_S);
1196 gmx_simd_store_r(f_buf + 0*GMX_SIMD_REAL_WIDTH, f_ix_S);
1197 gmx_simd_store_r(f_buf + 1*GMX_SIMD_REAL_WIDTH, f_iy_S);
1198 gmx_simd_store_r(f_buf + 2*GMX_SIMD_REAL_WIDTH, f_iz_S);
1199 gmx_simd_store_r(f_buf + 3*GMX_SIMD_REAL_WIDTH, f_kx_S);
1200 gmx_simd_store_r(f_buf + 4*GMX_SIMD_REAL_WIDTH, f_ky_S);
1201 gmx_simd_store_r(f_buf + 5*GMX_SIMD_REAL_WIDTH, f_kz_S);
1207 for (m = 0; m < DIM; m++)
1209 f[ai[s]][m] += f_buf[s + m*GMX_SIMD_REAL_WIDTH];
1210 f[aj[s]][m] -= f_buf[s + m*GMX_SIMD_REAL_WIDTH] + f_buf[s + (DIM+m)*GMX_SIMD_REAL_WIDTH];
1211 f[ak[s]][m] += f_buf[s + (DIM+m)*GMX_SIMD_REAL_WIDTH];
1216 while (s < GMX_SIMD_REAL_WIDTH && iu < nbonds);
1220 #endif /* GMX_SIMD_HAVE_REAL */
1222 real linear_angles(int nbonds,
1223 const t_iatom forceatoms[], const t_iparams forceparams[],
1224 const rvec x[], rvec f[], rvec fshift[],
1225 const t_pbc *pbc, const t_graph *g,
1226 real lambda, real *dvdlambda,
1227 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1228 int gmx_unused *global_atom_index)
1230 int i, m, ai, aj, ak, t1, t2, type;
1232 real L1, kA, kB, aA, aB, dr, dr2, va, vtot, a, b, klin;
1233 ivec jt, dt_ij, dt_kj;
1234 rvec r_ij, r_kj, r_ik, dx;
1238 for (i = 0; (i < nbonds); )
1240 type = forceatoms[i++];
1241 ai = forceatoms[i++];
1242 aj = forceatoms[i++];
1243 ak = forceatoms[i++];
1245 kA = forceparams[type].linangle.klinA;
1246 kB = forceparams[type].linangle.klinB;
1247 klin = L1*kA + lambda*kB;
1249 aA = forceparams[type].linangle.aA;
1250 aB = forceparams[type].linangle.aB;
1251 a = L1*aA+lambda*aB;
1254 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
1255 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
1256 rvec_sub(r_ij, r_kj, r_ik);
1259 for (m = 0; (m < DIM); m++)
1261 dr = -a * r_ij[m] - b * r_kj[m];
1266 f_j[m] = -(f_i[m]+f_k[m]);
1272 *dvdlambda += 0.5*(kB-kA)*dr2 + klin*(aB-aA)*iprod(dx, r_ik);
1278 copy_ivec(SHIFT_IVEC(g, aj), jt);
1280 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1281 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1282 t1 = IVEC2IS(dt_ij);
1283 t2 = IVEC2IS(dt_kj);
1285 rvec_inc(fshift[t1], f_i);
1286 rvec_inc(fshift[CENTRAL], f_j);
1287 rvec_inc(fshift[t2], f_k);
1292 real urey_bradley(int nbonds,
1293 const t_iatom forceatoms[], const t_iparams forceparams[],
1294 const rvec x[], rvec f[], rvec fshift[],
1295 const t_pbc *pbc, const t_graph *g,
1296 real lambda, real *dvdlambda,
1297 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1298 int gmx_unused *global_atom_index)
1300 int i, m, ai, aj, ak, t1, t2, type, ki;
1301 rvec r_ij, r_kj, r_ik;
1302 real cos_theta, cos_theta2, theta;
1303 real dVdt, va, vtot, dr, dr2, vbond, fbond, fik;
1304 real kthA, th0A, kUBA, r13A, kthB, th0B, kUBB, r13B;
1305 ivec jt, dt_ij, dt_kj, dt_ik;
1308 for (i = 0; (i < nbonds); )
1310 type = forceatoms[i++];
1311 ai = forceatoms[i++];
1312 aj = forceatoms[i++];
1313 ak = forceatoms[i++];
1314 th0A = forceparams[type].u_b.thetaA*DEG2RAD;
1315 kthA = forceparams[type].u_b.kthetaA;
1316 r13A = forceparams[type].u_b.r13A;
1317 kUBA = forceparams[type].u_b.kUBA;
1318 th0B = forceparams[type].u_b.thetaB*DEG2RAD;
1319 kthB = forceparams[type].u_b.kthetaB;
1320 r13B = forceparams[type].u_b.r13B;
1321 kUBB = forceparams[type].u_b.kUBB;
1323 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
1324 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
1326 *dvdlambda += harmonic(kthA, kthB, th0A, th0B, theta, lambda, &va, &dVdt); /* 21 */
1329 ki = pbc_rvec_sub(pbc, x[ai], x[ak], r_ik); /* 3 */
1330 dr2 = iprod(r_ik, r_ik); /* 5 */
1331 dr = dr2*gmx_invsqrt(dr2); /* 10 */
1333 *dvdlambda += harmonic(kUBA, kUBB, r13A, r13B, dr, lambda, &vbond, &fbond); /* 19 */
1335 cos_theta2 = sqr(cos_theta); /* 1 */
1343 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1344 sth = st*cos_theta; /* 1 */
1348 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1349 theta*RAD2DEG, va, dVdt);
1352 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1353 nrij2 = iprod(r_ij, r_ij);
1355 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
1356 cii = sth/nrij2; /* 10 */
1357 ckk = sth/nrkj2; /* 10 */
1359 for (m = 0; (m < DIM); m++) /* 39 */
1361 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
1362 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
1363 f_j[m] = -f_i[m]-f_k[m];
1370 copy_ivec(SHIFT_IVEC(g, aj), jt);
1372 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1373 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1374 t1 = IVEC2IS(dt_ij);
1375 t2 = IVEC2IS(dt_kj);
1377 rvec_inc(fshift[t1], f_i);
1378 rvec_inc(fshift[CENTRAL], f_j);
1379 rvec_inc(fshift[t2], f_k);
1381 /* Time for the bond calculations */
1387 vtot += vbond; /* 1*/
1388 fbond *= gmx_invsqrt(dr2); /* 6 */
1392 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ak), dt_ik);
1393 ki = IVEC2IS(dt_ik);
1395 for (m = 0; (m < DIM); m++) /* 15 */
1397 fik = fbond*r_ik[m];
1400 fshift[ki][m] += fik;
1401 fshift[CENTRAL][m] -= fik;
1407 real quartic_angles(int nbonds,
1408 const t_iatom forceatoms[], const t_iparams forceparams[],
1409 const rvec x[], rvec f[], rvec fshift[],
1410 const t_pbc *pbc, const t_graph *g,
1411 real gmx_unused lambda, real gmx_unused *dvdlambda,
1412 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1413 int gmx_unused *global_atom_index)
1415 int i, j, ai, aj, ak, t1, t2, type;
1417 real cos_theta, cos_theta2, theta, dt, dVdt, va, dtp, c, vtot;
1418 ivec jt, dt_ij, dt_kj;
1421 for (i = 0; (i < nbonds); )
1423 type = forceatoms[i++];
1424 ai = forceatoms[i++];
1425 aj = forceatoms[i++];
1426 ak = forceatoms[i++];
1428 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
1429 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
1431 dt = theta - forceparams[type].qangle.theta*DEG2RAD; /* 2 */
1434 va = forceparams[type].qangle.c[0];
1436 for (j = 1; j <= 4; j++)
1438 c = forceparams[type].qangle.c[j];
1447 cos_theta2 = sqr(cos_theta); /* 1 */
1456 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1457 sth = st*cos_theta; /* 1 */
1461 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1462 theta*RAD2DEG, va, dVdt);
1465 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1466 nrij2 = iprod(r_ij, r_ij);
1468 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
1469 cii = sth/nrij2; /* 10 */
1470 ckk = sth/nrkj2; /* 10 */
1472 for (m = 0; (m < DIM); m++) /* 39 */
1474 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
1475 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
1476 f_j[m] = -f_i[m]-f_k[m];
1483 copy_ivec(SHIFT_IVEC(g, aj), jt);
1485 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1486 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1487 t1 = IVEC2IS(dt_ij);
1488 t2 = IVEC2IS(dt_kj);
1490 rvec_inc(fshift[t1], f_i);
1491 rvec_inc(fshift[CENTRAL], f_j);
1492 rvec_inc(fshift[t2], f_k);
1498 real dih_angle(const rvec xi, const rvec xj, const rvec xk, const rvec xl,
1500 rvec r_ij, rvec r_kj, rvec r_kl, rvec m, rvec n,
1501 real *sign, int *t1, int *t2, int *t3)
1505 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
1506 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
1507 *t3 = pbc_rvec_sub(pbc, xk, xl, r_kl); /* 3 */
1509 cprod(r_ij, r_kj, m); /* 9 */
1510 cprod(r_kj, r_kl, n); /* 9 */
1511 phi = gmx_angle(m, n); /* 49 (assuming 25 for atan2) */
1512 ipr = iprod(r_ij, n); /* 5 */
1513 (*sign) = (ipr < 0.0) ? -1.0 : 1.0;
1514 phi = (*sign)*phi; /* 1 */
1520 #ifdef GMX_SIMD_HAVE_REAL
1522 /* As dih_angle above, but calculates 4 dihedral angles at once using SIMD,
1523 * also calculates the pre-factor required for the dihedral force update.
1524 * Note that bv and buf should be register aligned.
1526 static gmx_inline void
1527 dih_angle_simd(const rvec *x,
1528 const int *ai, const int *aj, const int *ak, const int *al,
1529 const pbc_simd_t *pbc,
1531 gmx_simd_real_t *phi_S,
1532 gmx_simd_real_t *mx_S, gmx_simd_real_t *my_S, gmx_simd_real_t *mz_S,
1533 gmx_simd_real_t *nx_S, gmx_simd_real_t *ny_S, gmx_simd_real_t *nz_S,
1534 gmx_simd_real_t *nrkj_m2_S,
1535 gmx_simd_real_t *nrkj_n2_S,
1540 gmx_simd_real_t rijx_S, rijy_S, rijz_S;
1541 gmx_simd_real_t rkjx_S, rkjy_S, rkjz_S;
1542 gmx_simd_real_t rklx_S, rkly_S, rklz_S;
1543 gmx_simd_real_t cx_S, cy_S, cz_S;
1544 gmx_simd_real_t cn_S;
1545 gmx_simd_real_t s_S;
1546 gmx_simd_real_t ipr_S;
1547 gmx_simd_real_t iprm_S, iprn_S;
1548 gmx_simd_real_t nrkj2_S, nrkj_1_S, nrkj_2_S, nrkj_S;
1549 gmx_simd_real_t toler_S;
1550 gmx_simd_real_t p_S, q_S;
1551 gmx_simd_real_t nrkj2_min_S;
1552 gmx_simd_real_t real_eps_S;
1554 /* Used to avoid division by zero.
1555 * We take into acount that we multiply the result by real_eps_S.
1557 nrkj2_min_S = gmx_simd_set1_r(GMX_REAL_MIN/(2*GMX_REAL_EPS));
1559 /* The value of the last significant bit (GMX_REAL_EPS is half of that) */
1560 real_eps_S = gmx_simd_set1_r(2*GMX_REAL_EPS);
1562 for (s = 0; s < GMX_SIMD_REAL_WIDTH; s++)
1564 /* If you can't use pbc_dx_simd below for PBC, e.g. because
1565 * you can't round in SIMD, use pbc_rvec_sub here.
1567 for (m = 0; m < DIM; m++)
1569 dr[s + (0*DIM + m)*GMX_SIMD_REAL_WIDTH] = x[ai[s]][m] - x[aj[s]][m];
1570 dr[s + (1*DIM + m)*GMX_SIMD_REAL_WIDTH] = x[ak[s]][m] - x[aj[s]][m];
1571 dr[s + (2*DIM + m)*GMX_SIMD_REAL_WIDTH] = x[ak[s]][m] - x[al[s]][m];
1575 rijx_S = gmx_simd_load_r(dr + 0*GMX_SIMD_REAL_WIDTH);
1576 rijy_S = gmx_simd_load_r(dr + 1*GMX_SIMD_REAL_WIDTH);
1577 rijz_S = gmx_simd_load_r(dr + 2*GMX_SIMD_REAL_WIDTH);
1578 rkjx_S = gmx_simd_load_r(dr + 3*GMX_SIMD_REAL_WIDTH);
1579 rkjy_S = gmx_simd_load_r(dr + 4*GMX_SIMD_REAL_WIDTH);
1580 rkjz_S = gmx_simd_load_r(dr + 5*GMX_SIMD_REAL_WIDTH);
1581 rklx_S = gmx_simd_load_r(dr + 6*GMX_SIMD_REAL_WIDTH);
1582 rkly_S = gmx_simd_load_r(dr + 7*GMX_SIMD_REAL_WIDTH);
1583 rklz_S = gmx_simd_load_r(dr + 8*GMX_SIMD_REAL_WIDTH);
1585 pbc_dx_simd(&rijx_S, &rijy_S, &rijz_S, pbc);
1586 pbc_dx_simd(&rkjx_S, &rkjy_S, &rkjz_S, pbc);
1587 pbc_dx_simd(&rklx_S, &rkly_S, &rklz_S, pbc);
1589 gmx_simd_cprod_r(rijx_S, rijy_S, rijz_S,
1590 rkjx_S, rkjy_S, rkjz_S,
1593 gmx_simd_cprod_r(rkjx_S, rkjy_S, rkjz_S,
1594 rklx_S, rkly_S, rklz_S,
1597 gmx_simd_cprod_r(*mx_S, *my_S, *mz_S,
1598 *nx_S, *ny_S, *nz_S,
1599 &cx_S, &cy_S, &cz_S);
1601 cn_S = gmx_simd_sqrt_r(gmx_simd_norm2_r(cx_S, cy_S, cz_S));
1603 s_S = gmx_simd_iprod_r(*mx_S, *my_S, *mz_S, *nx_S, *ny_S, *nz_S);
1605 /* Determine the dihedral angle, the sign might need correction */
1606 *phi_S = gmx_simd_atan2_r(cn_S, s_S);
1608 ipr_S = gmx_simd_iprod_r(rijx_S, rijy_S, rijz_S,
1609 *nx_S, *ny_S, *nz_S);
1611 iprm_S = gmx_simd_norm2_r(*mx_S, *my_S, *mz_S);
1612 iprn_S = gmx_simd_norm2_r(*nx_S, *ny_S, *nz_S);
1614 nrkj2_S = gmx_simd_norm2_r(rkjx_S, rkjy_S, rkjz_S);
1616 /* Avoid division by zero. When zero, the result is multiplied by 0
1617 * anyhow, so the 3 max below do not affect the final result.
1619 nrkj2_S = gmx_simd_max_r(nrkj2_S, nrkj2_min_S);
1620 nrkj_1_S = gmx_simd_invsqrt_r(nrkj2_S);
1621 nrkj_2_S = gmx_simd_mul_r(nrkj_1_S, nrkj_1_S);
1622 nrkj_S = gmx_simd_mul_r(nrkj2_S, nrkj_1_S);
1624 toler_S = gmx_simd_mul_r(nrkj2_S, real_eps_S);
1626 /* Here the plain-C code uses a conditional, but we can't do that in SIMD.
1627 * So we take a max with the tolerance instead. Since we multiply with
1628 * m or n later, the max does not affect the results.
1630 iprm_S = gmx_simd_max_r(iprm_S, toler_S);
1631 iprn_S = gmx_simd_max_r(iprn_S, toler_S);
1632 *nrkj_m2_S = gmx_simd_mul_r(nrkj_S, gmx_simd_inv_r(iprm_S));
1633 *nrkj_n2_S = gmx_simd_mul_r(nrkj_S, gmx_simd_inv_r(iprn_S));
1635 /* Set sign of phi_S with the sign of ipr_S; phi_S is currently positive */
1636 *phi_S = gmx_simd_xor_sign_r(*phi_S, ipr_S);
1637 p_S = gmx_simd_iprod_r(rijx_S, rijy_S, rijz_S,
1638 rkjx_S, rkjy_S, rkjz_S);
1639 p_S = gmx_simd_mul_r(p_S, nrkj_2_S);
1641 q_S = gmx_simd_iprod_r(rklx_S, rkly_S, rklz_S,
1642 rkjx_S, rkjy_S, rkjz_S);
1643 q_S = gmx_simd_mul_r(q_S, nrkj_2_S);
1645 gmx_simd_store_r(p, p_S);
1646 gmx_simd_store_r(q, q_S);
1649 #endif /* GMX_SIMD_HAVE_REAL */
1652 void do_dih_fup(int i, int j, int k, int l, real ddphi,
1653 rvec r_ij, rvec r_kj, rvec r_kl,
1654 rvec m, rvec n, rvec f[], rvec fshift[],
1655 const t_pbc *pbc, const t_graph *g,
1656 const rvec x[], int t1, int t2, int t3)
1659 rvec f_i, f_j, f_k, f_l;
1660 rvec uvec, vvec, svec, dx_jl;
1661 real iprm, iprn, nrkj, nrkj2, nrkj_1, nrkj_2;
1662 real a, b, p, q, toler;
1663 ivec jt, dt_ij, dt_kj, dt_lj;
1665 iprm = iprod(m, m); /* 5 */
1666 iprn = iprod(n, n); /* 5 */
1667 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1668 toler = nrkj2*GMX_REAL_EPS;
1669 if ((iprm > toler) && (iprn > toler))
1671 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1672 nrkj_2 = nrkj_1*nrkj_1; /* 1 */
1673 nrkj = nrkj2*nrkj_1; /* 1 */
1674 a = -ddphi*nrkj/iprm; /* 11 */
1675 svmul(a, m, f_i); /* 3 */
1676 b = ddphi*nrkj/iprn; /* 11 */
1677 svmul(b, n, f_l); /* 3 */
1678 p = iprod(r_ij, r_kj); /* 5 */
1679 p *= nrkj_2; /* 1 */
1680 q = iprod(r_kl, r_kj); /* 5 */
1681 q *= nrkj_2; /* 1 */
1682 svmul(p, f_i, uvec); /* 3 */
1683 svmul(q, f_l, vvec); /* 3 */
1684 rvec_sub(uvec, vvec, svec); /* 3 */
1685 rvec_sub(f_i, svec, f_j); /* 3 */
1686 rvec_add(f_l, svec, f_k); /* 3 */
1687 rvec_inc(f[i], f_i); /* 3 */
1688 rvec_dec(f[j], f_j); /* 3 */
1689 rvec_dec(f[k], f_k); /* 3 */
1690 rvec_inc(f[l], f_l); /* 3 */
1694 copy_ivec(SHIFT_IVEC(g, j), jt);
1695 ivec_sub(SHIFT_IVEC(g, i), jt, dt_ij);
1696 ivec_sub(SHIFT_IVEC(g, k), jt, dt_kj);
1697 ivec_sub(SHIFT_IVEC(g, l), jt, dt_lj);
1698 t1 = IVEC2IS(dt_ij);
1699 t2 = IVEC2IS(dt_kj);
1700 t3 = IVEC2IS(dt_lj);
1704 t3 = pbc_rvec_sub(pbc, x[l], x[j], dx_jl);
1711 rvec_inc(fshift[t1], f_i);
1712 rvec_dec(fshift[CENTRAL], f_j);
1713 rvec_dec(fshift[t2], f_k);
1714 rvec_inc(fshift[t3], f_l);
1719 /* As do_dih_fup above, but without shift forces */
1721 do_dih_fup_noshiftf(int i, int j, int k, int l, real ddphi,
1722 rvec r_ij, rvec r_kj, rvec r_kl,
1723 rvec m, rvec n, rvec f[])
1725 rvec f_i, f_j, f_k, f_l;
1726 rvec uvec, vvec, svec;
1727 real iprm, iprn, nrkj, nrkj2, nrkj_1, nrkj_2;
1728 real a, b, p, q, toler;
1730 iprm = iprod(m, m); /* 5 */
1731 iprn = iprod(n, n); /* 5 */
1732 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1733 toler = nrkj2*GMX_REAL_EPS;
1734 if ((iprm > toler) && (iprn > toler))
1736 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1737 nrkj_2 = nrkj_1*nrkj_1; /* 1 */
1738 nrkj = nrkj2*nrkj_1; /* 1 */
1739 a = -ddphi*nrkj/iprm; /* 11 */
1740 svmul(a, m, f_i); /* 3 */
1741 b = ddphi*nrkj/iprn; /* 11 */
1742 svmul(b, n, f_l); /* 3 */
1743 p = iprod(r_ij, r_kj); /* 5 */
1744 p *= nrkj_2; /* 1 */
1745 q = iprod(r_kl, r_kj); /* 5 */
1746 q *= nrkj_2; /* 1 */
1747 svmul(p, f_i, uvec); /* 3 */
1748 svmul(q, f_l, vvec); /* 3 */
1749 rvec_sub(uvec, vvec, svec); /* 3 */
1750 rvec_sub(f_i, svec, f_j); /* 3 */
1751 rvec_add(f_l, svec, f_k); /* 3 */
1752 rvec_inc(f[i], f_i); /* 3 */
1753 rvec_dec(f[j], f_j); /* 3 */
1754 rvec_dec(f[k], f_k); /* 3 */
1755 rvec_inc(f[l], f_l); /* 3 */
1759 /* As do_dih_fup_noshiftf above, but with pre-calculated pre-factors */
1760 static gmx_inline void
1761 do_dih_fup_noshiftf_precalc(int i, int j, int k, int l,
1763 real f_i_x, real f_i_y, real f_i_z,
1764 real mf_l_x, real mf_l_y, real mf_l_z,
1767 rvec f_i, f_j, f_k, f_l;
1768 rvec uvec, vvec, svec;
1776 svmul(p, f_i, uvec);
1777 svmul(q, f_l, vvec);
1778 rvec_sub(uvec, vvec, svec);
1779 rvec_sub(f_i, svec, f_j);
1780 rvec_add(f_l, svec, f_k);
1781 rvec_inc(f[i], f_i);
1782 rvec_dec(f[j], f_j);
1783 rvec_dec(f[k], f_k);
1784 rvec_inc(f[l], f_l);
1788 real dopdihs(real cpA, real cpB, real phiA, real phiB, int mult,
1789 real phi, real lambda, real *V, real *F)
1791 real v, dvdlambda, mdphi, v1, sdphi, ddphi;
1792 real L1 = 1.0 - lambda;
1793 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1794 real dph0 = (phiB - phiA)*DEG2RAD;
1795 real cp = L1*cpA + lambda*cpB;
1797 mdphi = mult*phi - ph0;
1799 ddphi = -cp*mult*sdphi;
1800 v1 = 1.0 + cos(mdphi);
1803 dvdlambda = (cpB - cpA)*v1 + cp*dph0*sdphi;
1810 /* That was 40 flops */
1814 dopdihs_noener(real cpA, real cpB, real phiA, real phiB, int mult,
1815 real phi, real lambda, real *F)
1817 real mdphi, sdphi, ddphi;
1818 real L1 = 1.0 - lambda;
1819 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1820 real cp = L1*cpA + lambda*cpB;
1822 mdphi = mult*phi - ph0;
1824 ddphi = -cp*mult*sdphi;
1828 /* That was 20 flops */
1832 dopdihs_mdphi(real cpA, real cpB, real phiA, real phiB, int mult,
1833 real phi, real lambda, real *cp, real *mdphi)
1835 real L1 = 1.0 - lambda;
1836 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1838 *cp = L1*cpA + lambda*cpB;
1840 *mdphi = mult*phi - ph0;
1843 static real dopdihs_min(real cpA, real cpB, real phiA, real phiB, int mult,
1844 real phi, real lambda, real *V, real *F)
1845 /* similar to dopdihs, except for a minus sign *
1846 * and a different treatment of mult/phi0 */
1848 real v, dvdlambda, mdphi, v1, sdphi, ddphi;
1849 real L1 = 1.0 - lambda;
1850 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1851 real dph0 = (phiB - phiA)*DEG2RAD;
1852 real cp = L1*cpA + lambda*cpB;
1854 mdphi = mult*(phi-ph0);
1856 ddphi = cp*mult*sdphi;
1857 v1 = 1.0-cos(mdphi);
1860 dvdlambda = (cpB-cpA)*v1 + cp*dph0*sdphi;
1867 /* That was 40 flops */
1870 real pdihs(int nbonds,
1871 const t_iatom forceatoms[], const t_iparams forceparams[],
1872 const rvec x[], rvec f[], rvec fshift[],
1873 const t_pbc *pbc, const t_graph *g,
1874 real lambda, real *dvdlambda,
1875 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1876 int gmx_unused *global_atom_index)
1878 int i, type, ai, aj, ak, al;
1880 rvec r_ij, r_kj, r_kl, m, n;
1881 real phi, sign, ddphi, vpd, vtot;
1885 for (i = 0; (i < nbonds); )
1887 type = forceatoms[i++];
1888 ai = forceatoms[i++];
1889 aj = forceatoms[i++];
1890 ak = forceatoms[i++];
1891 al = forceatoms[i++];
1893 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
1894 &sign, &t1, &t2, &t3); /* 84 */
1895 *dvdlambda += dopdihs(forceparams[type].pdihs.cpA,
1896 forceparams[type].pdihs.cpB,
1897 forceparams[type].pdihs.phiA,
1898 forceparams[type].pdihs.phiB,
1899 forceparams[type].pdihs.mult,
1900 phi, lambda, &vpd, &ddphi);
1903 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
1904 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
1907 fprintf(debug, "pdih: (%d,%d,%d,%d) phi=%g\n",
1908 ai, aj, ak, al, phi);
1915 void make_dp_periodic(real *dp) /* 1 flop? */
1917 /* dp cannot be outside (-pi,pi) */
1922 else if (*dp < -M_PI)
1929 /* As pdihs above, but without calculating energies and shift forces */
1931 pdihs_noener(int nbonds,
1932 const t_iatom forceatoms[], const t_iparams forceparams[],
1933 const rvec x[], rvec f[],
1934 const t_pbc gmx_unused *pbc, const t_graph gmx_unused *g,
1936 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1937 int gmx_unused *global_atom_index)
1939 int i, type, ai, aj, ak, al;
1941 rvec r_ij, r_kj, r_kl, m, n;
1942 real phi, sign, ddphi_tot, ddphi;
1944 for (i = 0; (i < nbonds); )
1946 ai = forceatoms[i+1];
1947 aj = forceatoms[i+2];
1948 ak = forceatoms[i+3];
1949 al = forceatoms[i+4];
1951 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
1952 &sign, &t1, &t2, &t3);
1956 /* Loop over dihedrals working on the same atoms,
1957 * so we avoid recalculating angles and force distributions.
1961 type = forceatoms[i];
1962 dopdihs_noener(forceparams[type].pdihs.cpA,
1963 forceparams[type].pdihs.cpB,
1964 forceparams[type].pdihs.phiA,
1965 forceparams[type].pdihs.phiB,
1966 forceparams[type].pdihs.mult,
1967 phi, lambda, &ddphi);
1972 while (i < nbonds &&
1973 forceatoms[i+1] == ai &&
1974 forceatoms[i+2] == aj &&
1975 forceatoms[i+3] == ak &&
1976 forceatoms[i+4] == al);
1978 do_dih_fup_noshiftf(ai, aj, ak, al, ddphi_tot, r_ij, r_kj, r_kl, m, n, f);
1983 #ifdef GMX_SIMD_HAVE_REAL
1985 /* As pdihs_noner above, but using SIMD to calculate many dihedrals at once */
1987 pdihs_noener_simd(int nbonds,
1988 const t_iatom forceatoms[], const t_iparams forceparams[],
1989 const rvec x[], rvec f[],
1990 const t_pbc *pbc, const t_graph gmx_unused *g,
1991 real gmx_unused lambda,
1992 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1993 int gmx_unused *global_atom_index)
1997 int type, ai[GMX_SIMD_REAL_WIDTH], aj[GMX_SIMD_REAL_WIDTH], ak[GMX_SIMD_REAL_WIDTH], al[GMX_SIMD_REAL_WIDTH];
1998 real dr_array[3*DIM*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *dr;
1999 real buf_array[7*GMX_SIMD_REAL_WIDTH+GMX_SIMD_REAL_WIDTH], *buf;
2000 real *cp, *phi0, *mult, *p, *q;
2001 gmx_simd_real_t phi0_S, phi_S;
2002 gmx_simd_real_t mx_S, my_S, mz_S;
2003 gmx_simd_real_t nx_S, ny_S, nz_S;
2004 gmx_simd_real_t nrkj_m2_S, nrkj_n2_S;
2005 gmx_simd_real_t cp_S, mdphi_S, mult_S;
2006 gmx_simd_real_t sin_S, cos_S;
2007 gmx_simd_real_t mddphi_S;
2008 gmx_simd_real_t sf_i_S, msf_l_S;
2009 pbc_simd_t pbc_simd;
2011 /* Ensure SIMD register alignment */
2012 dr = gmx_simd_align_r(dr_array);
2013 buf = gmx_simd_align_r(buf_array);
2015 /* Extract aligned pointer for parameters and variables */
2016 cp = buf + 0*GMX_SIMD_REAL_WIDTH;
2017 phi0 = buf + 1*GMX_SIMD_REAL_WIDTH;
2018 mult = buf + 2*GMX_SIMD_REAL_WIDTH;
2019 p = buf + 3*GMX_SIMD_REAL_WIDTH;
2020 q = buf + 4*GMX_SIMD_REAL_WIDTH;
2022 set_pbc_simd(pbc, &pbc_simd);
2024 /* nbonds is the number of dihedrals times nfa1, here we step GMX_SIMD_REAL_WIDTH dihs */
2025 for (i = 0; (i < nbonds); i += GMX_SIMD_REAL_WIDTH*nfa1)
2027 /* Collect atoms quadruplets for GMX_SIMD_REAL_WIDTH dihedrals.
2028 * iu indexes into forceatoms, we should not let iu go beyond nbonds.
2031 for (s = 0; s < GMX_SIMD_REAL_WIDTH; s++)
2033 type = forceatoms[iu];
2034 ai[s] = forceatoms[iu+1];
2035 aj[s] = forceatoms[iu+2];
2036 ak[s] = forceatoms[iu+3];
2037 al[s] = forceatoms[iu+4];
2039 cp[s] = forceparams[type].pdihs.cpA;
2040 phi0[s] = forceparams[type].pdihs.phiA*DEG2RAD;
2041 mult[s] = forceparams[type].pdihs.mult;
2043 /* At the end fill the arrays with identical entries */
2044 if (iu + nfa1 < nbonds)
2050 /* Caclulate GMX_SIMD_REAL_WIDTH dihedral angles at once */
2051 dih_angle_simd(x, ai, aj, ak, al, &pbc_simd,
2054 &mx_S, &my_S, &mz_S,
2055 &nx_S, &ny_S, &nz_S,
2060 cp_S = gmx_simd_load_r(cp);
2061 phi0_S = gmx_simd_load_r(phi0);
2062 mult_S = gmx_simd_load_r(mult);
2064 mdphi_S = gmx_simd_sub_r(gmx_simd_mul_r(mult_S, phi_S), phi0_S);
2066 /* Calculate GMX_SIMD_REAL_WIDTH sines at once */
2067 gmx_simd_sincos_r(mdphi_S, &sin_S, &cos_S);
2068 mddphi_S = gmx_simd_mul_r(gmx_simd_mul_r(cp_S, mult_S), sin_S);
2069 sf_i_S = gmx_simd_mul_r(mddphi_S, nrkj_m2_S);
2070 msf_l_S = gmx_simd_mul_r(mddphi_S, nrkj_n2_S);
2072 /* After this m?_S will contain f[i] */
2073 mx_S = gmx_simd_mul_r(sf_i_S, mx_S);
2074 my_S = gmx_simd_mul_r(sf_i_S, my_S);
2075 mz_S = gmx_simd_mul_r(sf_i_S, mz_S);
2077 /* After this m?_S will contain -f[l] */
2078 nx_S = gmx_simd_mul_r(msf_l_S, nx_S);
2079 ny_S = gmx_simd_mul_r(msf_l_S, ny_S);
2080 nz_S = gmx_simd_mul_r(msf_l_S, nz_S);
2082 gmx_simd_store_r(dr + 0*GMX_SIMD_REAL_WIDTH, mx_S);
2083 gmx_simd_store_r(dr + 1*GMX_SIMD_REAL_WIDTH, my_S);
2084 gmx_simd_store_r(dr + 2*GMX_SIMD_REAL_WIDTH, mz_S);
2085 gmx_simd_store_r(dr + 3*GMX_SIMD_REAL_WIDTH, nx_S);
2086 gmx_simd_store_r(dr + 4*GMX_SIMD_REAL_WIDTH, ny_S);
2087 gmx_simd_store_r(dr + 5*GMX_SIMD_REAL_WIDTH, nz_S);
2093 do_dih_fup_noshiftf_precalc(ai[s], aj[s], ak[s], al[s],
2095 dr[ XX *GMX_SIMD_REAL_WIDTH+s],
2096 dr[ YY *GMX_SIMD_REAL_WIDTH+s],
2097 dr[ ZZ *GMX_SIMD_REAL_WIDTH+s],
2098 dr[(DIM+XX)*GMX_SIMD_REAL_WIDTH+s],
2099 dr[(DIM+YY)*GMX_SIMD_REAL_WIDTH+s],
2100 dr[(DIM+ZZ)*GMX_SIMD_REAL_WIDTH+s],
2105 while (s < GMX_SIMD_REAL_WIDTH && iu < nbonds);
2109 #endif /* GMX_SIMD_HAVE_REAL */
2112 real idihs(int nbonds,
2113 const t_iatom forceatoms[], const t_iparams forceparams[],
2114 const rvec x[], rvec f[], rvec fshift[],
2115 const t_pbc *pbc, const t_graph *g,
2116 real lambda, real *dvdlambda,
2117 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2118 int gmx_unused *global_atom_index)
2120 int i, type, ai, aj, ak, al;
2122 real phi, phi0, dphi0, ddphi, sign, vtot;
2123 rvec r_ij, r_kj, r_kl, m, n;
2124 real L1, kk, dp, dp2, kA, kB, pA, pB, dvdl_term;
2129 for (i = 0; (i < nbonds); )
2131 type = forceatoms[i++];
2132 ai = forceatoms[i++];
2133 aj = forceatoms[i++];
2134 ak = forceatoms[i++];
2135 al = forceatoms[i++];
2137 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2138 &sign, &t1, &t2, &t3); /* 84 */
2140 /* phi can jump if phi0 is close to Pi/-Pi, which will cause huge
2141 * force changes if we just apply a normal harmonic.
2142 * Instead, we first calculate phi-phi0 and take it modulo (-Pi,Pi).
2143 * This means we will never have the periodicity problem, unless
2144 * the dihedral is Pi away from phiO, which is very unlikely due to
2147 kA = forceparams[type].harmonic.krA;
2148 kB = forceparams[type].harmonic.krB;
2149 pA = forceparams[type].harmonic.rA;
2150 pB = forceparams[type].harmonic.rB;
2152 kk = L1*kA + lambda*kB;
2153 phi0 = (L1*pA + lambda*pB)*DEG2RAD;
2154 dphi0 = (pB - pA)*DEG2RAD;
2158 make_dp_periodic(&dp);
2165 dvdl_term += 0.5*(kB - kA)*dp2 - kk*dphi0*dp;
2167 do_dih_fup(ai, aj, ak, al, -ddphi, r_ij, r_kj, r_kl, m, n,
2168 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2173 fprintf(debug, "idih: (%d,%d,%d,%d) phi=%g\n",
2174 ai, aj, ak, al, phi);
2179 *dvdlambda += dvdl_term;
2184 /*! \brief returns dx, rdist, and dpdl for functions posres() and fbposres()
2186 static void posres_dx(const rvec x, const rvec pos0A, const rvec pos0B,
2187 const rvec comA_sc, const rvec comB_sc,
2189 t_pbc *pbc, int refcoord_scaling, int npbcdim,
2190 rvec dx, rvec rdist, rvec dpdl)
2193 real posA, posB, L1, ref = 0.;
2198 for (m = 0; m < DIM; m++)
2204 switch (refcoord_scaling)
2208 rdist[m] = L1*posA + lambda*posB;
2209 dpdl[m] = posB - posA;
2212 /* Box relative coordinates are stored for dimensions with pbc */
2213 posA *= pbc->box[m][m];
2214 posB *= pbc->box[m][m];
2215 assert(npbcdim <= DIM);
2216 for (d = m+1; d < npbcdim; d++)
2218 posA += pos0A[d]*pbc->box[d][m];
2219 posB += pos0B[d]*pbc->box[d][m];
2221 ref = L1*posA + lambda*posB;
2223 dpdl[m] = posB - posA;
2226 ref = L1*comA_sc[m] + lambda*comB_sc[m];
2227 rdist[m] = L1*posA + lambda*posB;
2228 dpdl[m] = comB_sc[m] - comA_sc[m] + posB - posA;
2231 gmx_fatal(FARGS, "No such scaling method implemented");
2236 ref = L1*posA + lambda*posB;
2238 dpdl[m] = posB - posA;
2241 /* We do pbc_dx with ref+rdist,
2242 * since with only ref we can be up to half a box vector wrong.
2244 pos[m] = ref + rdist[m];
2249 pbc_dx(pbc, x, pos, dx);
2253 rvec_sub(x, pos, dx);
2257 /*! \brief Computes forces and potential for flat-bottom cylindrical restraints.
2258 * Returns the flat-bottom potential. */
2259 static real do_fbposres_cylinder(int fbdim, rvec fm, rvec dx, real rfb, real kk, gmx_bool bInvert)
2262 real dr, dr2, invdr, v, rfb2;
2268 for (d = 0; d < DIM; d++)
2277 ( (dr2 > rfb2 && bInvert == FALSE ) || (dr2 < rfb2 && bInvert == TRUE ) )
2282 v = 0.5*kk*sqr(dr - rfb);
2283 for (d = 0; d < DIM; d++)
2287 fm[d] = -kk*(dr-rfb)*dx[d]*invdr; /* Force pointing to the center */
2295 /*! \brief Adds forces of flat-bottomed positions restraints to f[]
2296 * and fixes vir_diag. Returns the flat-bottomed potential. */
2297 real fbposres(int nbonds,
2298 const t_iatom forceatoms[], const t_iparams forceparams[],
2299 const rvec x[], rvec f[], rvec vir_diag,
2301 int refcoord_scaling, int ePBC, rvec com)
2302 /* compute flat-bottomed positions restraints */
2304 int i, ai, m, d, type, npbcdim = 0, fbdim;
2305 const t_iparams *pr;
2307 real dr, dr2, rfb, rfb2, fact;
2308 rvec com_sc, rdist, dx, dpdl, fm;
2311 npbcdim = ePBC2npbcdim(ePBC);
2313 if (refcoord_scaling == erscCOM)
2316 for (m = 0; m < npbcdim; m++)
2318 assert(npbcdim <= DIM);
2319 for (d = m; d < npbcdim; d++)
2321 com_sc[m] += com[d]*pbc->box[d][m];
2327 for (i = 0; (i < nbonds); )
2329 type = forceatoms[i++];
2330 ai = forceatoms[i++];
2331 pr = &forceparams[type];
2333 /* same calculation as for normal posres, but with identical A and B states, and lambda==0 */
2334 posres_dx(x[ai], forceparams[type].fbposres.pos0, forceparams[type].fbposres.pos0,
2335 com_sc, com_sc, 0.0,
2336 pbc, refcoord_scaling, npbcdim,
2342 kk = pr->fbposres.k;
2343 rfb = pr->fbposres.r;
2346 /* with rfb<0, push particle out of the sphere/cylinder/layer */
2354 switch (pr->fbposres.geom)
2356 case efbposresSPHERE:
2357 /* spherical flat-bottom posres */
2360 ( (dr2 > rfb2 && bInvert == FALSE ) || (dr2 < rfb2 && bInvert == TRUE ) )
2364 v = 0.5*kk*sqr(dr - rfb);
2365 fact = -kk*(dr-rfb)/dr; /* Force pointing to the center pos0 */
2366 svmul(fact, dx, fm);
2369 case efbposresCYLINDERX:
2370 /* cylindrical flat-bottom posres in y-z plane. fm[XX] = 0. */
2372 v = do_fbposres_cylinder(fbdim, fm, dx, rfb, kk, bInvert);
2374 case efbposresCYLINDERY:
2375 /* cylindrical flat-bottom posres in x-z plane. fm[YY] = 0. */
2377 v = do_fbposres_cylinder(fbdim, fm, dx, rfb, kk, bInvert);
2379 case efbposresCYLINDER:
2380 /* equivalent to efbposresCYLINDERZ for backwards compatibility */
2381 case efbposresCYLINDERZ:
2382 /* cylindrical flat-bottom posres in x-y plane. fm[ZZ] = 0. */
2384 v = do_fbposres_cylinder(fbdim, fm, dx, rfb, kk, bInvert);
2386 case efbposresX: /* fbdim=XX */
2387 case efbposresY: /* fbdim=YY */
2388 case efbposresZ: /* fbdim=ZZ */
2389 /* 1D flat-bottom potential */
2390 fbdim = pr->fbposres.geom - efbposresX;
2392 if ( ( dr > rfb && bInvert == FALSE ) || ( 0 < dr && dr < rfb && bInvert == TRUE ) )
2394 v = 0.5*kk*sqr(dr - rfb);
2395 fm[fbdim] = -kk*(dr - rfb);
2397 else if ( (dr < (-rfb) && bInvert == FALSE ) || ( (-rfb) < dr && dr < 0 && bInvert == TRUE ))
2399 v = 0.5*kk*sqr(dr + rfb);
2400 fm[fbdim] = -kk*(dr + rfb);
2407 for (m = 0; (m < DIM); m++)
2410 /* Here we correct for the pbc_dx which included rdist */
2411 vir_diag[m] -= 0.5*(dx[m] + rdist[m])*fm[m];
2419 real posres(int nbonds,
2420 const t_iatom forceatoms[], const t_iparams forceparams[],
2421 const rvec x[], rvec f[], rvec vir_diag,
2423 real lambda, real *dvdlambda,
2424 int refcoord_scaling, int ePBC, rvec comA, rvec comB)
2426 int i, ai, m, d, type, npbcdim = 0;
2427 const t_iparams *pr;
2430 rvec comA_sc, comB_sc, rdist, dpdl, dx;
2431 gmx_bool bForceValid = TRUE;
2433 if ((f == NULL) || (vir_diag == NULL)) /* should both be null together! */
2435 bForceValid = FALSE;
2438 npbcdim = ePBC2npbcdim(ePBC);
2440 if (refcoord_scaling == erscCOM)
2442 clear_rvec(comA_sc);
2443 clear_rvec(comB_sc);
2444 for (m = 0; m < npbcdim; m++)
2446 assert(npbcdim <= DIM);
2447 for (d = m; d < npbcdim; d++)
2449 comA_sc[m] += comA[d]*pbc->box[d][m];
2450 comB_sc[m] += comB[d]*pbc->box[d][m];
2458 for (i = 0; (i < nbonds); )
2460 type = forceatoms[i++];
2461 ai = forceatoms[i++];
2462 pr = &forceparams[type];
2464 /* return dx, rdist, and dpdl */
2465 posres_dx(x[ai], forceparams[type].posres.pos0A, forceparams[type].posres.pos0B,
2466 comA_sc, comB_sc, lambda,
2467 pbc, refcoord_scaling, npbcdim,
2470 for (m = 0; (m < DIM); m++)
2472 kk = L1*pr->posres.fcA[m] + lambda*pr->posres.fcB[m];
2474 vtot += 0.5*kk*dx[m]*dx[m];
2476 0.5*(pr->posres.fcB[m] - pr->posres.fcA[m])*dx[m]*dx[m]
2479 /* Here we correct for the pbc_dx which included rdist */
2483 vir_diag[m] -= 0.5*(dx[m] + rdist[m])*fm;
2491 static real low_angres(int nbonds,
2492 const t_iatom forceatoms[], const t_iparams forceparams[],
2493 const rvec x[], rvec f[], rvec fshift[],
2494 const t_pbc *pbc, const t_graph *g,
2495 real lambda, real *dvdlambda,
2498 int i, m, type, ai, aj, ak, al;
2500 real phi, cos_phi, cos_phi2, vid, vtot, dVdphi;
2501 rvec r_ij, r_kl, f_i, f_k = {0, 0, 0};
2502 real st, sth, nrij2, nrkl2, c, cij, ckl;
2505 t2 = 0; /* avoid warning with gcc-3.3. It is never used uninitialized */
2508 ak = al = 0; /* to avoid warnings */
2509 for (i = 0; i < nbonds; )
2511 type = forceatoms[i++];
2512 ai = forceatoms[i++];
2513 aj = forceatoms[i++];
2514 t1 = pbc_rvec_sub(pbc, x[aj], x[ai], r_ij); /* 3 */
2517 ak = forceatoms[i++];
2518 al = forceatoms[i++];
2519 t2 = pbc_rvec_sub(pbc, x[al], x[ak], r_kl); /* 3 */
2528 cos_phi = cos_angle(r_ij, r_kl); /* 25 */
2529 phi = acos(cos_phi); /* 10 */
2531 *dvdlambda += dopdihs_min(forceparams[type].pdihs.cpA,
2532 forceparams[type].pdihs.cpB,
2533 forceparams[type].pdihs.phiA,
2534 forceparams[type].pdihs.phiB,
2535 forceparams[type].pdihs.mult,
2536 phi, lambda, &vid, &dVdphi); /* 40 */
2540 cos_phi2 = sqr(cos_phi); /* 1 */
2543 st = -dVdphi*gmx_invsqrt(1 - cos_phi2); /* 12 */
2544 sth = st*cos_phi; /* 1 */
2545 nrij2 = iprod(r_ij, r_ij); /* 5 */
2546 nrkl2 = iprod(r_kl, r_kl); /* 5 */
2548 c = st*gmx_invsqrt(nrij2*nrkl2); /* 11 */
2549 cij = sth/nrij2; /* 10 */
2550 ckl = sth/nrkl2; /* 10 */
2552 for (m = 0; m < DIM; m++) /* 18+18 */
2554 f_i[m] = (c*r_kl[m]-cij*r_ij[m]);
2559 f_k[m] = (c*r_ij[m]-ckl*r_kl[m]);
2567 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
2570 rvec_inc(fshift[t1], f_i);
2571 rvec_dec(fshift[CENTRAL], f_i);
2576 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, al), dt);
2579 rvec_inc(fshift[t2], f_k);
2580 rvec_dec(fshift[CENTRAL], f_k);
2585 return vtot; /* 184 / 157 (bZAxis) total */
2588 real angres(int nbonds,
2589 const t_iatom forceatoms[], const t_iparams forceparams[],
2590 const rvec x[], rvec f[], rvec fshift[],
2591 const t_pbc *pbc, const t_graph *g,
2592 real lambda, real *dvdlambda,
2593 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2594 int gmx_unused *global_atom_index)
2596 return low_angres(nbonds, forceatoms, forceparams, x, f, fshift, pbc, g,
2597 lambda, dvdlambda, FALSE);
2600 real angresz(int nbonds,
2601 const t_iatom forceatoms[], const t_iparams forceparams[],
2602 const rvec x[], rvec f[], rvec fshift[],
2603 const t_pbc *pbc, const t_graph *g,
2604 real lambda, real *dvdlambda,
2605 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2606 int gmx_unused *global_atom_index)
2608 return low_angres(nbonds, forceatoms, forceparams, x, f, fshift, pbc, g,
2609 lambda, dvdlambda, TRUE);
2612 real dihres(int nbonds,
2613 const t_iatom forceatoms[], const t_iparams forceparams[],
2614 const rvec x[], rvec f[], rvec fshift[],
2615 const t_pbc *pbc, const t_graph *g,
2616 real lambda, real *dvdlambda,
2617 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2618 int gmx_unused *global_atom_index)
2621 int ai, aj, ak, al, i, k, type, t1, t2, t3;
2622 real phi0A, phi0B, dphiA, dphiB, kfacA, kfacB, phi0, dphi, kfac;
2623 real phi, ddphi, ddp, ddp2, dp, sign, d2r, L1;
2624 rvec r_ij, r_kj, r_kl, m, n;
2631 for (i = 0; (i < nbonds); )
2633 type = forceatoms[i++];
2634 ai = forceatoms[i++];
2635 aj = forceatoms[i++];
2636 ak = forceatoms[i++];
2637 al = forceatoms[i++];
2639 phi0A = forceparams[type].dihres.phiA*d2r;
2640 dphiA = forceparams[type].dihres.dphiA*d2r;
2641 kfacA = forceparams[type].dihres.kfacA;
2643 phi0B = forceparams[type].dihres.phiB*d2r;
2644 dphiB = forceparams[type].dihres.dphiB*d2r;
2645 kfacB = forceparams[type].dihres.kfacB;
2647 phi0 = L1*phi0A + lambda*phi0B;
2648 dphi = L1*dphiA + lambda*dphiB;
2649 kfac = L1*kfacA + lambda*kfacB;
2651 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2652 &sign, &t1, &t2, &t3);
2657 fprintf(debug, "dihres[%d]: %d %d %d %d : phi=%f, dphi=%f, kfac=%f\n",
2658 k++, ai, aj, ak, al, phi0, dphi, kfac);
2660 /* phi can jump if phi0 is close to Pi/-Pi, which will cause huge
2661 * force changes if we just apply a normal harmonic.
2662 * Instead, we first calculate phi-phi0 and take it modulo (-Pi,Pi).
2663 * This means we will never have the periodicity problem, unless
2664 * the dihedral is Pi away from phiO, which is very unlikely due to
2668 make_dp_periodic(&dp);
2674 else if (dp < -dphi)
2686 vtot += 0.5*kfac*ddp2;
2689 *dvdlambda += 0.5*(kfacB - kfacA)*ddp2;
2690 /* lambda dependence from changing restraint distances */
2693 *dvdlambda -= kfac*ddp*((dphiB - dphiA)+(phi0B - phi0A));
2697 *dvdlambda += kfac*ddp*((dphiB - dphiA)-(phi0B - phi0A));
2699 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
2700 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2707 real unimplemented(int gmx_unused nbonds,
2708 const t_iatom gmx_unused forceatoms[], const t_iparams gmx_unused forceparams[],
2709 const rvec gmx_unused x[], rvec gmx_unused f[], rvec gmx_unused fshift[],
2710 const t_pbc gmx_unused *pbc, const t_graph gmx_unused *g,
2711 real gmx_unused lambda, real gmx_unused *dvdlambda,
2712 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2713 int gmx_unused *global_atom_index)
2715 gmx_impl("*** you are using a not implemented function");
2717 return 0.0; /* To make the compiler happy */
2720 real restrangles(int nbonds,
2721 const t_iatom forceatoms[], const t_iparams forceparams[],
2722 const rvec x[], rvec f[], rvec fshift[],
2723 const t_pbc *pbc, const t_graph *g,
2724 real gmx_unused lambda, real gmx_unused *dvdlambda,
2725 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2726 int gmx_unused *global_atom_index)
2728 int i, d, ai, aj, ak, type, m;
2731 ivec jt, dt_ij, dt_kj;
2733 real prefactor, ratio_ante, ratio_post;
2734 rvec delta_ante, delta_post, vec_temp;
2737 for (i = 0; (i < nbonds); )
2739 type = forceatoms[i++];
2740 ai = forceatoms[i++];
2741 aj = forceatoms[i++];
2742 ak = forceatoms[i++];
2744 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], vec_temp);
2745 pbc_rvec_sub(pbc, x[aj], x[ai], delta_ante);
2746 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], delta_post);
2749 /* This function computes factors needed for restricted angle potential.
2750 * The restricted angle potential is used in coarse-grained simulations to avoid singularities
2751 * when three particles align and the dihedral angle and dihedral potential
2752 * cannot be calculated. This potential is calculated using the formula:
2753 real restrangles(int nbonds,
2754 const t_iatom forceatoms[],const t_iparams forceparams[],
2755 const rvec x[],rvec f[],rvec fshift[],
2756 const t_pbc *pbc,const t_graph *g,
2757 real gmx_unused lambda,real gmx_unused *dvdlambda,
2758 const t_mdatoms gmx_unused *md,t_fcdata gmx_unused *fcd,
2759 int gmx_unused *global_atom_index)
2761 int i, d, ai, aj, ak, type, m;
2765 ivec jt,dt_ij,dt_kj;
2767 real prefactor, ratio_ante, ratio_post;
2768 rvec delta_ante, delta_post, vec_temp;
2771 for(i=0; (i<nbonds); )
2773 type = forceatoms[i++];
2774 ai = forceatoms[i++];
2775 aj = forceatoms[i++];
2776 ak = forceatoms[i++];
2778 * \f[V_{\rm ReB}(\theta_i) = \frac{1}{2} k_{\theta} \frac{(\cos\theta_i - \cos\theta_0)^2}
2779 * {\sin^2\theta_i}\f] ({eq:ReB} and ref \cite{MonicaGoga2013} from the manual).
2780 * For more explanations see comments file "restcbt.h". */
2782 compute_factors_restangles(type, forceparams, delta_ante, delta_post,
2783 &prefactor, &ratio_ante, &ratio_post, &v);
2785 /* Forces are computed per component */
2786 for (d = 0; d < DIM; d++)
2788 f_i[d] = prefactor * (ratio_ante * delta_ante[d] - delta_post[d]);
2789 f_j[d] = prefactor * ((ratio_post + 1.0) * delta_post[d] - (ratio_ante + 1.0) * delta_ante[d]);
2790 f_k[d] = prefactor * (delta_ante[d] - ratio_post * delta_post[d]);
2793 /* Computation of potential energy */
2799 for (m = 0; (m < DIM); m++)
2808 copy_ivec(SHIFT_IVEC(g, aj), jt);
2809 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
2810 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
2811 t1 = IVEC2IS(dt_ij);
2812 t2 = IVEC2IS(dt_kj);
2815 rvec_inc(fshift[t1], f_i);
2816 rvec_inc(fshift[CENTRAL], f_j);
2817 rvec_inc(fshift[t2], f_k);
2823 real restrdihs(int nbonds,
2824 const t_iatom forceatoms[], const t_iparams forceparams[],
2825 const rvec x[], rvec f[], rvec fshift[],
2826 const t_pbc *pbc, const t_graph *g,
2827 real gmx_unused lambda, real gmx_unused *dvlambda,
2828 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2829 int gmx_unused *global_atom_index)
2831 int i, d, type, ai, aj, ak, al;
2832 rvec f_i, f_j, f_k, f_l;
2834 ivec jt, dt_ij, dt_kj, dt_lj;
2837 rvec delta_ante, delta_crnt, delta_post, vec_temp;
2838 real factor_phi_ai_ante, factor_phi_ai_crnt, factor_phi_ai_post;
2839 real factor_phi_aj_ante, factor_phi_aj_crnt, factor_phi_aj_post;
2840 real factor_phi_ak_ante, factor_phi_ak_crnt, factor_phi_ak_post;
2841 real factor_phi_al_ante, factor_phi_al_crnt, factor_phi_al_post;
2846 for (i = 0; (i < nbonds); )
2848 type = forceatoms[i++];
2849 ai = forceatoms[i++];
2850 aj = forceatoms[i++];
2851 ak = forceatoms[i++];
2852 al = forceatoms[i++];
2854 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], vec_temp);
2855 pbc_rvec_sub(pbc, x[aj], x[ai], delta_ante);
2856 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], delta_crnt);
2857 pbc_rvec_sub(pbc, x[ak], x[al], vec_temp);
2858 pbc_rvec_sub(pbc, x[al], x[ak], delta_post);
2860 /* This function computes factors needed for restricted angle potential.
2861 * The restricted angle potential is used in coarse-grained simulations to avoid singularities
2862 * when three particles align and the dihedral angle and dihedral potential cannot be calculated.
2863 * This potential is calculated using the formula:
2864 * \f[V_{\rm ReB}(\theta_i) = \frac{1}{2} k_{\theta}
2865 * \frac{(\cos\theta_i - \cos\theta_0)^2}{\sin^2\theta_i}\f]
2866 * ({eq:ReB} and ref \cite{MonicaGoga2013} from the manual).
2867 * For more explanations see comments file "restcbt.h" */
2869 compute_factors_restrdihs(type, forceparams,
2870 delta_ante, delta_crnt, delta_post,
2871 &factor_phi_ai_ante, &factor_phi_ai_crnt, &factor_phi_ai_post,
2872 &factor_phi_aj_ante, &factor_phi_aj_crnt, &factor_phi_aj_post,
2873 &factor_phi_ak_ante, &factor_phi_ak_crnt, &factor_phi_ak_post,
2874 &factor_phi_al_ante, &factor_phi_al_crnt, &factor_phi_al_post,
2875 &prefactor_phi, &v);
2878 /* Computation of forces per component */
2879 for (d = 0; d < DIM; d++)
2881 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]);
2882 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]);
2883 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]);
2884 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]);
2886 /* Computation of the energy */
2892 /* Updating the forces */
2894 rvec_inc(f[ai], f_i);
2895 rvec_inc(f[aj], f_j);
2896 rvec_inc(f[ak], f_k);
2897 rvec_inc(f[al], f_l);
2900 /* Updating the fshift forces for the pressure coupling */
2903 copy_ivec(SHIFT_IVEC(g, aj), jt);
2904 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
2905 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
2906 ivec_sub(SHIFT_IVEC(g, al), jt, dt_lj);
2907 t1 = IVEC2IS(dt_ij);
2908 t2 = IVEC2IS(dt_kj);
2909 t3 = IVEC2IS(dt_lj);
2913 t3 = pbc_rvec_sub(pbc, x[al], x[aj], dx_jl);
2920 rvec_inc(fshift[t1], f_i);
2921 rvec_inc(fshift[CENTRAL], f_j);
2922 rvec_inc(fshift[t2], f_k);
2923 rvec_inc(fshift[t3], f_l);
2931 real cbtdihs(int nbonds,
2932 const t_iatom forceatoms[], const t_iparams forceparams[],
2933 const rvec x[], rvec f[], rvec fshift[],
2934 const t_pbc *pbc, const t_graph *g,
2935 real gmx_unused lambda, real gmx_unused *dvdlambda,
2936 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2937 int gmx_unused *global_atom_index)
2939 int type, ai, aj, ak, al, i, d;
2943 rvec f_i, f_j, f_k, f_l;
2944 ivec jt, dt_ij, dt_kj, dt_lj;
2946 rvec delta_ante, delta_crnt, delta_post;
2947 rvec f_phi_ai, f_phi_aj, f_phi_ak, f_phi_al;
2948 rvec f_theta_ante_ai, f_theta_ante_aj, f_theta_ante_ak;
2949 rvec f_theta_post_aj, f_theta_post_ak, f_theta_post_al;
2955 for (i = 0; (i < nbonds); )
2957 type = forceatoms[i++];
2958 ai = forceatoms[i++];
2959 aj = forceatoms[i++];
2960 ak = forceatoms[i++];
2961 al = forceatoms[i++];
2964 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], vec_temp);
2965 pbc_rvec_sub(pbc, x[aj], x[ai], delta_ante);
2966 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], vec_temp);
2967 pbc_rvec_sub(pbc, x[ak], x[aj], delta_crnt);
2968 pbc_rvec_sub(pbc, x[ak], x[al], vec_temp);
2969 pbc_rvec_sub(pbc, x[al], x[ak], delta_post);
2971 /* \brief Compute factors for CBT potential
2972 * The combined bending-torsion potential goes to zero in a very smooth manner, eliminating the numerical
2973 * instabilities, when three coarse-grained particles align and the dihedral angle and standard
2974 * dihedral potentials cannot be calculated. The CBT potential is calculated using the formula:
2975 * \f[V_{\rm CBT}(\theta_{i-1}, \theta_i, \phi_i) = k_{\phi} \sin^3\theta_{i-1} \sin^3\theta_{i}
2976 * \sum_{n=0}^4 { a_n \cos^n\phi_i}\f] ({eq:CBT} and ref \cite{MonicaGoga2013} from the manual).
2977 * It contains in its expression not only the dihedral angle \f$\phi\f$
2978 * but also \f[\theta_{i-1}\f] (theta_ante bellow) and \f[\theta_{i}\f] (theta_post bellow)
2979 * --- the adjacent bending angles.
2980 * For more explanations see comments file "restcbt.h". */
2982 compute_factors_cbtdihs(type, forceparams, delta_ante, delta_crnt, delta_post,
2983 f_phi_ai, f_phi_aj, f_phi_ak, f_phi_al,
2984 f_theta_ante_ai, f_theta_ante_aj, f_theta_ante_ak,
2985 f_theta_post_aj, f_theta_post_ak, f_theta_post_al,
2989 /* Acumulate the resuts per beads */
2990 for (d = 0; d < DIM; d++)
2992 f_i[d] = f_phi_ai[d] + f_theta_ante_ai[d];
2993 f_j[d] = f_phi_aj[d] + f_theta_ante_aj[d] + f_theta_post_aj[d];
2994 f_k[d] = f_phi_ak[d] + f_theta_ante_ak[d] + f_theta_post_ak[d];
2995 f_l[d] = f_phi_al[d] + f_theta_post_al[d];
2998 /* Compute the potential energy */
3003 /* Updating the forces */
3004 rvec_inc(f[ai], f_i);
3005 rvec_inc(f[aj], f_j);
3006 rvec_inc(f[ak], f_k);
3007 rvec_inc(f[al], f_l);
3010 /* Updating the fshift forces for the pressure coupling */
3013 copy_ivec(SHIFT_IVEC(g, aj), jt);
3014 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3015 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3016 ivec_sub(SHIFT_IVEC(g, al), jt, dt_lj);
3017 t1 = IVEC2IS(dt_ij);
3018 t2 = IVEC2IS(dt_kj);
3019 t3 = IVEC2IS(dt_lj);
3023 t3 = pbc_rvec_sub(pbc, x[al], x[aj], dx_jl);
3030 rvec_inc(fshift[t1], f_i);
3031 rvec_inc(fshift[CENTRAL], f_j);
3032 rvec_inc(fshift[t2], f_k);
3033 rvec_inc(fshift[t3], f_l);
3039 real rbdihs(int nbonds,
3040 const t_iatom forceatoms[], const t_iparams forceparams[],
3041 const rvec x[], rvec f[], rvec fshift[],
3042 const t_pbc *pbc, const t_graph *g,
3043 real lambda, real *dvdlambda,
3044 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3045 int gmx_unused *global_atom_index)
3047 const real c0 = 0.0, c1 = 1.0, c2 = 2.0, c3 = 3.0, c4 = 4.0, c5 = 5.0;
3048 int type, ai, aj, ak, al, i, j;
3050 rvec r_ij, r_kj, r_kl, m, n;
3051 real parmA[NR_RBDIHS];
3052 real parmB[NR_RBDIHS];
3053 real parm[NR_RBDIHS];
3054 real cos_phi, phi, rbp, rbpBA;
3055 real v, sign, ddphi, sin_phi;
3057 real L1 = 1.0-lambda;
3061 for (i = 0; (i < nbonds); )
3063 type = forceatoms[i++];
3064 ai = forceatoms[i++];
3065 aj = forceatoms[i++];
3066 ak = forceatoms[i++];
3067 al = forceatoms[i++];
3069 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
3070 &sign, &t1, &t2, &t3); /* 84 */
3072 /* Change to polymer convention */
3079 phi -= M_PI; /* 1 */
3083 /* Beware of accuracy loss, cannot use 1-sqrt(cos^2) ! */
3086 for (j = 0; (j < NR_RBDIHS); j++)
3088 parmA[j] = forceparams[type].rbdihs.rbcA[j];
3089 parmB[j] = forceparams[type].rbdihs.rbcB[j];
3090 parm[j] = L1*parmA[j]+lambda*parmB[j];
3092 /* Calculate cosine powers */
3093 /* Calculate the energy */
3094 /* Calculate the derivative */
3097 dvdl_term += (parmB[0]-parmA[0]);
3102 rbpBA = parmB[1]-parmA[1];
3103 ddphi += rbp*cosfac;
3106 dvdl_term += cosfac*rbpBA;
3108 rbpBA = parmB[2]-parmA[2];
3109 ddphi += c2*rbp*cosfac;
3112 dvdl_term += cosfac*rbpBA;
3114 rbpBA = parmB[3]-parmA[3];
3115 ddphi += c3*rbp*cosfac;
3118 dvdl_term += cosfac*rbpBA;
3120 rbpBA = parmB[4]-parmA[4];
3121 ddphi += c4*rbp*cosfac;
3124 dvdl_term += cosfac*rbpBA;
3126 rbpBA = parmB[5]-parmA[5];
3127 ddphi += c5*rbp*cosfac;
3130 dvdl_term += cosfac*rbpBA;
3132 ddphi = -ddphi*sin_phi; /* 11 */
3134 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
3135 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
3138 *dvdlambda += dvdl_term;
3143 int cmap_setup_grid_index(int ip, int grid_spacing, int *ipm1, int *ipp1, int *ipp2)
3149 ip = ip + grid_spacing - 1;
3151 else if (ip > grid_spacing)
3153 ip = ip - grid_spacing - 1;
3162 im1 = grid_spacing - 1;
3164 else if (ip == grid_spacing-2)
3168 else if (ip == grid_spacing-1)
3182 real cmap_dihs(int nbonds,
3183 const t_iatom forceatoms[], const t_iparams forceparams[],
3184 const gmx_cmap_t *cmap_grid,
3185 const rvec x[], rvec f[], rvec fshift[],
3186 const t_pbc *pbc, const t_graph *g,
3187 real gmx_unused lambda, real gmx_unused *dvdlambda,
3188 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3189 int gmx_unused *global_atom_index)
3191 int i, j, k, n, idx;
3192 int ai, aj, ak, al, am;
3193 int a1i, a1j, a1k, a1l, a2i, a2j, a2k, a2l;
3195 int t11, t21, t31, t12, t22, t32;
3196 int iphi1, ip1m1, ip1p1, ip1p2;
3197 int iphi2, ip2m1, ip2p1, ip2p2;
3199 int pos1, pos2, pos3, pos4;
3201 real ty[4], ty1[4], ty2[4], ty12[4], tc[16], tx[16];
3202 real phi1, cos_phi1, sin_phi1, sign1, xphi1;
3203 real phi2, cos_phi2, sin_phi2, sign2, xphi2;
3204 real dx, xx, tt, tu, e, df1, df2, vtot;
3205 real ra21, rb21, rg21, rg1, rgr1, ra2r1, rb2r1, rabr1;
3206 real ra22, rb22, rg22, rg2, rgr2, ra2r2, rb2r2, rabr2;
3207 real fg1, hg1, fga1, hgb1, gaa1, gbb1;
3208 real fg2, hg2, fga2, hgb2, gaa2, gbb2;
3211 rvec r1_ij, r1_kj, r1_kl, m1, n1;
3212 rvec r2_ij, r2_kj, r2_kl, m2, n2;
3213 rvec f1_i, f1_j, f1_k, f1_l;
3214 rvec f2_i, f2_j, f2_k, f2_l;
3215 rvec a1, b1, a2, b2;
3216 rvec f1, g1, h1, f2, g2, h2;
3217 rvec dtf1, dtg1, dth1, dtf2, dtg2, dth2;
3218 ivec jt1, dt1_ij, dt1_kj, dt1_lj;
3219 ivec jt2, dt2_ij, dt2_kj, dt2_lj;
3223 int loop_index[4][4] = {
3230 /* Total CMAP energy */
3233 for (n = 0; n < nbonds; )
3235 /* Five atoms are involved in the two torsions */
3236 type = forceatoms[n++];
3237 ai = forceatoms[n++];
3238 aj = forceatoms[n++];
3239 ak = forceatoms[n++];
3240 al = forceatoms[n++];
3241 am = forceatoms[n++];
3243 /* Which CMAP type is this */
3244 cmapA = forceparams[type].cmap.cmapA;
3245 cmapd = cmap_grid->cmapdata[cmapA].cmap;
3253 phi1 = dih_angle(x[a1i], x[a1j], x[a1k], x[a1l], pbc, r1_ij, r1_kj, r1_kl, m1, n1,
3254 &sign1, &t11, &t21, &t31); /* 84 */
3256 cos_phi1 = cos(phi1);
3258 a1[0] = r1_ij[1]*r1_kj[2]-r1_ij[2]*r1_kj[1];
3259 a1[1] = r1_ij[2]*r1_kj[0]-r1_ij[0]*r1_kj[2];
3260 a1[2] = r1_ij[0]*r1_kj[1]-r1_ij[1]*r1_kj[0]; /* 9 */
3262 b1[0] = r1_kl[1]*r1_kj[2]-r1_kl[2]*r1_kj[1];
3263 b1[1] = r1_kl[2]*r1_kj[0]-r1_kl[0]*r1_kj[2];
3264 b1[2] = r1_kl[0]*r1_kj[1]-r1_kl[1]*r1_kj[0]; /* 9 */
3266 pbc_rvec_sub(pbc, x[a1l], x[a1k], h1);
3268 ra21 = iprod(a1, a1); /* 5 */
3269 rb21 = iprod(b1, b1); /* 5 */
3270 rg21 = iprod(r1_kj, r1_kj); /* 5 */
3276 rabr1 = sqrt(ra2r1*rb2r1);
3278 sin_phi1 = rg1 * rabr1 * iprod(a1, h1) * (-1);
3280 if (cos_phi1 < -0.5 || cos_phi1 > 0.5)
3282 phi1 = asin(sin_phi1);
3292 phi1 = -M_PI - phi1;
3298 phi1 = acos(cos_phi1);
3306 xphi1 = phi1 + M_PI; /* 1 */
3308 /* Second torsion */
3314 phi2 = dih_angle(x[a2i], x[a2j], x[a2k], x[a2l], pbc, r2_ij, r2_kj, r2_kl, m2, n2,
3315 &sign2, &t12, &t22, &t32); /* 84 */
3317 cos_phi2 = cos(phi2);
3319 a2[0] = r2_ij[1]*r2_kj[2]-r2_ij[2]*r2_kj[1];
3320 a2[1] = r2_ij[2]*r2_kj[0]-r2_ij[0]*r2_kj[2];
3321 a2[2] = r2_ij[0]*r2_kj[1]-r2_ij[1]*r2_kj[0]; /* 9 */
3323 b2[0] = r2_kl[1]*r2_kj[2]-r2_kl[2]*r2_kj[1];
3324 b2[1] = r2_kl[2]*r2_kj[0]-r2_kl[0]*r2_kj[2];
3325 b2[2] = r2_kl[0]*r2_kj[1]-r2_kl[1]*r2_kj[0]; /* 9 */
3327 pbc_rvec_sub(pbc, x[a2l], x[a2k], h2);
3329 ra22 = iprod(a2, a2); /* 5 */
3330 rb22 = iprod(b2, b2); /* 5 */
3331 rg22 = iprod(r2_kj, r2_kj); /* 5 */
3337 rabr2 = sqrt(ra2r2*rb2r2);
3339 sin_phi2 = rg2 * rabr2 * iprod(a2, h2) * (-1);
3341 if (cos_phi2 < -0.5 || cos_phi2 > 0.5)
3343 phi2 = asin(sin_phi2);
3353 phi2 = -M_PI - phi2;
3359 phi2 = acos(cos_phi2);
3367 xphi2 = phi2 + M_PI; /* 1 */
3369 /* Range mangling */
3372 xphi1 = xphi1 + 2*M_PI;
3374 else if (xphi1 >= 2*M_PI)
3376 xphi1 = xphi1 - 2*M_PI;
3381 xphi2 = xphi2 + 2*M_PI;
3383 else if (xphi2 >= 2*M_PI)
3385 xphi2 = xphi2 - 2*M_PI;
3388 /* Number of grid points */
3389 dx = 2*M_PI / cmap_grid->grid_spacing;
3391 /* Where on the grid are we */
3392 iphi1 = static_cast<int>(xphi1/dx);
3393 iphi2 = static_cast<int>(xphi2/dx);
3395 iphi1 = cmap_setup_grid_index(iphi1, cmap_grid->grid_spacing, &ip1m1, &ip1p1, &ip1p2);
3396 iphi2 = cmap_setup_grid_index(iphi2, cmap_grid->grid_spacing, &ip2m1, &ip2p1, &ip2p2);
3398 pos1 = iphi1*cmap_grid->grid_spacing+iphi2;
3399 pos2 = ip1p1*cmap_grid->grid_spacing+iphi2;
3400 pos3 = ip1p1*cmap_grid->grid_spacing+ip2p1;
3401 pos4 = iphi1*cmap_grid->grid_spacing+ip2p1;
3403 ty[0] = cmapd[pos1*4];
3404 ty[1] = cmapd[pos2*4];
3405 ty[2] = cmapd[pos3*4];
3406 ty[3] = cmapd[pos4*4];
3408 ty1[0] = cmapd[pos1*4+1];
3409 ty1[1] = cmapd[pos2*4+1];
3410 ty1[2] = cmapd[pos3*4+1];
3411 ty1[3] = cmapd[pos4*4+1];
3413 ty2[0] = cmapd[pos1*4+2];
3414 ty2[1] = cmapd[pos2*4+2];
3415 ty2[2] = cmapd[pos3*4+2];
3416 ty2[3] = cmapd[pos4*4+2];
3418 ty12[0] = cmapd[pos1*4+3];
3419 ty12[1] = cmapd[pos2*4+3];
3420 ty12[2] = cmapd[pos3*4+3];
3421 ty12[3] = cmapd[pos4*4+3];
3423 /* Switch to degrees */
3424 dx = 360.0 / cmap_grid->grid_spacing;
3425 xphi1 = xphi1 * RAD2DEG;
3426 xphi2 = xphi2 * RAD2DEG;
3428 for (i = 0; i < 4; i++) /* 16 */
3431 tx[i+4] = ty1[i]*dx;
3432 tx[i+8] = ty2[i]*dx;
3433 tx[i+12] = ty12[i]*dx*dx;
3437 for (i = 0; i < 4; i++) /* 1056 */
3439 for (j = 0; j < 4; j++)
3442 for (k = 0; k < 16; k++)
3444 xx = xx + cmap_coeff_matrix[k*16+idx]*tx[k];
3452 tt = (xphi1-iphi1*dx)/dx;
3453 tu = (xphi2-iphi2*dx)/dx;
3459 for (i = 3; i >= 0; i--)
3461 l1 = loop_index[i][3];
3462 l2 = loop_index[i][2];
3463 l3 = loop_index[i][1];
3465 e = tt * e + ((tc[i*4+3]*tu+tc[i*4+2])*tu + tc[i*4+1])*tu+tc[i*4];
3466 df1 = tu * df1 + (3.0*tc[l1]*tt+2.0*tc[l2])*tt+tc[l3];
3467 df2 = tt * df2 + (3.0*tc[i*4+3]*tu+2.0*tc[i*4+2])*tu+tc[i*4+1];
3477 /* Do forces - first torsion */
3478 fg1 = iprod(r1_ij, r1_kj);
3479 hg1 = iprod(r1_kl, r1_kj);
3480 fga1 = fg1*ra2r1*rgr1;
3481 hgb1 = hg1*rb2r1*rgr1;
3485 for (i = 0; i < DIM; i++)
3487 dtf1[i] = gaa1 * a1[i];
3488 dtg1[i] = fga1 * a1[i] - hgb1 * b1[i];
3489 dth1[i] = gbb1 * b1[i];
3491 f1[i] = df1 * dtf1[i];
3492 g1[i] = df1 * dtg1[i];
3493 h1[i] = df1 * dth1[i];
3496 f1_j[i] = -f1[i] - g1[i];
3497 f1_k[i] = h1[i] + g1[i];
3500 f[a1i][i] = f[a1i][i] + f1_i[i];
3501 f[a1j][i] = f[a1j][i] + f1_j[i]; /* - f1[i] - g1[i] */
3502 f[a1k][i] = f[a1k][i] + f1_k[i]; /* h1[i] + g1[i] */
3503 f[a1l][i] = f[a1l][i] + f1_l[i]; /* h1[i] */
3506 /* Do forces - second torsion */
3507 fg2 = iprod(r2_ij, r2_kj);
3508 hg2 = iprod(r2_kl, r2_kj);
3509 fga2 = fg2*ra2r2*rgr2;
3510 hgb2 = hg2*rb2r2*rgr2;
3514 for (i = 0; i < DIM; i++)
3516 dtf2[i] = gaa2 * a2[i];
3517 dtg2[i] = fga2 * a2[i] - hgb2 * b2[i];
3518 dth2[i] = gbb2 * b2[i];
3520 f2[i] = df2 * dtf2[i];
3521 g2[i] = df2 * dtg2[i];
3522 h2[i] = df2 * dth2[i];
3525 f2_j[i] = -f2[i] - g2[i];
3526 f2_k[i] = h2[i] + g2[i];
3529 f[a2i][i] = f[a2i][i] + f2_i[i]; /* f2[i] */
3530 f[a2j][i] = f[a2j][i] + f2_j[i]; /* - f2[i] - g2[i] */
3531 f[a2k][i] = f[a2k][i] + f2_k[i]; /* h2[i] + g2[i] */
3532 f[a2l][i] = f[a2l][i] + f2_l[i]; /* - h2[i] */
3538 copy_ivec(SHIFT_IVEC(g, a1j), jt1);
3539 ivec_sub(SHIFT_IVEC(g, a1i), jt1, dt1_ij);
3540 ivec_sub(SHIFT_IVEC(g, a1k), jt1, dt1_kj);
3541 ivec_sub(SHIFT_IVEC(g, a1l), jt1, dt1_lj);
3542 t11 = IVEC2IS(dt1_ij);
3543 t21 = IVEC2IS(dt1_kj);
3544 t31 = IVEC2IS(dt1_lj);
3546 copy_ivec(SHIFT_IVEC(g, a2j), jt2);
3547 ivec_sub(SHIFT_IVEC(g, a2i), jt2, dt2_ij);
3548 ivec_sub(SHIFT_IVEC(g, a2k), jt2, dt2_kj);
3549 ivec_sub(SHIFT_IVEC(g, a2l), jt2, dt2_lj);
3550 t12 = IVEC2IS(dt2_ij);
3551 t22 = IVEC2IS(dt2_kj);
3552 t32 = IVEC2IS(dt2_lj);
3556 t31 = pbc_rvec_sub(pbc, x[a1l], x[a1j], h1);
3557 t32 = pbc_rvec_sub(pbc, x[a2l], x[a2j], h2);
3565 rvec_inc(fshift[t11], f1_i);
3566 rvec_inc(fshift[CENTRAL], f1_j);
3567 rvec_inc(fshift[t21], f1_k);
3568 rvec_inc(fshift[t31], f1_l);
3570 rvec_inc(fshift[t21], f2_i);
3571 rvec_inc(fshift[CENTRAL], f2_j);
3572 rvec_inc(fshift[t22], f2_k);
3573 rvec_inc(fshift[t32], f2_l);
3580 /***********************************************************
3582 * G R O M O S 9 6 F U N C T I O N S
3584 ***********************************************************/
3585 real g96harmonic(real kA, real kB, real xA, real xB, real x, real lambda,
3588 const real half = 0.5;
3589 real L1, kk, x0, dx, dx2;
3590 real v, f, dvdlambda;
3593 kk = L1*kA+lambda*kB;
3594 x0 = L1*xA+lambda*xB;
3601 dvdlambda = half*(kB-kA)*dx2 + (xA-xB)*kk*dx;
3608 /* That was 21 flops */
3611 real g96bonds(int nbonds,
3612 const t_iatom forceatoms[], const t_iparams forceparams[],
3613 const rvec x[], rvec f[], rvec fshift[],
3614 const t_pbc *pbc, const t_graph *g,
3615 real lambda, real *dvdlambda,
3616 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3617 int gmx_unused *global_atom_index)
3619 int i, m, ki, ai, aj, type;
3620 real dr2, fbond, vbond, fij, vtot;
3625 for (i = 0; (i < nbonds); )
3627 type = forceatoms[i++];
3628 ai = forceatoms[i++];
3629 aj = forceatoms[i++];
3631 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
3632 dr2 = iprod(dx, dx); /* 5 */
3634 *dvdlambda += g96harmonic(forceparams[type].harmonic.krA,
3635 forceparams[type].harmonic.krB,
3636 forceparams[type].harmonic.rA,
3637 forceparams[type].harmonic.rB,
3638 dr2, lambda, &vbond, &fbond);
3640 vtot += 0.5*vbond; /* 1*/
3644 fprintf(debug, "G96-BONDS: dr = %10g vbond = %10g fbond = %10g\n",
3645 sqrt(dr2), vbond, fbond);
3651 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
3654 for (m = 0; (m < DIM); m++) /* 15 */
3659 fshift[ki][m] += fij;
3660 fshift[CENTRAL][m] -= fij;
3666 real g96bond_angle(const rvec xi, const rvec xj, const rvec xk, const t_pbc *pbc,
3667 rvec r_ij, rvec r_kj,
3669 /* Return value is the angle between the bonds i-j and j-k */
3673 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
3674 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
3676 costh = cos_angle(r_ij, r_kj); /* 25 */
3681 real g96angles(int nbonds,
3682 const t_iatom forceatoms[], const t_iparams forceparams[],
3683 const rvec x[], rvec f[], rvec fshift[],
3684 const t_pbc *pbc, const t_graph *g,
3685 real lambda, real *dvdlambda,
3686 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3687 int gmx_unused *global_atom_index)
3689 int i, ai, aj, ak, type, m, t1, t2;
3691 real cos_theta, dVdt, va, vtot;
3692 real rij_1, rij_2, rkj_1, rkj_2, rijrkj_1;
3694 ivec jt, dt_ij, dt_kj;
3697 for (i = 0; (i < nbonds); )
3699 type = forceatoms[i++];
3700 ai = forceatoms[i++];
3701 aj = forceatoms[i++];
3702 ak = forceatoms[i++];
3704 cos_theta = g96bond_angle(x[ai], x[aj], x[ak], pbc, r_ij, r_kj, &t1, &t2);
3706 *dvdlambda += g96harmonic(forceparams[type].harmonic.krA,
3707 forceparams[type].harmonic.krB,
3708 forceparams[type].harmonic.rA,
3709 forceparams[type].harmonic.rB,
3710 cos_theta, lambda, &va, &dVdt);
3713 rij_1 = gmx_invsqrt(iprod(r_ij, r_ij));
3714 rkj_1 = gmx_invsqrt(iprod(r_kj, r_kj));
3715 rij_2 = rij_1*rij_1;
3716 rkj_2 = rkj_1*rkj_1;
3717 rijrkj_1 = rij_1*rkj_1; /* 23 */
3722 fprintf(debug, "G96ANGLES: costheta = %10g vth = %10g dV/dct = %10g\n",
3723 cos_theta, va, dVdt);
3726 for (m = 0; (m < DIM); m++) /* 42 */
3728 f_i[m] = dVdt*(r_kj[m]*rijrkj_1 - r_ij[m]*rij_2*cos_theta);
3729 f_k[m] = dVdt*(r_ij[m]*rijrkj_1 - r_kj[m]*rkj_2*cos_theta);
3730 f_j[m] = -f_i[m]-f_k[m];
3738 copy_ivec(SHIFT_IVEC(g, aj), jt);
3740 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3741 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3742 t1 = IVEC2IS(dt_ij);
3743 t2 = IVEC2IS(dt_kj);
3745 rvec_inc(fshift[t1], f_i);
3746 rvec_inc(fshift[CENTRAL], f_j);
3747 rvec_inc(fshift[t2], f_k); /* 9 */
3753 real cross_bond_bond(int nbonds,
3754 const t_iatom forceatoms[], const t_iparams forceparams[],
3755 const rvec x[], rvec f[], rvec fshift[],
3756 const t_pbc *pbc, const t_graph *g,
3757 real gmx_unused lambda, real gmx_unused *dvdlambda,
3758 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3759 int gmx_unused *global_atom_index)
3761 /* Potential from Lawrence and Skimmer, Chem. Phys. Lett. 372 (2003)
3764 int i, ai, aj, ak, type, m, t1, t2;
3766 real vtot, vrr, s1, s2, r1, r2, r1e, r2e, krr;
3768 ivec jt, dt_ij, dt_kj;
3771 for (i = 0; (i < nbonds); )
3773 type = forceatoms[i++];
3774 ai = forceatoms[i++];
3775 aj = forceatoms[i++];
3776 ak = forceatoms[i++];
3777 r1e = forceparams[type].cross_bb.r1e;
3778 r2e = forceparams[type].cross_bb.r2e;
3779 krr = forceparams[type].cross_bb.krr;
3781 /* Compute distance vectors ... */
3782 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
3783 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
3785 /* ... and their lengths */
3789 /* Deviations from ideality */
3793 /* Energy (can be negative!) */
3798 svmul(-krr*s2/r1, r_ij, f_i);
3799 svmul(-krr*s1/r2, r_kj, f_k);
3801 for (m = 0; (m < DIM); m++) /* 12 */
3803 f_j[m] = -f_i[m] - f_k[m];
3812 copy_ivec(SHIFT_IVEC(g, aj), jt);
3814 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3815 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3816 t1 = IVEC2IS(dt_ij);
3817 t2 = IVEC2IS(dt_kj);
3819 rvec_inc(fshift[t1], f_i);
3820 rvec_inc(fshift[CENTRAL], f_j);
3821 rvec_inc(fshift[t2], f_k); /* 9 */
3827 real cross_bond_angle(int nbonds,
3828 const t_iatom forceatoms[], const t_iparams forceparams[],
3829 const rvec x[], rvec f[], rvec fshift[],
3830 const t_pbc *pbc, const t_graph *g,
3831 real gmx_unused lambda, real gmx_unused *dvdlambda,
3832 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3833 int gmx_unused *global_atom_index)
3835 /* Potential from Lawrence and Skimmer, Chem. Phys. Lett. 372 (2003)
3838 int i, ai, aj, ak, type, m, t1, t2;
3839 rvec r_ij, r_kj, r_ik;
3840 real vtot, vrt, s1, s2, s3, r1, r2, r3, r1e, r2e, r3e, krt, k1, k2, k3;
3842 ivec jt, dt_ij, dt_kj;
3845 for (i = 0; (i < nbonds); )
3847 type = forceatoms[i++];
3848 ai = forceatoms[i++];
3849 aj = forceatoms[i++];
3850 ak = forceatoms[i++];
3851 r1e = forceparams[type].cross_ba.r1e;
3852 r2e = forceparams[type].cross_ba.r2e;
3853 r3e = forceparams[type].cross_ba.r3e;
3854 krt = forceparams[type].cross_ba.krt;
3856 /* Compute distance vectors ... */
3857 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
3858 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
3859 pbc_rvec_sub(pbc, x[ai], x[ak], r_ik);
3861 /* ... and their lengths */
3866 /* Deviations from ideality */
3871 /* Energy (can be negative!) */
3872 vrt = krt*s3*(s1+s2);
3878 k3 = -krt*(s1+s2)/r3;
3879 for (m = 0; (m < DIM); m++)
3881 f_i[m] = k1*r_ij[m] + k3*r_ik[m];
3882 f_k[m] = k2*r_kj[m] - k3*r_ik[m];
3883 f_j[m] = -f_i[m] - f_k[m];
3886 for (m = 0; (m < DIM); m++) /* 12 */
3896 copy_ivec(SHIFT_IVEC(g, aj), jt);
3898 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3899 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3900 t1 = IVEC2IS(dt_ij);
3901 t2 = IVEC2IS(dt_kj);
3903 rvec_inc(fshift[t1], f_i);
3904 rvec_inc(fshift[CENTRAL], f_j);
3905 rvec_inc(fshift[t2], f_k); /* 9 */
3911 static real bonded_tab(const char *type, int table_nr,
3912 const bondedtable_t *table, real kA, real kB, real r,
3913 real lambda, real *V, real *F)
3915 real k, tabscale, *VFtab, rt, eps, eps2, Yt, Ft, Geps, Heps2, Fp, VV, FF;
3919 k = (1.0 - lambda)*kA + lambda*kB;
3921 tabscale = table->scale;
3922 VFtab = table->data;
3925 n0 = static_cast<int>(rt);
3928 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",
3929 type, table_nr, r, n0, n0+1, table->n);
3936 Geps = VFtab[nnn+2]*eps;
3937 Heps2 = VFtab[nnn+3]*eps2;
3938 Fp = Ft + Geps + Heps2;
3940 FF = Fp + Geps + 2.0*Heps2;
3942 *F = -k*FF*tabscale;
3944 dvdlambda = (kB - kA)*VV;
3948 /* That was 22 flops */
3951 real tab_bonds(int nbonds,
3952 const t_iatom forceatoms[], const t_iparams forceparams[],
3953 const rvec x[], rvec f[], rvec fshift[],
3954 const t_pbc *pbc, const t_graph *g,
3955 real lambda, real *dvdlambda,
3956 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
3957 int gmx_unused *global_atom_index)
3959 int i, m, ki, ai, aj, type, table;
3960 real dr, dr2, fbond, vbond, fij, vtot;
3965 for (i = 0; (i < nbonds); )
3967 type = forceatoms[i++];
3968 ai = forceatoms[i++];
3969 aj = forceatoms[i++];
3971 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
3972 dr2 = iprod(dx, dx); /* 5 */
3973 dr = dr2*gmx_invsqrt(dr2); /* 10 */
3975 table = forceparams[type].tab.table;
3977 *dvdlambda += bonded_tab("bond", table,
3978 &fcd->bondtab[table],
3979 forceparams[type].tab.kA,
3980 forceparams[type].tab.kB,
3981 dr, lambda, &vbond, &fbond); /* 22 */
3989 vtot += vbond; /* 1*/
3990 fbond *= gmx_invsqrt(dr2); /* 6 */
3994 fprintf(debug, "TABBONDS: dr = %10g vbond = %10g fbond = %10g\n",
4000 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
4003 for (m = 0; (m < DIM); m++) /* 15 */
4008 fshift[ki][m] += fij;
4009 fshift[CENTRAL][m] -= fij;
4015 real tab_angles(int nbonds,
4016 const t_iatom forceatoms[], const t_iparams forceparams[],
4017 const rvec x[], rvec f[], rvec fshift[],
4018 const t_pbc *pbc, const t_graph *g,
4019 real lambda, real *dvdlambda,
4020 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
4021 int gmx_unused *global_atom_index)
4023 int i, ai, aj, ak, t1, t2, type, table;
4025 real cos_theta, cos_theta2, theta, dVdt, va, vtot;
4026 ivec jt, dt_ij, dt_kj;
4029 for (i = 0; (i < nbonds); )
4031 type = forceatoms[i++];
4032 ai = forceatoms[i++];
4033 aj = forceatoms[i++];
4034 ak = forceatoms[i++];
4036 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
4037 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
4039 table = forceparams[type].tab.table;
4041 *dvdlambda += bonded_tab("angle", table,
4042 &fcd->angletab[table],
4043 forceparams[type].tab.kA,
4044 forceparams[type].tab.kB,
4045 theta, lambda, &va, &dVdt); /* 22 */
4048 cos_theta2 = sqr(cos_theta); /* 1 */
4057 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
4058 sth = st*cos_theta; /* 1 */
4062 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
4063 theta*RAD2DEG, va, dVdt);
4066 nrkj2 = iprod(r_kj, r_kj); /* 5 */
4067 nrij2 = iprod(r_ij, r_ij);
4069 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
4070 cii = sth/nrij2; /* 10 */
4071 ckk = sth/nrkj2; /* 10 */
4073 for (m = 0; (m < DIM); m++) /* 39 */
4075 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
4076 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
4077 f_j[m] = -f_i[m]-f_k[m];
4084 copy_ivec(SHIFT_IVEC(g, aj), jt);
4086 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
4087 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
4088 t1 = IVEC2IS(dt_ij);
4089 t2 = IVEC2IS(dt_kj);
4091 rvec_inc(fshift[t1], f_i);
4092 rvec_inc(fshift[CENTRAL], f_j);
4093 rvec_inc(fshift[t2], f_k);
4099 real tab_dihs(int nbonds,
4100 const t_iatom forceatoms[], const t_iparams forceparams[],
4101 const rvec x[], rvec f[], rvec fshift[],
4102 const t_pbc *pbc, const t_graph *g,
4103 real lambda, real *dvdlambda,
4104 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
4105 int gmx_unused *global_atom_index)
4107 int i, type, ai, aj, ak, al, table;
4109 rvec r_ij, r_kj, r_kl, m, n;
4110 real phi, sign, ddphi, vpd, vtot;
4113 for (i = 0; (i < nbonds); )
4115 type = forceatoms[i++];
4116 ai = forceatoms[i++];
4117 aj = forceatoms[i++];
4118 ak = forceatoms[i++];
4119 al = forceatoms[i++];
4121 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
4122 &sign, &t1, &t2, &t3); /* 84 */
4124 table = forceparams[type].tab.table;
4126 /* Hopefully phi+M_PI never results in values < 0 */
4127 *dvdlambda += bonded_tab("dihedral", table,
4128 &fcd->dihtab[table],
4129 forceparams[type].tab.kA,
4130 forceparams[type].tab.kB,
4131 phi+M_PI, lambda, &vpd, &ddphi);
4134 do_dih_fup(ai, aj, ak, al, -ddphi, r_ij, r_kj, r_kl, m, n,
4135 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
4138 fprintf(debug, "pdih: (%d,%d,%d,%d) phi=%g\n",
4139 ai, aj, ak, al, phi);
4146 /* TODO This function could go away when idef is not a big bucket of
4148 gmx_bool ftype_is_bonded_potential(int ftype)
4151 (interaction_function[ftype].flags & IF_BOND) &&
4152 !(ftype == F_CONNBONDS || ftype == F_POSRES || ftype == F_FBPOSRES) &&
4153 (ftype < F_GB12 || ftype > F_GB14);
4156 static void zero_thread_forces(f_thread_t *f_t, int n,
4157 int nblock, int blocksize)
4159 int b, a0, a1, a, i, j;
4161 if (n > f_t->f_nalloc)
4163 f_t->f_nalloc = over_alloc_large(n);
4164 srenew(f_t->f, f_t->f_nalloc);
4167 if (f_t->red_mask != 0)
4169 for (b = 0; b < nblock; b++)
4171 if (f_t->red_mask && (1U<<b))
4174 a1 = std::min((b+1)*blocksize, n);
4175 for (a = a0; a < a1; a++)
4177 clear_rvec(f_t->f[a]);
4182 for (i = 0; i < SHIFTS; i++)
4184 clear_rvec(f_t->fshift[i]);
4186 for (i = 0; i < F_NRE; i++)
4190 for (i = 0; i < egNR; i++)
4192 for (j = 0; j < f_t->grpp.nener; j++)
4194 f_t->grpp.ener[i][j] = 0;
4197 for (i = 0; i < efptNR; i++)
4203 static void reduce_thread_force_buffer(int n, rvec *f,
4204 int nthreads, f_thread_t *f_t,
4205 int nblock, int block_size)
4207 /* The max thread number is arbitrary,
4208 * we used a fixed number to avoid memory management.
4209 * Using more than 16 threads is probably never useful performance wise.
4211 #define MAX_BONDED_THREADS 256
4214 if (nthreads > MAX_BONDED_THREADS)
4216 gmx_fatal(FARGS, "Can not reduce bonded forces on more than %d threads",
4217 MAX_BONDED_THREADS);
4220 /* This reduction can run on any number of threads,
4221 * independently of nthreads.
4223 #pragma omp parallel for num_threads(nthreads) schedule(static)
4224 for (b = 0; b < nblock; b++)
4226 rvec *fp[MAX_BONDED_THREADS];
4230 /* Determine which threads contribute to this block */
4232 for (ft = 1; ft < nthreads; ft++)
4234 if (f_t[ft].red_mask & (1U<<b))
4236 fp[nfb++] = f_t[ft].f;
4241 /* Reduce force buffers for threads that contribute */
4243 a1 = (b+1)*block_size;
4244 a1 = std::min(a1, n);
4245 for (a = a0; a < a1; a++)
4247 for (fb = 0; fb < nfb; fb++)
4249 rvec_inc(f[a], fp[fb][a]);
4256 static void reduce_thread_forces(int n, rvec *f, rvec *fshift,
4257 real *ener, gmx_grppairener_t *grpp, real *dvdl,
4258 int nthreads, f_thread_t *f_t,
4259 int nblock, int block_size,
4260 gmx_bool bCalcEnerVir,
4265 /* Reduce the bonded force buffer */
4266 reduce_thread_force_buffer(n, f, nthreads, f_t, nblock, block_size);
4269 /* When necessary, reduce energy and virial using one thread only */
4274 for (i = 0; i < SHIFTS; i++)
4276 for (t = 1; t < nthreads; t++)
4278 rvec_inc(fshift[i], f_t[t].fshift[i]);
4281 for (i = 0; i < F_NRE; i++)
4283 for (t = 1; t < nthreads; t++)
4285 ener[i] += f_t[t].ener[i];
4288 for (i = 0; i < egNR; i++)
4290 for (j = 0; j < f_t[1].grpp.nener; j++)
4292 for (t = 1; t < nthreads; t++)
4295 grpp->ener[i][j] += f_t[t].grpp.ener[i][j];
4301 for (i = 0; i < efptNR; i++)
4304 for (t = 1; t < nthreads; t++)
4306 dvdl[i] += f_t[t].dvdl[i];
4313 static real calc_one_bond(int thread,
4314 int ftype, const t_idef *idef,
4315 const rvec x[], rvec f[], rvec fshift[],
4317 const t_pbc *pbc, const t_graph *g,
4318 gmx_grppairener_t *grpp,
4320 real *lambda, real *dvdl,
4321 const t_mdatoms *md, t_fcdata *fcd,
4322 gmx_bool bCalcEnerVir,
4323 int *global_atom_index)
4325 int nat1, nbonds, efptFTYPE;
4330 if (IS_RESTRAINT_TYPE(ftype))
4332 efptFTYPE = efptRESTRAINT;
4336 efptFTYPE = efptBONDED;
4339 nat1 = interaction_function[ftype].nratoms + 1;
4340 nbonds = idef->il[ftype].nr/nat1;
4341 iatoms = idef->il[ftype].iatoms;
4343 nb0 = idef->il_thread_division[ftype*(idef->nthreads+1)+thread];
4344 nbn = idef->il_thread_division[ftype*(idef->nthreads+1)+thread+1] - nb0;
4346 if (!IS_LISTED_LJ_C(ftype))
4348 if (ftype == F_CMAP)
4350 v = cmap_dihs(nbn, iatoms+nb0,
4351 idef->iparams, &idef->cmap_grid,
4353 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
4354 md, fcd, global_atom_index);
4356 #ifdef GMX_SIMD_HAVE_REAL
4357 else if (ftype == F_ANGLES &&
4358 !bCalcEnerVir && fr->efep == efepNO)
4360 /* No energies, shift forces, dvdl */
4361 angles_noener_simd(nbn, idef->il[ftype].iatoms+nb0,
4364 pbc, g, lambda[efptFTYPE], md, fcd,
4369 else if (ftype == F_PDIHS &&
4370 !bCalcEnerVir && fr->efep == efepNO)
4372 /* No energies, shift forces, dvdl */
4373 #ifdef GMX_SIMD_HAVE_REAL
4378 (nbn, idef->il[ftype].iatoms+nb0,
4381 pbc, g, lambda[efptFTYPE], md, fcd,
4387 v = interaction_function[ftype].ifunc(nbn, iatoms+nb0,
4390 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
4391 md, fcd, global_atom_index);
4396 v = do_nonbonded_listed(ftype, nbn, iatoms+nb0, idef->iparams, x, f, fshift,
4397 pbc, g, lambda, dvdl, md, fr, grpp, global_atom_index);
4402 inc_nrnb(nrnb, interaction_function[ftype].nrnb_ind, nbonds);
4408 void calc_bonds(const gmx_multisim_t *ms,
4410 const rvec x[], history_t *hist,
4411 rvec f[], t_forcerec *fr,
4412 const t_pbc *pbc, const t_graph *g,
4413 gmx_enerdata_t *enerd, t_nrnb *nrnb,
4415 const t_mdatoms *md,
4416 t_fcdata *fcd, int *global_atom_index,
4419 gmx_bool bCalcEnerVir;
4421 real dvdl[efptNR]; /* The dummy array is to have a place to store the dhdl at other values
4422 of lambda, which will be thrown away in the end*/
4423 const t_pbc *pbc_null;
4426 assert(fr->nthreads == idef->nthreads);
4428 bCalcEnerVir = (force_flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY));
4430 for (i = 0; i < efptNR; i++)
4446 p_graph(debug, "Bondage is fun", g);
4450 /* Do pre force calculation stuff which might require communication */
4451 if (idef->il[F_ORIRES].nr)
4453 enerd->term[F_ORIRESDEV] =
4454 calc_orires_dev(ms, idef->il[F_ORIRES].nr,
4455 idef->il[F_ORIRES].iatoms,
4456 idef->iparams, md, x,
4457 pbc_null, fcd, hist);
4459 if (idef->il[F_DISRES].nr)
4461 calc_disres_R_6(idef->il[F_DISRES].nr,
4462 idef->il[F_DISRES].iatoms,
4463 idef->iparams, x, pbc_null,
4466 if (fcd->disres.nsystems > 1)
4468 gmx_sum_sim(2*fcd->disres.nres, fcd->disres.Rt_6, ms);
4473 #pragma omp parallel for num_threads(fr->nthreads) schedule(static)
4474 for (thread = 0; thread < fr->nthreads; thread++)
4481 gmx_grppairener_t *grpp;
4486 fshift = fr->fshift;
4488 grpp = &enerd->grpp;
4493 zero_thread_forces(&fr->f_t[thread], fr->natoms_force,
4494 fr->red_nblock, 1<<fr->red_ashift);
4496 ft = fr->f_t[thread].f;
4497 fshift = fr->f_t[thread].fshift;
4498 epot = fr->f_t[thread].ener;
4499 grpp = &fr->f_t[thread].grpp;
4500 dvdlt = fr->f_t[thread].dvdl;
4502 /* Loop over all bonded force types to calculate the bonded forces */
4503 for (ftype = 0; (ftype < F_NRE); ftype++)
4505 if (idef->il[ftype].nr > 0 && ftype_is_bonded_potential(ftype))
4507 v = calc_one_bond(thread, ftype, idef, x,
4508 ft, fshift, fr, pbc_null, g, grpp,
4509 nrnb, lambda, dvdlt,
4510 md, fcd, bCalcEnerVir,
4516 if (fr->nthreads > 1)
4518 reduce_thread_forces(fr->natoms_force, f, fr->fshift,
4519 enerd->term, &enerd->grpp, dvdl,
4520 fr->nthreads, fr->f_t,
4521 fr->red_nblock, 1<<fr->red_ashift,
4523 force_flags & GMX_FORCE_DHDL);
4525 if (force_flags & GMX_FORCE_DHDL)
4527 for (i = 0; i < efptNR; i++)
4529 enerd->dvdl_nonlin[i] += dvdl[i];
4533 /* Copy the sum of violations for the distance restraints from fcd */
4536 enerd->term[F_DISRESVIOL] = fcd->disres.sumviol;
4541 void calc_bonds_lambda(const t_idef *idef,
4544 const t_pbc *pbc, const t_graph *g,
4545 gmx_grppairener_t *grpp, real *epot, t_nrnb *nrnb,
4547 const t_mdatoms *md,
4549 int *global_atom_index)
4551 int ftype, nr_nonperturbed, nr;
4553 real dvdl_dum[efptNR];
4555 const t_pbc *pbc_null;
4567 /* Copy the whole idef, so we can modify the contents locally */
4569 idef_fe.nthreads = 1;
4570 snew(idef_fe.il_thread_division, F_NRE*(idef_fe.nthreads+1));
4572 /* We already have the forces, so we use temp buffers here */
4573 snew(f, fr->natoms_force);
4574 snew(fshift, SHIFTS);
4576 /* Loop over all bonded force types to calculate the bonded energies */
4577 for (ftype = 0; (ftype < F_NRE); ftype++)
4579 if (ftype_is_bonded_potential(ftype))
4581 /* Set the work range of thread 0 to the perturbed bondeds only */
4582 nr_nonperturbed = idef->il[ftype].nr_nonperturbed;
4583 nr = idef->il[ftype].nr;
4584 idef_fe.il_thread_division[ftype*2+0] = nr_nonperturbed;
4585 idef_fe.il_thread_division[ftype*2+1] = nr;
4587 /* This is only to get the flop count correct */
4588 idef_fe.il[ftype].nr = nr - nr_nonperturbed;
4590 if (nr - nr_nonperturbed > 0)
4592 v = calc_one_bond(0, ftype, &idef_fe,
4593 x, f, fshift, fr, pbc_null, g,
4594 grpp, nrnb, lambda, dvdl_dum,
4605 sfree(idef_fe.il_thread_division);