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
51 #include "gmx_fatal.h"
57 #include "nonbonded.h"
59 #if !defined GMX_DOUBLE && defined GMX_X86_SSE2
60 #include "gmx_x86_simd_single.h"
61 #define SSE_PROPER_DIHEDRALS
64 /* Find a better place for this? */
65 const int cmap_coeff_matrix[] = {
66 1, 0, -3, 2, 0, 0, 0, 0, -3, 0, 9, -6, 2, 0, -6, 4,
67 0, 0, 0, 0, 0, 0, 0, 0, 3, 0, -9, 6, -2, 0, 6, -4,
68 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 9, -6, 0, 0, -6, 4,
69 0, 0, 3, -2, 0, 0, 0, 0, 0, 0, -9, 6, 0, 0, 6, -4,
70 0, 0, 0, 0, 1, 0, -3, 2, -2, 0, 6, -4, 1, 0, -3, 2,
71 0, 0, 0, 0, 0, 0, 0, 0, -1, 0, 3, -2, 1, 0, -3, 2,
72 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -3, 2, 0, 0, 3, -2,
73 0, 0, 0, 0, 0, 0, 3, -2, 0, 0, -6, 4, 0, 0, 3, -2,
74 0, 1, -2, 1, 0, 0, 0, 0, 0, -3, 6, -3, 0, 2, -4, 2,
75 0, 0, 0, 0, 0, 0, 0, 0, 0, 3, -6, 3, 0, -2, 4, -2,
76 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -3, 3, 0, 0, 2, -2,
77 0, 0, -1, 1, 0, 0, 0, 0, 0, 0, 3, -3, 0, 0, -2, 2,
78 0, 0, 0, 0, 0, 1, -2, 1, 0, -2, 4, -2, 0, 1, -2, 1,
79 0, 0, 0, 0, 0, 0, 0, 0, 0, -1, 2, -1, 0, 1, -2, 1,
80 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, -1, 0, 0, -1, 1,
81 0, 0, 0, 0, 0, 0, -1, 1, 0, 0, 2, -2, 0, 0, -1, 1
86 int glatnr(int *global_atom_index, int i)
90 if (global_atom_index == NULL)
96 atnr = global_atom_index[i] + 1;
102 static int pbc_rvec_sub(const t_pbc *pbc, const rvec xi, const rvec xj, rvec dx)
106 return pbc_dx_aiuc(pbc, xi, xj, dx);
110 rvec_sub(xi, xj, dx);
116 * Morse potential bond by Frank Everdij
118 * Three parameters needed:
120 * b0 = equilibrium distance in nm
121 * be = beta in nm^-1 (actually, it's nu_e*Sqrt(2*pi*pi*mu/D_e))
122 * cb = well depth in kJ/mol
124 * Note: the potential is referenced to be +cb at infinite separation
125 * and zero at the equilibrium distance!
128 real morse_bonds(int nbonds,
129 const t_iatom forceatoms[], const t_iparams forceparams[],
130 const rvec x[], rvec f[], rvec fshift[],
131 const t_pbc *pbc, const t_graph *g,
132 real lambda, real *dvdlambda,
133 const t_mdatoms *md, t_fcdata *fcd,
134 int *global_atom_index)
136 const real one = 1.0;
137 const real two = 2.0;
138 real dr, dr2, temp, omtemp, cbomtemp, fbond, vbond, fij, vtot;
139 real b0, be, cb, b0A, beA, cbA, b0B, beB, cbB, L1;
141 int i, m, ki, type, ai, aj;
145 for (i = 0; (i < nbonds); )
147 type = forceatoms[i++];
148 ai = forceatoms[i++];
149 aj = forceatoms[i++];
151 b0A = forceparams[type].morse.b0A;
152 beA = forceparams[type].morse.betaA;
153 cbA = forceparams[type].morse.cbA;
155 b0B = forceparams[type].morse.b0B;
156 beB = forceparams[type].morse.betaB;
157 cbB = forceparams[type].morse.cbB;
159 L1 = one-lambda; /* 1 */
160 b0 = L1*b0A + lambda*b0B; /* 3 */
161 be = L1*beA + lambda*beB; /* 3 */
162 cb = L1*cbA + lambda*cbB; /* 3 */
164 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
165 dr2 = iprod(dx, dx); /* 5 */
166 dr = dr2*gmx_invsqrt(dr2); /* 10 */
167 temp = exp(-be*(dr-b0)); /* 12 */
171 /* bonds are constrainted. This may _not_ include bond constraints if they are lambda dependent */
172 *dvdlambda += cbB-cbA;
176 omtemp = one-temp; /* 1 */
177 cbomtemp = cb*omtemp; /* 1 */
178 vbond = cbomtemp*omtemp; /* 1 */
179 fbond = -two*be*temp*cbomtemp*gmx_invsqrt(dr2); /* 9 */
180 vtot += vbond; /* 1 */
182 *dvdlambda += (cbB - cbA) * omtemp * omtemp - (2-2*omtemp)*omtemp * cb * ((b0B-b0A)*be - (beB-beA)*(dr-b0)); /* 15 */
186 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
190 for (m = 0; (m < DIM); m++) /* 15 */
195 fshift[ki][m] += fij;
196 fshift[CENTRAL][m] -= fij;
202 real cubic_bonds(int nbonds,
203 const t_iatom forceatoms[], const t_iparams forceparams[],
204 const rvec x[], rvec f[], rvec fshift[],
205 const t_pbc *pbc, const t_graph *g,
206 real lambda, real *dvdlambda,
207 const t_mdatoms *md, t_fcdata *fcd,
208 int *global_atom_index)
210 const real three = 3.0;
211 const real two = 2.0;
213 real dr, dr2, dist, kdist, kdist2, fbond, vbond, fij, vtot;
215 int i, m, ki, type, ai, aj;
219 for (i = 0; (i < nbonds); )
221 type = forceatoms[i++];
222 ai = forceatoms[i++];
223 aj = forceatoms[i++];
225 b0 = forceparams[type].cubic.b0;
226 kb = forceparams[type].cubic.kb;
227 kcub = forceparams[type].cubic.kcub;
229 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
230 dr2 = iprod(dx, dx); /* 5 */
237 dr = dr2*gmx_invsqrt(dr2); /* 10 */
242 vbond = kdist2 + kcub*kdist2*dist;
243 fbond = -(two*kdist + three*kdist2*kcub)/dr;
245 vtot += vbond; /* 21 */
249 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
252 for (m = 0; (m < DIM); m++) /* 15 */
257 fshift[ki][m] += fij;
258 fshift[CENTRAL][m] -= fij;
264 real FENE_bonds(int nbonds,
265 const t_iatom forceatoms[], const t_iparams forceparams[],
266 const rvec x[], rvec f[], rvec fshift[],
267 const t_pbc *pbc, const t_graph *g,
268 real lambda, real *dvdlambda,
269 const t_mdatoms *md, t_fcdata *fcd,
270 int *global_atom_index)
272 const real half = 0.5;
273 const real one = 1.0;
275 real dr, dr2, bm2, omdr2obm2, fbond, vbond, fij, vtot;
277 int i, m, ki, type, ai, aj;
281 for (i = 0; (i < nbonds); )
283 type = forceatoms[i++];
284 ai = forceatoms[i++];
285 aj = forceatoms[i++];
287 bm = forceparams[type].fene.bm;
288 kb = forceparams[type].fene.kb;
290 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
291 dr2 = iprod(dx, dx); /* 5 */
303 "r^2 (%f) >= bm^2 (%f) in FENE bond between atoms %d and %d",
305 glatnr(global_atom_index, ai),
306 glatnr(global_atom_index, aj));
309 omdr2obm2 = one - dr2/bm2;
311 vbond = -half*kb*bm2*log(omdr2obm2);
312 fbond = -kb/omdr2obm2;
314 vtot += vbond; /* 35 */
318 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
321 for (m = 0; (m < DIM); m++) /* 15 */
326 fshift[ki][m] += fij;
327 fshift[CENTRAL][m] -= fij;
333 real harmonic(real kA, real kB, real xA, real xB, real x, real lambda,
336 const real half = 0.5;
337 real L1, kk, x0, dx, dx2;
338 real v, f, dvdlambda;
341 kk = L1*kA+lambda*kB;
342 x0 = L1*xA+lambda*xB;
349 dvdlambda = half*(kB-kA)*dx2 + (xA-xB)*kk*dx;
356 /* That was 19 flops */
360 real bonds(int nbonds,
361 const t_iatom forceatoms[], const t_iparams forceparams[],
362 const rvec x[], rvec f[], rvec fshift[],
363 const t_pbc *pbc, const t_graph *g,
364 real lambda, real *dvdlambda,
365 const t_mdatoms *md, t_fcdata *fcd,
366 int *global_atom_index)
368 int i, m, ki, ai, aj, type;
369 real dr, dr2, fbond, vbond, fij, vtot;
374 for (i = 0; (i < nbonds); )
376 type = forceatoms[i++];
377 ai = forceatoms[i++];
378 aj = forceatoms[i++];
380 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
381 dr2 = iprod(dx, dx); /* 5 */
382 dr = dr2*gmx_invsqrt(dr2); /* 10 */
384 *dvdlambda += harmonic(forceparams[type].harmonic.krA,
385 forceparams[type].harmonic.krB,
386 forceparams[type].harmonic.rA,
387 forceparams[type].harmonic.rB,
388 dr, lambda, &vbond, &fbond); /* 19 */
396 vtot += vbond; /* 1*/
397 fbond *= gmx_invsqrt(dr2); /* 6 */
401 fprintf(debug, "BONDS: dr = %10g vbond = %10g fbond = %10g\n",
407 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
410 for (m = 0; (m < DIM); m++) /* 15 */
415 fshift[ki][m] += fij;
416 fshift[CENTRAL][m] -= fij;
422 real restraint_bonds(int nbonds,
423 const t_iatom forceatoms[], const t_iparams forceparams[],
424 const rvec x[], rvec f[], rvec fshift[],
425 const t_pbc *pbc, const t_graph *g,
426 real lambda, real *dvdlambda,
427 const t_mdatoms *md, t_fcdata *fcd,
428 int *global_atom_index)
430 int i, m, ki, ai, aj, type;
431 real dr, dr2, fbond, vbond, fij, vtot;
433 real low, dlow, up1, dup1, up2, dup2, k, dk;
441 for (i = 0; (i < nbonds); )
443 type = forceatoms[i++];
444 ai = forceatoms[i++];
445 aj = forceatoms[i++];
447 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
448 dr2 = iprod(dx, dx); /* 5 */
449 dr = dr2*gmx_invsqrt(dr2); /* 10 */
451 low = L1*forceparams[type].restraint.lowA + lambda*forceparams[type].restraint.lowB;
452 dlow = -forceparams[type].restraint.lowA + forceparams[type].restraint.lowB;
453 up1 = L1*forceparams[type].restraint.up1A + lambda*forceparams[type].restraint.up1B;
454 dup1 = -forceparams[type].restraint.up1A + forceparams[type].restraint.up1B;
455 up2 = L1*forceparams[type].restraint.up2A + lambda*forceparams[type].restraint.up2B;
456 dup2 = -forceparams[type].restraint.up2A + forceparams[type].restraint.up2B;
457 k = L1*forceparams[type].restraint.kA + lambda*forceparams[type].restraint.kB;
458 dk = -forceparams[type].restraint.kA + forceparams[type].restraint.kB;
467 *dvdlambda += 0.5*dk*drh2 - k*dlow*drh;
480 *dvdlambda += 0.5*dk*drh2 - k*dup1*drh;
485 vbond = k*(up2 - up1)*(0.5*(up2 - up1) + drh);
486 fbond = -k*(up2 - up1);
487 *dvdlambda += dk*(up2 - up1)*(0.5*(up2 - up1) + drh)
488 + k*(dup2 - dup1)*(up2 - up1 + drh)
489 - k*(up2 - up1)*dup2;
497 vtot += vbond; /* 1*/
498 fbond *= gmx_invsqrt(dr2); /* 6 */
502 fprintf(debug, "BONDS: dr = %10g vbond = %10g fbond = %10g\n",
508 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
511 for (m = 0; (m < DIM); m++) /* 15 */
516 fshift[ki][m] += fij;
517 fshift[CENTRAL][m] -= fij;
524 real polarize(int nbonds,
525 const t_iatom forceatoms[], const t_iparams forceparams[],
526 const rvec x[], rvec f[], rvec fshift[],
527 const t_pbc *pbc, const t_graph *g,
528 real lambda, real *dvdlambda,
529 const t_mdatoms *md, t_fcdata *fcd,
530 int *global_atom_index)
532 int i, m, ki, ai, aj, type;
533 real dr, dr2, fbond, vbond, fij, vtot, ksh;
538 for (i = 0; (i < nbonds); )
540 type = forceatoms[i++];
541 ai = forceatoms[i++];
542 aj = forceatoms[i++];
543 ksh = sqr(md->chargeA[aj])*ONE_4PI_EPS0/forceparams[type].polarize.alpha;
546 fprintf(debug, "POL: local ai = %d aj = %d ksh = %.3f\n", ai, aj, ksh);
549 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
550 dr2 = iprod(dx, dx); /* 5 */
551 dr = dr2*gmx_invsqrt(dr2); /* 10 */
553 *dvdlambda += harmonic(ksh, ksh, 0, 0, dr, lambda, &vbond, &fbond); /* 19 */
560 vtot += vbond; /* 1*/
561 fbond *= gmx_invsqrt(dr2); /* 6 */
565 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
568 for (m = 0; (m < DIM); m++) /* 15 */
573 fshift[ki][m] += fij;
574 fshift[CENTRAL][m] -= fij;
580 real anharm_polarize(int nbonds,
581 const t_iatom forceatoms[], const t_iparams forceparams[],
582 const rvec x[], rvec f[], rvec fshift[],
583 const t_pbc *pbc, const t_graph *g,
584 real lambda, real *dvdlambda,
585 const t_mdatoms *md, t_fcdata *fcd,
586 int *global_atom_index)
588 int i, m, ki, ai, aj, type;
589 real dr, dr2, fbond, vbond, fij, vtot, ksh, khyp, drcut, ddr, ddr3;
594 for (i = 0; (i < nbonds); )
596 type = forceatoms[i++];
597 ai = forceatoms[i++];
598 aj = forceatoms[i++];
599 ksh = sqr(md->chargeA[aj])*ONE_4PI_EPS0/forceparams[type].anharm_polarize.alpha; /* 7*/
600 khyp = forceparams[type].anharm_polarize.khyp;
601 drcut = forceparams[type].anharm_polarize.drcut;
604 fprintf(debug, "POL: local ai = %d aj = %d ksh = %.3f\n", ai, aj, ksh);
607 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
608 dr2 = iprod(dx, dx); /* 5 */
609 dr = dr2*gmx_invsqrt(dr2); /* 10 */
611 *dvdlambda += harmonic(ksh, ksh, 0, 0, dr, lambda, &vbond, &fbond); /* 19 */
622 vbond += khyp*ddr*ddr3;
623 fbond -= 4*khyp*ddr3;
625 fbond *= gmx_invsqrt(dr2); /* 6 */
626 vtot += vbond; /* 1*/
630 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
633 for (m = 0; (m < DIM); m++) /* 15 */
638 fshift[ki][m] += fij;
639 fshift[CENTRAL][m] -= fij;
645 real water_pol(int nbonds,
646 const t_iatom forceatoms[], const t_iparams forceparams[],
647 const rvec x[], rvec f[], rvec fshift[],
648 const t_pbc *pbc, const t_graph *g,
649 real lambda, real *dvdlambda,
650 const t_mdatoms *md, t_fcdata *fcd,
651 int *global_atom_index)
653 /* This routine implements anisotropic polarizibility for water, through
654 * a shell connected to a dummy with spring constant that differ in the
655 * three spatial dimensions in the molecular frame.
657 int i, m, aO, aH1, aH2, aD, aS, type, type0;
658 rvec dOH1, dOH2, dHH, dOD, dDS, nW, kk, dx, kdx, proj;
662 real vtot, fij, r_HH, r_OD, r_nW, tx, ty, tz, qS;
667 type0 = forceatoms[0];
669 qS = md->chargeA[aS];
670 kk[XX] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_x;
671 kk[YY] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_y;
672 kk[ZZ] = sqr(qS)*ONE_4PI_EPS0/forceparams[type0].wpol.al_z;
673 r_HH = 1.0/forceparams[type0].wpol.rHH;
674 r_OD = 1.0/forceparams[type0].wpol.rOD;
677 fprintf(debug, "WPOL: qS = %10.5f aS = %5d\n", qS, aS);
678 fprintf(debug, "WPOL: kk = %10.3f %10.3f %10.3f\n",
679 kk[XX], kk[YY], kk[ZZ]);
680 fprintf(debug, "WPOL: rOH = %10.3f rHH = %10.3f rOD = %10.3f\n",
681 forceparams[type0].wpol.rOH,
682 forceparams[type0].wpol.rHH,
683 forceparams[type0].wpol.rOD);
685 for (i = 0; (i < nbonds); i += 6)
687 type = forceatoms[i];
690 gmx_fatal(FARGS, "Sorry, type = %d, type0 = %d, file = %s, line = %d",
691 type, type0, __FILE__, __LINE__);
693 aO = forceatoms[i+1];
694 aH1 = forceatoms[i+2];
695 aH2 = forceatoms[i+3];
696 aD = forceatoms[i+4];
697 aS = forceatoms[i+5];
699 /* Compute vectors describing the water frame */
700 rvec_sub(x[aH1], x[aO], dOH1);
701 rvec_sub(x[aH2], x[aO], dOH2);
702 rvec_sub(x[aH2], x[aH1], dHH);
703 rvec_sub(x[aD], x[aO], dOD);
704 rvec_sub(x[aS], x[aD], dDS);
705 cprod(dOH1, dOH2, nW);
707 /* Compute inverse length of normal vector
708 * (this one could be precomputed, but I'm too lazy now)
710 r_nW = gmx_invsqrt(iprod(nW, nW));
711 /* This is for precision, but does not make a big difference,
714 r_OD = gmx_invsqrt(iprod(dOD, dOD));
716 /* Normalize the vectors in the water frame */
718 svmul(r_HH, dHH, dHH);
719 svmul(r_OD, dOD, dOD);
721 /* Compute displacement of shell along components of the vector */
722 dx[ZZ] = iprod(dDS, dOD);
723 /* Compute projection on the XY plane: dDS - dx[ZZ]*dOD */
724 for (m = 0; (m < DIM); m++)
726 proj[m] = dDS[m]-dx[ZZ]*dOD[m];
729 /*dx[XX] = iprod(dDS,nW);
730 dx[YY] = iprod(dDS,dHH);*/
731 dx[XX] = iprod(proj, nW);
732 for (m = 0; (m < DIM); m++)
734 proj[m] -= dx[XX]*nW[m];
736 dx[YY] = iprod(proj, dHH);
741 fprintf(debug, "WPOL: dx2=%10g dy2=%10g dz2=%10g sum=%10g dDS^2=%10g\n",
742 sqr(dx[XX]), sqr(dx[YY]), sqr(dx[ZZ]), iprod(dx, dx), iprod(dDS, dDS));
743 fprintf(debug, "WPOL: dHH=(%10g,%10g,%10g)\n", dHH[XX], dHH[YY], dHH[ZZ]);
744 fprintf(debug, "WPOL: dOD=(%10g,%10g,%10g), 1/r_OD = %10g\n",
745 dOD[XX], dOD[YY], dOD[ZZ], 1/r_OD);
746 fprintf(debug, "WPOL: nW =(%10g,%10g,%10g), 1/r_nW = %10g\n",
747 nW[XX], nW[YY], nW[ZZ], 1/r_nW);
748 fprintf(debug, "WPOL: dx =%10g, dy =%10g, dz =%10g\n",
749 dx[XX], dx[YY], dx[ZZ]);
750 fprintf(debug, "WPOL: dDSx=%10g, dDSy=%10g, dDSz=%10g\n",
751 dDS[XX], dDS[YY], dDS[ZZ]);
754 /* Now compute the forces and energy */
755 kdx[XX] = kk[XX]*dx[XX];
756 kdx[YY] = kk[YY]*dx[YY];
757 kdx[ZZ] = kk[ZZ]*dx[ZZ];
758 vtot += iprod(dx, kdx);
759 for (m = 0; (m < DIM); m++)
761 /* This is a tensor operation but written out for speed */
775 fprintf(debug, "WPOL: vwpol=%g\n", 0.5*iprod(dx, kdx));
776 fprintf(debug, "WPOL: df = (%10g, %10g, %10g)\n", df[XX], df[YY], df[ZZ]);
784 static real do_1_thole(const rvec xi, const rvec xj, rvec fi, rvec fj,
785 const t_pbc *pbc, real qq,
786 rvec fshift[], real afac)
789 real r12sq, r12_1, r12n, r12bar, v0, v1, fscal, ebar, fff;
792 t = pbc_rvec_sub(pbc, xi, xj, r12); /* 3 */
794 r12sq = iprod(r12, r12); /* 5 */
795 r12_1 = gmx_invsqrt(r12sq); /* 5 */
796 r12bar = afac/r12_1; /* 5 */
797 v0 = qq*ONE_4PI_EPS0*r12_1; /* 2 */
798 ebar = exp(-r12bar); /* 5 */
799 v1 = (1-(1+0.5*r12bar)*ebar); /* 4 */
800 fscal = ((v0*r12_1)*v1 - v0*0.5*afac*ebar*(r12bar+1))*r12_1; /* 9 */
803 fprintf(debug, "THOLE: v0 = %.3f v1 = %.3f r12= % .3f r12bar = %.3f fscal = %.3f ebar = %.3f\n", v0, v1, 1/r12_1, r12bar, fscal, ebar);
806 for (m = 0; (m < DIM); m++)
812 fshift[CENTRAL][m] -= fff;
815 return v0*v1; /* 1 */
819 real thole_pol(int nbonds,
820 const t_iatom forceatoms[], const t_iparams forceparams[],
821 const rvec x[], rvec f[], rvec fshift[],
822 const t_pbc *pbc, const t_graph *g,
823 real lambda, real *dvdlambda,
824 const t_mdatoms *md, t_fcdata *fcd,
825 int *global_atom_index)
827 /* Interaction between two pairs of particles with opposite charge */
828 int i, type, a1, da1, a2, da2;
829 real q1, q2, qq, a, al1, al2, afac;
832 for (i = 0; (i < nbonds); )
834 type = forceatoms[i++];
835 a1 = forceatoms[i++];
836 da1 = forceatoms[i++];
837 a2 = forceatoms[i++];
838 da2 = forceatoms[i++];
839 q1 = md->chargeA[da1];
840 q2 = md->chargeA[da2];
841 a = forceparams[type].thole.a;
842 al1 = forceparams[type].thole.alpha1;
843 al2 = forceparams[type].thole.alpha2;
845 afac = a*pow(al1*al2, -1.0/6.0);
846 V += do_1_thole(x[a1], x[a2], f[a1], f[a2], pbc, qq, fshift, afac);
847 V += do_1_thole(x[da1], x[a2], f[da1], f[a2], pbc, -qq, fshift, afac);
848 V += do_1_thole(x[a1], x[da2], f[a1], f[da2], pbc, -qq, fshift, afac);
849 V += do_1_thole(x[da1], x[da2], f[da1], f[da2], pbc, qq, fshift, afac);
855 real bond_angle(const rvec xi, const rvec xj, const rvec xk, const t_pbc *pbc,
856 rvec r_ij, rvec r_kj, real *costh,
858 /* Return value is the angle between the bonds i-j and j-k */
863 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
864 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
866 *costh = cos_angle(r_ij, r_kj); /* 25 */
867 th = acos(*costh); /* 10 */
872 real angles(int nbonds,
873 const t_iatom forceatoms[], const t_iparams forceparams[],
874 const rvec x[], rvec f[], rvec fshift[],
875 const t_pbc *pbc, const t_graph *g,
876 real lambda, real *dvdlambda,
877 const t_mdatoms *md, t_fcdata *fcd,
878 int *global_atom_index)
880 int i, ai, aj, ak, t1, t2, type;
882 real cos_theta, cos_theta2, theta, dVdt, va, vtot;
883 ivec jt, dt_ij, dt_kj;
886 for (i = 0; i < nbonds; )
888 type = forceatoms[i++];
889 ai = forceatoms[i++];
890 aj = forceatoms[i++];
891 ak = forceatoms[i++];
893 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
894 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
896 *dvdlambda += harmonic(forceparams[type].harmonic.krA,
897 forceparams[type].harmonic.krB,
898 forceparams[type].harmonic.rA*DEG2RAD,
899 forceparams[type].harmonic.rB*DEG2RAD,
900 theta, lambda, &va, &dVdt); /* 21 */
903 cos_theta2 = sqr(cos_theta);
913 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
914 sth = st*cos_theta; /* 1 */
918 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
919 theta*RAD2DEG, va, dVdt);
922 nrij2 = iprod(r_ij, r_ij); /* 5 */
923 nrkj2 = iprod(r_kj, r_kj); /* 5 */
925 nrij_1 = gmx_invsqrt(nrij2); /* 10 */
926 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
928 cik = st*nrij_1*nrkj_1; /* 2 */
929 cii = sth*nrij_1*nrij_1; /* 2 */
930 ckk = sth*nrkj_1*nrkj_1; /* 2 */
932 for (m = 0; m < DIM; m++)
934 f_i[m] = -(cik*r_kj[m] - cii*r_ij[m]);
935 f_k[m] = -(cik*r_ij[m] - ckk*r_kj[m]);
936 f_j[m] = -f_i[m] - f_k[m];
943 copy_ivec(SHIFT_IVEC(g, aj), jt);
945 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
946 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
950 rvec_inc(fshift[t1], f_i);
951 rvec_inc(fshift[CENTRAL], f_j);
952 rvec_inc(fshift[t2], f_k);
959 real linear_angles(int nbonds,
960 const t_iatom forceatoms[], const t_iparams forceparams[],
961 const rvec x[], rvec f[], rvec fshift[],
962 const t_pbc *pbc, const t_graph *g,
963 real lambda, real *dvdlambda,
964 const t_mdatoms *md, t_fcdata *fcd,
965 int *global_atom_index)
967 int i, m, ai, aj, ak, t1, t2, type;
969 real L1, kA, kB, aA, aB, dr, dr2, va, vtot, a, b, klin;
970 ivec jt, dt_ij, dt_kj;
971 rvec r_ij, r_kj, r_ik, dx;
975 for (i = 0; (i < nbonds); )
977 type = forceatoms[i++];
978 ai = forceatoms[i++];
979 aj = forceatoms[i++];
980 ak = forceatoms[i++];
982 kA = forceparams[type].linangle.klinA;
983 kB = forceparams[type].linangle.klinB;
984 klin = L1*kA + lambda*kB;
986 aA = forceparams[type].linangle.aA;
987 aB = forceparams[type].linangle.aB;
991 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
992 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
993 rvec_sub(r_ij, r_kj, r_ik);
996 for (m = 0; (m < DIM); m++)
998 dr = -a * r_ij[m] - b * r_kj[m];
1003 f_j[m] = -(f_i[m]+f_k[m]);
1009 *dvdlambda += 0.5*(kB-kA)*dr2 + klin*(aB-aA)*iprod(dx, r_ik);
1015 copy_ivec(SHIFT_IVEC(g, aj), jt);
1017 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1018 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1019 t1 = IVEC2IS(dt_ij);
1020 t2 = IVEC2IS(dt_kj);
1022 rvec_inc(fshift[t1], f_i);
1023 rvec_inc(fshift[CENTRAL], f_j);
1024 rvec_inc(fshift[t2], f_k);
1029 real urey_bradley(int nbonds,
1030 const t_iatom forceatoms[], const t_iparams forceparams[],
1031 const rvec x[], rvec f[], rvec fshift[],
1032 const t_pbc *pbc, const t_graph *g,
1033 real lambda, real *dvdlambda,
1034 const t_mdatoms *md, t_fcdata *fcd,
1035 int *global_atom_index)
1037 int i, m, ai, aj, ak, t1, t2, type, ki;
1038 rvec r_ij, r_kj, r_ik;
1039 real cos_theta, cos_theta2, theta;
1040 real dVdt, va, vtot, dr, dr2, vbond, fbond, fik;
1041 real kthA, th0A, kUBA, r13A, kthB, th0B, kUBB, r13B;
1042 ivec jt, dt_ij, dt_kj, dt_ik;
1045 for (i = 0; (i < nbonds); )
1047 type = forceatoms[i++];
1048 ai = forceatoms[i++];
1049 aj = forceatoms[i++];
1050 ak = forceatoms[i++];
1051 th0A = forceparams[type].u_b.thetaA*DEG2RAD;
1052 kthA = forceparams[type].u_b.kthetaA;
1053 r13A = forceparams[type].u_b.r13A;
1054 kUBA = forceparams[type].u_b.kUBA;
1055 th0B = forceparams[type].u_b.thetaB*DEG2RAD;
1056 kthB = forceparams[type].u_b.kthetaB;
1057 r13B = forceparams[type].u_b.r13B;
1058 kUBB = forceparams[type].u_b.kUBB;
1060 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
1061 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
1063 *dvdlambda += harmonic(kthA, kthB, th0A, th0B, theta, lambda, &va, &dVdt); /* 21 */
1066 ki = pbc_rvec_sub(pbc, x[ai], x[ak], r_ik); /* 3 */
1067 dr2 = iprod(r_ik, r_ik); /* 5 */
1068 dr = dr2*gmx_invsqrt(dr2); /* 10 */
1070 *dvdlambda += harmonic(kUBA, kUBB, r13A, r13B, dr, lambda, &vbond, &fbond); /* 19 */
1072 cos_theta2 = sqr(cos_theta); /* 1 */
1080 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1081 sth = st*cos_theta; /* 1 */
1085 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1086 theta*RAD2DEG, va, dVdt);
1089 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1090 nrij2 = iprod(r_ij, r_ij);
1092 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
1093 cii = sth/nrij2; /* 10 */
1094 ckk = sth/nrkj2; /* 10 */
1096 for (m = 0; (m < DIM); m++) /* 39 */
1098 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
1099 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
1100 f_j[m] = -f_i[m]-f_k[m];
1107 copy_ivec(SHIFT_IVEC(g, aj), jt);
1109 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1110 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1111 t1 = IVEC2IS(dt_ij);
1112 t2 = IVEC2IS(dt_kj);
1114 rvec_inc(fshift[t1], f_i);
1115 rvec_inc(fshift[CENTRAL], f_j);
1116 rvec_inc(fshift[t2], f_k);
1118 /* Time for the bond calculations */
1124 vtot += vbond; /* 1*/
1125 fbond *= gmx_invsqrt(dr2); /* 6 */
1129 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, ak), dt_ik);
1130 ki = IVEC2IS(dt_ik);
1132 for (m = 0; (m < DIM); m++) /* 15 */
1134 fik = fbond*r_ik[m];
1137 fshift[ki][m] += fik;
1138 fshift[CENTRAL][m] -= fik;
1144 real quartic_angles(int nbonds,
1145 const t_iatom forceatoms[], const t_iparams forceparams[],
1146 const rvec x[], rvec f[], rvec fshift[],
1147 const t_pbc *pbc, const t_graph *g,
1148 real lambda, real *dvdlambda,
1149 const t_mdatoms *md, t_fcdata *fcd,
1150 int *global_atom_index)
1152 int i, j, ai, aj, ak, t1, t2, type;
1154 real cos_theta, cos_theta2, theta, dt, dVdt, va, dtp, c, vtot;
1155 ivec jt, dt_ij, dt_kj;
1158 for (i = 0; (i < nbonds); )
1160 type = forceatoms[i++];
1161 ai = forceatoms[i++];
1162 aj = forceatoms[i++];
1163 ak = forceatoms[i++];
1165 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
1166 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
1168 dt = theta - forceparams[type].qangle.theta*DEG2RAD; /* 2 */
1171 va = forceparams[type].qangle.c[0];
1173 for (j = 1; j <= 4; j++)
1175 c = forceparams[type].qangle.c[j];
1184 cos_theta2 = sqr(cos_theta); /* 1 */
1193 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
1194 sth = st*cos_theta; /* 1 */
1198 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
1199 theta*RAD2DEG, va, dVdt);
1202 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1203 nrij2 = iprod(r_ij, r_ij);
1205 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
1206 cii = sth/nrij2; /* 10 */
1207 ckk = sth/nrkj2; /* 10 */
1209 for (m = 0; (m < DIM); m++) /* 39 */
1211 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
1212 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
1213 f_j[m] = -f_i[m]-f_k[m];
1220 copy_ivec(SHIFT_IVEC(g, aj), jt);
1222 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
1223 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
1224 t1 = IVEC2IS(dt_ij);
1225 t2 = IVEC2IS(dt_kj);
1227 rvec_inc(fshift[t1], f_i);
1228 rvec_inc(fshift[CENTRAL], f_j);
1229 rvec_inc(fshift[t2], f_k);
1235 real dih_angle(const rvec xi, const rvec xj, const rvec xk, const rvec xl,
1237 rvec r_ij, rvec r_kj, rvec r_kl, rvec m, rvec n,
1238 real *sign, int *t1, int *t2, int *t3)
1242 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
1243 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
1244 *t3 = pbc_rvec_sub(pbc, xk, xl, r_kl); /* 3 */
1246 cprod(r_ij, r_kj, m); /* 9 */
1247 cprod(r_kj, r_kl, n); /* 9 */
1248 phi = gmx_angle(m, n); /* 49 (assuming 25 for atan2) */
1249 ipr = iprod(r_ij, n); /* 5 */
1250 (*sign) = (ipr < 0.0) ? -1.0 : 1.0;
1251 phi = (*sign)*phi; /* 1 */
1257 #ifdef SSE_PROPER_DIHEDRALS
1259 /* x86 SIMD inner-product of 4 float vectors */
1260 #define GMX_MM_IPROD_PS(ax, ay, az, bx, by, bz) \
1261 _mm_add_ps(_mm_add_ps(_mm_mul_ps(ax, bx), _mm_mul_ps(ay, by)), _mm_mul_ps(az, bz))
1263 /* x86 SIMD norm^2 of 4 float vectors */
1264 #define GMX_MM_NORM2_PS(ax, ay, az) GMX_MM_IPROD_PS(ax, ay, az, ax, ay, az)
1266 /* x86 SIMD cross-product of 4 float vectors */
1267 #define GMX_MM_CPROD_PS(ax, ay, az, bx, by, bz, cx, cy, cz) \
1269 cx = _mm_sub_ps(_mm_mul_ps(ay, bz), _mm_mul_ps(az, by)); \
1270 cy = _mm_sub_ps(_mm_mul_ps(az, bx), _mm_mul_ps(ax, bz)); \
1271 cz = _mm_sub_ps(_mm_mul_ps(ax, by), _mm_mul_ps(ay, bx)); \
1274 /* load 4 rvec's into 3 x86 SIMD float registers */
1275 #define load_rvec4(r0, r1, r2, r3, rx_SSE, ry_SSE, rz_SSE) \
1278 rx_SSE = _mm_load_ps(r0); \
1279 ry_SSE = _mm_load_ps(r1); \
1280 rz_SSE = _mm_load_ps(r2); \
1281 tmp = _mm_load_ps(r3); \
1282 _MM_TRANSPOSE4_PS(rx_SSE, ry_SSE, rz_SSE, tmp); \
1285 #define store_rvec4(rx_SSE, ry_SSE, rz_SSE, r0, r1, r2, r3) \
1287 __m128 tmp = _mm_setzero_ps(); \
1288 _MM_TRANSPOSE4_PS(rx_SSE, ry_SSE, rz_SSE, tmp); \
1289 _mm_store_ps(r0, rx_SSE); \
1290 _mm_store_ps(r1, ry_SSE); \
1291 _mm_store_ps(r2, rz_SSE); \
1292 _mm_store_ps(r3, tmp ); \
1295 /* An rvec in a structure which can be allocated 16-byte aligned */
1301 /* As dih_angle above, but calculates 4 dihedral angles at once using SSE,
1302 * also calculates the pre-factor required for the dihedral force update.
1303 * Note that bv and buf should be 16-byte aligned.
1306 dih_angle_sse(const rvec *x,
1307 int ai[4], int aj[4], int ak[4], int al[4],
1309 int t1[4], int t2[4], int t3[4],
1314 __m128 rijx_SSE, rijy_SSE, rijz_SSE;
1315 __m128 rkjx_SSE, rkjy_SSE, rkjz_SSE;
1316 __m128 rklx_SSE, rkly_SSE, rklz_SSE;
1317 __m128 mx_SSE, my_SSE, mz_SSE;
1318 __m128 nx_SSE, ny_SSE, nz_SSE;
1319 __m128 cx_SSE, cy_SSE, cz_SSE;
1325 __m128 iprm_SSE, iprn_SSE;
1326 __m128 nrkj2_SSE, nrkj_1_SSE, nrkj_2_SSE, nrkj_SSE;
1327 __m128 nrkj_m2_SSE, nrkj_n2_SSE;
1328 __m128 p_SSE, q_SSE;
1329 __m128 fmin_SSE = _mm_set1_ps(GMX_FLOAT_MIN);
1331 for (s = 0; s < 4; s++)
1333 t1[s] = pbc_rvec_sub(pbc, x[ai[s]], x[aj[s]], bv[0+s].v);
1334 t2[s] = pbc_rvec_sub(pbc, x[ak[s]], x[aj[s]], bv[4+s].v);
1335 t3[s] = pbc_rvec_sub(pbc, x[ak[s]], x[al[s]], bv[8+s].v);
1338 load_rvec4(bv[0].v, bv[1].v, bv[2].v, bv[3].v, rijx_SSE, rijy_SSE, rijz_SSE);
1339 load_rvec4(bv[4].v, bv[5].v, bv[6].v, bv[7].v, rkjx_SSE, rkjy_SSE, rkjz_SSE);
1340 load_rvec4(bv[8].v, bv[9].v, bv[10].v, bv[11].v, rklx_SSE, rkly_SSE, rklz_SSE);
1342 GMX_MM_CPROD_PS(rijx_SSE, rijy_SSE, rijz_SSE,
1343 rkjx_SSE, rkjy_SSE, rkjz_SSE,
1344 mx_SSE, my_SSE, mz_SSE);
1346 GMX_MM_CPROD_PS(rkjx_SSE, rkjy_SSE, rkjz_SSE,
1347 rklx_SSE, rkly_SSE, rklz_SSE,
1348 nx_SSE, ny_SSE, nz_SSE);
1350 GMX_MM_CPROD_PS(mx_SSE, my_SSE, mz_SSE,
1351 nx_SSE, ny_SSE, nz_SSE,
1352 cx_SSE, cy_SSE, cz_SSE);
1354 cn_SSE = gmx_mm_sqrt_ps(GMX_MM_NORM2_PS(cx_SSE, cy_SSE, cz_SSE));
1356 s_SSE = GMX_MM_IPROD_PS(mx_SSE, my_SSE, mz_SSE, nx_SSE, ny_SSE, nz_SSE);
1358 phi_SSE = gmx_mm_atan2_ps(cn_SSE, s_SSE);
1359 _mm_store_ps(buf+16, phi_SSE);
1361 ipr_SSE = GMX_MM_IPROD_PS(rijx_SSE, rijy_SSE, rijz_SSE,
1362 nx_SSE, ny_SSE, nz_SSE);
1364 signs = _mm_movemask_ps(ipr_SSE);
1366 for (s = 0; s < 4; s++)
1370 buf[16+s] = -buf[16+s];
1374 iprm_SSE = GMX_MM_NORM2_PS(mx_SSE, my_SSE, mz_SSE);
1375 iprn_SSE = GMX_MM_NORM2_PS(nx_SSE, ny_SSE, nz_SSE);
1377 /* store_rvec4 messes with the input, don't use it after this! */
1378 store_rvec4(mx_SSE, my_SSE, mz_SSE, bv[0].v, bv[1].v, bv[2].v, bv[3].v);
1379 store_rvec4(nx_SSE, ny_SSE, nz_SSE, bv[4].v, bv[5].v, bv[6].v, bv[7].v);
1381 nrkj2_SSE = GMX_MM_NORM2_PS(rkjx_SSE, rkjy_SSE, rkjz_SSE);
1383 /* Avoid division by zero. When zero, the result is multiplied by 0
1384 * anyhow, so the 3 max below do not affect the final result.
1386 nrkj2_SSE = _mm_max_ps(nrkj2_SSE, fmin_SSE);
1387 nrkj_1_SSE = gmx_mm_invsqrt_ps(nrkj2_SSE);
1388 nrkj_2_SSE = _mm_mul_ps(nrkj_1_SSE, nrkj_1_SSE);
1389 nrkj_SSE = _mm_mul_ps(nrkj2_SSE, nrkj_1_SSE);
1391 iprm_SSE = _mm_max_ps(iprm_SSE, fmin_SSE);
1392 iprn_SSE = _mm_max_ps(iprn_SSE, fmin_SSE);
1393 nrkj_m2_SSE = _mm_mul_ps(nrkj_SSE, gmx_mm_inv_ps(iprm_SSE));
1394 nrkj_n2_SSE = _mm_mul_ps(nrkj_SSE, gmx_mm_inv_ps(iprn_SSE));
1396 _mm_store_ps(buf+0, nrkj_m2_SSE);
1397 _mm_store_ps(buf+4, nrkj_n2_SSE);
1399 p_SSE = GMX_MM_IPROD_PS(rijx_SSE, rijy_SSE, rijz_SSE,
1400 rkjx_SSE, rkjy_SSE, rkjz_SSE);
1401 p_SSE = _mm_mul_ps(p_SSE, nrkj_2_SSE);
1403 q_SSE = GMX_MM_IPROD_PS(rklx_SSE, rkly_SSE, rklz_SSE,
1404 rkjx_SSE, rkjy_SSE, rkjz_SSE);
1405 q_SSE = _mm_mul_ps(q_SSE, nrkj_2_SSE);
1407 _mm_store_ps(buf+8, p_SSE);
1408 _mm_store_ps(buf+12, q_SSE);
1411 #endif /* SSE_PROPER_DIHEDRALS */
1414 void do_dih_fup(int i, int j, int k, int l, real ddphi,
1415 rvec r_ij, rvec r_kj, rvec r_kl,
1416 rvec m, rvec n, rvec f[], rvec fshift[],
1417 const t_pbc *pbc, const t_graph *g,
1418 const rvec x[], int t1, int t2, int t3)
1421 rvec f_i, f_j, f_k, f_l;
1422 rvec uvec, vvec, svec, dx_jl;
1423 real iprm, iprn, nrkj, nrkj2, nrkj_1, nrkj_2;
1424 real a, b, p, q, toler;
1425 ivec jt, dt_ij, dt_kj, dt_lj;
1427 iprm = iprod(m, m); /* 5 */
1428 iprn = iprod(n, n); /* 5 */
1429 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1430 toler = nrkj2*GMX_REAL_EPS;
1431 if ((iprm > toler) && (iprn > toler))
1433 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1434 nrkj_2 = nrkj_1*nrkj_1; /* 1 */
1435 nrkj = nrkj2*nrkj_1; /* 1 */
1436 a = -ddphi*nrkj/iprm; /* 11 */
1437 svmul(a, m, f_i); /* 3 */
1438 b = ddphi*nrkj/iprn; /* 11 */
1439 svmul(b, n, f_l); /* 3 */
1440 p = iprod(r_ij, r_kj); /* 5 */
1441 p *= nrkj_2; /* 1 */
1442 q = iprod(r_kl, r_kj); /* 5 */
1443 q *= nrkj_2; /* 1 */
1444 svmul(p, f_i, uvec); /* 3 */
1445 svmul(q, f_l, vvec); /* 3 */
1446 rvec_sub(uvec, vvec, svec); /* 3 */
1447 rvec_sub(f_i, svec, f_j); /* 3 */
1448 rvec_add(f_l, svec, f_k); /* 3 */
1449 rvec_inc(f[i], f_i); /* 3 */
1450 rvec_dec(f[j], f_j); /* 3 */
1451 rvec_dec(f[k], f_k); /* 3 */
1452 rvec_inc(f[l], f_l); /* 3 */
1456 copy_ivec(SHIFT_IVEC(g, j), jt);
1457 ivec_sub(SHIFT_IVEC(g, i), jt, dt_ij);
1458 ivec_sub(SHIFT_IVEC(g, k), jt, dt_kj);
1459 ivec_sub(SHIFT_IVEC(g, l), jt, dt_lj);
1460 t1 = IVEC2IS(dt_ij);
1461 t2 = IVEC2IS(dt_kj);
1462 t3 = IVEC2IS(dt_lj);
1466 t3 = pbc_rvec_sub(pbc, x[l], x[j], dx_jl);
1473 rvec_inc(fshift[t1], f_i);
1474 rvec_dec(fshift[CENTRAL], f_j);
1475 rvec_dec(fshift[t2], f_k);
1476 rvec_inc(fshift[t3], f_l);
1481 /* As do_dih_fup above, but without shift forces */
1483 do_dih_fup_noshiftf(int i, int j, int k, int l, real ddphi,
1484 rvec r_ij, rvec r_kj, rvec r_kl,
1485 rvec m, rvec n, rvec f[])
1487 rvec f_i, f_j, f_k, f_l;
1488 rvec uvec, vvec, svec, dx_jl;
1489 real iprm, iprn, nrkj, nrkj2, nrkj_1, nrkj_2;
1490 real a, b, p, q, toler;
1491 ivec jt, dt_ij, dt_kj, dt_lj;
1493 iprm = iprod(m, m); /* 5 */
1494 iprn = iprod(n, n); /* 5 */
1495 nrkj2 = iprod(r_kj, r_kj); /* 5 */
1496 toler = nrkj2*GMX_REAL_EPS;
1497 if ((iprm > toler) && (iprn > toler))
1499 nrkj_1 = gmx_invsqrt(nrkj2); /* 10 */
1500 nrkj_2 = nrkj_1*nrkj_1; /* 1 */
1501 nrkj = nrkj2*nrkj_1; /* 1 */
1502 a = -ddphi*nrkj/iprm; /* 11 */
1503 svmul(a, m, f_i); /* 3 */
1504 b = ddphi*nrkj/iprn; /* 11 */
1505 svmul(b, n, f_l); /* 3 */
1506 p = iprod(r_ij, r_kj); /* 5 */
1507 p *= nrkj_2; /* 1 */
1508 q = iprod(r_kl, r_kj); /* 5 */
1509 q *= nrkj_2; /* 1 */
1510 svmul(p, f_i, uvec); /* 3 */
1511 svmul(q, f_l, vvec); /* 3 */
1512 rvec_sub(uvec, vvec, svec); /* 3 */
1513 rvec_sub(f_i, svec, f_j); /* 3 */
1514 rvec_add(f_l, svec, f_k); /* 3 */
1515 rvec_inc(f[i], f_i); /* 3 */
1516 rvec_dec(f[j], f_j); /* 3 */
1517 rvec_dec(f[k], f_k); /* 3 */
1518 rvec_inc(f[l], f_l); /* 3 */
1522 /* As do_dih_fup_noshiftf above, but with pre-calculated pre-factors */
1524 do_dih_fup_noshiftf_precalc(int i, int j, int k, int l, real ddphi,
1525 real nrkj_m2, real nrkj_n2,
1527 rvec m, rvec n, rvec f[])
1529 rvec f_i, f_j, f_k, f_l;
1530 rvec uvec, vvec, svec, dx_jl;
1532 ivec jt, dt_ij, dt_kj, dt_lj;
1538 svmul(p, f_i, uvec);
1539 svmul(q, f_l, vvec);
1540 rvec_sub(uvec, vvec, svec);
1541 rvec_sub(f_i, svec, f_j);
1542 rvec_add(f_l, svec, f_k);
1543 rvec_inc(f[i], f_i);
1544 rvec_dec(f[j], f_j);
1545 rvec_dec(f[k], f_k);
1546 rvec_inc(f[l], f_l);
1550 real dopdihs(real cpA, real cpB, real phiA, real phiB, int mult,
1551 real phi, real lambda, real *V, real *F)
1553 real v, dvdlambda, mdphi, v1, sdphi, ddphi;
1554 real L1 = 1.0 - lambda;
1555 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1556 real dph0 = (phiB - phiA)*DEG2RAD;
1557 real cp = L1*cpA + lambda*cpB;
1559 mdphi = mult*phi - ph0;
1561 ddphi = -cp*mult*sdphi;
1562 v1 = 1.0 + cos(mdphi);
1565 dvdlambda = (cpB - cpA)*v1 + cp*dph0*sdphi;
1572 /* That was 40 flops */
1576 dopdihs_noener(real cpA, real cpB, real phiA, real phiB, int mult,
1577 real phi, real lambda, real *F)
1579 real mdphi, sdphi, ddphi;
1580 real L1 = 1.0 - lambda;
1581 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1582 real cp = L1*cpA + lambda*cpB;
1584 mdphi = mult*phi - ph0;
1586 ddphi = -cp*mult*sdphi;
1590 /* That was 20 flops */
1594 dopdihs_mdphi(real cpA, real cpB, real phiA, real phiB, int mult,
1595 real phi, real lambda, real *cp, real *mdphi)
1597 real L1 = 1.0 - lambda;
1598 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1600 *cp = L1*cpA + lambda*cpB;
1602 *mdphi = mult*phi - ph0;
1605 static real dopdihs_min(real cpA, real cpB, real phiA, real phiB, int mult,
1606 real phi, real lambda, real *V, real *F)
1607 /* similar to dopdihs, except for a minus sign *
1608 * and a different treatment of mult/phi0 */
1610 real v, dvdlambda, mdphi, v1, sdphi, ddphi;
1611 real L1 = 1.0 - lambda;
1612 real ph0 = (L1*phiA + lambda*phiB)*DEG2RAD;
1613 real dph0 = (phiB - phiA)*DEG2RAD;
1614 real cp = L1*cpA + lambda*cpB;
1616 mdphi = mult*(phi-ph0);
1618 ddphi = cp*mult*sdphi;
1619 v1 = 1.0-cos(mdphi);
1622 dvdlambda = (cpB-cpA)*v1 + cp*dph0*sdphi;
1629 /* That was 40 flops */
1632 real pdihs(int nbonds,
1633 const t_iatom forceatoms[], const t_iparams forceparams[],
1634 const rvec x[], rvec f[], rvec fshift[],
1635 const t_pbc *pbc, const t_graph *g,
1636 real lambda, real *dvdlambda,
1637 const t_mdatoms *md, t_fcdata *fcd,
1638 int *global_atom_index)
1640 int i, type, ai, aj, ak, al;
1642 rvec r_ij, r_kj, r_kl, m, n;
1643 real phi, sign, ddphi, vpd, vtot;
1647 for (i = 0; (i < nbonds); )
1649 type = forceatoms[i++];
1650 ai = forceatoms[i++];
1651 aj = forceatoms[i++];
1652 ak = forceatoms[i++];
1653 al = forceatoms[i++];
1655 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
1656 &sign, &t1, &t2, &t3); /* 84 */
1657 *dvdlambda += dopdihs(forceparams[type].pdihs.cpA,
1658 forceparams[type].pdihs.cpB,
1659 forceparams[type].pdihs.phiA,
1660 forceparams[type].pdihs.phiB,
1661 forceparams[type].pdihs.mult,
1662 phi, lambda, &vpd, &ddphi);
1665 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
1666 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
1669 fprintf(debug, "pdih: (%d,%d,%d,%d) phi=%g\n",
1670 ai, aj, ak, al, phi);
1677 void make_dp_periodic(real *dp) /* 1 flop? */
1679 /* dp cannot be outside (-pi,pi) */
1684 else if (*dp < -M_PI)
1691 /* As pdihs above, but without calculating energies and shift forces */
1693 pdihs_noener(int nbonds,
1694 const t_iatom forceatoms[], const t_iparams forceparams[],
1695 const rvec x[], rvec f[],
1696 const t_pbc *pbc, const t_graph *g,
1698 const t_mdatoms *md, t_fcdata *fcd,
1699 int *global_atom_index)
1701 int i, type, ai, aj, ak, al;
1703 rvec r_ij, r_kj, r_kl, m, n;
1704 real phi, sign, ddphi_tot, ddphi;
1706 for (i = 0; (i < nbonds); )
1708 ai = forceatoms[i+1];
1709 aj = forceatoms[i+2];
1710 ak = forceatoms[i+3];
1711 al = forceatoms[i+4];
1713 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
1714 &sign, &t1, &t2, &t3);
1718 /* Loop over dihedrals working on the same atoms,
1719 * so we avoid recalculating angles and force distributions.
1723 type = forceatoms[i];
1724 dopdihs_noener(forceparams[type].pdihs.cpA,
1725 forceparams[type].pdihs.cpB,
1726 forceparams[type].pdihs.phiA,
1727 forceparams[type].pdihs.phiB,
1728 forceparams[type].pdihs.mult,
1729 phi, lambda, &ddphi);
1734 while (i < nbonds &&
1735 forceatoms[i+1] == ai &&
1736 forceatoms[i+2] == aj &&
1737 forceatoms[i+3] == ak &&
1738 forceatoms[i+4] == al);
1740 do_dih_fup_noshiftf(ai, aj, ak, al, ddphi_tot, r_ij, r_kj, r_kl, m, n, f);
1745 #ifdef SSE_PROPER_DIHEDRALS
1747 /* As pdihs_noner above, but using SSE to calculate 4 dihedrals at once */
1749 pdihs_noener_sse(int nbonds,
1750 const t_iatom forceatoms[], const t_iparams forceparams[],
1751 const rvec x[], rvec f[],
1752 const t_pbc *pbc, const t_graph *g,
1754 const t_mdatoms *md, t_fcdata *fcd,
1755 int *global_atom_index)
1758 int type, ai[4], aj[4], ak[4], al[4];
1759 int t1[4], t2[4], t3[4];
1761 real cp[4], mdphi[4];
1763 rvec_sse_t rs_array[13], *rs;
1764 real buf_array[24], *buf;
1765 __m128 mdphi_SSE, sin_SSE, cos_SSE;
1767 /* Ensure 16-byte alignment */
1768 rs = (rvec_sse_t *)(((size_t)(rs_array +1)) & (~((size_t)15)));
1769 buf = (float *)(((size_t)(buf_array+3)) & (~((size_t)15)));
1771 for (i = 0; (i < nbonds); i += 20)
1773 /* Collect atoms quadruplets for 4 dihedrals */
1775 for (s = 0; s < 4; s++)
1777 ai[s] = forceatoms[i4+1];
1778 aj[s] = forceatoms[i4+2];
1779 ak[s] = forceatoms[i4+3];
1780 al[s] = forceatoms[i4+4];
1781 /* At the end fill the arrays with identical entries */
1782 if (i4 + 5 < nbonds)
1788 /* Caclulate 4 dihedral angles at once */
1789 dih_angle_sse(x, ai, aj, ak, al, pbc, t1, t2, t3, rs, buf);
1792 for (s = 0; s < 4; s++)
1796 /* Calculate the coefficient and angle deviation */
1797 type = forceatoms[i4];
1798 dopdihs_mdphi(forceparams[type].pdihs.cpA,
1799 forceparams[type].pdihs.cpB,
1800 forceparams[type].pdihs.phiA,
1801 forceparams[type].pdihs.phiB,
1802 forceparams[type].pdihs.mult,
1803 buf[16+s], lambda, &cp[s], &buf[16+s]);
1804 mult[s] = forceparams[type].pdihs.mult;
1813 /* Calculate 4 sines at once */
1814 mdphi_SSE = _mm_load_ps(buf+16);
1815 gmx_mm_sincos_ps(mdphi_SSE, &sin_SSE, &cos_SSE);
1816 _mm_store_ps(buf+16, sin_SSE);
1822 ddphi = -cp[s]*mult[s]*buf[16+s];
1824 do_dih_fup_noshiftf_precalc(ai[s], aj[s], ak[s], al[s], ddphi,
1825 buf[ 0+s], buf[ 4+s],
1826 buf[ 8+s], buf[12+s],
1827 rs[0+s].v, rs[4+s].v,
1832 while (s < 4 && i4 < nbonds);
1836 #endif /* SSE_PROPER_DIHEDRALS */
1839 real idihs(int nbonds,
1840 const t_iatom forceatoms[], const t_iparams forceparams[],
1841 const rvec x[], rvec f[], rvec fshift[],
1842 const t_pbc *pbc, const t_graph *g,
1843 real lambda, real *dvdlambda,
1844 const t_mdatoms *md, t_fcdata *fcd,
1845 int *global_atom_index)
1847 int i, type, ai, aj, ak, al;
1849 real phi, phi0, dphi0, ddphi, sign, vtot;
1850 rvec r_ij, r_kj, r_kl, m, n;
1851 real L1, kk, dp, dp2, kA, kB, pA, pB, dvdl_term;
1856 for (i = 0; (i < nbonds); )
1858 type = forceatoms[i++];
1859 ai = forceatoms[i++];
1860 aj = forceatoms[i++];
1861 ak = forceatoms[i++];
1862 al = forceatoms[i++];
1864 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
1865 &sign, &t1, &t2, &t3); /* 84 */
1867 /* phi can jump if phi0 is close to Pi/-Pi, which will cause huge
1868 * force changes if we just apply a normal harmonic.
1869 * Instead, we first calculate phi-phi0 and take it modulo (-Pi,Pi).
1870 * This means we will never have the periodicity problem, unless
1871 * the dihedral is Pi away from phiO, which is very unlikely due to
1874 kA = forceparams[type].harmonic.krA;
1875 kB = forceparams[type].harmonic.krB;
1876 pA = forceparams[type].harmonic.rA;
1877 pB = forceparams[type].harmonic.rB;
1879 kk = L1*kA + lambda*kB;
1880 phi0 = (L1*pA + lambda*pB)*DEG2RAD;
1881 dphi0 = (pB - pA)*DEG2RAD;
1885 make_dp_periodic(&dp);
1892 dvdl_term += 0.5*(kB - kA)*dp2 - kk*dphi0*dp;
1894 do_dih_fup(ai, aj, ak, al, (real)(-ddphi), r_ij, r_kj, r_kl, m, n,
1895 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
1900 fprintf(debug, "idih: (%d,%d,%d,%d) phi=%g\n",
1901 ai, aj, ak, al, phi);
1906 *dvdlambda += dvdl_term;
1911 /*! \brief returns dx, rdist, and dpdl for functions posres() and fbposres()
1913 static void posres_dx(const rvec x, const rvec pos0A, const rvec pos0B,
1914 const rvec comA_sc, const rvec comB_sc,
1916 t_pbc *pbc, int refcoord_scaling, int npbcdim,
1917 rvec dx, rvec rdist, rvec dpdl)
1920 real posA, posB, L1, ref = 0.;
1925 for (m = 0; m < DIM; m++)
1931 switch (refcoord_scaling)
1935 rdist[m] = L1*posA + lambda*posB;
1936 dpdl[m] = posB - posA;
1939 /* Box relative coordinates are stored for dimensions with pbc */
1940 posA *= pbc->box[m][m];
1941 posB *= pbc->box[m][m];
1942 for (d = m+1; d < npbcdim; d++)
1944 posA += pos0A[d]*pbc->box[d][m];
1945 posB += pos0B[d]*pbc->box[d][m];
1947 ref = L1*posA + lambda*posB;
1949 dpdl[m] = posB - posA;
1952 ref = L1*comA_sc[m] + lambda*comB_sc[m];
1953 rdist[m] = L1*posA + lambda*posB;
1954 dpdl[m] = comB_sc[m] - comA_sc[m] + posB - posA;
1957 gmx_fatal(FARGS, "No such scaling method implemented");
1962 ref = L1*posA + lambda*posB;
1964 dpdl[m] = posB - posA;
1967 /* We do pbc_dx with ref+rdist,
1968 * since with only ref we can be up to half a box vector wrong.
1970 pos[m] = ref + rdist[m];
1975 pbc_dx(pbc, x, pos, dx);
1979 rvec_sub(x, pos, dx);
1983 /*! \brief Adds forces of flat-bottomed positions restraints to f[]
1984 * and fixes vir_diag. Returns the flat-bottomed potential. */
1985 real fbposres(int nbonds,
1986 const t_iatom forceatoms[], const t_iparams forceparams[],
1987 const rvec x[], rvec f[], rvec vir_diag,
1989 int refcoord_scaling, int ePBC, rvec com)
1990 /* compute flat-bottomed positions restraints */
1992 int i, ai, m, d, type, npbcdim = 0, fbdim;
1993 const t_iparams *pr;
1995 real ref = 0, dr, dr2, rpot, rfb, rfb2, fact, invdr;
1996 rvec com_sc, rdist, pos, dx, dpdl, fm;
1999 npbcdim = ePBC2npbcdim(ePBC);
2001 if (refcoord_scaling == erscCOM)
2004 for (m = 0; m < npbcdim; m++)
2006 for (d = m; d < npbcdim; d++)
2008 com_sc[m] += com[d]*pbc->box[d][m];
2014 for (i = 0; (i < nbonds); )
2016 type = forceatoms[i++];
2017 ai = forceatoms[i++];
2018 pr = &forceparams[type];
2020 /* same calculation as for normal posres, but with identical A and B states, and lambda==0 */
2021 posres_dx(x[ai], forceparams[type].fbposres.pos0, forceparams[type].fbposres.pos0,
2022 com_sc, com_sc, 0.0,
2023 pbc, refcoord_scaling, npbcdim,
2029 kk = pr->fbposres.k;
2030 rfb = pr->fbposres.r;
2033 /* with rfb<0, push particle out of the sphere/cylinder/layer */
2041 switch (pr->fbposres.geom)
2043 case efbposresSPHERE:
2044 /* spherical flat-bottom posres */
2047 ( (dr2 > rfb2 && bInvert == FALSE ) || (dr2 < rfb2 && bInvert == TRUE ) )
2051 v = 0.5*kk*sqr(dr - rfb);
2052 fact = -kk*(dr-rfb)/dr; /* Force pointing to the center pos0 */
2053 svmul(fact, dx, fm);
2056 case efbposresCYLINDER:
2057 /* cylidrical flat-bottom posres in x-y plane. fm[ZZ] = 0. */
2058 dr2 = sqr(dx[XX])+sqr(dx[YY]);
2060 ( (dr2 > rfb2 && bInvert == FALSE ) || (dr2 < rfb2 && bInvert == TRUE ) )
2065 v = 0.5*kk*sqr(dr - rfb);
2066 fm[XX] = -kk*(dr-rfb)*dx[XX]*invdr; /* Force pointing to the center */
2067 fm[YY] = -kk*(dr-rfb)*dx[YY]*invdr;
2070 case efbposresX: /* fbdim=XX */
2071 case efbposresY: /* fbdim=YY */
2072 case efbposresZ: /* fbdim=ZZ */
2073 /* 1D flat-bottom potential */
2074 fbdim = pr->fbposres.geom - efbposresX;
2076 if ( ( dr > rfb && bInvert == FALSE ) || ( 0 < dr && dr < rfb && bInvert == TRUE ) )
2078 v = 0.5*kk*sqr(dr - rfb);
2079 fm[fbdim] = -kk*(dr - rfb);
2081 else if ( (dr < (-rfb) && bInvert == FALSE ) || ( (-rfb) < dr && dr < 0 && bInvert == TRUE ))
2083 v = 0.5*kk*sqr(dr + rfb);
2084 fm[fbdim] = -kk*(dr + rfb);
2091 for (m = 0; (m < DIM); m++)
2094 /* Here we correct for the pbc_dx which included rdist */
2095 vir_diag[m] -= 0.5*(dx[m] + rdist[m])*fm[m];
2103 real posres(int nbonds,
2104 const t_iatom forceatoms[], const t_iparams forceparams[],
2105 const rvec x[], rvec f[], rvec vir_diag,
2107 real lambda, real *dvdlambda,
2108 int refcoord_scaling, int ePBC, rvec comA, rvec comB)
2110 int i, ai, m, d, type, ki, npbcdim = 0;
2111 const t_iparams *pr;
2114 real posA, posB, ref = 0;
2115 rvec comA_sc, comB_sc, rdist, dpdl, pos, dx;
2116 gmx_bool bForceValid = TRUE;
2118 if ((f == NULL) || (vir_diag == NULL)) /* should both be null together! */
2120 bForceValid = FALSE;
2123 npbcdim = ePBC2npbcdim(ePBC);
2125 if (refcoord_scaling == erscCOM)
2127 clear_rvec(comA_sc);
2128 clear_rvec(comB_sc);
2129 for (m = 0; m < npbcdim; m++)
2131 for (d = m; d < npbcdim; d++)
2133 comA_sc[m] += comA[d]*pbc->box[d][m];
2134 comB_sc[m] += comB[d]*pbc->box[d][m];
2142 for (i = 0; (i < nbonds); )
2144 type = forceatoms[i++];
2145 ai = forceatoms[i++];
2146 pr = &forceparams[type];
2148 /* return dx, rdist, and dpdl */
2149 posres_dx(x[ai], forceparams[type].posres.pos0A, forceparams[type].posres.pos0B,
2150 comA_sc, comB_sc, lambda,
2151 pbc, refcoord_scaling, npbcdim,
2154 for (m = 0; (m < DIM); m++)
2156 kk = L1*pr->posres.fcA[m] + lambda*pr->posres.fcB[m];
2158 vtot += 0.5*kk*dx[m]*dx[m];
2160 0.5*(pr->posres.fcB[m] - pr->posres.fcA[m])*dx[m]*dx[m]
2163 /* Here we correct for the pbc_dx which included rdist */
2167 vir_diag[m] -= 0.5*(dx[m] + rdist[m])*fm;
2175 static real low_angres(int nbonds,
2176 const t_iatom forceatoms[], const t_iparams forceparams[],
2177 const rvec x[], rvec f[], rvec fshift[],
2178 const t_pbc *pbc, const t_graph *g,
2179 real lambda, real *dvdlambda,
2182 int i, m, type, ai, aj, ak, al;
2184 real phi, cos_phi, cos_phi2, vid, vtot, dVdphi;
2185 rvec r_ij, r_kl, f_i, f_k = {0, 0, 0};
2186 real st, sth, nrij2, nrkl2, c, cij, ckl;
2189 t2 = 0; /* avoid warning with gcc-3.3. It is never used uninitialized */
2192 ak = al = 0; /* to avoid warnings */
2193 for (i = 0; i < nbonds; )
2195 type = forceatoms[i++];
2196 ai = forceatoms[i++];
2197 aj = forceatoms[i++];
2198 t1 = pbc_rvec_sub(pbc, x[aj], x[ai], r_ij); /* 3 */
2201 ak = forceatoms[i++];
2202 al = forceatoms[i++];
2203 t2 = pbc_rvec_sub(pbc, x[al], x[ak], r_kl); /* 3 */
2212 cos_phi = cos_angle(r_ij, r_kl); /* 25 */
2213 phi = acos(cos_phi); /* 10 */
2215 *dvdlambda += dopdihs_min(forceparams[type].pdihs.cpA,
2216 forceparams[type].pdihs.cpB,
2217 forceparams[type].pdihs.phiA,
2218 forceparams[type].pdihs.phiB,
2219 forceparams[type].pdihs.mult,
2220 phi, lambda, &vid, &dVdphi); /* 40 */
2224 cos_phi2 = sqr(cos_phi); /* 1 */
2227 st = -dVdphi*gmx_invsqrt(1 - cos_phi2); /* 12 */
2228 sth = st*cos_phi; /* 1 */
2229 nrij2 = iprod(r_ij, r_ij); /* 5 */
2230 nrkl2 = iprod(r_kl, r_kl); /* 5 */
2232 c = st*gmx_invsqrt(nrij2*nrkl2); /* 11 */
2233 cij = sth/nrij2; /* 10 */
2234 ckl = sth/nrkl2; /* 10 */
2236 for (m = 0; m < DIM; m++) /* 18+18 */
2238 f_i[m] = (c*r_kl[m]-cij*r_ij[m]);
2243 f_k[m] = (c*r_ij[m]-ckl*r_kl[m]);
2251 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
2254 rvec_inc(fshift[t1], f_i);
2255 rvec_dec(fshift[CENTRAL], f_i);
2260 ivec_sub(SHIFT_IVEC(g, ak), SHIFT_IVEC(g, al), dt);
2263 rvec_inc(fshift[t2], f_k);
2264 rvec_dec(fshift[CENTRAL], f_k);
2269 return vtot; /* 184 / 157 (bZAxis) total */
2272 real angres(int nbonds,
2273 const t_iatom forceatoms[], const t_iparams forceparams[],
2274 const rvec x[], rvec f[], rvec fshift[],
2275 const t_pbc *pbc, const t_graph *g,
2276 real lambda, real *dvdlambda,
2277 const t_mdatoms *md, t_fcdata *fcd,
2278 int *global_atom_index)
2280 return low_angres(nbonds, forceatoms, forceparams, x, f, fshift, pbc, g,
2281 lambda, dvdlambda, FALSE);
2284 real angresz(int nbonds,
2285 const t_iatom forceatoms[], const t_iparams forceparams[],
2286 const rvec x[], rvec f[], rvec fshift[],
2287 const t_pbc *pbc, const t_graph *g,
2288 real lambda, real *dvdlambda,
2289 const t_mdatoms *md, t_fcdata *fcd,
2290 int *global_atom_index)
2292 return low_angres(nbonds, forceatoms, forceparams, x, f, fshift, pbc, g,
2293 lambda, dvdlambda, TRUE);
2296 real dihres(int nbonds,
2297 const t_iatom forceatoms[], const t_iparams forceparams[],
2298 const rvec x[], rvec f[], rvec fshift[],
2299 const t_pbc *pbc, const t_graph *g,
2300 real lambda, real *dvdlambda,
2301 const t_mdatoms *md, t_fcdata *fcd,
2302 int *global_atom_index)
2305 int ai, aj, ak, al, i, k, type, t1, t2, t3;
2306 real phi0A, phi0B, dphiA, dphiB, kfacA, kfacB, phi0, dphi, kfac;
2307 real phi, ddphi, ddp, ddp2, dp, sign, d2r, fc, L1;
2308 rvec r_ij, r_kj, r_kl, m, n;
2315 for (i = 0; (i < nbonds); )
2317 type = forceatoms[i++];
2318 ai = forceatoms[i++];
2319 aj = forceatoms[i++];
2320 ak = forceatoms[i++];
2321 al = forceatoms[i++];
2323 phi0A = forceparams[type].dihres.phiA*d2r;
2324 dphiA = forceparams[type].dihres.dphiA*d2r;
2325 kfacA = forceparams[type].dihres.kfacA;
2327 phi0B = forceparams[type].dihres.phiB*d2r;
2328 dphiB = forceparams[type].dihres.dphiB*d2r;
2329 kfacB = forceparams[type].dihres.kfacB;
2331 phi0 = L1*phi0A + lambda*phi0B;
2332 dphi = L1*dphiA + lambda*dphiB;
2333 kfac = L1*kfacA + lambda*kfacB;
2335 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2336 &sign, &t1, &t2, &t3);
2341 fprintf(debug, "dihres[%d]: %d %d %d %d : phi=%f, dphi=%f, kfac=%f\n",
2342 k++, ai, aj, ak, al, phi0, dphi, kfac);
2344 /* phi can jump if phi0 is close to Pi/-Pi, which will cause huge
2345 * force changes if we just apply a normal harmonic.
2346 * Instead, we first calculate phi-phi0 and take it modulo (-Pi,Pi).
2347 * This means we will never have the periodicity problem, unless
2348 * the dihedral is Pi away from phiO, which is very unlikely due to
2352 make_dp_periodic(&dp);
2358 else if (dp < -dphi)
2370 vtot += 0.5*kfac*ddp2;
2373 *dvdlambda += 0.5*(kfacB - kfacA)*ddp2;
2374 /* lambda dependence from changing restraint distances */
2377 *dvdlambda -= kfac*ddp*((dphiB - dphiA)+(phi0B - phi0A));
2381 *dvdlambda += kfac*ddp*((dphiB - dphiA)-(phi0B - phi0A));
2383 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
2384 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2391 real unimplemented(int nbonds,
2392 const t_iatom forceatoms[], const t_iparams forceparams[],
2393 const rvec x[], rvec f[], rvec fshift[],
2394 const t_pbc *pbc, const t_graph *g,
2395 real lambda, real *dvdlambda,
2396 const t_mdatoms *md, t_fcdata *fcd,
2397 int *global_atom_index)
2399 gmx_impl("*** you are using a not implemented function");
2401 return 0.0; /* To make the compiler happy */
2404 real rbdihs(int nbonds,
2405 const t_iatom forceatoms[], const t_iparams forceparams[],
2406 const rvec x[], rvec f[], rvec fshift[],
2407 const t_pbc *pbc, const t_graph *g,
2408 real lambda, real *dvdlambda,
2409 const t_mdatoms *md, t_fcdata *fcd,
2410 int *global_atom_index)
2412 const real c0 = 0.0, c1 = 1.0, c2 = 2.0, c3 = 3.0, c4 = 4.0, c5 = 5.0;
2413 int type, ai, aj, ak, al, i, j;
2415 rvec r_ij, r_kj, r_kl, m, n;
2416 real parmA[NR_RBDIHS];
2417 real parmB[NR_RBDIHS];
2418 real parm[NR_RBDIHS];
2419 real cos_phi, phi, rbp, rbpBA;
2420 real v, sign, ddphi, sin_phi;
2422 real L1 = 1.0-lambda;
2426 for (i = 0; (i < nbonds); )
2428 type = forceatoms[i++];
2429 ai = forceatoms[i++];
2430 aj = forceatoms[i++];
2431 ak = forceatoms[i++];
2432 al = forceatoms[i++];
2434 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
2435 &sign, &t1, &t2, &t3); /* 84 */
2437 /* Change to polymer convention */
2444 phi -= M_PI; /* 1 */
2448 /* Beware of accuracy loss, cannot use 1-sqrt(cos^2) ! */
2451 for (j = 0; (j < NR_RBDIHS); j++)
2453 parmA[j] = forceparams[type].rbdihs.rbcA[j];
2454 parmB[j] = forceparams[type].rbdihs.rbcB[j];
2455 parm[j] = L1*parmA[j]+lambda*parmB[j];
2457 /* Calculate cosine powers */
2458 /* Calculate the energy */
2459 /* Calculate the derivative */
2462 dvdl_term += (parmB[0]-parmA[0]);
2467 rbpBA = parmB[1]-parmA[1];
2468 ddphi += rbp*cosfac;
2471 dvdl_term += cosfac*rbpBA;
2473 rbpBA = parmB[2]-parmA[2];
2474 ddphi += c2*rbp*cosfac;
2477 dvdl_term += cosfac*rbpBA;
2479 rbpBA = parmB[3]-parmA[3];
2480 ddphi += c3*rbp*cosfac;
2483 dvdl_term += cosfac*rbpBA;
2485 rbpBA = parmB[4]-parmA[4];
2486 ddphi += c4*rbp*cosfac;
2489 dvdl_term += cosfac*rbpBA;
2491 rbpBA = parmB[5]-parmA[5];
2492 ddphi += c5*rbp*cosfac;
2495 dvdl_term += cosfac*rbpBA;
2497 ddphi = -ddphi*sin_phi; /* 11 */
2499 do_dih_fup(ai, aj, ak, al, ddphi, r_ij, r_kj, r_kl, m, n,
2500 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
2503 *dvdlambda += dvdl_term;
2508 int cmap_setup_grid_index(int ip, int grid_spacing, int *ipm1, int *ipp1, int *ipp2)
2514 ip = ip + grid_spacing - 1;
2516 else if (ip > grid_spacing)
2518 ip = ip - grid_spacing - 1;
2527 im1 = grid_spacing - 1;
2529 else if (ip == grid_spacing-2)
2533 else if (ip == grid_spacing-1)
2547 real cmap_dihs(int nbonds,
2548 const t_iatom forceatoms[], const t_iparams forceparams[],
2549 const gmx_cmap_t *cmap_grid,
2550 const rvec x[], rvec f[], rvec fshift[],
2551 const t_pbc *pbc, const t_graph *g,
2552 real lambda, real *dvdlambda,
2553 const t_mdatoms *md, t_fcdata *fcd,
2554 int *global_atom_index)
2556 int i, j, k, n, idx;
2557 int ai, aj, ak, al, am;
2558 int a1i, a1j, a1k, a1l, a2i, a2j, a2k, a2l;
2560 int t11, t21, t31, t12, t22, t32;
2561 int iphi1, ip1m1, ip1p1, ip1p2;
2562 int iphi2, ip2m1, ip2p1, ip2p2;
2564 int pos1, pos2, pos3, pos4, tmp;
2566 real ty[4], ty1[4], ty2[4], ty12[4], tc[16], tx[16];
2567 real phi1, psi1, cos_phi1, sin_phi1, sign1, xphi1;
2568 real phi2, psi2, cos_phi2, sin_phi2, sign2, xphi2;
2569 real dx, xx, tt, tu, e, df1, df2, ddf1, ddf2, ddf12, vtot;
2570 real ra21, rb21, rg21, rg1, rgr1, ra2r1, rb2r1, rabr1;
2571 real ra22, rb22, rg22, rg2, rgr2, ra2r2, rb2r2, rabr2;
2572 real fg1, hg1, fga1, hgb1, gaa1, gbb1;
2573 real fg2, hg2, fga2, hgb2, gaa2, gbb2;
2576 rvec r1_ij, r1_kj, r1_kl, m1, n1;
2577 rvec r2_ij, r2_kj, r2_kl, m2, n2;
2578 rvec f1_i, f1_j, f1_k, f1_l;
2579 rvec f2_i, f2_j, f2_k, f2_l;
2580 rvec a1, b1, a2, b2;
2581 rvec f1, g1, h1, f2, g2, h2;
2582 rvec dtf1, dtg1, dth1, dtf2, dtg2, dth2;
2583 ivec jt1, dt1_ij, dt1_kj, dt1_lj;
2584 ivec jt2, dt2_ij, dt2_kj, dt2_lj;
2588 int loop_index[4][4] = {
2595 /* Total CMAP energy */
2598 for (n = 0; n < nbonds; )
2600 /* Five atoms are involved in the two torsions */
2601 type = forceatoms[n++];
2602 ai = forceatoms[n++];
2603 aj = forceatoms[n++];
2604 ak = forceatoms[n++];
2605 al = forceatoms[n++];
2606 am = forceatoms[n++];
2608 /* Which CMAP type is this */
2609 cmapA = forceparams[type].cmap.cmapA;
2610 cmapd = cmap_grid->cmapdata[cmapA].cmap;
2618 phi1 = dih_angle(x[a1i], x[a1j], x[a1k], x[a1l], pbc, r1_ij, r1_kj, r1_kl, m1, n1,
2619 &sign1, &t11, &t21, &t31); /* 84 */
2621 cos_phi1 = cos(phi1);
2623 a1[0] = r1_ij[1]*r1_kj[2]-r1_ij[2]*r1_kj[1];
2624 a1[1] = r1_ij[2]*r1_kj[0]-r1_ij[0]*r1_kj[2];
2625 a1[2] = r1_ij[0]*r1_kj[1]-r1_ij[1]*r1_kj[0]; /* 9 */
2627 b1[0] = r1_kl[1]*r1_kj[2]-r1_kl[2]*r1_kj[1];
2628 b1[1] = r1_kl[2]*r1_kj[0]-r1_kl[0]*r1_kj[2];
2629 b1[2] = r1_kl[0]*r1_kj[1]-r1_kl[1]*r1_kj[0]; /* 9 */
2631 tmp = pbc_rvec_sub(pbc, x[a1l], x[a1k], h1);
2633 ra21 = iprod(a1, a1); /* 5 */
2634 rb21 = iprod(b1, b1); /* 5 */
2635 rg21 = iprod(r1_kj, r1_kj); /* 5 */
2641 rabr1 = sqrt(ra2r1*rb2r1);
2643 sin_phi1 = rg1 * rabr1 * iprod(a1, h1) * (-1);
2645 if (cos_phi1 < -0.5 || cos_phi1 > 0.5)
2647 phi1 = asin(sin_phi1);
2657 phi1 = -M_PI - phi1;
2663 phi1 = acos(cos_phi1);
2671 xphi1 = phi1 + M_PI; /* 1 */
2673 /* Second torsion */
2679 phi2 = dih_angle(x[a2i], x[a2j], x[a2k], x[a2l], pbc, r2_ij, r2_kj, r2_kl, m2, n2,
2680 &sign2, &t12, &t22, &t32); /* 84 */
2682 cos_phi2 = cos(phi2);
2684 a2[0] = r2_ij[1]*r2_kj[2]-r2_ij[2]*r2_kj[1];
2685 a2[1] = r2_ij[2]*r2_kj[0]-r2_ij[0]*r2_kj[2];
2686 a2[2] = r2_ij[0]*r2_kj[1]-r2_ij[1]*r2_kj[0]; /* 9 */
2688 b2[0] = r2_kl[1]*r2_kj[2]-r2_kl[2]*r2_kj[1];
2689 b2[1] = r2_kl[2]*r2_kj[0]-r2_kl[0]*r2_kj[2];
2690 b2[2] = r2_kl[0]*r2_kj[1]-r2_kl[1]*r2_kj[0]; /* 9 */
2692 tmp = pbc_rvec_sub(pbc, x[a2l], x[a2k], h2);
2694 ra22 = iprod(a2, a2); /* 5 */
2695 rb22 = iprod(b2, b2); /* 5 */
2696 rg22 = iprod(r2_kj, r2_kj); /* 5 */
2702 rabr2 = sqrt(ra2r2*rb2r2);
2704 sin_phi2 = rg2 * rabr2 * iprod(a2, h2) * (-1);
2706 if (cos_phi2 < -0.5 || cos_phi2 > 0.5)
2708 phi2 = asin(sin_phi2);
2718 phi2 = -M_PI - phi2;
2724 phi2 = acos(cos_phi2);
2732 xphi2 = phi2 + M_PI; /* 1 */
2734 /* Range mangling */
2737 xphi1 = xphi1 + 2*M_PI;
2739 else if (xphi1 >= 2*M_PI)
2741 xphi1 = xphi1 - 2*M_PI;
2746 xphi2 = xphi2 + 2*M_PI;
2748 else if (xphi2 >= 2*M_PI)
2750 xphi2 = xphi2 - 2*M_PI;
2753 /* Number of grid points */
2754 dx = 2*M_PI / cmap_grid->grid_spacing;
2756 /* Where on the grid are we */
2757 iphi1 = (int)(xphi1/dx);
2758 iphi2 = (int)(xphi2/dx);
2760 iphi1 = cmap_setup_grid_index(iphi1, cmap_grid->grid_spacing, &ip1m1, &ip1p1, &ip1p2);
2761 iphi2 = cmap_setup_grid_index(iphi2, cmap_grid->grid_spacing, &ip2m1, &ip2p1, &ip2p2);
2763 pos1 = iphi1*cmap_grid->grid_spacing+iphi2;
2764 pos2 = ip1p1*cmap_grid->grid_spacing+iphi2;
2765 pos3 = ip1p1*cmap_grid->grid_spacing+ip2p1;
2766 pos4 = iphi1*cmap_grid->grid_spacing+ip2p1;
2768 ty[0] = cmapd[pos1*4];
2769 ty[1] = cmapd[pos2*4];
2770 ty[2] = cmapd[pos3*4];
2771 ty[3] = cmapd[pos4*4];
2773 ty1[0] = cmapd[pos1*4+1];
2774 ty1[1] = cmapd[pos2*4+1];
2775 ty1[2] = cmapd[pos3*4+1];
2776 ty1[3] = cmapd[pos4*4+1];
2778 ty2[0] = cmapd[pos1*4+2];
2779 ty2[1] = cmapd[pos2*4+2];
2780 ty2[2] = cmapd[pos3*4+2];
2781 ty2[3] = cmapd[pos4*4+2];
2783 ty12[0] = cmapd[pos1*4+3];
2784 ty12[1] = cmapd[pos2*4+3];
2785 ty12[2] = cmapd[pos3*4+3];
2786 ty12[3] = cmapd[pos4*4+3];
2788 /* Switch to degrees */
2789 dx = 360.0 / cmap_grid->grid_spacing;
2790 xphi1 = xphi1 * RAD2DEG;
2791 xphi2 = xphi2 * RAD2DEG;
2793 for (i = 0; i < 4; i++) /* 16 */
2796 tx[i+4] = ty1[i]*dx;
2797 tx[i+8] = ty2[i]*dx;
2798 tx[i+12] = ty12[i]*dx*dx;
2802 for (i = 0; i < 4; i++) /* 1056 */
2804 for (j = 0; j < 4; j++)
2807 for (k = 0; k < 16; k++)
2809 xx = xx + cmap_coeff_matrix[k*16+idx]*tx[k];
2817 tt = (xphi1-iphi1*dx)/dx;
2818 tu = (xphi2-iphi2*dx)/dx;
2827 for (i = 3; i >= 0; i--)
2829 l1 = loop_index[i][3];
2830 l2 = loop_index[i][2];
2831 l3 = loop_index[i][1];
2833 e = tt * e + ((tc[i*4+3]*tu+tc[i*4+2])*tu + tc[i*4+1])*tu+tc[i*4];
2834 df1 = tu * df1 + (3.0*tc[l1]*tt+2.0*tc[l2])*tt+tc[l3];
2835 df2 = tt * df2 + (3.0*tc[i*4+3]*tu+2.0*tc[i*4+2])*tu+tc[i*4+1];
2836 ddf1 = tu * ddf1 + 2.0*3.0*tc[l1]*tt+2.0*tc[l2];
2837 ddf2 = tt * ddf2 + 2.0*3.0*tc[4*i+3]*tu+2.0*tc[4*i+2];
2840 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) +
2841 3.0*tu*tu*(tc[7]+2.0*tc[11]*tt+3.0*tc[15]*tt*tt);
2846 ddf1 = ddf1 * fac * fac;
2847 ddf2 = ddf2 * fac * fac;
2848 ddf12 = ddf12 * fac * fac;
2853 /* Do forces - first torsion */
2854 fg1 = iprod(r1_ij, r1_kj);
2855 hg1 = iprod(r1_kl, r1_kj);
2856 fga1 = fg1*ra2r1*rgr1;
2857 hgb1 = hg1*rb2r1*rgr1;
2861 for (i = 0; i < DIM; i++)
2863 dtf1[i] = gaa1 * a1[i];
2864 dtg1[i] = fga1 * a1[i] - hgb1 * b1[i];
2865 dth1[i] = gbb1 * b1[i];
2867 f1[i] = df1 * dtf1[i];
2868 g1[i] = df1 * dtg1[i];
2869 h1[i] = df1 * dth1[i];
2872 f1_j[i] = -f1[i] - g1[i];
2873 f1_k[i] = h1[i] + g1[i];
2876 f[a1i][i] = f[a1i][i] + f1_i[i];
2877 f[a1j][i] = f[a1j][i] + f1_j[i]; /* - f1[i] - g1[i] */
2878 f[a1k][i] = f[a1k][i] + f1_k[i]; /* h1[i] + g1[i] */
2879 f[a1l][i] = f[a1l][i] + f1_l[i]; /* h1[i] */
2882 /* Do forces - second torsion */
2883 fg2 = iprod(r2_ij, r2_kj);
2884 hg2 = iprod(r2_kl, r2_kj);
2885 fga2 = fg2*ra2r2*rgr2;
2886 hgb2 = hg2*rb2r2*rgr2;
2890 for (i = 0; i < DIM; i++)
2892 dtf2[i] = gaa2 * a2[i];
2893 dtg2[i] = fga2 * a2[i] - hgb2 * b2[i];
2894 dth2[i] = gbb2 * b2[i];
2896 f2[i] = df2 * dtf2[i];
2897 g2[i] = df2 * dtg2[i];
2898 h2[i] = df2 * dth2[i];
2901 f2_j[i] = -f2[i] - g2[i];
2902 f2_k[i] = h2[i] + g2[i];
2905 f[a2i][i] = f[a2i][i] + f2_i[i]; /* f2[i] */
2906 f[a2j][i] = f[a2j][i] + f2_j[i]; /* - f2[i] - g2[i] */
2907 f[a2k][i] = f[a2k][i] + f2_k[i]; /* h2[i] + g2[i] */
2908 f[a2l][i] = f[a2l][i] + f2_l[i]; /* - h2[i] */
2914 copy_ivec(SHIFT_IVEC(g, a1j), jt1);
2915 ivec_sub(SHIFT_IVEC(g, a1i), jt1, dt1_ij);
2916 ivec_sub(SHIFT_IVEC(g, a1k), jt1, dt1_kj);
2917 ivec_sub(SHIFT_IVEC(g, a1l), jt1, dt1_lj);
2918 t11 = IVEC2IS(dt1_ij);
2919 t21 = IVEC2IS(dt1_kj);
2920 t31 = IVEC2IS(dt1_lj);
2922 copy_ivec(SHIFT_IVEC(g, a2j), jt2);
2923 ivec_sub(SHIFT_IVEC(g, a2i), jt2, dt2_ij);
2924 ivec_sub(SHIFT_IVEC(g, a2k), jt2, dt2_kj);
2925 ivec_sub(SHIFT_IVEC(g, a2l), jt2, dt2_lj);
2926 t12 = IVEC2IS(dt2_ij);
2927 t22 = IVEC2IS(dt2_kj);
2928 t32 = IVEC2IS(dt2_lj);
2932 t31 = pbc_rvec_sub(pbc, x[a1l], x[a1j], h1);
2933 t32 = pbc_rvec_sub(pbc, x[a2l], x[a2j], h2);
2941 rvec_inc(fshift[t11], f1_i);
2942 rvec_inc(fshift[CENTRAL], f1_j);
2943 rvec_inc(fshift[t21], f1_k);
2944 rvec_inc(fshift[t31], f1_l);
2946 rvec_inc(fshift[t21], f2_i);
2947 rvec_inc(fshift[CENTRAL], f2_j);
2948 rvec_inc(fshift[t22], f2_k);
2949 rvec_inc(fshift[t32], f2_l);
2956 /***********************************************************
2958 * G R O M O S 9 6 F U N C T I O N S
2960 ***********************************************************/
2961 real g96harmonic(real kA, real kB, real xA, real xB, real x, real lambda,
2964 const real half = 0.5;
2965 real L1, kk, x0, dx, dx2;
2966 real v, f, dvdlambda;
2969 kk = L1*kA+lambda*kB;
2970 x0 = L1*xA+lambda*xB;
2977 dvdlambda = half*(kB-kA)*dx2 + (xA-xB)*kk*dx;
2984 /* That was 21 flops */
2987 real g96bonds(int nbonds,
2988 const t_iatom forceatoms[], const t_iparams forceparams[],
2989 const rvec x[], rvec f[], rvec fshift[],
2990 const t_pbc *pbc, const t_graph *g,
2991 real lambda, real *dvdlambda,
2992 const t_mdatoms *md, t_fcdata *fcd,
2993 int *global_atom_index)
2995 int i, m, ki, ai, aj, type;
2996 real dr2, fbond, vbond, fij, vtot;
3001 for (i = 0; (i < nbonds); )
3003 type = forceatoms[i++];
3004 ai = forceatoms[i++];
3005 aj = forceatoms[i++];
3007 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
3008 dr2 = iprod(dx, dx); /* 5 */
3010 *dvdlambda += g96harmonic(forceparams[type].harmonic.krA,
3011 forceparams[type].harmonic.krB,
3012 forceparams[type].harmonic.rA,
3013 forceparams[type].harmonic.rB,
3014 dr2, lambda, &vbond, &fbond);
3016 vtot += 0.5*vbond; /* 1*/
3020 fprintf(debug, "G96-BONDS: dr = %10g vbond = %10g fbond = %10g\n",
3021 sqrt(dr2), vbond, fbond);
3027 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
3030 for (m = 0; (m < DIM); m++) /* 15 */
3035 fshift[ki][m] += fij;
3036 fshift[CENTRAL][m] -= fij;
3042 real g96bond_angle(const rvec xi, const rvec xj, const rvec xk, const t_pbc *pbc,
3043 rvec r_ij, rvec r_kj,
3045 /* Return value is the angle between the bonds i-j and j-k */
3049 *t1 = pbc_rvec_sub(pbc, xi, xj, r_ij); /* 3 */
3050 *t2 = pbc_rvec_sub(pbc, xk, xj, r_kj); /* 3 */
3052 costh = cos_angle(r_ij, r_kj); /* 25 */
3057 real g96angles(int nbonds,
3058 const t_iatom forceatoms[], const t_iparams forceparams[],
3059 const rvec x[], rvec f[], rvec fshift[],
3060 const t_pbc *pbc, const t_graph *g,
3061 real lambda, real *dvdlambda,
3062 const t_mdatoms *md, t_fcdata *fcd,
3063 int *global_atom_index)
3065 int i, ai, aj, ak, type, m, t1, t2;
3067 real cos_theta, dVdt, va, vtot;
3068 real rij_1, rij_2, rkj_1, rkj_2, rijrkj_1;
3070 ivec jt, dt_ij, dt_kj;
3073 for (i = 0; (i < nbonds); )
3075 type = forceatoms[i++];
3076 ai = forceatoms[i++];
3077 aj = forceatoms[i++];
3078 ak = forceatoms[i++];
3080 cos_theta = g96bond_angle(x[ai], x[aj], x[ak], pbc, r_ij, r_kj, &t1, &t2);
3082 *dvdlambda += g96harmonic(forceparams[type].harmonic.krA,
3083 forceparams[type].harmonic.krB,
3084 forceparams[type].harmonic.rA,
3085 forceparams[type].harmonic.rB,
3086 cos_theta, lambda, &va, &dVdt);
3089 rij_1 = gmx_invsqrt(iprod(r_ij, r_ij));
3090 rkj_1 = gmx_invsqrt(iprod(r_kj, r_kj));
3091 rij_2 = rij_1*rij_1;
3092 rkj_2 = rkj_1*rkj_1;
3093 rijrkj_1 = rij_1*rkj_1; /* 23 */
3098 fprintf(debug, "G96ANGLES: costheta = %10g vth = %10g dV/dct = %10g\n",
3099 cos_theta, va, dVdt);
3102 for (m = 0; (m < DIM); m++) /* 42 */
3104 f_i[m] = dVdt*(r_kj[m]*rijrkj_1 - r_ij[m]*rij_2*cos_theta);
3105 f_k[m] = dVdt*(r_ij[m]*rijrkj_1 - r_kj[m]*rkj_2*cos_theta);
3106 f_j[m] = -f_i[m]-f_k[m];
3114 copy_ivec(SHIFT_IVEC(g, aj), jt);
3116 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3117 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3118 t1 = IVEC2IS(dt_ij);
3119 t2 = IVEC2IS(dt_kj);
3121 rvec_inc(fshift[t1], f_i);
3122 rvec_inc(fshift[CENTRAL], f_j);
3123 rvec_inc(fshift[t2], f_k); /* 9 */
3129 real cross_bond_bond(int nbonds,
3130 const t_iatom forceatoms[], const t_iparams forceparams[],
3131 const rvec x[], rvec f[], rvec fshift[],
3132 const t_pbc *pbc, const t_graph *g,
3133 real lambda, real *dvdlambda,
3134 const t_mdatoms *md, t_fcdata *fcd,
3135 int *global_atom_index)
3137 /* Potential from Lawrence and Skimmer, Chem. Phys. Lett. 372 (2003)
3140 int i, ai, aj, ak, type, m, t1, t2;
3142 real vtot, vrr, s1, s2, r1, r2, r1e, r2e, krr;
3144 ivec jt, dt_ij, dt_kj;
3147 for (i = 0; (i < nbonds); )
3149 type = forceatoms[i++];
3150 ai = forceatoms[i++];
3151 aj = forceatoms[i++];
3152 ak = forceatoms[i++];
3153 r1e = forceparams[type].cross_bb.r1e;
3154 r2e = forceparams[type].cross_bb.r2e;
3155 krr = forceparams[type].cross_bb.krr;
3157 /* Compute distance vectors ... */
3158 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
3159 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
3161 /* ... and their lengths */
3165 /* Deviations from ideality */
3169 /* Energy (can be negative!) */
3174 svmul(-krr*s2/r1, r_ij, f_i);
3175 svmul(-krr*s1/r2, r_kj, f_k);
3177 for (m = 0; (m < DIM); m++) /* 12 */
3179 f_j[m] = -f_i[m] - f_k[m];
3188 copy_ivec(SHIFT_IVEC(g, aj), jt);
3190 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3191 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3192 t1 = IVEC2IS(dt_ij);
3193 t2 = IVEC2IS(dt_kj);
3195 rvec_inc(fshift[t1], f_i);
3196 rvec_inc(fshift[CENTRAL], f_j);
3197 rvec_inc(fshift[t2], f_k); /* 9 */
3203 real cross_bond_angle(int nbonds,
3204 const t_iatom forceatoms[], const t_iparams forceparams[],
3205 const rvec x[], rvec f[], rvec fshift[],
3206 const t_pbc *pbc, const t_graph *g,
3207 real lambda, real *dvdlambda,
3208 const t_mdatoms *md, t_fcdata *fcd,
3209 int *global_atom_index)
3211 /* Potential from Lawrence and Skimmer, Chem. Phys. Lett. 372 (2003)
3214 int i, ai, aj, ak, type, m, t1, t2, t3;
3215 rvec r_ij, r_kj, r_ik;
3216 real vtot, vrt, s1, s2, s3, r1, r2, r3, r1e, r2e, r3e, krt, k1, k2, k3;
3218 ivec jt, dt_ij, dt_kj;
3221 for (i = 0; (i < nbonds); )
3223 type = forceatoms[i++];
3224 ai = forceatoms[i++];
3225 aj = forceatoms[i++];
3226 ak = forceatoms[i++];
3227 r1e = forceparams[type].cross_ba.r1e;
3228 r2e = forceparams[type].cross_ba.r2e;
3229 r3e = forceparams[type].cross_ba.r3e;
3230 krt = forceparams[type].cross_ba.krt;
3232 /* Compute distance vectors ... */
3233 t1 = pbc_rvec_sub(pbc, x[ai], x[aj], r_ij);
3234 t2 = pbc_rvec_sub(pbc, x[ak], x[aj], r_kj);
3235 t3 = pbc_rvec_sub(pbc, x[ai], x[ak], r_ik);
3237 /* ... and their lengths */
3242 /* Deviations from ideality */
3247 /* Energy (can be negative!) */
3248 vrt = krt*s3*(s1+s2);
3254 k3 = -krt*(s1+s2)/r3;
3255 for (m = 0; (m < DIM); m++)
3257 f_i[m] = k1*r_ij[m] + k3*r_ik[m];
3258 f_k[m] = k2*r_kj[m] - k3*r_ik[m];
3259 f_j[m] = -f_i[m] - f_k[m];
3262 for (m = 0; (m < DIM); m++) /* 12 */
3272 copy_ivec(SHIFT_IVEC(g, aj), jt);
3274 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3275 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3276 t1 = IVEC2IS(dt_ij);
3277 t2 = IVEC2IS(dt_kj);
3279 rvec_inc(fshift[t1], f_i);
3280 rvec_inc(fshift[CENTRAL], f_j);
3281 rvec_inc(fshift[t2], f_k); /* 9 */
3287 static real bonded_tab(const char *type, int table_nr,
3288 const bondedtable_t *table, real kA, real kB, real r,
3289 real lambda, real *V, real *F)
3291 real k, tabscale, *VFtab, rt, eps, eps2, Yt, Ft, Geps, Heps2, Fp, VV, FF;
3293 real v, f, dvdlambda;
3295 k = (1.0 - lambda)*kA + lambda*kB;
3297 tabscale = table->scale;
3298 VFtab = table->data;
3304 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",
3305 type, table_nr, r, n0, n0+1, table->n);
3312 Geps = VFtab[nnn+2]*eps;
3313 Heps2 = VFtab[nnn+3]*eps2;
3314 Fp = Ft + Geps + Heps2;
3316 FF = Fp + Geps + 2.0*Heps2;
3318 *F = -k*FF*tabscale;
3320 dvdlambda = (kB - kA)*VV;
3324 /* That was 22 flops */
3327 real tab_bonds(int nbonds,
3328 const t_iatom forceatoms[], const t_iparams forceparams[],
3329 const rvec x[], rvec f[], rvec fshift[],
3330 const t_pbc *pbc, const t_graph *g,
3331 real lambda, real *dvdlambda,
3332 const t_mdatoms *md, t_fcdata *fcd,
3333 int *global_atom_index)
3335 int i, m, ki, ai, aj, type, table;
3336 real dr, dr2, fbond, vbond, fij, vtot;
3341 for (i = 0; (i < nbonds); )
3343 type = forceatoms[i++];
3344 ai = forceatoms[i++];
3345 aj = forceatoms[i++];
3347 ki = pbc_rvec_sub(pbc, x[ai], x[aj], dx); /* 3 */
3348 dr2 = iprod(dx, dx); /* 5 */
3349 dr = dr2*gmx_invsqrt(dr2); /* 10 */
3351 table = forceparams[type].tab.table;
3353 *dvdlambda += bonded_tab("bond", table,
3354 &fcd->bondtab[table],
3355 forceparams[type].tab.kA,
3356 forceparams[type].tab.kB,
3357 dr, lambda, &vbond, &fbond); /* 22 */
3365 vtot += vbond; /* 1*/
3366 fbond *= gmx_invsqrt(dr2); /* 6 */
3370 fprintf(debug, "TABBONDS: dr = %10g vbond = %10g fbond = %10g\n",
3376 ivec_sub(SHIFT_IVEC(g, ai), SHIFT_IVEC(g, aj), dt);
3379 for (m = 0; (m < DIM); m++) /* 15 */
3384 fshift[ki][m] += fij;
3385 fshift[CENTRAL][m] -= fij;
3391 real tab_angles(int nbonds,
3392 const t_iatom forceatoms[], const t_iparams forceparams[],
3393 const rvec x[], rvec f[], rvec fshift[],
3394 const t_pbc *pbc, const t_graph *g,
3395 real lambda, real *dvdlambda,
3396 const t_mdatoms *md, t_fcdata *fcd,
3397 int *global_atom_index)
3399 int i, ai, aj, ak, t1, t2, type, table;
3401 real cos_theta, cos_theta2, theta, dVdt, va, vtot;
3402 ivec jt, dt_ij, dt_kj;
3405 for (i = 0; (i < nbonds); )
3407 type = forceatoms[i++];
3408 ai = forceatoms[i++];
3409 aj = forceatoms[i++];
3410 ak = forceatoms[i++];
3412 theta = bond_angle(x[ai], x[aj], x[ak], pbc,
3413 r_ij, r_kj, &cos_theta, &t1, &t2); /* 41 */
3415 table = forceparams[type].tab.table;
3417 *dvdlambda += bonded_tab("angle", table,
3418 &fcd->angletab[table],
3419 forceparams[type].tab.kA,
3420 forceparams[type].tab.kB,
3421 theta, lambda, &va, &dVdt); /* 22 */
3424 cos_theta2 = sqr(cos_theta); /* 1 */
3433 st = dVdt*gmx_invsqrt(1 - cos_theta2); /* 12 */
3434 sth = st*cos_theta; /* 1 */
3438 fprintf(debug, "ANGLES: theta = %10g vth = %10g dV/dtheta = %10g\n",
3439 theta*RAD2DEG, va, dVdt);
3442 nrkj2 = iprod(r_kj, r_kj); /* 5 */
3443 nrij2 = iprod(r_ij, r_ij);
3445 cik = st*gmx_invsqrt(nrkj2*nrij2); /* 12 */
3446 cii = sth/nrij2; /* 10 */
3447 ckk = sth/nrkj2; /* 10 */
3449 for (m = 0; (m < DIM); m++) /* 39 */
3451 f_i[m] = -(cik*r_kj[m]-cii*r_ij[m]);
3452 f_k[m] = -(cik*r_ij[m]-ckk*r_kj[m]);
3453 f_j[m] = -f_i[m]-f_k[m];
3460 copy_ivec(SHIFT_IVEC(g, aj), jt);
3462 ivec_sub(SHIFT_IVEC(g, ai), jt, dt_ij);
3463 ivec_sub(SHIFT_IVEC(g, ak), jt, dt_kj);
3464 t1 = IVEC2IS(dt_ij);
3465 t2 = IVEC2IS(dt_kj);
3467 rvec_inc(fshift[t1], f_i);
3468 rvec_inc(fshift[CENTRAL], f_j);
3469 rvec_inc(fshift[t2], f_k);
3475 real tab_dihs(int nbonds,
3476 const t_iatom forceatoms[], const t_iparams forceparams[],
3477 const rvec x[], rvec f[], rvec fshift[],
3478 const t_pbc *pbc, const t_graph *g,
3479 real lambda, real *dvdlambda,
3480 const t_mdatoms *md, t_fcdata *fcd,
3481 int *global_atom_index)
3483 int i, type, ai, aj, ak, al, table;
3485 rvec r_ij, r_kj, r_kl, m, n;
3486 real phi, sign, ddphi, vpd, vtot;
3489 for (i = 0; (i < nbonds); )
3491 type = forceatoms[i++];
3492 ai = forceatoms[i++];
3493 aj = forceatoms[i++];
3494 ak = forceatoms[i++];
3495 al = forceatoms[i++];
3497 phi = dih_angle(x[ai], x[aj], x[ak], x[al], pbc, r_ij, r_kj, r_kl, m, n,
3498 &sign, &t1, &t2, &t3); /* 84 */
3500 table = forceparams[type].tab.table;
3502 /* Hopefully phi+M_PI never results in values < 0 */
3503 *dvdlambda += bonded_tab("dihedral", table,
3504 &fcd->dihtab[table],
3505 forceparams[type].tab.kA,
3506 forceparams[type].tab.kB,
3507 phi+M_PI, lambda, &vpd, &ddphi);
3510 do_dih_fup(ai, aj, ak, al, -ddphi, r_ij, r_kj, r_kl, m, n,
3511 f, fshift, pbc, g, x, t1, t2, t3); /* 112 */
3514 fprintf(debug, "pdih: (%d,%d,%d,%d) phi=%g\n",
3515 ai, aj, ak, al, phi);
3523 calc_bonded_reduction_mask(const t_idef *idef,
3528 int ftype, nb, nat1, nb0, nb1, i, a;
3532 for (ftype = 0; ftype < F_NRE; ftype++)
3534 if (interaction_function[ftype].flags & IF_BOND &&
3535 !(ftype == F_CONNBONDS || ftype == F_POSRES) &&
3536 (ftype<F_GB12 || ftype>F_GB14))
3538 nb = idef->il[ftype].nr;
3541 nat1 = interaction_function[ftype].nratoms + 1;
3543 /* Divide this interaction equally over the threads.
3544 * This is not stored: should match division in calc_bonds.
3546 nb0 = (((nb/nat1)* t )/nt)*nat1;
3547 nb1 = (((nb/nat1)*(t+1))/nt)*nat1;
3549 for (i = nb0; i < nb1; i += nat1)
3551 for (a = 1; a < nat1; a++)
3553 mask |= (1U << (idef->il[ftype].iatoms[i+a]>>shift));
3563 void init_bonded_thread_force_reduction(t_forcerec *fr,
3566 #define MAX_BLOCK_BITS 32
3570 if (fr->nthreads <= 1)
3577 /* We divide the force array in a maximum of 32 blocks.
3578 * Minimum force block reduction size is 2^6=64.
3581 while (fr->natoms_force > (int)(MAX_BLOCK_BITS*(1U<<fr->red_ashift)))
3587 fprintf(debug, "bonded force buffer block atom shift %d bits\n",
3591 /* Determine to which blocks each thread's bonded force calculation
3592 * contributes. Store this is a mask for each thread.
3594 #pragma omp parallel for num_threads(fr->nthreads) schedule(static)
3595 for (t = 1; t < fr->nthreads; t++)
3597 fr->f_t[t].red_mask =
3598 calc_bonded_reduction_mask(idef, fr->red_ashift, t, fr->nthreads);
3601 /* Determine the maximum number of blocks we need to reduce over */
3604 for (t = 0; t < fr->nthreads; t++)
3607 for (b = 0; b < MAX_BLOCK_BITS; b++)
3609 if (fr->f_t[t].red_mask & (1U<<b))
3611 fr->red_nblock = max(fr->red_nblock, b+1);
3617 fprintf(debug, "thread %d flags %x count %d\n",
3618 t, fr->f_t[t].red_mask, c);
3624 fprintf(debug, "Number of blocks to reduce: %d of size %d\n",
3625 fr->red_nblock, 1<<fr->red_ashift);
3626 fprintf(debug, "Reduction density %.2f density/#thread %.2f\n",
3627 ctot*(1<<fr->red_ashift)/(double)fr->natoms_force,
3628 ctot*(1<<fr->red_ashift)/(double)(fr->natoms_force*fr->nthreads));
3632 static void zero_thread_forces(f_thread_t *f_t, int n,
3633 int nblock, int blocksize)
3635 int b, a0, a1, a, i, j;
3637 if (n > f_t->f_nalloc)
3639 f_t->f_nalloc = over_alloc_large(n);
3640 srenew(f_t->f, f_t->f_nalloc);
3643 if (f_t->red_mask != 0)
3645 for (b = 0; b < nblock; b++)
3647 if (f_t->red_mask && (1U<<b))
3650 a1 = min((b+1)*blocksize, n);
3651 for (a = a0; a < a1; a++)
3653 clear_rvec(f_t->f[a]);
3658 for (i = 0; i < SHIFTS; i++)
3660 clear_rvec(f_t->fshift[i]);
3662 for (i = 0; i < F_NRE; i++)
3666 for (i = 0; i < egNR; i++)
3668 for (j = 0; j < f_t->grpp.nener; j++)
3670 f_t->grpp.ener[i][j] = 0;
3673 for (i = 0; i < efptNR; i++)
3679 static void reduce_thread_force_buffer(int n, rvec *f,
3680 int nthreads, f_thread_t *f_t,
3681 int nblock, int block_size)
3683 /* The max thread number is arbitrary,
3684 * we used a fixed number to avoid memory management.
3685 * Using more than 16 threads is probably never useful performance wise.
3687 #define MAX_BONDED_THREADS 256
3690 if (nthreads > MAX_BONDED_THREADS)
3692 gmx_fatal(FARGS, "Can not reduce bonded forces on more than %d threads",
3693 MAX_BONDED_THREADS);
3696 /* This reduction can run on any number of threads,
3697 * independently of nthreads.
3699 #pragma omp parallel for num_threads(nthreads) schedule(static)
3700 for (b = 0; b < nblock; b++)
3702 rvec *fp[MAX_BONDED_THREADS];
3706 /* Determine which threads contribute to this block */
3708 for (ft = 1; ft < nthreads; ft++)
3710 if (f_t[ft].red_mask & (1U<<b))
3712 fp[nfb++] = f_t[ft].f;
3717 /* Reduce force buffers for threads that contribute */
3719 a1 = (b+1)*block_size;
3721 for (a = a0; a < a1; a++)
3723 for (fb = 0; fb < nfb; fb++)
3725 rvec_inc(f[a], fp[fb][a]);
3732 static void reduce_thread_forces(int n, rvec *f, rvec *fshift,
3733 real *ener, gmx_grppairener_t *grpp, real *dvdl,
3734 int nthreads, f_thread_t *f_t,
3735 int nblock, int block_size,
3736 gmx_bool bCalcEnerVir,
3741 /* Reduce the bonded force buffer */
3742 reduce_thread_force_buffer(n, f, nthreads, f_t, nblock, block_size);
3745 /* When necessary, reduce energy and virial using one thread only */
3750 for (i = 0; i < SHIFTS; i++)
3752 for (t = 1; t < nthreads; t++)
3754 rvec_inc(fshift[i], f_t[t].fshift[i]);
3757 for (i = 0; i < F_NRE; i++)
3759 for (t = 1; t < nthreads; t++)
3761 ener[i] += f_t[t].ener[i];
3764 for (i = 0; i < egNR; i++)
3766 for (j = 0; j < f_t[1].grpp.nener; j++)
3768 for (t = 1; t < nthreads; t++)
3771 grpp->ener[i][j] += f_t[t].grpp.ener[i][j];
3777 for (i = 0; i < efptNR; i++)
3780 for (t = 1; t < nthreads; t++)
3782 dvdl[i] += f_t[t].dvdl[i];
3789 static real calc_one_bond(FILE *fplog, int thread,
3790 int ftype, const t_idef *idef,
3791 rvec x[], rvec f[], rvec fshift[],
3793 const t_pbc *pbc, const t_graph *g,
3794 gmx_enerdata_t *enerd, gmx_grppairener_t *grpp,
3796 real *lambda, real *dvdl,
3797 const t_mdatoms *md, t_fcdata *fcd,
3798 gmx_bool bCalcEnerVir,
3799 int *global_atom_index, gmx_bool bPrintSepPot)
3801 int ind, nat1, nbonds, efptFTYPE;
3806 if (IS_RESTRAINT_TYPE(ftype))
3808 efptFTYPE = efptRESTRAINT;
3812 efptFTYPE = efptBONDED;
3815 if (interaction_function[ftype].flags & IF_BOND &&
3816 !(ftype == F_CONNBONDS || ftype == F_POSRES))
3818 ind = interaction_function[ftype].nrnb_ind;
3819 nat1 = interaction_function[ftype].nratoms + 1;
3820 nbonds = idef->il[ftype].nr/nat1;
3821 iatoms = idef->il[ftype].iatoms;
3823 nb0 = ((nbonds* thread )/(fr->nthreads))*nat1;
3824 nbn = ((nbonds*(thread+1))/(fr->nthreads))*nat1 - nb0;
3826 if (!IS_LISTED_LJ_C(ftype))
3828 if (ftype == F_CMAP)
3830 v = cmap_dihs(nbn, iatoms+nb0,
3831 idef->iparams, &idef->cmap_grid,
3832 (const rvec*)x, f, fshift,
3833 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
3834 md, fcd, global_atom_index);
3836 else if (ftype == F_PDIHS &&
3837 !bCalcEnerVir && fr->efep == efepNO)
3839 /* No energies, shift forces, dvdl */
3840 #ifndef SSE_PROPER_DIHEDRALS
3845 (nbn, idef->il[ftype].iatoms+nb0,
3848 pbc, g, lambda[efptFTYPE], md, fcd,
3854 v = interaction_function[ftype].ifunc(nbn, iatoms+nb0,
3856 (const rvec*)x, f, fshift,
3857 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]),
3858 md, fcd, global_atom_index);
3862 fprintf(fplog, " %-23s #%4d V %12.5e dVdl %12.5e\n",
3863 interaction_function[ftype].longname,
3864 nbonds/nat1, v, lambda[efptFTYPE]);
3869 v = do_nonbonded_listed(ftype, nbn, iatoms+nb0, idef->iparams, (const rvec*)x, f, fshift,
3870 pbc, g, lambda, dvdl, md, fr, grpp, global_atom_index);
3872 enerd->dvdl_nonlin[efptCOUL] += dvdl[efptCOUL];
3873 enerd->dvdl_nonlin[efptVDW] += dvdl[efptVDW];
3877 fprintf(fplog, " %-5s + %-15s #%4d dVdl %12.5e\n",
3878 interaction_function[ftype].longname,
3879 interaction_function[F_LJ14].longname, nbonds/nat1, dvdl[efptVDW]);
3880 fprintf(fplog, " %-5s + %-15s #%4d dVdl %12.5e\n",
3881 interaction_function[ftype].longname,
3882 interaction_function[F_COUL14].longname, nbonds/nat1, dvdl[efptCOUL]);
3885 if (ind != -1 && thread == 0)
3887 inc_nrnb(nrnb, ind, nbonds);
3894 /* WARNING! THIS FUNCTION MUST EXACTLY TRACK THE calc
3895 function, or horrible things will happen when doing free energy
3896 calculations! In a good coding world, this would not be a
3897 different function, but for speed reasons, it needs to be made a
3898 separate function. TODO for 5.0 - figure out a way to reorganize
3899 to reduce duplication.
3902 static real calc_one_bond_foreign(FILE *fplog, int ftype, const t_idef *idef,
3903 rvec x[], rvec f[], t_forcerec *fr,
3904 const t_pbc *pbc, const t_graph *g,
3905 gmx_grppairener_t *grpp, t_nrnb *nrnb,
3906 real *lambda, real *dvdl,
3907 const t_mdatoms *md, t_fcdata *fcd,
3908 int *global_atom_index, gmx_bool bPrintSepPot)
3910 int ind, nat1, nbonds, efptFTYPE, nbonds_np;
3914 if (IS_RESTRAINT_TYPE(ftype))
3916 efptFTYPE = efptRESTRAINT;
3920 efptFTYPE = efptBONDED;
3923 if (ftype < F_GB12 || ftype > F_GB14)
3925 if (interaction_function[ftype].flags & IF_BOND &&
3926 !(ftype == F_CONNBONDS || ftype == F_POSRES || ftype == F_FBPOSRES))
3928 ind = interaction_function[ftype].nrnb_ind;
3929 nat1 = interaction_function[ftype].nratoms+1;
3930 nbonds_np = idef->il[ftype].nr_nonperturbed;
3931 nbonds = idef->il[ftype].nr - nbonds_np;
3932 iatoms = idef->il[ftype].iatoms + nbonds_np;
3935 if (!IS_LISTED_LJ_C(ftype))
3937 if (ftype == F_CMAP)
3939 v = cmap_dihs(nbonds, iatoms,
3940 idef->iparams, &idef->cmap_grid,
3941 (const rvec*)x, f, fr->fshift,
3942 pbc, g, lambda[efptFTYPE], &(dvdl[efptFTYPE]), md, fcd,
3947 v = interaction_function[ftype].ifunc(nbonds, iatoms,
3949 (const rvec*)x, f, fr->fshift,
3950 pbc, g, lambda[efptFTYPE], &dvdl[efptFTYPE],
3951 md, fcd, global_atom_index);
3956 v = do_nonbonded_listed(ftype, nbonds, iatoms,
3958 (const rvec*)x, f, fr->fshift,
3959 pbc, g, lambda, dvdl,
3960 md, fr, grpp, global_atom_index);
3964 inc_nrnb(nrnb, ind, nbonds/nat1);
3972 void calc_bonds(FILE *fplog, const gmx_multisim_t *ms,
3974 rvec x[], history_t *hist,
3975 rvec f[], t_forcerec *fr,
3976 const t_pbc *pbc, const t_graph *g,
3977 gmx_enerdata_t *enerd, t_nrnb *nrnb,
3979 const t_mdatoms *md,
3980 t_fcdata *fcd, int *global_atom_index,
3981 t_atomtypes *atype, gmx_genborn_t *born,
3983 gmx_bool bPrintSepPot, gmx_large_int_t step)
3985 gmx_bool bCalcEnerVir;
3987 real v, dvdl[efptNR], dvdl_dum[efptNR]; /* The dummy array is to have a place to store the dhdl at other values
3988 of lambda, which will be thrown away in the end*/
3989 const t_pbc *pbc_null;
3993 bCalcEnerVir = (force_flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY));
3995 for (i = 0; i < efptNR; i++)
4009 fprintf(fplog, "Step %s: bonded V and dVdl for this node\n",
4010 gmx_step_str(step, buf));
4016 p_graph(debug, "Bondage is fun", g);
4020 /* Do pre force calculation stuff which might require communication */
4021 if (idef->il[F_ORIRES].nr)
4023 enerd->term[F_ORIRESDEV] =
4024 calc_orires_dev(ms, idef->il[F_ORIRES].nr,
4025 idef->il[F_ORIRES].iatoms,
4026 idef->iparams, md, (const rvec*)x,
4027 pbc_null, fcd, hist);
4029 if (idef->il[F_DISRES].nr)
4031 calc_disres_R_6(ms, idef->il[F_DISRES].nr,
4032 idef->il[F_DISRES].iatoms,
4033 idef->iparams, (const rvec*)x, pbc_null,
4037 #pragma omp parallel for num_threads(fr->nthreads) schedule(static)
4038 for (thread = 0; thread < fr->nthreads; thread++)
4040 int ftype, nbonds, ind, nat1;
4045 gmx_grppairener_t *grpp;
4051 fshift = fr->fshift;
4053 grpp = &enerd->grpp;
4058 zero_thread_forces(&fr->f_t[thread], fr->natoms_force,
4059 fr->red_nblock, 1<<fr->red_ashift);
4061 ft = fr->f_t[thread].f;
4062 fshift = fr->f_t[thread].fshift;
4063 epot = fr->f_t[thread].ener;
4064 grpp = &fr->f_t[thread].grpp;
4065 dvdlt = fr->f_t[thread].dvdl;
4067 /* Loop over all bonded force types to calculate the bonded forces */
4068 for (ftype = 0; (ftype < F_NRE); ftype++)
4070 if (idef->il[ftype].nr > 0 &&
4071 (interaction_function[ftype].flags & IF_BOND) &&
4072 (ftype < F_GB12 || ftype > F_GB14) &&
4073 !(ftype == F_CONNBONDS || ftype == F_POSRES))
4075 v = calc_one_bond(fplog, thread, ftype, idef, x,
4076 ft, fshift, fr, pbc_null, g, enerd, grpp,
4077 nrnb, lambda, dvdlt,
4078 md, fcd, bCalcEnerVir,
4079 global_atom_index, bPrintSepPot);
4084 if (fr->nthreads > 1)
4086 reduce_thread_forces(fr->natoms_force, f, fr->fshift,
4087 enerd->term, &enerd->grpp, dvdl,
4088 fr->nthreads, fr->f_t,
4089 fr->red_nblock, 1<<fr->red_ashift,
4091 force_flags & GMX_FORCE_DHDL);
4093 if (force_flags & GMX_FORCE_DHDL)
4095 for (i = 0; i < efptNR; i++)
4097 enerd->dvdl_nonlin[i] += dvdl[i];
4101 /* Copy the sum of violations for the distance restraints from fcd */
4104 enerd->term[F_DISRESVIOL] = fcd->disres.sumviol;
4109 void calc_bonds_lambda(FILE *fplog,
4113 const t_pbc *pbc, const t_graph *g,
4114 gmx_grppairener_t *grpp, real *epot, t_nrnb *nrnb,
4116 const t_mdatoms *md,
4118 int *global_atom_index)
4120 int i, ftype, nbonds_np, nbonds, ind, nat;
4122 real dvdl_dum[efptNR];
4123 rvec *f, *fshift_orig;
4124 const t_pbc *pbc_null;
4136 snew(f, fr->natoms_force);
4137 /* We want to preserve the fshift array in forcerec */
4138 fshift_orig = fr->fshift;
4139 snew(fr->fshift, SHIFTS);
4141 /* Loop over all bonded force types to calculate the bonded forces */
4142 for (ftype = 0; (ftype < F_NRE); ftype++)
4144 v = calc_one_bond_foreign(fplog, ftype, idef, x,
4145 f, fr, pbc_null, g, grpp, nrnb, lambda, dvdl_dum,
4146 md, fcd, global_atom_index, FALSE);
4151 fr->fshift = fshift_orig;