1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROningen Mixture of Alchemy and Childrens' Stories
52 #include "gmx_fatal.h"
58 #include "nonbonded.h"
60 /* Include the SIMD macro file and then check for support */
61 #include "gromacs/simd/macros.h"
62 #if defined GMX_HAVE_SIMD_MACROS && defined GMX_SIMD_HAVE_TRIGONOMETRIC
64 #include "gromacs/simd/vector_operations.h"
67 /* Find a better place for this? */
68 const int cmap_coeff_matrix[] = {
69 1, 0, -3, 2, 0, 0, 0, 0, -3, 0, 9, -6, 2, 0, -6, 4,
70 0, 0, 0, 0, 0, 0, 0, 0, 3, 0, -9, 6, -2, 0, 6, -4,
71 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 9, -6, 0, 0, -6, 4,
72 0, 0, 3, -2, 0, 0, 0, 0, 0, 0, -9, 6, 0, 0, 6, -4,
73 0, 0, 0, 0, 1, 0, -3, 2, -2, 0, 6, -4, 1, 0, -3, 2,
74 0, 0, 0, 0, 0, 0, 0, 0, -1, 0, 3, -2, 1, 0, -3, 2,
75 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -3, 2, 0, 0, 3, -2,
76 0, 0, 0, 0, 0, 0, 3, -2, 0, 0, -6, 4, 0, 0, 3, -2,
77 0, 1, -2, 1, 0, 0, 0, 0, 0, -3, 6, -3, 0, 2, -4, 2,
78 0, 0, 0, 0, 0, 0, 0, 0, 0, 3, -6, 3, 0, -2, 4, -2,
79 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -3, 3, 0, 0, 2, -2,
80 0, 0, -1, 1, 0, 0, 0, 0, 0, 0, 3, -3, 0, 0, -2, 2,
81 0, 0, 0, 0, 0, 1, -2, 1, 0, -2, 4, -2, 0, 1, -2, 1,
82 0, 0, 0, 0, 0, 0, 0, 0, 0, -1, 2, -1, 0, 1, -2, 1,
83 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, -1, 0, 0, -1, 1,
84 0, 0, 0, 0, 0, 0, -1, 1, 0, 0, 2, -2, 0, 0, -1, 1
89 int glatnr(int *global_atom_index, int i)
93 if (global_atom_index == NULL)
99 atnr = global_atom_index[i] + 1;
105 static int pbc_rvec_sub(const t_pbc *pbc, const rvec xi, const rvec xj, rvec dx)
109 return pbc_dx_aiuc(pbc, xi, xj, dx);
113 rvec_sub(xi, xj, dx);
120 /* SIMD PBC data structure, containing 1/boxdiag and the box vectors */
133 /* Set the SIMD pbc data from a normal t_pbc struct */
134 static void set_pbc_simd(const t_pbc *pbc, pbc_simd_t *pbc_simd)
139 /* Setting inv_bdiag to 0 effectively turns off PBC */
140 clear_rvec(inv_bdiag);
143 for (d = 0; d < pbc->ndim_ePBC; d++)
145 inv_bdiag[d] = 1.0/pbc->box[d][d];
149 pbc_simd->inv_bzz = gmx_set1_pr(inv_bdiag[ZZ]);
150 pbc_simd->inv_byy = gmx_set1_pr(inv_bdiag[YY]);
151 pbc_simd->inv_bxx = gmx_set1_pr(inv_bdiag[XX]);
155 pbc_simd->bzx = gmx_set1_pr(pbc->box[ZZ][XX]);
156 pbc_simd->bzy = gmx_set1_pr(pbc->box[ZZ][YY]);
157 pbc_simd->bzz = gmx_set1_pr(pbc->box[ZZ][ZZ]);
158 pbc_simd->byx = gmx_set1_pr(pbc->box[YY][XX]);
159 pbc_simd->byy = gmx_set1_pr(pbc->box[YY][YY]);
160 pbc_simd->bxx = gmx_set1_pr(pbc->box[XX][XX]);
164 pbc_simd->bzx = gmx_setzero_pr();
165 pbc_simd->bzy = gmx_setzero_pr();
166 pbc_simd->bzz = gmx_setzero_pr();
167 pbc_simd->byx = gmx_setzero_pr();
168 pbc_simd->byy = gmx_setzero_pr();
169 pbc_simd->bxx = gmx_setzero_pr();
173 /* Correct distance vector *dx,*dy,*dz for PBC using SIMD */
174 static gmx_inline void
175 pbc_dx_simd(gmx_mm_pr *dx, gmx_mm_pr *dy, gmx_mm_pr *dz,
176 const pbc_simd_t *pbc)
180 sh = gmx_round_pr(gmx_mul_pr(*dz, pbc->inv_bzz));
181 *dx = gmx_nmsub_pr(sh, pbc->bzx, *dx);
182 *dy = gmx_nmsub_pr(sh, pbc->bzy, *dy);
183 *dz = gmx_nmsub_pr(sh, pbc->bzz, *dz);
185 sh = gmx_round_pr(gmx_mul_pr(*dy, pbc->inv_byy));
186 *dx = gmx_nmsub_pr(sh, pbc->byx, *dx);
187 *dy = gmx_nmsub_pr(sh, pbc->byy, *dy);
189 sh = gmx_round_pr(gmx_mul_pr(*dx, pbc->inv_bxx));
190 *dx = gmx_nmsub_pr(sh, pbc->bxx, *dx);
193 #endif /* SIMD_BONDEDS */
196 * Morse potential bond by Frank Everdij
198 * Three parameters needed:
200 * b0 = equilibrium distance in nm
201 * be = beta in nm^-1 (actually, it's nu_e*Sqrt(2*pi*pi*mu/D_e))
202 * cb = well depth in kJ/mol
204 * Note: the potential is referenced to be +cb at infinite separation
205 * and zero at the equilibrium distance!
208 real morse_bonds(int nbonds,
209 const t_iatom forceatoms[], const t_iparams forceparams[],
210 const rvec x[], rvec f[], rvec fshift[],
211 const t_pbc *pbc, const t_graph *g,
212 real lambda, real *dvdlambda,
213 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
214 int gmx_unused *global_atom_index)
216 const real one = 1.0;
217 const real two = 2.0;
218 real dr, dr2, temp, omtemp, cbomtemp, fbond, vbond, fij, vtot;
219 real b0, be, cb, b0A, beA, cbA, b0B, beB, cbB, L1;
221 int i, m, ki, type, ai, aj;
225 for (i = 0; (i < nbonds); )
227 type = forceatoms[i++];
228 ai = forceatoms[i++];
229 aj = forceatoms[i++];
231 b0A = forceparams[type].morse.b0A;
232 beA = forceparams[type].morse.betaA;
233 cbA = forceparams[type].morse.cbA;
235 b0B = forceparams[type].morse.b0B;
236 beB = forceparams[type].morse.betaB;
237 cbB = forceparams[type].morse.cbB;
239 L1 = one-lambda; /* 1 */
240 b0 = L1*b0A + lambda*b0B; /* 3 */
241 be = L1*beA + lambda*beB; /* 3 */
242 cb = L1*cbA + lambda*cbB; /* 3 */
244 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
245 dr2 = iprod(dx, dx); /* 5 */
246 dr = dr2*gmx_invsqrt(dr2); /* 10 */
247 temp = exp(-be*(dr-b0)); /* 12 */
251 /* bonds are constrainted. This may _not_ include bond constraints if they are lambda dependent */
252 *dvdlambda += cbB-cbA;
256 omtemp = one-temp; /* 1 */
257 cbomtemp = cb*omtemp; /* 1 */
258 vbond = cbomtemp*omtemp; /* 1 */
259 fbond = -two*be*temp*cbomtemp*gmx_invsqrt(dr2); /* 9 */
260 vtot += vbond; /* 1 */
262 *dvdlambda += (cbB - cbA) * omtemp * omtemp - (2-2*omtemp)*omtemp * cb * ((b0B-b0A)*be - (beB-beA)*(dr-b0)); /* 15 */
266 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
270 for (m = 0; (m < DIM); m++) /* 15 */
275 fshift[ki][m] += fij;
276 fshift[CENTRAL][m] -= fij;
282 real cubic_bonds(int nbonds,
283 const t_iatom forceatoms[], const t_iparams forceparams[],
284 const rvec x[], rvec f[], rvec fshift[],
285 const t_pbc *pbc, const t_graph *g,
286 real gmx_unused lambda, real gmx_unused *dvdlambda,
287 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
288 int gmx_unused *global_atom_index)
290 const real three = 3.0;
291 const real two = 2.0;
293 real dr, dr2, dist, kdist, kdist2, fbond, vbond, fij, vtot;
295 int i, m, ki, type, ai, aj;
299 for (i = 0; (i < nbonds); )
301 type = forceatoms[i++];
302 ai = forceatoms[i++];
303 aj = forceatoms[i++];
305 b0 = forceparams[type].cubic.b0;
306 kb = forceparams[type].cubic.kb;
307 kcub = forceparams[type].cubic.kcub;
309 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
310 dr2 = iprod(dx, dx); /* 5 */
317 dr = dr2*gmx_invsqrt(dr2); /* 10 */
322 vbond = kdist2 + kcub*kdist2*dist;
323 fbond = -(two*kdist + three*kdist2*kcub)/dr;
325 vtot += vbond; /* 21 */
329 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
332 for (m = 0; (m < DIM); m++) /* 15 */
337 fshift[ki][m] += fij;
338 fshift[CENTRAL][m] -= fij;
344 real FENE_bonds(int nbonds,
345 const t_iatom forceatoms[], const t_iparams forceparams[],
346 const rvec x[], rvec f[], rvec fshift[],
347 const t_pbc *pbc, const t_graph *g,
348 real gmx_unused lambda, real gmx_unused *dvdlambda,
349 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
350 int *global_atom_index)
352 const real half = 0.5;
353 const real one = 1.0;
355 real dr, dr2, bm2, omdr2obm2, fbond, vbond, fij, vtot;
357 int i, m, ki, type, ai, aj;
361 for (i = 0; (i < nbonds); )
363 type = forceatoms[i++];
364 ai = forceatoms[i++];
365 aj = forceatoms[i++];
367 bm = forceparams[type].fene.bm;
368 kb = forceparams[type].fene.kb;
370 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
371 dr2 = iprod(dx, dx); /* 5 */
383 "r^2 (%f) >= bm^2 (%f) in FENE bond between atoms %d and %d",
385 glatnr(global_atom_index, ai),
386 glatnr(global_atom_index, aj));
389 omdr2obm2 = one - dr2/bm2;
391 vbond = -half*kb*bm2*log(omdr2obm2);
392 fbond = -kb/omdr2obm2;
394 vtot += vbond; /* 35 */
398 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
401 for (m = 0; (m < DIM); m++) /* 15 */
406 fshift[ki][m] += fij;
407 fshift[CENTRAL][m] -= fij;
413 real harmonic(real kA, real kB, real xA, real xB, real x, real lambda,
416 const real half = 0.5;
417 real L1, kk, x0, dx, dx2;
418 real v, f, dvdlambda;
421 kk = L1*kA+lambda*kB;
422 x0 = L1*xA+lambda*xB;
429 dvdlambda = half*(kB-kA)*dx2 + (xA-xB)*kk*dx;
436 /* That was 19 flops */
440 real bonds(int nbonds,
441 const t_iatom forceatoms[], const t_iparams forceparams[],
442 const rvec x[], rvec f[], rvec fshift[],
443 const t_pbc *pbc, const t_graph *g,
444 real lambda, real *dvdlambda,
445 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
446 int gmx_unused *global_atom_index)
448 int i, m, ki, ai, aj, type;
449 real dr, dr2, fbond, vbond, fij, vtot;
454 for (i = 0; (i < nbonds); )
456 type = forceatoms[i++];
457 ai = forceatoms[i++];
458 aj = forceatoms[i++];
460 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
461 dr2 = iprod(dx, dx); /* 5 */
462 dr = dr2*gmx_invsqrt(dr2); /* 10 */
464 *dvdlambda += harmonic(forceparams[type].harmonic.krA,
465 forceparams[type].harmonic.krB,
466 forceparams[type].harmonic.rA,
467 forceparams[type].harmonic.rB,
468 dr, lambda, &vbond, &fbond); /* 19 */
476 vtot += vbond; /* 1*/
477 fbond *= gmx_invsqrt(dr2); /* 6 */
481 fprintf(debug, "BONDS: dr = %10g vbond = %10g fbond = %10g\n",
487 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
490 for (m = 0; (m < DIM); m++) /* 15 */
495 fshift[ki][m] += fij;
496 fshift[CENTRAL][m] -= fij;
502 real restraint_bonds(int nbonds,
503 const t_iatom forceatoms[], const t_iparams forceparams[],
504 const rvec x[], rvec f[], rvec fshift[],
505 const t_pbc *pbc, const t_graph *g,
506 real lambda, real *dvdlambda,
507 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
508 int gmx_unused *global_atom_index)
510 int i, m, ki, ai, aj, type;
511 real dr, dr2, fbond, vbond, fij, vtot;
513 real low, dlow, up1, dup1, up2, dup2, k, dk;
521 for (i = 0; (i < nbonds); )
523 type = forceatoms[i++];
524 ai = forceatoms[i++];
525 aj = forceatoms[i++];
527 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
528 dr2 = iprod(dx, dx); /* 5 */
529 dr = dr2*gmx_invsqrt(dr2); /* 10 */
531 low = L1*forceparams[type].restraint.lowA + lambda*forceparams[type].restraint.lowB;
532 dlow = -forceparams[type].restraint.lowA + forceparams[type].restraint.lowB;
533 up1 = L1*forceparams[type].restraint.up1A + lambda*forceparams[type].restraint.up1B;
534 dup1 = -forceparams[type].restraint.up1A + forceparams[type].restraint.up1B;
535 up2 = L1*forceparams[type].restraint.up2A + lambda*forceparams[type].restraint.up2B;
536 dup2 = -forceparams[type].restraint.up2A + forceparams[type].restraint.up2B;
537 k = L1*forceparams[type].restraint.kA + lambda*forceparams[type].restraint.kB;
538 dk = -forceparams[type].restraint.kA + forceparams[type].restraint.kB;
547 *dvdlambda += 0.5*dk*drh2 - k*dlow*drh;
560 *dvdlambda += 0.5*dk*drh2 - k*dup1*drh;
565 vbond = k*(up2 - up1)*(0.5*(up2 - up1) + drh);
566 fbond = -k*(up2 - up1);
567 *dvdlambda += dk*(up2 - up1)*(0.5*(up2 - up1) + drh)
568 + k*(dup2 - dup1)*(up2 - up1 + drh)
569 - k*(up2 - up1)*dup2;
577 vtot += vbond; /* 1*/
578 fbond *= gmx_invsqrt(dr2); /* 6 */
582 fprintf(debug, "BONDS: dr = %10g vbond = %10g fbond = %10g\n",
588 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
591 for (m = 0; (m < DIM); m++) /* 15 */
596 fshift[ki][m] += fij;
597 fshift[CENTRAL][m] -= fij;
604 real polarize(int nbonds,
605 const t_iatom forceatoms[], const t_iparams forceparams[],
606 const rvec x[], rvec f[], rvec fshift[],
607 const t_pbc *pbc, const t_graph *g,
608 real lambda, real *dvdlambda,
609 const t_mdatoms *md, t_fcdata gmx_unused *fcd,
610 int gmx_unused *global_atom_index)
612 int i, m, ki, ai, aj, type;
613 real dr, dr2, fbond, vbond, fij, vtot, ksh;
618 for (i = 0; (i < nbonds); )
620 type = forceatoms[i++];
621 ai = forceatoms[i++];
622 aj = forceatoms[i++];
623 ksh = sqr(md->chargeA[aj])*ONE_4PI_EPS0/forceparams[type].polarize.alpha;
626 fprintf(debug, "POL: local ai = %d aj = %d ksh = %.3f\n", ai, aj, ksh);
629 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
630 dr2 = iprod(dx, dx); /* 5 */
631 dr = dr2*gmx_invsqrt(dr2); /* 10 */
633 *dvdlambda += harmonic(ksh, ksh, 0, 0, dr, lambda, &vbond, &fbond); /* 19 */
640 vtot += vbond; /* 1*/
641 fbond *= gmx_invsqrt(dr2); /* 6 */
645 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
648 for (m = 0; (m < DIM); m++) /* 15 */
653 fshift[ki][m] += fij;
654 fshift[CENTRAL][m] -= fij;
660 real anharm_polarize(int nbonds,
661 const t_iatom forceatoms[], const t_iparams forceparams[],
662 const rvec x[], rvec f[], rvec fshift[],
663 const t_pbc *pbc, const t_graph *g,
664 real lambda, real *dvdlambda,
665 const t_mdatoms *md, t_fcdata gmx_unused *fcd,
666 int gmx_unused *global_atom_index)
668 int i, m, ki, ai, aj, type;
669 real dr, dr2, fbond, vbond, fij, vtot, ksh, khyp, drcut, ddr, ddr3;
674 for (i = 0; (i < nbonds); )
676 type = forceatoms[i++];
677 ai = forceatoms[i++];
678 aj = forceatoms[i++];
679 ksh = sqr(md->chargeA[aj])*ONE_4PI_EPS0/forceparams[type].anharm_polarize.alpha; /* 7*/
680 khyp = forceparams[type].anharm_polarize.khyp;
681 drcut = forceparams[type].anharm_polarize.drcut;
684 fprintf(debug, "POL: local ai = %d aj = %d ksh = %.3f\n", ai, aj, ksh);
687 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
688 dr2 = iprod(dx, dx); /* 5 */
689 dr = dr2*gmx_invsqrt(dr2); /* 10 */
691 *dvdlambda += harmonic(ksh, ksh, 0, 0, dr, lambda, &vbond, &fbond); /* 19 */
702 vbond += khyp*ddr*ddr3;
703 fbond -= 4*khyp*ddr3;
705 fbond *= gmx_invsqrt(dr2); /* 6 */
706 vtot += vbond; /* 1*/
710 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
713 for (m = 0; (m < DIM); m++) /* 15 */
718 fshift[ki][m] += fij;
719 fshift[CENTRAL][m] -= fij;
725 real water_pol(int nbonds,
726 const t_iatom forceatoms[], const t_iparams forceparams[],
727 const rvec x[], rvec f[], rvec gmx_unused fshift[],
728 const t_pbc gmx_unused *pbc, const t_graph gmx_unused *g,
729 real gmx_unused lambda, real gmx_unused *dvdlambda,
730 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
731 int gmx_unused *global_atom_index)
733 /* This routine implements anisotropic polarizibility for water, through
734 * a shell connected to a dummy with spring constant that differ in the
735 * three spatial dimensions in the molecular frame.
737 int i, m, aO, aH1, aH2, aD, aS, type, type0;
738 rvec dOH1, dOH2, dHH, dOD, dDS, nW, kk, dx, kdx, proj;
742 real vtot, fij, r_HH, r_OD, r_nW, tx, ty, tz, qS;
747 type0 = forceatoms[0];
749 qS = md->chargeA[aS];
750 kk[XX] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_x;
751 kk[YY] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_y;
752 kk[ZZ] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_z;
753 r_HH = 1.0/forceparams[type0].wpol.rHH;
754 r_OD = 1.0/forceparams[type0].wpol.rOD;
757 fprintf(debug, "WPOL: qS = %10.5f aS = %5d\n", qS, aS);
758 fprintf(debug, "WPOL: kk = %10.3f %10.3f %10.3f\n",
759 kk[XX], kk[YY], kk[ZZ]);
760 fprintf(debug, "WPOL: rOH = %10.3f rHH = %10.3f rOD = %10.3f\n",
761 forceparams[type0].wpol.rOH,
762 forceparams[type0].wpol.rHH,
763 forceparams[type0].wpol.rOD);
765 for (i = 0; (i < nbonds); i += 6)
767 type = forceatoms[i];
770 gmx_fatal(FARGS, "Sorry, type = %d, type0 = %d, file = %s, line = %d",
771 type, type0, __FILE__, __LINE__);
773 aO = forceatoms[i+1];
774 aH1 = forceatoms[i+2];
775 aH2 = forceatoms[i+3];
776 aD = forceatoms[i+4];
777 aS = forceatoms[i+5];
779 /* Compute vectors describing the water frame */
780 rvec_sub(x[aH1], x[aO], dOH1);
781 rvec_sub(x[aH2], x[aO], dOH2);
782 rvec_sub(x[aH2], x[aH1], dHH);
783 rvec_sub(x[aD], x[aO], dOD);
784 rvec_sub(x[aS], x[aD], dDS);
785 cprod(dOH1, dOH2, nW);
787 /* Compute inverse length of normal vector
788 * (this one could be precomputed, but I'm too lazy now)
790 r_nW = gmx_invsqrt(iprod(nW, nW));
791 /* This is for precision, but does not make a big difference,
794 r_OD = gmx_invsqrt(iprod(dOD, dOD));
796 /* Normalize the vectors in the water frame */
798 svmul(r_HH, dHH, dHH);
799 svmul(r_OD, dOD, dOD);
801 /* Compute displacement of shell along components of the vector */
802 dx[ZZ] = iprod(dDS, dOD);
803 /* Compute projection on the XY plane: dDS - dx[ZZ]*dOD */
804 for (m = 0; (m < DIM); m++)
806 proj[m] = dDS[m]-dx[ZZ]*dOD[m];
809 /*dx[XX] = iprod(dDS,nW);
810 dx[YY] = iprod(dDS,dHH);*/
811 dx[XX] = iprod(proj, nW);
812 for (m = 0; (m < DIM); m++)
814 proj[m] -= dx[XX]*nW[m];
816 dx[YY] = iprod(proj, dHH);
821 fprintf(debug, "WPOL: dx2=%10g dy2=%10g dz2=%10g sum=%10g dDS^2=%10g\n",
822 sqr(dx[XX]), sqr(dx[YY]), sqr(dx[ZZ]), iprod(dx, dx), iprod(dDS, dDS));
823 fprintf(debug, "WPOL: dHH=(%10g,%10g,%10g)\n", dHH[XX], dHH[YY], dHH[ZZ]);
824 fprintf(debug, "WPOL: dOD=(%10g,%10g,%10g), 1/r_OD = %10g\n",
825 dOD[XX], dOD[YY], dOD[ZZ], 1/r_OD);
826 fprintf(debug, "WPOL: nW =(%10g,%10g,%10g), 1/r_nW = %10g\n",
827 nW[XX], nW[YY], nW[ZZ], 1/r_nW);
828 fprintf(debug, "WPOL: dx =%10g, dy =%10g, dz =%10g\n",
829 dx[XX], dx[YY], dx[ZZ]);
830 fprintf(debug, "WPOL: dDSx=%10g, dDSy=%10g, dDSz=%10g\n",
831 dDS[XX], dDS[YY], dDS[ZZ]);
834 /* Now compute the forces and energy */
835 kdx[XX] = kk[XX]*dx[XX];
836 kdx[YY] = kk[YY]*dx[YY];
837 kdx[ZZ] = kk[ZZ]*dx[ZZ];
838 vtot += iprod(dx, kdx);
839 for (m = 0; (m < DIM); m++)
841 /* This is a tensor operation but written out for speed */
855 fprintf(debug, "WPOL: vwpol=%g\n", 0.5*iprod(dx, kdx));
856 fprintf(debug, "WPOL: df = (%10g, %10g, %10g)\n", df[XX], df[YY], df[ZZ]);
864 static real do_1_thole(const rvec xi, const rvec xj, rvec fi, rvec fj,
865 const t_pbc *pbc, real qq,
866 rvec fshift[], real afac)
869 real r12sq, r12_1, r12n, r12bar, v0, v1, fscal, ebar, fff;
872 t = pbc_rvec_sub(pbc, xi, xj, r12); /* 3 */
874 r12sq = iprod(r12, r12); /* 5 */
875 r12_1 = gmx_invsqrt(r12sq); /* 5 */
876 r12bar = afac/r12_1; /* 5 */
877 v0 = qq*ONE_4PI_EPS0*r12_1; /* 2 */
878 ebar = exp(-r12bar); /* 5 */
879 v1 = (1-(1+0.5*r12bar)*ebar); /* 4 */
880 fscal = ((v0*r12_1)*v1 - v0*0.5*afac*ebar*(r12bar+1))*r12_1; /* 9 */
883 fprintf(debug, "THOLE: v0 = %.3f v1 = %.3f r12= % .3f r12bar = %.3f fscal = %.3f ebar = %.3f\n", v0, v1, 1/r12_1, r12bar, fscal, ebar);
886 for (m = 0; (m < DIM); m++)
892 fshift[CENTRAL][m] -= fff;
895 return v0*v1; /* 1 */
899 real thole_pol(int nbonds,
900 const t_iatom forceatoms[], const t_iparams forceparams[],
901 const rvec x[], rvec f[], rvec fshift[],
902 const t_pbc *pbc, const t_graph gmx_unused *g,
903 real gmx_unused lambda, real gmx_unused *dvdlambda,
904 const t_mdatoms *md, t_fcdata gmx_unused *fcd,
905 int gmx_unused *global_atom_index)
907 /* Interaction between two pairs of particles with opposite charge */
908 int i, type, a1, da1, a2, da2;
909 real q1, q2, qq, a, al1, al2, afac;
912 for (i = 0; (i < nbonds); )
914 type = forceatoms[i++];
915 a1 = forceatoms[i++];
916 da1 = forceatoms[i++];
917 a2 = forceatoms[i++];
918 da2 = forceatoms[i++];
919 q1 = md->chargeA[da1];
920 q2 = md->chargeA[da2];
921 a = forceparams[type].thole.a;
922 al1 = forceparams[type].thole.alpha1;
923 al2 = forceparams[type].thole.alpha2;
925 afac = a*pow(al1*al2, -1.0/6.0);
926 V += do_1_thole(x[a1], x[a2], f[a1], f[a2], pbc, qq, fshift, afac);
927 V += do_1_thole(x[da1], x[a2], f[da1], f[a2], pbc, -qq, fshift, afac);
928 V += do_1_thole(x[a1], x[da2], f[a1], f[da2], pbc, -qq, fshift, afac);
929 V += do_1_thole(x[da1], x[da2], f[da1], f[da2], pbc, qq, fshift, afac);
935 real bond_angle(const rvec xi, const rvec xj, const rvec xk, const t_pbc *pbc,
936 rvec r_ij, rvec r_kj, real *costh,
938 /* Return value is the angle between the bonds i-j and j-k */
943 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
944 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
946 *costh = cos_angle(r_ij, r_kj); /* 25 */
947 th = acos(*costh); /* 10 */
952 real angles(int nbonds,
953 const t_iatom forceatoms[], const t_iparams forceparams[],
954 const rvec x[], rvec f[], rvec fshift[],
955 const t_pbc *pbc, const t_graph *g,
956 real lambda, real *dvdlambda,
957 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
958 int gmx_unused *global_atom_index)
960 int i, ai, aj, ak, t1, t2, type;
962 real cos_theta, cos_theta2, theta, dVdt, va, vtot;
963 ivec jt, dt_ij, dt_kj;
966 for (i = 0; i < nbonds; )
968 type = forceatoms[i++];
969 ai = forceatoms[i++];
970 aj = forceatoms[i++];
971 ak = forceatoms[i++];
973 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
974 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
976 *dvdlambda += harmonic(forceparams[type].harmonic.krA,
977 forceparams[type].harmonic.krB,
978 forceparams[type].harmonic.rA*DEG2RAD,
979 forceparams[type].harmonic.rB*DEG2RAD,
980 theta, lambda, &va, &dVdt); /* 21 */
983 cos_theta2 = sqr(cos_theta);
993 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
994 sth = st*cos_theta; /* 1 */
998 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
999 theta*RAD2DEG, va, dVdt);
1002 nrij2 = iprod(r_ij, r_ij); /* 5 */
1003 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1005 nrij_1 = gmx_invsqrt(nrij2); /* 10 */
1006 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1008 cik = st*nrij_1*nrkj_1; /* 2 */
1009 cii = sth*nrij_1*nrij_1; /* 2 */
1010 ckk = sth*nrkj_1*nrkj_1; /* 2 */
1012 for (m = 0; m < DIM; m++)
1014 f_i[m] = -(cik*r_kj[m] - cii*r_ij[m]);
1015 f_k[m] = -(cik*r_ij[m] - ckk*r_kj[m]);
1016 f_j[m] = -f_i[m] - f_k[m];
1023 copy_ivec(SHIFT_IVEC(g, aj), jt);
1025 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1026 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1027 t1 = IVEC2IS(dt_ij);
1028 t2 = IVEC2IS(dt_kj);
1030 rvec_inc(fshift[t1], f_i);
1031 rvec_inc(fshift[CENTRAL], f_j);
1032 rvec_inc(fshift[t2], f_k);
1041 /* As angles, but using SIMD to calculate many dihedrals at once.
1042 * This routines does not calculate energies and shift forces.
1044 static gmx_inline void
1045 angles_noener_simd(int nbonds,
1046 const t_iatom forceatoms[], const t_iparams forceparams[],
1047 const rvec x[], rvec f[],
1048 const t_pbc *pbc, const t_graph gmx_unused *g,
1049 real gmx_unused lambda,
1050 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1051 int gmx_unused *global_atom_index)
1053 #define UNROLL GMX_SIMD_WIDTH_HERE
1056 int type, ai[UNROLL], aj[UNROLL], ak[UNROLL];
1057 real coeff_array[2*UNROLL+UNROLL], *coeff;
1058 real dr_array[2*DIM*UNROLL+UNROLL], *dr;
1059 real f_buf_array[6*UNROLL+UNROLL], *f_buf;
1060 gmx_mm_pr k_S, theta0_S;
1061 gmx_mm_pr rijx_S, rijy_S, rijz_S;
1062 gmx_mm_pr rkjx_S, rkjy_S, rkjz_S;
1064 gmx_mm_pr min_one_plus_eps_S;
1065 gmx_mm_pr rij_rkj_S;
1066 gmx_mm_pr nrij2_S, nrij_1_S;
1067 gmx_mm_pr nrkj2_S, nrkj_1_S;
1068 gmx_mm_pr cos_S, invsin_S;
1070 gmx_mm_pr st_S, sth_S;
1071 gmx_mm_pr cik_S, cii_S, ckk_S;
1072 gmx_mm_pr f_ix_S, f_iy_S, f_iz_S;
1073 gmx_mm_pr f_kx_S, f_ky_S, f_kz_S;
1074 pbc_simd_t pbc_simd;
1076 /* Ensure register memory alignment */
1077 coeff = gmx_simd_align_real(coeff_array);
1078 dr = gmx_simd_align_real(dr_array);
1079 f_buf = gmx_simd_align_real(f_buf_array);
1081 set_pbc_simd(pbc, &pbc_simd);
1083 one_S = gmx_set1_pr(1.0);
1085 /* The smallest number > -1 */
1086 min_one_plus_eps_S = gmx_set1_pr(-1.0 + 2*GMX_REAL_EPS);
1088 /* nbonds is the number of angles times nfa1, here we step UNROLL angles */
1089 for (i = 0; (i < nbonds); i += UNROLL*nfa1)
1091 /* Collect atoms for UNROLL angles.
1092 * iu indexes into forceatoms, we should not let iu go beyond nbonds.
1095 for (s = 0; s < UNROLL; s++)
1097 type = forceatoms[iu];
1098 ai[s] = forceatoms[iu+1];
1099 aj[s] = forceatoms[iu+2];
1100 ak[s] = forceatoms[iu+3];
1102 coeff[s] = forceparams[type].harmonic.krA;
1103 coeff[UNROLL+s] = forceparams[type].harmonic.rA*DEG2RAD;
1105 /* If you can't use pbc_dx_simd below for PBC, e.g. because
1106 * you can't round in SIMD, use pbc_rvec_sub here.
1108 /* Store the non PBC corrected distances packed and aligned */
1109 for (m = 0; m < DIM; m++)
1111 dr[s + m *UNROLL] = x[ai[s]][m] - x[aj[s]][m];
1112 dr[s + (DIM+m)*UNROLL] = x[ak[s]][m] - x[aj[s]][m];
1115 /* At the end fill the arrays with identical entries */
1116 if (iu + nfa1 < nbonds)
1122 k_S = gmx_load_pr(coeff);
1123 theta0_S = gmx_load_pr(coeff+UNROLL);
1125 rijx_S = gmx_load_pr(dr + 0*UNROLL);
1126 rijy_S = gmx_load_pr(dr + 1*UNROLL);
1127 rijz_S = gmx_load_pr(dr + 2*UNROLL);
1128 rkjx_S = gmx_load_pr(dr + 3*UNROLL);
1129 rkjy_S = gmx_load_pr(dr + 4*UNROLL);
1130 rkjz_S = gmx_load_pr(dr + 5*UNROLL);
1132 pbc_dx_simd(&rijx_S, &rijy_S, &rijz_S, &pbc_simd);
1133 pbc_dx_simd(&rkjx_S, &rkjy_S, &rkjz_S, &pbc_simd);
1135 rij_rkj_S = gmx_iprod_pr(rijx_S, rijy_S, rijz_S,
1136 rkjx_S, rkjy_S, rkjz_S);
1138 nrij2_S = gmx_norm2_pr(rijx_S, rijy_S, rijz_S);
1139 nrkj2_S = gmx_norm2_pr(rkjx_S, rkjy_S, rkjz_S);
1141 nrij_1_S = gmx_invsqrt_pr(nrij2_S);
1142 nrkj_1_S = gmx_invsqrt_pr(nrkj2_S);
1144 cos_S = gmx_mul_pr(rij_rkj_S, gmx_mul_pr(nrij_1_S, nrkj_1_S));
1146 /* To allow for 180 degrees, we take the max of cos and -1 + 1bit,
1147 * so we can safely get the 1/sin from 1/sqrt(1 - cos^2).
1148 * This also ensures that rounding errors would cause the argument
1149 * of gmx_acos_pr to be < -1.
1150 * Note that we do not take precautions for cos(0)=1, so the outer
1151 * atoms in an angle should not be on top of each other.
1153 cos_S = gmx_max_pr(cos_S, min_one_plus_eps_S);
1155 theta_S = gmx_acos_pr(cos_S);
1157 invsin_S = gmx_invsqrt_pr(gmx_sub_pr(one_S, gmx_mul_pr(cos_S, cos_S)));
1159 st_S = gmx_mul_pr(gmx_mul_pr(k_S, gmx_sub_pr(theta0_S, theta_S)),
1161 sth_S = gmx_mul_pr(st_S, cos_S);
1163 cik_S = gmx_mul_pr(st_S, gmx_mul_pr(nrij_1_S, nrkj_1_S));
1164 cii_S = gmx_mul_pr(sth_S, gmx_mul_pr(nrij_1_S, nrij_1_S));
1165 ckk_S = gmx_mul_pr(sth_S, gmx_mul_pr(nrkj_1_S, nrkj_1_S));
1167 f_ix_S = gmx_mul_pr(cii_S, rijx_S);
1168 f_ix_S = gmx_nmsub_pr(cik_S, rkjx_S, f_ix_S);
1169 f_iy_S = gmx_mul_pr(cii_S, rijy_S);
1170 f_iy_S = gmx_nmsub_pr(cik_S, rkjy_S, f_iy_S);
1171 f_iz_S = gmx_mul_pr(cii_S, rijz_S);
1172 f_iz_S = gmx_nmsub_pr(cik_S, rkjz_S, f_iz_S);
1173 f_kx_S = gmx_mul_pr(ckk_S, rkjx_S);
1174 f_kx_S = gmx_nmsub_pr(cik_S, rijx_S, f_kx_S);
1175 f_ky_S = gmx_mul_pr(ckk_S, rkjy_S);
1176 f_ky_S = gmx_nmsub_pr(cik_S, rijy_S, f_ky_S);
1177 f_kz_S = gmx_mul_pr(ckk_S, rkjz_S);
1178 f_kz_S = gmx_nmsub_pr(cik_S, rijz_S, f_kz_S);
1180 gmx_store_pr(f_buf + 0*UNROLL, f_ix_S);
1181 gmx_store_pr(f_buf + 1*UNROLL, f_iy_S);
1182 gmx_store_pr(f_buf + 2*UNROLL, f_iz_S);
1183 gmx_store_pr(f_buf + 3*UNROLL, f_kx_S);
1184 gmx_store_pr(f_buf + 4*UNROLL, f_ky_S);
1185 gmx_store_pr(f_buf + 5*UNROLL, f_kz_S);
1191 for (m = 0; m < DIM; m++)
1193 f[ai[s]][m] += f_buf[s + m*UNROLL];
1194 f[aj[s]][m] -= f_buf[s + m*UNROLL] + f_buf[s + (DIM+m)*UNROLL];
1195 f[ak[s]][m] += f_buf[s + (DIM+m)*UNROLL];
1200 while (s < UNROLL && iu < nbonds);
1205 #endif /* SIMD_BONDEDS */
1207 real linear_angles(int nbonds,
1208 const t_iatom forceatoms[], const t_iparams forceparams[],
1209 const rvec x[], rvec f[], rvec fshift[],
1210 const t_pbc *pbc, const t_graph *g,
1211 real lambda, real *dvdlambda,
1212 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1213 int gmx_unused *global_atom_index)
1215 int i, m, ai, aj, ak, t1, t2, type;
1217 real L1, kA, kB, aA, aB, dr, dr2, va, vtot, a, b, klin;
1218 ivec jt, dt_ij, dt_kj;
1219 rvec r_ij, r_kj, r_ik, dx;
1223 for (i = 0; (i < nbonds); )
1225 type = forceatoms[i++];
1226 ai = forceatoms[i++];
1227 aj = forceatoms[i++];
1228 ak = forceatoms[i++];
1230 kA = forceparams[type].linangle.klinA;
1231 kB = forceparams[type].linangle.klinB;
1232 klin = L1*kA + lambda*kB;
1234 aA = forceparams[type].linangle.aA;
1235 aB = forceparams[type].linangle.aB;
1236 a = L1*aA+lambda*aB;
1239 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
1240 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
1241 rvec_sub(r_ij, r_kj, r_ik);
1244 for (m = 0; (m < DIM); m++)
1246 dr = -a * r_ij[m] - b * r_kj[m];
1251 f_j[m] = -(f_i[m]+f_k[m]);
1257 *dvdlambda += 0.5*(kB-kA)*dr2 + klin*(aB-aA)*iprod(dx, r_ik);
1263 copy_ivec(SHIFT_IVEC(g, aj), jt);
1265 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1266 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1267 t1 = IVEC2IS(dt_ij);
1268 t2 = IVEC2IS(dt_kj);
1270 rvec_inc(fshift[t1], f_i);
1271 rvec_inc(fshift[CENTRAL], f_j);
1272 rvec_inc(fshift[t2], f_k);
1277 real urey_bradley(int nbonds,
1278 const t_iatom forceatoms[], const t_iparams forceparams[],
1279 const rvec x[], rvec f[], rvec fshift[],
1280 const t_pbc *pbc, const t_graph *g,
1281 real lambda, real *dvdlambda,
1282 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1283 int gmx_unused *global_atom_index)
1285 int i, m, ai, aj, ak, t1, t2, type, ki;
1286 rvec r_ij, r_kj, r_ik;
1287 real cos_theta, cos_theta2, theta;
1288 real dVdt, va, vtot, dr, dr2, vbond, fbond, fik;
1289 real kthA, th0A, kUBA, r13A, kthB, th0B, kUBB, r13B;
1290 ivec jt, dt_ij, dt_kj, dt_ik;
1293 for (i = 0; (i < nbonds); )
1295 type = forceatoms[i++];
1296 ai = forceatoms[i++];
1297 aj = forceatoms[i++];
1298 ak = forceatoms[i++];
1299 th0A = forceparams[type].u_b.thetaA*DEG2RAD;
1300 kthA = forceparams[type].u_b.kthetaA;
1301 r13A = forceparams[type].u_b.r13A;
1302 kUBA = forceparams[type].u_b.kUBA;
1303 th0B = forceparams[type].u_b.thetaB*DEG2RAD;
1304 kthB = forceparams[type].u_b.kthetaB;
1305 r13B = forceparams[type].u_b.r13B;
1306 kUBB = forceparams[type].u_b.kUBB;
1308 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
1309 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
1311 *dvdlambda += harmonic(kthA, kthB, th0A, th0B, theta, lambda, &va, &dVdt); /* 21 */
1314 ki = pbc_rvec_sub(pbc, x[ai], x[ak], r_ik); /* 3 */
1315 dr2 = iprod(r_ik, r_ik); /* 5 */
1316 dr = dr2*gmx_invsqrt(dr2); /* 10 */
1318 *dvdlambda += harmonic(kUBA, kUBB, r13A, r13B, dr, lambda, &vbond, &fbond); /* 19 */
1320 cos_theta2 = sqr(cos_theta); /* 1 */
1328 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1329 sth = st*cos_theta; /* 1 */
1333 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1334 theta*RAD2DEG, va, dVdt);
1337 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1338 nrij2 = iprod(r_ij, r_ij);
1340 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
1341 cii = sth/nrij2; /* 10 */
1342 ckk = sth/nrkj2; /* 10 */
1344 for (m = 0; (m < DIM); m++) /* 39 */
1346 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
1347 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
1348 f_j[m] = -f_i[m]-f_k[m];
1355 copy_ivec(SHIFT_IVEC(g, aj), jt);
1357 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1358 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1359 t1 = IVEC2IS(dt_ij);
1360 t2 = IVEC2IS(dt_kj);
1362 rvec_inc(fshift[t1], f_i);
1363 rvec_inc(fshift[CENTRAL], f_j);
1364 rvec_inc(fshift[t2], f_k);
1366 /* Time for the bond calculations */
1372 vtot += vbond; /* 1*/
1373 fbond *= gmx_invsqrt(dr2); /* 6 */
1377 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ak), dt_ik);
1378 ki = IVEC2IS(dt_ik);
1380 for (m = 0; (m < DIM); m++) /* 15 */
1382 fik = fbond*r_ik[m];
1385 fshift[ki][m] += fik;
1386 fshift[CENTRAL][m] -= fik;
1392 real quartic_angles(int nbonds,
1393 const t_iatom forceatoms[], const t_iparams forceparams[],
1394 const rvec x[], rvec f[], rvec fshift[],
1395 const t_pbc *pbc, const t_graph *g,
1396 real gmx_unused lambda, real gmx_unused *dvdlambda,
1397 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1398 int gmx_unused *global_atom_index)
1400 int i, j, ai, aj, ak, t1, t2, type;
1402 real cos_theta, cos_theta2, theta, dt, dVdt, va, dtp, c, vtot;
1403 ivec jt, dt_ij, dt_kj;
1406 for (i = 0; (i < nbonds); )
1408 type = forceatoms[i++];
1409 ai = forceatoms[i++];
1410 aj = forceatoms[i++];
1411 ak = forceatoms[i++];
1413 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
1414 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
1416 dt = theta - forceparams[type].qangle.theta*DEG2RAD; /* 2 */
1419 va = forceparams[type].qangle.c[0];
1421 for (j = 1; j <= 4; j++)
1423 c = forceparams[type].qangle.c[j];
1432 cos_theta2 = sqr(cos_theta); /* 1 */
1441 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1442 sth = st*cos_theta; /* 1 */
1446 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1447 theta*RAD2DEG, va, dVdt);
1450 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1451 nrij2 = iprod(r_ij, r_ij);
1453 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
1454 cii = sth/nrij2; /* 10 */
1455 ckk = sth/nrkj2; /* 10 */
1457 for (m = 0; (m < DIM); m++) /* 39 */
1459 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
1460 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
1461 f_j[m] = -f_i[m]-f_k[m];
1468 copy_ivec(SHIFT_IVEC(g, aj), jt);
1470 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1471 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1472 t1 = IVEC2IS(dt_ij);
1473 t2 = IVEC2IS(dt_kj);
1475 rvec_inc(fshift[t1], f_i);
1476 rvec_inc(fshift[CENTRAL], f_j);
1477 rvec_inc(fshift[t2], f_k);
1483 real dih_angle(const rvec xi, const rvec xj, const rvec xk, const rvec xl,
1485 rvec r_ij, rvec r_kj, rvec r_kl, rvec m, rvec n,
1486 real *sign, int *t1, int *t2, int *t3)
1490 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
1491 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
1492 *t3 = pbc_rvec_sub(pbc, xk, xl, r_kl); /* 3 */
1494 cprod(r_ij, r_kj, m); /* 9 */
1495 cprod(r_kj, r_kl, n); /* 9 */
1496 phi = gmx_angle(m, n); /* 49 (assuming 25 for atan2) */
1497 ipr = iprod(r_ij, n); /* 5 */
1498 (*sign) = (ipr < 0.0) ? -1.0 : 1.0;
1499 phi = (*sign)*phi; /* 1 */
1507 /* As dih_angle above, but calculates 4 dihedral angles at once using SIMD,
1508 * also calculates the pre-factor required for the dihedral force update.
1509 * Note that bv and buf should be register aligned.
1511 static gmx_inline void
1512 dih_angle_simd(const rvec *x,
1513 const int *ai, const int *aj, const int *ak, const int *al,
1514 const pbc_simd_t *pbc,
1517 gmx_mm_pr *mx_S, gmx_mm_pr *my_S, gmx_mm_pr *mz_S,
1518 gmx_mm_pr *nx_S, gmx_mm_pr *ny_S, gmx_mm_pr *nz_S,
1519 gmx_mm_pr *nrkj_m2_S,
1520 gmx_mm_pr *nrkj_n2_S,
1524 #define UNROLL GMX_SIMD_WIDTH_HERE
1526 gmx_mm_pr rijx_S, rijy_S, rijz_S;
1527 gmx_mm_pr rkjx_S, rkjy_S, rkjz_S;
1528 gmx_mm_pr rklx_S, rkly_S, rklz_S;
1529 gmx_mm_pr cx_S, cy_S, cz_S;
1533 gmx_mm_pr iprm_S, iprn_S;
1534 gmx_mm_pr nrkj2_S, nrkj_1_S, nrkj_2_S, nrkj_S;
1537 gmx_mm_pr nrkj2_min_S;
1538 gmx_mm_pr real_eps_S;
1540 /* Used to avoid division by zero.
1541 * We take into acount that we multiply the result by real_eps_S.
1543 nrkj2_min_S = gmx_set1_pr(GMX_REAL_MIN/(2*GMX_REAL_EPS));
1545 /* The value of the last significant bit (GMX_REAL_EPS is half of that) */
1546 real_eps_S = gmx_set1_pr(2*GMX_REAL_EPS);
1548 for (s = 0; s < UNROLL; s++)
1550 /* If you can't use pbc_dx_simd below for PBC, e.g. because
1551 * you can't round in SIMD, use pbc_rvec_sub here.
1553 for (m = 0; m < DIM; m++)
1555 dr[s + (0*DIM + m)*UNROLL] = x[ai[s]][m] - x[aj[s]][m];
1556 dr[s + (1*DIM + m)*UNROLL] = x[ak[s]][m] - x[aj[s]][m];
1557 dr[s + (2*DIM + m)*UNROLL] = x[ak[s]][m] - x[al[s]][m];
1561 rijx_S = gmx_load_pr(dr + 0*UNROLL);
1562 rijy_S = gmx_load_pr(dr + 1*UNROLL);
1563 rijz_S = gmx_load_pr(dr + 2*UNROLL);
1564 rkjx_S = gmx_load_pr(dr + 3*UNROLL);
1565 rkjy_S = gmx_load_pr(dr + 4*UNROLL);
1566 rkjz_S = gmx_load_pr(dr + 5*UNROLL);
1567 rklx_S = gmx_load_pr(dr + 6*UNROLL);
1568 rkly_S = gmx_load_pr(dr + 7*UNROLL);
1569 rklz_S = gmx_load_pr(dr + 8*UNROLL);
1571 pbc_dx_simd(&rijx_S, &rijy_S, &rijz_S, pbc);
1572 pbc_dx_simd(&rkjx_S, &rkjy_S, &rkjz_S, pbc);
1573 pbc_dx_simd(&rklx_S, &rkly_S, &rklz_S, pbc);
1575 gmx_cprod_pr(rijx_S, rijy_S, rijz_S,
1576 rkjx_S, rkjy_S, rkjz_S,
1579 gmx_cprod_pr(rkjx_S, rkjy_S, rkjz_S,
1580 rklx_S, rkly_S, rklz_S,
1583 gmx_cprod_pr(*mx_S, *my_S, *mz_S,
1584 *nx_S, *ny_S, *nz_S,
1585 &cx_S, &cy_S, &cz_S);
1587 cn_S = gmx_sqrt_pr(gmx_norm2_pr(cx_S, cy_S, cz_S));
1589 s_S = gmx_iprod_pr(*mx_S, *my_S, *mz_S, *nx_S, *ny_S, *nz_S);
1591 /* Determine the dihedral angle, the sign might need correction */
1592 *phi_S = gmx_atan2_pr(cn_S, s_S);
1594 ipr_S = gmx_iprod_pr(rijx_S, rijy_S, rijz_S,
1595 *nx_S, *ny_S, *nz_S);
1597 iprm_S = gmx_norm2_pr(*mx_S, *my_S, *mz_S);
1598 iprn_S = gmx_norm2_pr(*nx_S, *ny_S, *nz_S);
1600 nrkj2_S = gmx_norm2_pr(rkjx_S, rkjy_S, rkjz_S);
1602 /* Avoid division by zero. When zero, the result is multiplied by 0
1603 * anyhow, so the 3 max below do not affect the final result.
1605 nrkj2_S = gmx_max_pr(nrkj2_S, nrkj2_min_S);
1606 nrkj_1_S = gmx_invsqrt_pr(nrkj2_S);
1607 nrkj_2_S = gmx_mul_pr(nrkj_1_S, nrkj_1_S);
1608 nrkj_S = gmx_mul_pr(nrkj2_S, nrkj_1_S);
1610 toler_S = gmx_mul_pr(nrkj2_S, real_eps_S);
1612 /* Here the plain-C code uses a conditional, but we can't do that in SIMD.
1613 * So we take a max with the tolerance instead. Since we multiply with
1614 * m or n later, the max does not affect the results.
1616 iprm_S = gmx_max_pr(iprm_S, toler_S);
1617 iprn_S = gmx_max_pr(iprn_S, toler_S);
1618 *nrkj_m2_S = gmx_mul_pr(nrkj_S, gmx_inv_pr(iprm_S));
1619 *nrkj_n2_S = gmx_mul_pr(nrkj_S, gmx_inv_pr(iprn_S));
1621 /* Set sign of phi_S with the sign of ipr_S; phi_S is currently positive */
1622 *phi_S = gmx_cpsgn_nonneg_pr(ipr_S, *phi_S);
1624 p_S = gmx_iprod_pr(rijx_S, rijy_S, rijz_S,
1625 rkjx_S, rkjy_S, rkjz_S);
1626 p_S = gmx_mul_pr(p_S, nrkj_2_S);
1628 q_S = gmx_iprod_pr(rklx_S, rkly_S, rklz_S,
1629 rkjx_S, rkjy_S, rkjz_S);
1630 q_S = gmx_mul_pr(q_S, nrkj_2_S);
1632 gmx_store_pr(p, p_S);
1633 gmx_store_pr(q, q_S);
1637 #endif /* SIMD_BONDEDS */
1640 void do_dih_fup(int i, int j, int k, int l, real ddphi,
1641 rvec r_ij, rvec r_kj, rvec r_kl,
1642 rvec m, rvec n, rvec f[], rvec fshift[],
1643 const t_pbc *pbc, const t_graph *g,
1644 const rvec x[], int t1, int t2, int t3)
1647 rvec f_i, f_j, f_k, f_l;
1648 rvec uvec, vvec, svec, dx_jl;
1649 real iprm, iprn, nrkj, nrkj2, nrkj_1, nrkj_2;
1650 real a, b, p, q, toler;
1651 ivec jt, dt_ij, dt_kj, dt_lj;
1653 iprm = iprod(m, m); /* 5 */
1654 iprn = iprod(n, n); /* 5 */
1655 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1656 toler = nrkj2*GMX_REAL_EPS;
1657 if ((iprm > toler) && (iprn > toler))
1659 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1660 nrkj_2 = nrkj_1*nrkj_1; /* 1 */
1661 nrkj = nrkj2*nrkj_1; /* 1 */
1662 a = -ddphi*nrkj/iprm; /* 11 */
1663 svmul(a, m, f_i); /* 3 */
1664 b = ddphi*nrkj/iprn; /* 11 */
1665 svmul(b, n, f_l); /* 3 */
1666 p = iprod(r_ij, r_kj); /* 5 */
1667 p *= nrkj_2; /* 1 */
1668 q = iprod(r_kl, r_kj); /* 5 */
1669 q *= nrkj_2; /* 1 */
1670 svmul(p, f_i, uvec); /* 3 */
1671 svmul(q, f_l, vvec); /* 3 */
1672 rvec_sub(uvec, vvec, svec); /* 3 */
1673 rvec_sub(f_i, svec, f_j); /* 3 */
1674 rvec_add(f_l, svec, f_k); /* 3 */
1675 rvec_inc(f[i], f_i); /* 3 */
1676 rvec_dec(f[j], f_j); /* 3 */
1677 rvec_dec(f[k], f_k); /* 3 */
1678 rvec_inc(f[l], f_l); /* 3 */
1682 copy_ivec(SHIFT_IVEC(g, j), jt);
1683 ivec_sub(SHIFT_IVEC(g, i), jt, dt_ij);
1684 ivec_sub(SHIFT_IVEC(g, k), jt, dt_kj);
1685 ivec_sub(SHIFT_IVEC(g, l), jt, dt_lj);
1686 t1 = IVEC2IS(dt_ij);
1687 t2 = IVEC2IS(dt_kj);
1688 t3 = IVEC2IS(dt_lj);
1692 t3 = pbc_rvec_sub(pbc, x[l], x[j], dx_jl);
1699 rvec_inc(fshift[t1], f_i);
1700 rvec_dec(fshift[CENTRAL], f_j);
1701 rvec_dec(fshift[t2], f_k);
1702 rvec_inc(fshift[t3], f_l);
1707 /* As do_dih_fup above, but without shift forces */
1709 do_dih_fup_noshiftf(int i, int j, int k, int l, real ddphi,
1710 rvec r_ij, rvec r_kj, rvec r_kl,
1711 rvec m, rvec n, rvec f[])
1713 rvec f_i, f_j, f_k, f_l;
1714 rvec uvec, vvec, svec, dx_jl;
1715 real iprm, iprn, nrkj, nrkj2, nrkj_1, nrkj_2;
1716 real a, b, p, q, toler;
1717 ivec jt, dt_ij, dt_kj, dt_lj;
1719 iprm = iprod(m, m); /* 5 */
1720 iprn = iprod(n, n); /* 5 */
1721 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1722 toler = nrkj2*GMX_REAL_EPS;
1723 if ((iprm > toler) && (iprn > toler))
1725 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1726 nrkj_2 = nrkj_1*nrkj_1; /* 1 */
1727 nrkj = nrkj2*nrkj_1; /* 1 */
1728 a = -ddphi*nrkj/iprm; /* 11 */
1729 svmul(a, m, f_i); /* 3 */
1730 b = ddphi*nrkj/iprn; /* 11 */
1731 svmul(b, n, f_l); /* 3 */
1732 p = iprod(r_ij, r_kj); /* 5 */
1733 p *= nrkj_2; /* 1 */
1734 q = iprod(r_kl, r_kj); /* 5 */
1735 q *= nrkj_2; /* 1 */
1736 svmul(p, f_i, uvec); /* 3 */
1737 svmul(q, f_l, vvec); /* 3 */
1738 rvec_sub(uvec, vvec, svec); /* 3 */
1739 rvec_sub(f_i, svec, f_j); /* 3 */
1740 rvec_add(f_l, svec, f_k); /* 3 */
1741 rvec_inc(f[i], f_i); /* 3 */
1742 rvec_dec(f[j], f_j); /* 3 */
1743 rvec_dec(f[k], f_k); /* 3 */
1744 rvec_inc(f[l], f_l); /* 3 */
1748 /* As do_dih_fup_noshiftf above, but with pre-calculated pre-factors */
1749 static gmx_inline void
1750 do_dih_fup_noshiftf_precalc(int i, int j, int k, int l,
1752 real f_i_x, real f_i_y, real f_i_z,
1753 real mf_l_x, real mf_l_y, real mf_l_z,
1756 rvec f_i, f_j, f_k, f_l;
1757 rvec uvec, vvec, svec;
1765 svmul(p, f_i, uvec);
1766 svmul(q, f_l, vvec);
1767 rvec_sub(uvec, vvec, svec);
1768 rvec_sub(f_i, svec, f_j);
1769 rvec_add(f_l, svec, f_k);
1770 rvec_inc(f[i], f_i);
1771 rvec_dec(f[j], f_j);
1772 rvec_dec(f[k], f_k);
1773 rvec_inc(f[l], f_l);
1777 real dopdihs(real cpA, real cpB, real phiA, real phiB, int mult,
1778 real phi, real lambda, real *V, real *F)
1780 real v, dvdlambda, mdphi, v1, sdphi, ddphi;
1781 real L1 = 1.0 - lambda;
1782 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1783 real dph0 = (phiB - phiA)*DEG2RAD;
1784 real cp = L1*cpA + lambda*cpB;
1786 mdphi = mult*phi - ph0;
1788 ddphi = -cp*mult*sdphi;
1789 v1 = 1.0 + cos(mdphi);
1792 dvdlambda = (cpB - cpA)*v1 + cp*dph0*sdphi;
1799 /* That was 40 flops */
1803 dopdihs_noener(real cpA, real cpB, real phiA, real phiB, int mult,
1804 real phi, real lambda, real *F)
1806 real mdphi, sdphi, ddphi;
1807 real L1 = 1.0 - lambda;
1808 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1809 real cp = L1*cpA + lambda*cpB;
1811 mdphi = mult*phi - ph0;
1813 ddphi = -cp*mult*sdphi;
1817 /* That was 20 flops */
1821 dopdihs_mdphi(real cpA, real cpB, real phiA, real phiB, int mult,
1822 real phi, real lambda, real *cp, real *mdphi)
1824 real L1 = 1.0 - lambda;
1825 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1827 *cp = L1*cpA + lambda*cpB;
1829 *mdphi = mult*phi - ph0;
1832 static real dopdihs_min(real cpA, real cpB, real phiA, real phiB, int mult,
1833 real phi, real lambda, real *V, real *F)
1834 /* similar to dopdihs, except for a minus sign *
1835 * and a different treatment of mult/phi0 */
1837 real v, dvdlambda, mdphi, v1, sdphi, ddphi;
1838 real L1 = 1.0 - lambda;
1839 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1840 real dph0 = (phiB - phiA)*DEG2RAD;
1841 real cp = L1*cpA + lambda*cpB;
1843 mdphi = mult*(phi-ph0);
1845 ddphi = cp*mult*sdphi;
1846 v1 = 1.0-cos(mdphi);
1849 dvdlambda = (cpB-cpA)*v1 + cp*dph0*sdphi;
1856 /* That was 40 flops */
1859 real pdihs(int nbonds,
1860 const t_iatom forceatoms[], const t_iparams forceparams[],
1861 const rvec x[], rvec f[], rvec fshift[],
1862 const t_pbc *pbc, const t_graph *g,
1863 real lambda, real *dvdlambda,
1864 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1865 int gmx_unused *global_atom_index)
1867 int i, type, ai, aj, ak, al;
1869 rvec r_ij, r_kj, r_kl, m, n;
1870 real phi, sign, ddphi, vpd, vtot;
1874 for (i = 0; (i < nbonds); )
1876 type = forceatoms[i++];
1877 ai = forceatoms[i++];
1878 aj = forceatoms[i++];
1879 ak = forceatoms[i++];
1880 al = forceatoms[i++];
1882 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
1883 &sign, &t1, &t2, &t3); /* 84 */
1884 *dvdlambda += dopdihs(forceparams[type].pdihs.cpA,
1885 forceparams[type].pdihs.cpB,
1886 forceparams[type].pdihs.phiA,
1887 forceparams[type].pdihs.phiB,
1888 forceparams[type].pdihs.mult,
1889 phi, lambda, &vpd, &ddphi);
1892 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
1893 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
1896 fprintf(debug, "pdih: (%d,%d,%d,%d) phi=%g\n",
1897 ai, aj, ak, al, phi);
1904 void make_dp_periodic(real *dp) /* 1 flop? */
1906 /* dp cannot be outside (-pi,pi) */
1911 else if (*dp < -M_PI)
1918 /* As pdihs above, but without calculating energies and shift forces */
1920 pdihs_noener(int nbonds,
1921 const t_iatom forceatoms[], const t_iparams forceparams[],
1922 const rvec x[], rvec f[],
1923 const t_pbc gmx_unused *pbc, const t_graph gmx_unused *g,
1925 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1926 int gmx_unused *global_atom_index)
1928 int i, type, ai, aj, ak, al;
1930 rvec r_ij, r_kj, r_kl, m, n;
1931 real phi, sign, ddphi_tot, ddphi;
1933 for (i = 0; (i < nbonds); )
1935 ai = forceatoms[i+1];
1936 aj = forceatoms[i+2];
1937 ak = forceatoms[i+3];
1938 al = forceatoms[i+4];
1940 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
1941 &sign, &t1, &t2, &t3);
1945 /* Loop over dihedrals working on the same atoms,
1946 * so we avoid recalculating angles and force distributions.
1950 type = forceatoms[i];
1951 dopdihs_noener(forceparams[type].pdihs.cpA,
1952 forceparams[type].pdihs.cpB,
1953 forceparams[type].pdihs.phiA,
1954 forceparams[type].pdihs.phiB,
1955 forceparams[type].pdihs.mult,
1956 phi, lambda, &ddphi);
1961 while (i < nbonds &&
1962 forceatoms[i+1] == ai &&
1963 forceatoms[i+2] == aj &&
1964 forceatoms[i+3] == ak &&
1965 forceatoms[i+4] == al);
1967 do_dih_fup_noshiftf(ai, aj, ak, al, ddphi_tot, r_ij, r_kj, r_kl, m, n, f);
1974 /* As pdihs_noner above, but using SIMD to calculate many dihedrals at once */
1976 pdihs_noener_simd(int nbonds,
1977 const t_iatom forceatoms[], const t_iparams forceparams[],
1978 const rvec x[], rvec f[],
1979 const t_pbc *pbc, const t_graph gmx_unused *g,
1980 real gmx_unused lambda,
1981 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
1982 int gmx_unused *global_atom_index)
1984 #define UNROLL GMX_SIMD_WIDTH_HERE
1987 int type, ai[UNROLL], aj[UNROLL], ak[UNROLL], al[UNROLL];
1988 int t1[UNROLL], t2[UNROLL], t3[UNROLL];
1990 real dr_array[3*DIM*UNROLL+UNROLL], *dr;
1991 real buf_array[7*UNROLL+UNROLL], *buf;
1992 real *cp, *phi0, *mult, *phi, *p, *q, *sf_i, *msf_l;
1993 gmx_mm_pr phi0_S, phi_S;
1994 gmx_mm_pr mx_S, my_S, mz_S;
1995 gmx_mm_pr nx_S, ny_S, nz_S;
1996 gmx_mm_pr nrkj_m2_S, nrkj_n2_S;
1997 gmx_mm_pr cp_S, mdphi_S, mult_S;
1998 gmx_mm_pr sin_S, cos_S;
2000 gmx_mm_pr sf_i_S, msf_l_S;
2001 pbc_simd_t pbc_simd;
2003 /* Ensure SIMD register alignment */
2004 dr = gmx_simd_align_real(dr_array);
2005 buf = gmx_simd_align_real(buf_array);
2007 /* Extract aligned pointer for parameters and variables */
2008 cp = buf + 0*UNROLL;
2009 phi0 = buf + 1*UNROLL;
2010 mult = buf + 2*UNROLL;
2013 sf_i = buf + 5*UNROLL;
2014 msf_l = buf + 6*UNROLL;
2016 set_pbc_simd(pbc, &pbc_simd);
2018 /* nbonds is the number of dihedrals times nfa1, here we step UNROLL dihs */
2019 for (i = 0; (i < nbonds); i += UNROLL*nfa1)
2021 /* Collect atoms quadruplets for UNROLL dihedrals.
2022 * iu indexes into forceatoms, we should not let iu go beyond nbonds.
2025 for (s = 0; s < UNROLL; s++)
2027 type = forceatoms[iu];
2028 ai[s] = forceatoms[iu+1];
2029 aj[s] = forceatoms[iu+2];
2030 ak[s] = forceatoms[iu+3];
2031 al[s] = forceatoms[iu+4];
2033 cp[s] = forceparams[type].pdihs.cpA;
2034 phi0[s] = forceparams[type].pdihs.phiA*DEG2RAD;
2035 mult[s] = forceparams[type].pdihs.mult;
2037 /* At the end fill the arrays with identical entries */
2038 if (iu + nfa1 < nbonds)
2044 /* Caclulate UNROLL dihedral angles at once */
2045 dih_angle_simd(x, ai, aj, ak, al, &pbc_simd,
2048 &mx_S, &my_S, &mz_S,
2049 &nx_S, &ny_S, &nz_S,
2054 cp_S = gmx_load_pr(cp);
2055 phi0_S = gmx_load_pr(phi0);
2056 mult_S = gmx_load_pr(mult);
2058 mdphi_S = gmx_sub_pr(gmx_mul_pr(mult_S, phi_S), phi0_S);
2060 /* Calculate UNROLL sines at once */
2061 gmx_sincos_pr(mdphi_S, &sin_S, &cos_S);
2062 mddphi_S = gmx_mul_pr(gmx_mul_pr(cp_S, mult_S), sin_S);
2063 sf_i_S = gmx_mul_pr(mddphi_S, nrkj_m2_S);
2064 msf_l_S = gmx_mul_pr(mddphi_S, nrkj_n2_S);
2066 /* After this m?_S will contain f[i] */
2067 mx_S = gmx_mul_pr(sf_i_S, mx_S);
2068 my_S = gmx_mul_pr(sf_i_S, my_S);
2069 mz_S = gmx_mul_pr(sf_i_S, mz_S);
2071 /* After this m?_S will contain -f[l] */
2072 nx_S = gmx_mul_pr(msf_l_S, nx_S);
2073 ny_S = gmx_mul_pr(msf_l_S, ny_S);
2074 nz_S = gmx_mul_pr(msf_l_S, nz_S);
2076 gmx_store_pr(dr + 0*UNROLL, mx_S);
2077 gmx_store_pr(dr + 1*UNROLL, my_S);
2078 gmx_store_pr(dr + 2*UNROLL, mz_S);
2079 gmx_store_pr(dr + 3*UNROLL, nx_S);
2080 gmx_store_pr(dr + 4*UNROLL, ny_S);
2081 gmx_store_pr(dr + 5*UNROLL, nz_S);
2087 do_dih_fup_noshiftf_precalc(ai[s], aj[s], ak[s], al[s],
2092 dr[(DIM+XX)*UNROLL+s],
2093 dr[(DIM+YY)*UNROLL+s],
2094 dr[(DIM+ZZ)*UNROLL+s],
2099 while (s < UNROLL && iu < nbonds);
2104 #endif /* SIMD_BONDEDS */
2107 real idihs(int nbonds,
2108 const t_iatom forceatoms[], const t_iparams forceparams[],
2109 const rvec x[], rvec f[], rvec fshift[],
2110 const t_pbc *pbc, const t_graph *g,
2111 real lambda, real *dvdlambda,
2112 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2113 int gmx_unused *global_atom_index)
2115 int i, type, ai, aj, ak, al;
2117 real phi, phi0, dphi0, ddphi, sign, vtot;
2118 rvec r_ij, r_kj, r_kl, m, n;
2119 real L1, kk, dp, dp2, kA, kB, pA, pB, dvdl_term;
2124 for (i = 0; (i < nbonds); )
2126 type = forceatoms[i++];
2127 ai = forceatoms[i++];
2128 aj = forceatoms[i++];
2129 ak = forceatoms[i++];
2130 al = forceatoms[i++];
2132 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2133 &sign, &t1, &t2, &t3); /* 84 */
2135 /* phi can jump if phi0 is close to Pi/-Pi, which will cause huge
2136 * force changes if we just apply a normal harmonic.
2137 * Instead, we first calculate phi-phi0 and take it modulo (-Pi,Pi).
2138 * This means we will never have the periodicity problem, unless
2139 * the dihedral is Pi away from phiO, which is very unlikely due to
2142 kA = forceparams[type].harmonic.krA;
2143 kB = forceparams[type].harmonic.krB;
2144 pA = forceparams[type].harmonic.rA;
2145 pB = forceparams[type].harmonic.rB;
2147 kk = L1*kA + lambda*kB;
2148 phi0 = (L1*pA + lambda*pB)*DEG2RAD;
2149 dphi0 = (pB - pA)*DEG2RAD;
2153 make_dp_periodic(&dp);
2160 dvdl_term += 0.5*(kB - kA)*dp2 - kk*dphi0*dp;
2162 do_dih_fup(ai, aj, ak, al, (real)(-ddphi), r_ij, r_kj, r_kl, m, n,
2163 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2168 fprintf(debug, "idih: (%d,%d,%d,%d) phi=%g\n",
2169 ai, aj, ak, al, phi);
2174 *dvdlambda += dvdl_term;
2179 /*! \brief returns dx, rdist, and dpdl for functions posres() and fbposres()
2181 static void posres_dx(const rvec x, const rvec pos0A, const rvec pos0B,
2182 const rvec comA_sc, const rvec comB_sc,
2184 t_pbc *pbc, int refcoord_scaling, int npbcdim,
2185 rvec dx, rvec rdist, rvec dpdl)
2188 real posA, posB, L1, ref = 0.;
2193 for (m = 0; m < DIM; m++)
2199 switch (refcoord_scaling)
2203 rdist[m] = L1*posA + lambda*posB;
2204 dpdl[m] = posB - posA;
2207 /* Box relative coordinates are stored for dimensions with pbc */
2208 posA *= pbc->box[m][m];
2209 posB *= pbc->box[m][m];
2210 for (d = m+1; d < npbcdim; d++)
2212 posA += pos0A[d]*pbc->box[d][m];
2213 posB += pos0B[d]*pbc->box[d][m];
2215 ref = L1*posA + lambda*posB;
2217 dpdl[m] = posB - posA;
2220 ref = L1*comA_sc[m] + lambda*comB_sc[m];
2221 rdist[m] = L1*posA + lambda*posB;
2222 dpdl[m] = comB_sc[m] - comA_sc[m] + posB - posA;
2225 gmx_fatal(FARGS, "No such scaling method implemented");
2230 ref = L1*posA + lambda*posB;
2232 dpdl[m] = posB - posA;
2235 /* We do pbc_dx with ref+rdist,
2236 * since with only ref we can be up to half a box vector wrong.
2238 pos[m] = ref + rdist[m];
2243 pbc_dx(pbc, x, pos, dx);
2247 rvec_sub(x, pos, dx);
2251 /*! \brief Adds forces of flat-bottomed positions restraints to f[]
2252 * and fixes vir_diag. Returns the flat-bottomed potential. */
2253 real fbposres(int nbonds,
2254 const t_iatom forceatoms[], const t_iparams forceparams[],
2255 const rvec x[], rvec f[], rvec vir_diag,
2257 int refcoord_scaling, int ePBC, rvec com)
2258 /* compute flat-bottomed positions restraints */
2260 int i, ai, m, d, type, npbcdim = 0, fbdim;
2261 const t_iparams *pr;
2263 real ref = 0, dr, dr2, rpot, rfb, rfb2, fact, invdr;
2264 rvec com_sc, rdist, pos, dx, dpdl, fm;
2267 npbcdim = ePBC2npbcdim(ePBC);
2269 if (refcoord_scaling == erscCOM)
2272 for (m = 0; m < npbcdim; m++)
2274 for (d = m; d < npbcdim; d++)
2276 com_sc[m] += com[d]*pbc->box[d][m];
2282 for (i = 0; (i < nbonds); )
2284 type = forceatoms[i++];
2285 ai = forceatoms[i++];
2286 pr = &forceparams[type];
2288 /* same calculation as for normal posres, but with identical A and B states, and lambda==0 */
2289 posres_dx(x[ai], forceparams[type].fbposres.pos0, forceparams[type].fbposres.pos0,
2290 com_sc, com_sc, 0.0,
2291 pbc, refcoord_scaling, npbcdim,
2297 kk = pr->fbposres.k;
2298 rfb = pr->fbposres.r;
2301 /* with rfb<0, push particle out of the sphere/cylinder/layer */
2309 switch (pr->fbposres.geom)
2311 case efbposresSPHERE:
2312 /* spherical flat-bottom posres */
2315 ( (dr2 > rfb2 && bInvert == FALSE ) || (dr2 < rfb2 && bInvert == TRUE ) )
2319 v = 0.5*kk*sqr(dr - rfb);
2320 fact = -kk*(dr-rfb)/dr; /* Force pointing to the center pos0 */
2321 svmul(fact, dx, fm);
2324 case efbposresCYLINDER:
2325 /* cylidrical flat-bottom posres in x-y plane. fm[ZZ] = 0. */
2326 dr2 = sqr(dx[XX])+sqr(dx[YY]);
2328 ( (dr2 > rfb2 && bInvert == FALSE ) || (dr2 < rfb2 && bInvert == TRUE ) )
2333 v = 0.5*kk*sqr(dr - rfb);
2334 fm[XX] = -kk*(dr-rfb)*dx[XX]*invdr; /* Force pointing to the center */
2335 fm[YY] = -kk*(dr-rfb)*dx[YY]*invdr;
2338 case efbposresX: /* fbdim=XX */
2339 case efbposresY: /* fbdim=YY */
2340 case efbposresZ: /* fbdim=ZZ */
2341 /* 1D flat-bottom potential */
2342 fbdim = pr->fbposres.geom - efbposresX;
2344 if ( ( dr > rfb && bInvert == FALSE ) || ( 0 < dr && dr < rfb && bInvert == TRUE ) )
2346 v = 0.5*kk*sqr(dr - rfb);
2347 fm[fbdim] = -kk*(dr - rfb);
2349 else if ( (dr < (-rfb) && bInvert == FALSE ) || ( (-rfb) < dr && dr < 0 && bInvert == TRUE ))
2351 v = 0.5*kk*sqr(dr + rfb);
2352 fm[fbdim] = -kk*(dr + rfb);
2359 for (m = 0; (m < DIM); m++)
2362 /* Here we correct for the pbc_dx which included rdist */
2363 vir_diag[m] -= 0.5*(dx[m] + rdist[m])*fm[m];
2371 real posres(int nbonds,
2372 const t_iatom forceatoms[], const t_iparams forceparams[],
2373 const rvec x[], rvec f[], rvec vir_diag,
2375 real lambda, real *dvdlambda,
2376 int refcoord_scaling, int ePBC, rvec comA, rvec comB)
2378 int i, ai, m, d, type, ki, npbcdim = 0;
2379 const t_iparams *pr;
2382 real posA, posB, ref = 0;
2383 rvec comA_sc, comB_sc, rdist, dpdl, pos, dx;
2384 gmx_bool bForceValid = TRUE;
2386 if ((f == NULL) || (vir_diag == NULL)) /* should both be null together! */
2388 bForceValid = FALSE;
2391 npbcdim = ePBC2npbcdim(ePBC);
2393 if (refcoord_scaling == erscCOM)
2395 clear_rvec(comA_sc);
2396 clear_rvec(comB_sc);
2397 for (m = 0; m < npbcdim; m++)
2399 for (d = m; d < npbcdim; d++)
2401 comA_sc[m] += comA[d]*pbc->box[d][m];
2402 comB_sc[m] += comB[d]*pbc->box[d][m];
2410 for (i = 0; (i < nbonds); )
2412 type = forceatoms[i++];
2413 ai = forceatoms[i++];
2414 pr = &forceparams[type];
2416 /* return dx, rdist, and dpdl */
2417 posres_dx(x[ai], forceparams[type].posres.pos0A, forceparams[type].posres.pos0B,
2418 comA_sc, comB_sc, lambda,
2419 pbc, refcoord_scaling, npbcdim,
2422 for (m = 0; (m < DIM); m++)
2424 kk = L1*pr->posres.fcA[m] + lambda*pr->posres.fcB[m];
2426 vtot += 0.5*kk*dx[m]*dx[m];
2428 0.5*(pr->posres.fcB[m] - pr->posres.fcA[m])*dx[m]*dx[m]
2431 /* Here we correct for the pbc_dx which included rdist */
2435 vir_diag[m] -= 0.5*(dx[m] + rdist[m])*fm;
2443 static real low_angres(int nbonds,
2444 const t_iatom forceatoms[], const t_iparams forceparams[],
2445 const rvec x[], rvec f[], rvec fshift[],
2446 const t_pbc *pbc, const t_graph *g,
2447 real lambda, real *dvdlambda,
2450 int i, m, type, ai, aj, ak, al;
2452 real phi, cos_phi, cos_phi2, vid, vtot, dVdphi;
2453 rvec r_ij, r_kl, f_i, f_k = {0, 0, 0};
2454 real st, sth, nrij2, nrkl2, c, cij, ckl;
2457 t2 = 0; /* avoid warning with gcc-3.3. It is never used uninitialized */
2460 ak = al = 0; /* to avoid warnings */
2461 for (i = 0; i < nbonds; )
2463 type = forceatoms[i++];
2464 ai = forceatoms[i++];
2465 aj = forceatoms[i++];
2466 t1 = pbc_rvec_sub(pbc, x[aj], x[ai], r_ij); /* 3 */
2469 ak = forceatoms[i++];
2470 al = forceatoms[i++];
2471 t2 = pbc_rvec_sub(pbc, x[al], x[ak], r_kl); /* 3 */
2480 cos_phi = cos_angle(r_ij, r_kl); /* 25 */
2481 phi = acos(cos_phi); /* 10 */
2483 *dvdlambda += dopdihs_min(forceparams[type].pdihs.cpA,
2484 forceparams[type].pdihs.cpB,
2485 forceparams[type].pdihs.phiA,
2486 forceparams[type].pdihs.phiB,
2487 forceparams[type].pdihs.mult,
2488 phi, lambda, &vid, &dVdphi); /* 40 */
2492 cos_phi2 = sqr(cos_phi); /* 1 */
2495 st = -dVdphi*gmx_invsqrt(1 - cos_phi2); /* 12 */
2496 sth = st*cos_phi; /* 1 */
2497 nrij2 = iprod(r_ij, r_ij); /* 5 */
2498 nrkl2 = iprod(r_kl, r_kl); /* 5 */
2500 c = st*gmx_invsqrt(nrij2*nrkl2); /* 11 */
2501 cij = sth/nrij2; /* 10 */
2502 ckl = sth/nrkl2; /* 10 */
2504 for (m = 0; m < DIM; m++) /* 18+18 */
2506 f_i[m] = (c*r_kl[m]-cij*r_ij[m]);
2511 f_k[m] = (c*r_ij[m]-ckl*r_kl[m]);
2519 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
2522 rvec_inc(fshift[t1], f_i);
2523 rvec_dec(fshift[CENTRAL], f_i);
2528 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, al), dt);
2531 rvec_inc(fshift[t2], f_k);
2532 rvec_dec(fshift[CENTRAL], f_k);
2537 return vtot; /* 184 / 157 (bZAxis) total */
2540 real angres(int nbonds,
2541 const t_iatom forceatoms[], const t_iparams forceparams[],
2542 const rvec x[], rvec f[], rvec fshift[],
2543 const t_pbc *pbc, const t_graph *g,
2544 real lambda, real *dvdlambda,
2545 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2546 int gmx_unused *global_atom_index)
2548 return low_angres(nbonds, forceatoms, forceparams, x, f, fshift, pbc, g,
2549 lambda, dvdlambda, FALSE);
2552 real angresz(int nbonds,
2553 const t_iatom forceatoms[], const t_iparams forceparams[],
2554 const rvec x[], rvec f[], rvec fshift[],
2555 const t_pbc *pbc, const t_graph *g,
2556 real lambda, real *dvdlambda,
2557 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2558 int gmx_unused *global_atom_index)
2560 return low_angres(nbonds, forceatoms, forceparams, x, f, fshift, pbc, g,
2561 lambda, dvdlambda, TRUE);
2564 real dihres(int nbonds,
2565 const t_iatom forceatoms[], const t_iparams forceparams[],
2566 const rvec x[], rvec f[], rvec fshift[],
2567 const t_pbc *pbc, const t_graph *g,
2568 real lambda, real *dvdlambda,
2569 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2570 int gmx_unused *global_atom_index)
2573 int ai, aj, ak, al, i, k, type, t1, t2, t3;
2574 real phi0A, phi0B, dphiA, dphiB, kfacA, kfacB, phi0, dphi, kfac;
2575 real phi, ddphi, ddp, ddp2, dp, sign, d2r, fc, L1;
2576 rvec r_ij, r_kj, r_kl, m, n;
2583 for (i = 0; (i < nbonds); )
2585 type = forceatoms[i++];
2586 ai = forceatoms[i++];
2587 aj = forceatoms[i++];
2588 ak = forceatoms[i++];
2589 al = forceatoms[i++];
2591 phi0A = forceparams[type].dihres.phiA*d2r;
2592 dphiA = forceparams[type].dihres.dphiA*d2r;
2593 kfacA = forceparams[type].dihres.kfacA;
2595 phi0B = forceparams[type].dihres.phiB*d2r;
2596 dphiB = forceparams[type].dihres.dphiB*d2r;
2597 kfacB = forceparams[type].dihres.kfacB;
2599 phi0 = L1*phi0A + lambda*phi0B;
2600 dphi = L1*dphiA + lambda*dphiB;
2601 kfac = L1*kfacA + lambda*kfacB;
2603 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2604 &sign, &t1, &t2, &t3);
2609 fprintf(debug, "dihres[%d]: %d %d %d %d : phi=%f, dphi=%f, kfac=%f\n",
2610 k++, ai, aj, ak, al, phi0, dphi, kfac);
2612 /* phi can jump if phi0 is close to Pi/-Pi, which will cause huge
2613 * force changes if we just apply a normal harmonic.
2614 * Instead, we first calculate phi-phi0 and take it modulo (-Pi,Pi).
2615 * This means we will never have the periodicity problem, unless
2616 * the dihedral is Pi away from phiO, which is very unlikely due to
2620 make_dp_periodic(&dp);
2626 else if (dp < -dphi)
2638 vtot += 0.5*kfac*ddp2;
2641 *dvdlambda += 0.5*(kfacB - kfacA)*ddp2;
2642 /* lambda dependence from changing restraint distances */
2645 *dvdlambda -= kfac*ddp*((dphiB - dphiA)+(phi0B - phi0A));
2649 *dvdlambda += kfac*ddp*((dphiB - dphiA)-(phi0B - phi0A));
2651 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
2652 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2659 real unimplemented(int gmx_unused nbonds,
2660 const t_iatom gmx_unused forceatoms[], const t_iparams gmx_unused forceparams[],
2661 const rvec gmx_unused x[], rvec gmx_unused f[], rvec gmx_unused fshift[],
2662 const t_pbc gmx_unused *pbc, const t_graph gmx_unused *g,
2663 real gmx_unused lambda, real gmx_unused *dvdlambda,
2664 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2665 int gmx_unused *global_atom_index)
2667 gmx_impl("*** you are using a not implemented function");
2669 return 0.0; /* To make the compiler happy */
2672 real rbdihs(int nbonds,
2673 const t_iatom forceatoms[], const t_iparams forceparams[],
2674 const rvec x[], rvec f[], rvec fshift[],
2675 const t_pbc *pbc, const t_graph *g,
2676 real lambda, real *dvdlambda,
2677 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2678 int gmx_unused *global_atom_index)
2680 const real c0 = 0.0, c1 = 1.0, c2 = 2.0, c3 = 3.0, c4 = 4.0, c5 = 5.0;
2681 int type, ai, aj, ak, al, i, j;
2683 rvec r_ij, r_kj, r_kl, m, n;
2684 real parmA[NR_RBDIHS];
2685 real parmB[NR_RBDIHS];
2686 real parm[NR_RBDIHS];
2687 real cos_phi, phi, rbp, rbpBA;
2688 real v, sign, ddphi, sin_phi;
2690 real L1 = 1.0-lambda;
2694 for (i = 0; (i < nbonds); )
2696 type = forceatoms[i++];
2697 ai = forceatoms[i++];
2698 aj = forceatoms[i++];
2699 ak = forceatoms[i++];
2700 al = forceatoms[i++];
2702 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2703 &sign, &t1, &t2, &t3); /* 84 */
2705 /* Change to polymer convention */
2712 phi -= M_PI; /* 1 */
2716 /* Beware of accuracy loss, cannot use 1-sqrt(cos^2) ! */
2719 for (j = 0; (j < NR_RBDIHS); j++)
2721 parmA[j] = forceparams[type].rbdihs.rbcA[j];
2722 parmB[j] = forceparams[type].rbdihs.rbcB[j];
2723 parm[j] = L1*parmA[j]+lambda*parmB[j];
2725 /* Calculate cosine powers */
2726 /* Calculate the energy */
2727 /* Calculate the derivative */
2730 dvdl_term += (parmB[0]-parmA[0]);
2735 rbpBA = parmB[1]-parmA[1];
2736 ddphi += rbp*cosfac;
2739 dvdl_term += cosfac*rbpBA;
2741 rbpBA = parmB[2]-parmA[2];
2742 ddphi += c2*rbp*cosfac;
2745 dvdl_term += cosfac*rbpBA;
2747 rbpBA = parmB[3]-parmA[3];
2748 ddphi += c3*rbp*cosfac;
2751 dvdl_term += cosfac*rbpBA;
2753 rbpBA = parmB[4]-parmA[4];
2754 ddphi += c4*rbp*cosfac;
2757 dvdl_term += cosfac*rbpBA;
2759 rbpBA = parmB[5]-parmA[5];
2760 ddphi += c5*rbp*cosfac;
2763 dvdl_term += cosfac*rbpBA;
2765 ddphi = -ddphi*sin_phi; /* 11 */
2767 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
2768 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2771 *dvdlambda += dvdl_term;
2776 int cmap_setup_grid_index(int ip, int grid_spacing, int *ipm1, int *ipp1, int *ipp2)
2782 ip = ip + grid_spacing - 1;
2784 else if (ip > grid_spacing)
2786 ip = ip - grid_spacing - 1;
2795 im1 = grid_spacing - 1;
2797 else if (ip == grid_spacing-2)
2801 else if (ip == grid_spacing-1)
2815 real cmap_dihs(int nbonds,
2816 const t_iatom forceatoms[], const t_iparams forceparams[],
2817 const gmx_cmap_t *cmap_grid,
2818 const rvec x[], rvec f[], rvec fshift[],
2819 const t_pbc *pbc, const t_graph *g,
2820 real gmx_unused lambda, real gmx_unused *dvdlambda,
2821 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
2822 int gmx_unused *global_atom_index)
2824 int i, j, k, n, idx;
2825 int ai, aj, ak, al, am;
2826 int a1i, a1j, a1k, a1l, a2i, a2j, a2k, a2l;
2828 int t11, t21, t31, t12, t22, t32;
2829 int iphi1, ip1m1, ip1p1, ip1p2;
2830 int iphi2, ip2m1, ip2p1, ip2p2;
2832 int pos1, pos2, pos3, pos4, tmp;
2834 real ty[4], ty1[4], ty2[4], ty12[4], tc[16], tx[16];
2835 real phi1, psi1, cos_phi1, sin_phi1, sign1, xphi1;
2836 real phi2, psi2, cos_phi2, sin_phi2, sign2, xphi2;
2837 real dx, xx, tt, tu, e, df1, df2, ddf1, ddf2, ddf12, vtot;
2838 real ra21, rb21, rg21, rg1, rgr1, ra2r1, rb2r1, rabr1;
2839 real ra22, rb22, rg22, rg2, rgr2, ra2r2, rb2r2, rabr2;
2840 real fg1, hg1, fga1, hgb1, gaa1, gbb1;
2841 real fg2, hg2, fga2, hgb2, gaa2, gbb2;
2844 rvec r1_ij, r1_kj, r1_kl, m1, n1;
2845 rvec r2_ij, r2_kj, r2_kl, m2, n2;
2846 rvec f1_i, f1_j, f1_k, f1_l;
2847 rvec f2_i, f2_j, f2_k, f2_l;
2848 rvec a1, b1, a2, b2;
2849 rvec f1, g1, h1, f2, g2, h2;
2850 rvec dtf1, dtg1, dth1, dtf2, dtg2, dth2;
2851 ivec jt1, dt1_ij, dt1_kj, dt1_lj;
2852 ivec jt2, dt2_ij, dt2_kj, dt2_lj;
2856 int loop_index[4][4] = {
2863 /* Total CMAP energy */
2866 for (n = 0; n < nbonds; )
2868 /* Five atoms are involved in the two torsions */
2869 type = forceatoms[n++];
2870 ai = forceatoms[n++];
2871 aj = forceatoms[n++];
2872 ak = forceatoms[n++];
2873 al = forceatoms[n++];
2874 am = forceatoms[n++];
2876 /* Which CMAP type is this */
2877 cmapA = forceparams[type].cmap.cmapA;
2878 cmapd = cmap_grid->cmapdata[cmapA].cmap;
2886 phi1 = dih_angle(x[a1i], x[a1j], x[a1k], x[a1l], pbc, r1_ij, r1_kj, r1_kl, m1, n1,
2887 &sign1, &t11, &t21, &t31); /* 84 */
2889 cos_phi1 = cos(phi1);
2891 a1[0] = r1_ij[1]*r1_kj[2]-r1_ij[2]*r1_kj[1];
2892 a1[1] = r1_ij[2]*r1_kj[0]-r1_ij[0]*r1_kj[2];
2893 a1[2] = r1_ij[0]*r1_kj[1]-r1_ij[1]*r1_kj[0]; /* 9 */
2895 b1[0] = r1_kl[1]*r1_kj[2]-r1_kl[2]*r1_kj[1];
2896 b1[1] = r1_kl[2]*r1_kj[0]-r1_kl[0]*r1_kj[2];
2897 b1[2] = r1_kl[0]*r1_kj[1]-r1_kl[1]*r1_kj[0]; /* 9 */
2899 tmp = pbc_rvec_sub(pbc, x[a1l], x[a1k], h1);
2901 ra21 = iprod(a1, a1); /* 5 */
2902 rb21 = iprod(b1, b1); /* 5 */
2903 rg21 = iprod(r1_kj, r1_kj); /* 5 */
2909 rabr1 = sqrt(ra2r1*rb2r1);
2911 sin_phi1 = rg1 * rabr1 * iprod(a1, h1) * (-1);
2913 if (cos_phi1 < -0.5 || cos_phi1 > 0.5)
2915 phi1 = asin(sin_phi1);
2925 phi1 = -M_PI - phi1;
2931 phi1 = acos(cos_phi1);
2939 xphi1 = phi1 + M_PI; /* 1 */
2941 /* Second torsion */
2947 phi2 = dih_angle(x[a2i], x[a2j], x[a2k], x[a2l], pbc, r2_ij, r2_kj, r2_kl, m2, n2,
2948 &sign2, &t12, &t22, &t32); /* 84 */
2950 cos_phi2 = cos(phi2);
2952 a2[0] = r2_ij[1]*r2_kj[2]-r2_ij[2]*r2_kj[1];
2953 a2[1] = r2_ij[2]*r2_kj[0]-r2_ij[0]*r2_kj[2];
2954 a2[2] = r2_ij[0]*r2_kj[1]-r2_ij[1]*r2_kj[0]; /* 9 */
2956 b2[0] = r2_kl[1]*r2_kj[2]-r2_kl[2]*r2_kj[1];
2957 b2[1] = r2_kl[2]*r2_kj[0]-r2_kl[0]*r2_kj[2];
2958 b2[2] = r2_kl[0]*r2_kj[1]-r2_kl[1]*r2_kj[0]; /* 9 */
2960 tmp = pbc_rvec_sub(pbc, x[a2l], x[a2k], h2);
2962 ra22 = iprod(a2, a2); /* 5 */
2963 rb22 = iprod(b2, b2); /* 5 */
2964 rg22 = iprod(r2_kj, r2_kj); /* 5 */
2970 rabr2 = sqrt(ra2r2*rb2r2);
2972 sin_phi2 = rg2 * rabr2 * iprod(a2, h2) * (-1);
2974 if (cos_phi2 < -0.5 || cos_phi2 > 0.5)
2976 phi2 = asin(sin_phi2);
2986 phi2 = -M_PI - phi2;
2992 phi2 = acos(cos_phi2);
3000 xphi2 = phi2 + M_PI; /* 1 */
3002 /* Range mangling */
3005 xphi1 = xphi1 + 2*M_PI;
3007 else if (xphi1 >= 2*M_PI)
3009 xphi1 = xphi1 - 2*M_PI;
3014 xphi2 = xphi2 + 2*M_PI;
3016 else if (xphi2 >= 2*M_PI)
3018 xphi2 = xphi2 - 2*M_PI;
3021 /* Number of grid points */
3022 dx = 2*M_PI / cmap_grid->grid_spacing;
3024 /* Where on the grid are we */
3025 iphi1 = (int)(xphi1/dx);
3026 iphi2 = (int)(xphi2/dx);
3028 iphi1 = cmap_setup_grid_index(iphi1, cmap_grid->grid_spacing, &ip1m1, &ip1p1, &ip1p2);
3029 iphi2 = cmap_setup_grid_index(iphi2, cmap_grid->grid_spacing, &ip2m1, &ip2p1, &ip2p2);
3031 pos1 = iphi1*cmap_grid->grid_spacing+iphi2;
3032 pos2 = ip1p1*cmap_grid->grid_spacing+iphi2;
3033 pos3 = ip1p1*cmap_grid->grid_spacing+ip2p1;
3034 pos4 = iphi1*cmap_grid->grid_spacing+ip2p1;
3036 ty[0] = cmapd[pos1*4];
3037 ty[1] = cmapd[pos2*4];
3038 ty[2] = cmapd[pos3*4];
3039 ty[3] = cmapd[pos4*4];
3041 ty1[0] = cmapd[pos1*4+1];
3042 ty1[1] = cmapd[pos2*4+1];
3043 ty1[2] = cmapd[pos3*4+1];
3044 ty1[3] = cmapd[pos4*4+1];
3046 ty2[0] = cmapd[pos1*4+2];
3047 ty2[1] = cmapd[pos2*4+2];
3048 ty2[2] = cmapd[pos3*4+2];
3049 ty2[3] = cmapd[pos4*4+2];
3051 ty12[0] = cmapd[pos1*4+3];
3052 ty12[1] = cmapd[pos2*4+3];
3053 ty12[2] = cmapd[pos3*4+3];
3054 ty12[3] = cmapd[pos4*4+3];
3056 /* Switch to degrees */
3057 dx = 360.0 / cmap_grid->grid_spacing;
3058 xphi1 = xphi1 * RAD2DEG;
3059 xphi2 = xphi2 * RAD2DEG;
3061 for (i = 0; i < 4; i++) /* 16 */
3064 tx[i+4] = ty1[i]*dx;
3065 tx[i+8] = ty2[i]*dx;
3066 tx[i+12] = ty12[i]*dx*dx;
3070 for (i = 0; i < 4; i++) /* 1056 */
3072 for (j = 0; j < 4; j++)
3075 for (k = 0; k < 16; k++)
3077 xx = xx + cmap_coeff_matrix[k*16+idx]*tx[k];
3085 tt = (xphi1-iphi1*dx)/dx;
3086 tu = (xphi2-iphi2*dx)/dx;
3095 for (i = 3; i >= 0; i--)
3097 l1 = loop_index[i][3];
3098 l2 = loop_index[i][2];
3099 l3 = loop_index[i][1];
3101 e = tt * e + ((tc[i*4+3]*tu+tc[i*4+2])*tu + tc[i*4+1])*tu+tc[i*4];
3102 df1 = tu * df1 + (3.0*tc[l1]*tt+2.0*tc[l2])*tt+tc[l3];
3103 df2 = tt * df2 + (3.0*tc[i*4+3]*tu+2.0*tc[i*4+2])*tu+tc[i*4+1];
3104 ddf1 = tu * ddf1 + 2.0*3.0*tc[l1]*tt+2.0*tc[l2];
3105 ddf2 = tt * ddf2 + 2.0*3.0*tc[4*i+3]*tu+2.0*tc[4*i+2];
3108 ddf12 = tc[5] + 2.0*tc[9]*tt + 3.0*tc[13]*tt*tt + 2.0*tu*(tc[6]+2.0*tc[10]*tt+3.0*tc[14]*tt*tt) +
3109 3.0*tu*tu*(tc[7]+2.0*tc[11]*tt+3.0*tc[15]*tt*tt);
3114 ddf1 = ddf1 * fac * fac;
3115 ddf2 = ddf2 * fac * fac;
3116 ddf12 = ddf12 * fac * fac;
3121 /* Do forces - first torsion */
3122 fg1 = iprod(r1_ij, r1_kj);
3123 hg1 = iprod(r1_kl, r1_kj);
3124 fga1 = fg1*ra2r1*rgr1;
3125 hgb1 = hg1*rb2r1*rgr1;
3129 for (i = 0; i < DIM; i++)
3131 dtf1[i] = gaa1 * a1[i];
3132 dtg1[i] = fga1 * a1[i] - hgb1 * b1[i];
3133 dth1[i] = gbb1 * b1[i];
3135 f1[i] = df1 * dtf1[i];
3136 g1[i] = df1 * dtg1[i];
3137 h1[i] = df1 * dth1[i];
3140 f1_j[i] = -f1[i] - g1[i];
3141 f1_k[i] = h1[i] + g1[i];
3144 f[a1i][i] = f[a1i][i] + f1_i[i];
3145 f[a1j][i] = f[a1j][i] + f1_j[i]; /* - f1[i] - g1[i] */
3146 f[a1k][i] = f[a1k][i] + f1_k[i]; /* h1[i] + g1[i] */
3147 f[a1l][i] = f[a1l][i] + f1_l[i]; /* h1[i] */
3150 /* Do forces - second torsion */
3151 fg2 = iprod(r2_ij, r2_kj);
3152 hg2 = iprod(r2_kl, r2_kj);
3153 fga2 = fg2*ra2r2*rgr2;
3154 hgb2 = hg2*rb2r2*rgr2;
3158 for (i = 0; i < DIM; i++)
3160 dtf2[i] = gaa2 * a2[i];
3161 dtg2[i] = fga2 * a2[i] - hgb2 * b2[i];
3162 dth2[i] = gbb2 * b2[i];
3164 f2[i] = df2 * dtf2[i];
3165 g2[i] = df2 * dtg2[i];
3166 h2[i] = df2 * dth2[i];
3169 f2_j[i] = -f2[i] - g2[i];
3170 f2_k[i] = h2[i] + g2[i];
3173 f[a2i][i] = f[a2i][i] + f2_i[i]; /* f2[i] */
3174 f[a2j][i] = f[a2j][i] + f2_j[i]; /* - f2[i] - g2[i] */
3175 f[a2k][i] = f[a2k][i] + f2_k[i]; /* h2[i] + g2[i] */
3176 f[a2l][i] = f[a2l][i] + f2_l[i]; /* - h2[i] */
3182 copy_ivec(SHIFT_IVEC(g, a1j), jt1);
3183 ivec_sub(SHIFT_IVEC(g, a1i), jt1, dt1_ij);
3184 ivec_sub(SHIFT_IVEC(g, a1k), jt1, dt1_kj);
3185 ivec_sub(SHIFT_IVEC(g, a1l), jt1, dt1_lj);
3186 t11 = IVEC2IS(dt1_ij);
3187 t21 = IVEC2IS(dt1_kj);
3188 t31 = IVEC2IS(dt1_lj);
3190 copy_ivec(SHIFT_IVEC(g, a2j), jt2);
3191 ivec_sub(SHIFT_IVEC(g, a2i), jt2, dt2_ij);
3192 ivec_sub(SHIFT_IVEC(g, a2k), jt2, dt2_kj);
3193 ivec_sub(SHIFT_IVEC(g, a2l), jt2, dt2_lj);
3194 t12 = IVEC2IS(dt2_ij);
3195 t22 = IVEC2IS(dt2_kj);
3196 t32 = IVEC2IS(dt2_lj);
3200 t31 = pbc_rvec_sub(pbc, x[a1l], x[a1j], h1);
3201 t32 = pbc_rvec_sub(pbc, x[a2l], x[a2j], h2);
3209 rvec_inc(fshift[t11], f1_i);
3210 rvec_inc(fshift[CENTRAL], f1_j);
3211 rvec_inc(fshift[t21], f1_k);
3212 rvec_inc(fshift[t31], f1_l);
3214 rvec_inc(fshift[t21], f2_i);
3215 rvec_inc(fshift[CENTRAL], f2_j);
3216 rvec_inc(fshift[t22], f2_k);
3217 rvec_inc(fshift[t32], f2_l);
3224 /***********************************************************
3226 * G R O M O S 9 6 F U N C T I O N S
3228 ***********************************************************/
3229 real g96harmonic(real kA, real kB, real xA, real xB, real x, real lambda,
3232 const real half = 0.5;
3233 real L1, kk, x0, dx, dx2;
3234 real v, f, dvdlambda;
3237 kk = L1*kA+lambda*kB;
3238 x0 = L1*xA+lambda*xB;
3245 dvdlambda = half*(kB-kA)*dx2 + (xA-xB)*kk*dx;
3252 /* That was 21 flops */
3255 real g96bonds(int nbonds,
3256 const t_iatom forceatoms[], const t_iparams forceparams[],
3257 const rvec x[], rvec f[], rvec fshift[],
3258 const t_pbc *pbc, const t_graph *g,
3259 real lambda, real *dvdlambda,
3260 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3261 int gmx_unused *global_atom_index)
3263 int i, m, ki, ai, aj, type;
3264 real dr2, fbond, vbond, fij, vtot;
3269 for (i = 0; (i < nbonds); )
3271 type = forceatoms[i++];
3272 ai = forceatoms[i++];
3273 aj = forceatoms[i++];
3275 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
3276 dr2 = iprod(dx, dx); /* 5 */
3278 *dvdlambda += g96harmonic(forceparams[type].harmonic.krA,
3279 forceparams[type].harmonic.krB,
3280 forceparams[type].harmonic.rA,
3281 forceparams[type].harmonic.rB,
3282 dr2, lambda, &vbond, &fbond);
3284 vtot += 0.5*vbond; /* 1*/
3288 fprintf(debug, "G96-BONDS: dr = %10g vbond = %10g fbond = %10g\n",
3289 sqrt(dr2), vbond, fbond);
3295 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
3298 for (m = 0; (m < DIM); m++) /* 15 */
3303 fshift[ki][m] += fij;
3304 fshift[CENTRAL][m] -= fij;
3310 real g96bond_angle(const rvec xi, const rvec xj, const rvec xk, const t_pbc *pbc,
3311 rvec r_ij, rvec r_kj,
3313 /* Return value is the angle between the bonds i-j and j-k */
3317 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
3318 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
3320 costh = cos_angle(r_ij, r_kj); /* 25 */
3325 real g96angles(int nbonds,
3326 const t_iatom forceatoms[], const t_iparams forceparams[],
3327 const rvec x[], rvec f[], rvec fshift[],
3328 const t_pbc *pbc, const t_graph *g,
3329 real lambda, real *dvdlambda,
3330 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3331 int gmx_unused *global_atom_index)
3333 int i, ai, aj, ak, type, m, t1, t2;
3335 real cos_theta, dVdt, va, vtot;
3336 real rij_1, rij_2, rkj_1, rkj_2, rijrkj_1;
3338 ivec jt, dt_ij, dt_kj;
3341 for (i = 0; (i < nbonds); )
3343 type = forceatoms[i++];
3344 ai = forceatoms[i++];
3345 aj = forceatoms[i++];
3346 ak = forceatoms[i++];
3348 cos_theta = g96bond_angle(x[ai], x[aj], x[ak], pbc, r_ij, r_kj, &t1, &t2);
3350 *dvdlambda += g96harmonic(forceparams[type].harmonic.krA,
3351 forceparams[type].harmonic.krB,
3352 forceparams[type].harmonic.rA,
3353 forceparams[type].harmonic.rB,
3354 cos_theta, lambda, &va, &dVdt);
3357 rij_1 = gmx_invsqrt(iprod(r_ij, r_ij));
3358 rkj_1 = gmx_invsqrt(iprod(r_kj, r_kj));
3359 rij_2 = rij_1*rij_1;
3360 rkj_2 = rkj_1*rkj_1;
3361 rijrkj_1 = rij_1*rkj_1; /* 23 */
3366 fprintf(debug, "G96ANGLES: costheta = %10g vth = %10g dV/dct = %10g\n",
3367 cos_theta, va, dVdt);
3370 for (m = 0; (m < DIM); m++) /* 42 */
3372 f_i[m] = dVdt*(r_kj[m]*rijrkj_1 - r_ij[m]*rij_2*cos_theta);
3373 f_k[m] = dVdt*(r_ij[m]*rijrkj_1 - r_kj[m]*rkj_2*cos_theta);
3374 f_j[m] = -f_i[m]-f_k[m];
3382 copy_ivec(SHIFT_IVEC(g, aj), jt);
3384 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3385 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3386 t1 = IVEC2IS(dt_ij);
3387 t2 = IVEC2IS(dt_kj);
3389 rvec_inc(fshift[t1], f_i);
3390 rvec_inc(fshift[CENTRAL], f_j);
3391 rvec_inc(fshift[t2], f_k); /* 9 */
3397 real cross_bond_bond(int nbonds,
3398 const t_iatom forceatoms[], const t_iparams forceparams[],
3399 const rvec x[], rvec f[], rvec fshift[],
3400 const t_pbc *pbc, const t_graph *g,
3401 real gmx_unused lambda, real gmx_unused *dvdlambda,
3402 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3403 int gmx_unused *global_atom_index)
3405 /* Potential from Lawrence and Skimmer, Chem. Phys. Lett. 372 (2003)
3408 int i, ai, aj, ak, type, m, t1, t2;
3410 real vtot, vrr, s1, s2, r1, r2, r1e, r2e, krr;
3412 ivec jt, dt_ij, dt_kj;
3415 for (i = 0; (i < nbonds); )
3417 type = forceatoms[i++];
3418 ai = forceatoms[i++];
3419 aj = forceatoms[i++];
3420 ak = forceatoms[i++];
3421 r1e = forceparams[type].cross_bb.r1e;
3422 r2e = forceparams[type].cross_bb.r2e;
3423 krr = forceparams[type].cross_bb.krr;
3425 /* Compute distance vectors ... */
3426 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
3427 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
3429 /* ... and their lengths */
3433 /* Deviations from ideality */
3437 /* Energy (can be negative!) */
3442 svmul(-krr*s2/r1, r_ij, f_i);
3443 svmul(-krr*s1/r2, r_kj, f_k);
3445 for (m = 0; (m < DIM); m++) /* 12 */
3447 f_j[m] = -f_i[m] - f_k[m];
3456 copy_ivec(SHIFT_IVEC(g, aj), jt);
3458 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3459 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3460 t1 = IVEC2IS(dt_ij);
3461 t2 = IVEC2IS(dt_kj);
3463 rvec_inc(fshift[t1], f_i);
3464 rvec_inc(fshift[CENTRAL], f_j);
3465 rvec_inc(fshift[t2], f_k); /* 9 */
3471 real cross_bond_angle(int nbonds,
3472 const t_iatom forceatoms[], const t_iparams forceparams[],
3473 const rvec x[], rvec f[], rvec fshift[],
3474 const t_pbc *pbc, const t_graph *g,
3475 real gmx_unused lambda, real gmx_unused *dvdlambda,
3476 const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
3477 int gmx_unused *global_atom_index)
3479 /* Potential from Lawrence and Skimmer, Chem. Phys. Lett. 372 (2003)
3482 int i, ai, aj, ak, type, m, t1, t2, t3;
3483 rvec r_ij, r_kj, r_ik;
3484 real vtot, vrt, s1, s2, s3, r1, r2, r3, r1e, r2e, r3e, krt, k1, k2, k3;
3486 ivec jt, dt_ij, dt_kj;
3489 for (i = 0; (i < nbonds); )
3491 type = forceatoms[i++];
3492 ai = forceatoms[i++];
3493 aj = forceatoms[i++];
3494 ak = forceatoms[i++];
3495 r1e = forceparams[type].cross_ba.r1e;
3496 r2e = forceparams[type].cross_ba.r2e;
3497 r3e = forceparams[type].cross_ba.r3e;
3498 krt = forceparams[type].cross_ba.krt;
3500 /* Compute distance vectors ... */
3501 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
3502 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
3503 t3 = pbc_rvec_sub(pbc, x[ai], x[ak], r_ik);
3505 /* ... and their lengths */
3510 /* Deviations from ideality */
3515 /* Energy (can be negative!) */
3516 vrt = krt*s3*(s1+s2);
3522 k3 = -krt*(s1+s2)/r3;
3523 for (m = 0; (m < DIM); m++)
3525 f_i[m] = k1*r_ij[m] + k3*r_ik[m];
3526 f_k[m] = k2*r_kj[m] - k3*r_ik[m];
3527 f_j[m] = -f_i[m] - f_k[m];
3530 for (m = 0; (m < DIM); m++) /* 12 */
3540 copy_ivec(SHIFT_IVEC(g, aj), jt);
3542 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3543 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3544 t1 = IVEC2IS(dt_ij);
3545 t2 = IVEC2IS(dt_kj);
3547 rvec_inc(fshift[t1], f_i);
3548 rvec_inc(fshift[CENTRAL], f_j);
3549 rvec_inc(fshift[t2], f_k); /* 9 */
3555 static real bonded_tab(const char *type, int table_nr,
3556 const bondedtable_t *table, real kA, real kB, real r,
3557 real lambda, real *V, real *F)
3559 real k, tabscale, *VFtab, rt, eps, eps2, Yt, Ft, Geps, Heps2, Fp, VV, FF;
3561 real v, f, dvdlambda;
3563 k = (1.0 - lambda)*kA + lambda*kB;
3565 tabscale = table->scale;
3566 VFtab = table->data;
3572 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",
3573 type, table_nr, r, n0, n0+1, table->n);
3580 Geps = VFtab[nnn+2]*eps;
3581 Heps2 = VFtab[nnn+3]*eps2;
3582 Fp = Ft + Geps + Heps2;
3584 FF = Fp + Geps + 2.0*Heps2;
3586 *F = -k*FF*tabscale;
3588 dvdlambda = (kB - kA)*VV;
3592 /* That was 22 flops */
3595 real tab_bonds(int nbonds,
3596 const t_iatom forceatoms[], const t_iparams forceparams[],
3597 const rvec x[], rvec f[], rvec fshift[],
3598 const t_pbc *pbc, const t_graph *g,
3599 real lambda, real *dvdlambda,
3600 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
3601 int gmx_unused *global_atom_index)
3603 int i, m, ki, ai, aj, type, table;
3604 real dr, dr2, fbond, vbond, fij, vtot;
3609 for (i = 0; (i < nbonds); )
3611 type = forceatoms[i++];
3612 ai = forceatoms[i++];
3613 aj = forceatoms[i++];
3615 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
3616 dr2 = iprod(dx, dx); /* 5 */
3617 dr = dr2*gmx_invsqrt(dr2); /* 10 */
3619 table = forceparams[type].tab.table;
3621 *dvdlambda += bonded_tab("bond", table,
3622 &fcd->bondtab[table],
3623 forceparams[type].tab.kA,
3624 forceparams[type].tab.kB,
3625 dr, lambda, &vbond, &fbond); /* 22 */
3633 vtot += vbond; /* 1*/
3634 fbond *= gmx_invsqrt(dr2); /* 6 */
3638 fprintf(debug, "TABBONDS: dr = %10g vbond = %10g fbond = %10g\n",
3644 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
3647 for (m = 0; (m < DIM); m++) /* 15 */
3652 fshift[ki][m] += fij;
3653 fshift[CENTRAL][m] -= fij;
3659 real tab_angles(int nbonds,
3660 const t_iatom forceatoms[], const t_iparams forceparams[],
3661 const rvec x[], rvec f[], rvec fshift[],
3662 const t_pbc *pbc, const t_graph *g,
3663 real lambda, real *dvdlambda,
3664 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
3665 int gmx_unused *global_atom_index)
3667 int i, ai, aj, ak, t1, t2, type, table;
3669 real cos_theta, cos_theta2, theta, dVdt, va, vtot;
3670 ivec jt, dt_ij, dt_kj;
3673 for (i = 0; (i < nbonds); )
3675 type = forceatoms[i++];
3676 ai = forceatoms[i++];
3677 aj = forceatoms[i++];
3678 ak = forceatoms[i++];
3680 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
3681 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
3683 table = forceparams[type].tab.table;
3685 *dvdlambda += bonded_tab("angle", table,
3686 &fcd->angletab[table],
3687 forceparams[type].tab.kA,
3688 forceparams[type].tab.kB,
3689 theta, lambda, &va, &dVdt); /* 22 */
3692 cos_theta2 = sqr(cos_theta); /* 1 */
3701 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
3702 sth = st*cos_theta; /* 1 */
3706 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
3707 theta*RAD2DEG, va, dVdt);
3710 nrkj2 = iprod(r_kj, r_kj); /* 5 */
3711 nrij2 = iprod(r_ij, r_ij);
3713 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
3714 cii = sth/nrij2; /* 10 */
3715 ckk = sth/nrkj2; /* 10 */
3717 for (m = 0; (m < DIM); m++) /* 39 */
3719 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
3720 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
3721 f_j[m] = -f_i[m]-f_k[m];
3728 copy_ivec(SHIFT_IVEC(g, aj), jt);
3730 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3731 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3732 t1 = IVEC2IS(dt_ij);
3733 t2 = IVEC2IS(dt_kj);
3735 rvec_inc(fshift[t1], f_i);
3736 rvec_inc(fshift[CENTRAL], f_j);
3737 rvec_inc(fshift[t2], f_k);
3743 real tab_dihs(int nbonds,
3744 const t_iatom forceatoms[], const t_iparams forceparams[],
3745 const rvec x[], rvec f[], rvec fshift[],
3746 const t_pbc *pbc, const t_graph *g,
3747 real lambda, real *dvdlambda,
3748 const t_mdatoms gmx_unused *md, t_fcdata *fcd,
3749 int gmx_unused *global_atom_index)
3751 int i, type, ai, aj, ak, al, table;
3753 rvec r_ij, r_kj, r_kl, m, n;
3754 real phi, sign, ddphi, vpd, vtot;
3757 for (i = 0; (i < nbonds); )
3759 type = forceatoms[i++];
3760 ai = forceatoms[i++];
3761 aj = forceatoms[i++];
3762 ak = forceatoms[i++];
3763 al = forceatoms[i++];
3765 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
3766 &sign, &t1, &t2, &t3); /* 84 */
3768 table = forceparams[type].tab.table;
3770 /* Hopefully phi+M_PI never results in values < 0 */
3771 *dvdlambda += bonded_tab("dihedral", table,
3772 &fcd->dihtab[table],
3773 forceparams[type].tab.kA,
3774 forceparams[type].tab.kB,
3775 phi+M_PI, lambda, &vpd, &ddphi);
3778 do_dih_fup(ai, aj, ak, al, -ddphi, r_ij, r_kj, r_kl, m, n,
3779 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
3782 fprintf(debug, "pdih: (%d,%d,%d,%d) phi=%g\n",
3783 ai, aj, ak, al, phi);
3790 /* Return if this is a potential calculated in bondfree.c,
3791 * i.e. an interaction that actually calculates a potential and
3792 * works on multiple atoms (not e.g. a connection or a position restraint).
3794 static gmx_inline gmx_bool ftype_is_bonded_potential(int ftype)
3797 (interaction_function[ftype].flags & IF_BOND) &&
3798 !(ftype == F_CONNBONDS || ftype == F_POSRES) &&
3799 (ftype < F_GB12 || ftype > F_GB14);
3802 static void divide_bondeds_over_threads(t_idef *idef, int nthreads)
3809 idef->nthreads = nthreads;
3811 if (F_NRE*(nthreads+1) > idef->il_thread_division_nalloc)
3813 idef->il_thread_division_nalloc = F_NRE*(nthreads+1);
3814 snew(idef->il_thread_division, idef->il_thread_division_nalloc);
3817 for (ftype = 0; ftype < F_NRE; ftype++)
3819 if (ftype_is_bonded_potential(ftype))
3821 nat1 = interaction_function[ftype].nratoms + 1;
3823 for (t = 0; t <= nthreads; t++)
3825 /* Divide the interactions equally over the threads.
3826 * When the different types of bonded interactions
3827 * are distributed roughly equally over the threads,
3828 * this should lead to well localized output into
3829 * the force buffer on each thread.
3830 * If this is not the case, a more advanced scheme
3831 * (not implemented yet) will do better.
3833 il_nr_thread = (((idef->il[ftype].nr/nat1)*t)/nthreads)*nat1;
3835 /* Ensure that distance restraint pairs with the same label
3836 * end up on the same thread.
3837 * This is slighlty tricky code, since the next for iteration
3838 * may have an initial il_nr_thread lower than the final value
3839 * in the previous iteration, but this will anyhow be increased
3840 * to the approriate value again by this while loop.
3842 while (ftype == F_DISRES &&
3844 il_nr_thread < idef->il[ftype].nr &&
3845 idef->iparams[idef->il[ftype].iatoms[il_nr_thread]].disres.label ==
3846 idef->iparams[idef->il[ftype].iatoms[il_nr_thread-nat1]].disres.label)
3848 il_nr_thread += nat1;
3851 idef->il_thread_division[ftype*(nthreads+1)+t] = il_nr_thread;
3858 calc_bonded_reduction_mask(const t_idef *idef,
3863 int ftype, nb, nat1, nb0, nb1, i, a;
3867 for (ftype = 0; ftype < F_NRE; ftype++)
3869 if (ftype_is_bonded_potential(ftype))
3871 nb = idef->il[ftype].nr;
3874 nat1 = interaction_function[ftype].nratoms + 1;
3876 /* Divide this interaction equally over the threads.
3877 * This is not stored: should match division in calc_bonds.
3879 nb0 = idef->il_thread_division[ftype*(nt+1)+t];
3880 nb1 = idef->il_thread_division[ftype*(nt+1)+t+1];
3882 for (i = nb0; i < nb1; i += nat1)
3884 for (a = 1; a < nat1; a++)
3886 mask |= (1U << (idef->il[ftype].iatoms[i+a]>>shift));
3896 void setup_bonded_threading(t_forcerec *fr, t_idef *idef)
3898 #define MAX_BLOCK_BITS 32
3902 assert(fr->nthreads >= 1);
3904 /* Divide the bonded interaction over the threads */
3905 divide_bondeds_over_threads(idef, fr->nthreads);
3907 if (fr->nthreads == 1)
3914 /* We divide the force array in a maximum of 32 blocks.
3915 * Minimum force block reduction size is 2^6=64.
3918 while (fr->natoms_force > (int)(MAX_BLOCK_BITS*(1U<<fr->red_ashift)))
3924 fprintf(debug, "bonded force buffer block atom shift %d bits\n",
3928 /* Determine to which blocks each thread's bonded force calculation
3929 * contributes. Store this is a mask for each thread.
3931 #pragma omp parallel for num_threads(fr->nthreads) schedule(static)
3932 for (t = 1; t < fr->nthreads; t++)
3934 fr->f_t[t].red_mask =
3935 calc_bonded_reduction_mask(idef, fr->red_ashift, t, fr->nthreads);
3938 /* Determine the maximum number of blocks we need to reduce over */
3941 for (t = 0; t < fr->nthreads; t++)
3944 for (b = 0; b < MAX_BLOCK_BITS; b++)
3946 if (fr->f_t[t].red_mask & (1U<<b))
3948 fr->red_nblock = max(fr->red_nblock, b+1);
3954 fprintf(debug, "thread %d flags %x count %d\n",
3955 t, fr->f_t[t].red_mask, c);
3961 fprintf(debug, "Number of blocks to reduce: %d of size %d\n",
3962 fr->red_nblock, 1<<fr->red_ashift);
3963 fprintf(debug, "Reduction density %.2f density/#thread %.2f\n",
3964 ctot*(1<<fr->red_ashift)/(double)fr->natoms_force,
3965 ctot*(1<<fr->red_ashift)/(double)(fr->natoms_force*fr->nthreads));
3969 static void zero_thread_forces(f_thread_t *f_t, int n,
3970 int nblock, int blocksize)
3972 int b, a0, a1, a, i, j;
3974 if (n > f_t->f_nalloc)
3976 f_t->f_nalloc = over_alloc_large(n);
3977 srenew(f_t->f, f_t->f_nalloc);
3980 if (f_t->red_mask != 0)
3982 for (b = 0; b < nblock; b++)
3984 if (f_t->red_mask && (1U<<b))
3987 a1 = min((b+1)*blocksize, n);
3988 for (a = a0; a < a1; a++)
3990 clear_rvec(f_t->f[a]);
3995 for (i = 0; i < SHIFTS; i++)
3997 clear_rvec(f_t->fshift[i]);
3999 for (i = 0; i < F_NRE; i++)
4003 for (i = 0; i < egNR; i++)
4005 for (j = 0; j < f_t->grpp.nener; j++)
4007 f_t->grpp.ener[i][j] = 0;
4010 for (i = 0; i < efptNR; i++)
4016 static void reduce_thread_force_buffer(int n, rvec *f,
4017 int nthreads, f_thread_t *f_t,
4018 int nblock, int block_size)
4020 /* The max thread number is arbitrary,
4021 * we used a fixed number to avoid memory management.
4022 * Using more than 16 threads is probably never useful performance wise.
4024 #define MAX_BONDED_THREADS 256
4027 if (nthreads > MAX_BONDED_THREADS)
4029 gmx_fatal(FARGS, "Can not reduce bonded forces on more than %d threads",
4030 MAX_BONDED_THREADS);
4033 /* This reduction can run on any number of threads,
4034 * independently of nthreads.
4036 #pragma omp parallel for num_threads(nthreads) schedule(static)
4037 for (b = 0; b < nblock; b++)
4039 rvec *fp[MAX_BONDED_THREADS];
4043 /* Determine which threads contribute to this block */
4045 for (ft = 1; ft < nthreads; ft++)
4047 if (f_t[ft].red_mask & (1U<<b))
4049 fp[nfb++] = f_t[ft].f;
4054 /* Reduce force buffers for threads that contribute */
4056 a1 = (b+1)*block_size;
4058 for (a = a0; a < a1; a++)
4060 for (fb = 0; fb < nfb; fb++)
4062 rvec_inc(f[a], fp[fb][a]);
4069 static void reduce_thread_forces(int n, rvec *f, rvec *fshift,
4070 real *ener, gmx_grppairener_t *grpp, real *dvdl,
4071 int nthreads, f_thread_t *f_t,
4072 int nblock, int block_size,
4073 gmx_bool bCalcEnerVir,
4078 /* Reduce the bonded force buffer */
4079 reduce_thread_force_buffer(n, f, nthreads, f_t, nblock, block_size);
4082 /* When necessary, reduce energy and virial using one thread only */
4087 for (i = 0; i < SHIFTS; i++)
4089 for (t = 1; t < nthreads; t++)
4091 rvec_inc(fshift[i], f_t[t].fshift[i]);
4094 for (i = 0; i < F_NRE; i++)
4096 for (t = 1; t < nthreads; t++)
4098 ener[i] += f_t[t].ener[i];
4101 for (i = 0; i < egNR; i++)
4103 for (j = 0; j < f_t[1].grpp.nener; j++)
4105 for (t = 1; t < nthreads; t++)
4108 grpp->ener[i][j] += f_t[t].grpp.ener[i][j];
4114 for (i = 0; i < efptNR; i++)
4117 for (t = 1; t < nthreads; t++)
4119 dvdl[i] += f_t[t].dvdl[i];
4126 static real calc_one_bond(FILE *fplog, int thread,
4127 int ftype, const t_idef *idef,
4128 rvec x[], rvec f[], rvec fshift[],
4130 const t_pbc *pbc, const t_graph *g,
4131 gmx_grppairener_t *grpp,
4133 real *lambda, real *dvdl,
4134 const t_mdatoms *md, t_fcdata *fcd,
4135 gmx_bool bCalcEnerVir,
4136 int *global_atom_index, gmx_bool bPrintSepPot)
4138 int nat1, nbonds, efptFTYPE;
4143 if (IS_RESTRAINT_TYPE(ftype))
4145 efptFTYPE = efptRESTRAINT;
4149 efptFTYPE = efptBONDED;
4152 nat1 = interaction_function[ftype].nratoms + 1;
4153 nbonds = idef->il[ftype].nr/nat1;
4154 iatoms = idef->il[ftype].iatoms;
4156 nb0 = idef->il_thread_division[ftype*(idef->nthreads+1)+thread];
4157 nbn = idef->il_thread_division[ftype*(idef->nthreads+1)+thread+1] - nb0;
4159 if (!IS_LISTED_LJ_C(ftype))
4161 if (ftype == F_CMAP)
4163 v = cmap_dihs(nbn, iatoms+nb0,
4164 idef->iparams, &idef->cmap_grid,
4165 (const rvec*)x, f, fshift,
4166 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
4167 md, fcd, global_atom_index);
4170 else if (ftype == F_ANGLES &&
4171 !bCalcEnerVir && fr->efep == efepNO)
4173 /* No energies, shift forces, dvdl */
4174 angles_noener_simd(nbn, idef->il[ftype].iatoms+nb0,
4177 pbc, g, lambda[efptFTYPE], md, fcd,
4182 else if (ftype == F_PDIHS &&
4183 !bCalcEnerVir && fr->efep == efepNO)
4185 /* No energies, shift forces, dvdl */
4186 #ifndef SIMD_BONDEDS
4191 (nbn, idef->il[ftype].iatoms+nb0,
4194 pbc, g, lambda[efptFTYPE], md, fcd,
4200 v = interaction_function[ftype].ifunc(nbn, iatoms+nb0,
4202 (const rvec*)x, f, fshift,
4203 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
4204 md, fcd, global_atom_index);
4208 fprintf(fplog, " %-23s #%4d V %12.5e dVdl %12.5e\n",
4209 interaction_function[ftype].longname,
4210 nbonds, v, lambda[efptFTYPE]);
4215 v = do_nonbonded_listed(ftype, nbn, iatoms+nb0, idef->iparams, (const rvec*)x, f, fshift,
4216 pbc, g, lambda, dvdl, md, fr, grpp, global_atom_index);
4220 fprintf(fplog, " %-5s + %-15s #%4d dVdl %12.5e\n",
4221 interaction_function[ftype].longname,
4222 interaction_function[F_LJ14].longname, nbonds, dvdl[efptVDW]);
4223 fprintf(fplog, " %-5s + %-15s #%4d dVdl %12.5e\n",
4224 interaction_function[ftype].longname,
4225 interaction_function[F_COUL14].longname, nbonds, dvdl[efptCOUL]);
4231 inc_nrnb(nrnb, interaction_function[ftype].nrnb_ind, nbonds);
4237 void calc_bonds(FILE *fplog, const gmx_multisim_t *ms,
4239 rvec x[], history_t *hist,
4240 rvec f[], t_forcerec *fr,
4241 const t_pbc *pbc, const t_graph *g,
4242 gmx_enerdata_t *enerd, t_nrnb *nrnb,
4244 const t_mdatoms *md,
4245 t_fcdata *fcd, int *global_atom_index,
4246 t_atomtypes gmx_unused *atype, gmx_genborn_t gmx_unused *born,
4248 gmx_bool bPrintSepPot, gmx_int64_t step)
4250 gmx_bool bCalcEnerVir;
4252 real v, dvdl[efptNR], dvdl_dum[efptNR]; /* The dummy array is to have a place to store the dhdl at other values
4253 of lambda, which will be thrown away in the end*/
4254 const t_pbc *pbc_null;
4258 assert(fr->nthreads == idef->nthreads);
4260 bCalcEnerVir = (force_flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY));
4262 for (i = 0; i < efptNR; i++)
4276 fprintf(fplog, "Step %s: bonded V and dVdl for this node\n",
4277 gmx_step_str(step, buf));
4283 p_graph(debug, "Bondage is fun", g);
4287 /* Do pre force calculation stuff which might require communication */
4288 if (idef->il[F_ORIRES].nr)
4290 enerd->term[F_ORIRESDEV] =
4291 calc_orires_dev(ms, idef->il[F_ORIRES].nr,
4292 idef->il[F_ORIRES].iatoms,
4293 idef->iparams, md, (const rvec*)x,
4294 pbc_null, fcd, hist);
4296 if (idef->il[F_DISRES].nr)
4298 calc_disres_R_6(idef->il[F_DISRES].nr,
4299 idef->il[F_DISRES].iatoms,
4300 idef->iparams, (const rvec*)x, pbc_null,
4303 if (fcd->disres.nsystems > 1)
4305 gmx_sum_sim(2*fcd->disres.nres, fcd->disres.Rt_6, ms);
4310 #pragma omp parallel for num_threads(fr->nthreads) schedule(static)
4311 for (thread = 0; thread < fr->nthreads; thread++)
4318 gmx_grppairener_t *grpp;
4323 fshift = fr->fshift;
4325 grpp = &enerd->grpp;
4330 zero_thread_forces(&fr->f_t[thread], fr->natoms_force,
4331 fr->red_nblock, 1<<fr->red_ashift);
4333 ft = fr->f_t[thread].f;
4334 fshift = fr->f_t[thread].fshift;
4335 epot = fr->f_t[thread].ener;
4336 grpp = &fr->f_t[thread].grpp;
4337 dvdlt = fr->f_t[thread].dvdl;
4339 /* Loop over all bonded force types to calculate the bonded forces */
4340 for (ftype = 0; (ftype < F_NRE); ftype++)
4342 if (idef->il[ftype].nr > 0 && ftype_is_bonded_potential(ftype))
4344 v = calc_one_bond(fplog, thread, ftype, idef, x,
4345 ft, fshift, fr, pbc_null, g, grpp,
4346 nrnb, lambda, dvdlt,
4347 md, fcd, bCalcEnerVir,
4348 global_atom_index, bPrintSepPot);
4353 if (fr->nthreads > 1)
4355 reduce_thread_forces(fr->natoms_force, f, fr->fshift,
4356 enerd->term, &enerd->grpp, dvdl,
4357 fr->nthreads, fr->f_t,
4358 fr->red_nblock, 1<<fr->red_ashift,
4360 force_flags & GMX_FORCE_DHDL);
4362 if (force_flags & GMX_FORCE_DHDL)
4364 for (i = 0; i < efptNR; i++)
4366 enerd->dvdl_nonlin[i] += dvdl[i];
4370 /* Copy the sum of violations for the distance restraints from fcd */
4373 enerd->term[F_DISRESVIOL] = fcd->disres.sumviol;
4378 void calc_bonds_lambda(FILE *fplog,
4382 const t_pbc *pbc, const t_graph *g,
4383 gmx_grppairener_t *grpp, real *epot, t_nrnb *nrnb,
4385 const t_mdatoms *md,
4387 int *global_atom_index)
4389 int i, ftype, nr_nonperturbed, nr;
4391 real dvdl_dum[efptNR];
4393 const t_pbc *pbc_null;
4405 /* Copy the whole idef, so we can modify the contents locally */
4407 idef_fe.nthreads = 1;
4408 snew(idef_fe.il_thread_division, F_NRE*(idef_fe.nthreads+1));
4410 /* We already have the forces, so we use temp buffers here */
4411 snew(f, fr->natoms_force);
4412 snew(fshift, SHIFTS);
4414 /* Loop over all bonded force types to calculate the bonded energies */
4415 for (ftype = 0; (ftype < F_NRE); ftype++)
4417 if (ftype_is_bonded_potential(ftype))
4419 /* Set the work range of thread 0 to the perturbed bondeds only */
4420 nr_nonperturbed = idef->il[ftype].nr_nonperturbed;
4421 nr = idef->il[ftype].nr;
4422 idef_fe.il_thread_division[ftype*2+0] = nr_nonperturbed;
4423 idef_fe.il_thread_division[ftype*2+1] = nr;
4425 /* This is only to get the flop count correct */
4426 idef_fe.il[ftype].nr = nr - nr_nonperturbed;
4428 if (nr - nr_nonperturbed > 0)
4430 v = calc_one_bond(fplog, 0, ftype, &idef_fe,
4431 x, f, fshift, fr, pbc_null, g,
4432 grpp, nrnb, lambda, dvdl_dum,
4434 global_atom_index, FALSE);
4443 sfree(idef_fe.il_thread_division);