3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * GROningen Mixture of Alchemy and Childrens' Stories
53 real calc_ewaldcoeff(real rc, real dtol)
55 real x = 5, low, high;
64 while (gmx_erfc(x*rc) > dtol);
66 n = i+60; /* search tolerance is 2^-60 */
69 for (i = 0; i < n; i++)
72 if (gmx_erfc(x*rc) > dtol)
86 real ewald_LRcorrection(int start, int end,
87 t_commrec *cr, int thread, t_forcerec *fr,
88 real *chargeA, real *chargeB,
89 gmx_bool calc_excl_corr,
90 t_blocka *excl, rvec x[],
91 matrix box, rvec mu_tot[],
92 int ewald_geometry, real epsilon_surface,
94 real lambda, real *dvdlambda)
96 int i, i1, i2, j, k, m, iv, jv, q;
98 double q2sumA, q2sumB, Vexcl, dvdl_excl; /* Necessary for precision */
100 real v, vc, qiA, qiB, dr, dr2, rinv, fscal, enercorr;
101 real Vself[2], Vdipole[2], rinv2, ewc = fr->ewaldcoeff, ewcdr;
102 rvec df, dx, mutot[2], dipcorrA, dipcorrB;
104 real vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
105 real L1, dipole_coeff, qqA, qqB, qqL, vr0;
108 real tabscale = fr->tabscale;
109 real eps, eps2, VV, FF, F, Y, Geps, Heps2, Fp, fijC, r1t;
110 real *VFtab = fr->coulvdwtab;
113 gmx_bool bFreeEnergy = (chargeB != NULL);
114 gmx_bool bMolPBC = fr->bMolPBC;
116 one_4pi_eps = ONE_4PI_EPS0/fr->epsilon_r;
117 vr0 = ewc*M_2_SQRTPI;
128 /* Note that we have to transform back to gromacs units, since
129 * mu_tot contains the dipole in debye units (for output).
131 for (i = 0; (i < DIM); i++)
133 mutot[0][i] = mu_tot[0][i]*DEBYE2ENM;
134 mutot[1][i] = mu_tot[1][i]*DEBYE2ENM;
139 switch (ewald_geometry)
142 if (epsilon_surface != 0)
145 2*M_PI*ONE_4PI_EPS0/((2*epsilon_surface + fr->epsilon_r)*vol);
146 for (i = 0; (i < DIM); i++)
148 dipcorrA[i] = 2*dipole_coeff*mutot[0][i];
149 dipcorrB[i] = 2*dipole_coeff*mutot[1][i];
154 dipole_coeff = 2*M_PI*one_4pi_eps/vol;
155 dipcorrA[ZZ] = 2*dipole_coeff*mutot[0][ZZ];
156 dipcorrB[ZZ] = 2*dipole_coeff*mutot[1][ZZ];
159 gmx_incons("Unsupported Ewald geometry");
164 fprintf(debug, "dipcorr = %8.3f %8.3f %8.3f\n",
165 dipcorrA[XX], dipcorrA[YY], dipcorrA[ZZ]);
166 fprintf(debug, "mutot = %8.3f %8.3f %8.3f\n",
167 mutot[0][XX], mutot[0][YY], mutot[0][ZZ]);
171 if ((calc_excl_corr || dipole_coeff != 0) && !bFreeEnergy)
173 for (i = start; (i < end); i++)
175 /* Initiate local variables (for this i-particle) to 0 */
176 qiA = chargeA[i]*one_4pi_eps;
181 i2 = excl->index[i+1];
183 /* Loop over excluded neighbours */
184 for (j = i1; (j < i2); j++)
188 * First we must test whether k <> i, and then, because the
189 * exclusions are all listed twice i->k and k->i we must select
190 * just one of the two.
191 * As a minor optimization we only compute forces when the charges
196 qqA = qiA*chargeA[k];
199 rvec_sub(x[i], x[k], dx);
202 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
203 for (m = DIM-1; (m >= 0); m--)
205 if (dx[m] > 0.5*box[m][m])
207 rvec_dec(dx, box[m]);
209 else if (dx[m] < -0.5*box[m][m])
211 rvec_inc(dx, box[m]);
216 /* Distance between two excluded particles may be zero in the
221 rinv = gmx_invsqrt(dr2);
234 Geps = eps*VFtab[nnn+2];
235 Heps2 = eps2*VFtab[nnn+3];
238 FF = Fp+Geps+2.0*Heps2;
243 fscal = vc*rinv2+fijC*tabscale*rinv;
244 /* End of tabulated interaction part */
247 /* This is the code you would want instead if not using
251 vc = qqA*gmx_erf(ewcdr)*rinv;
254 /* Relative accuracy at R_ERF_R_INACC of 3e-10 */
255 #define R_ERF_R_INACC 0.006
257 /* Relative accuracy at R_ERF_R_INACC of 2e-5 */
258 #define R_ERF_R_INACC 0.1
260 if (ewcdr > R_ERF_R_INACC)
262 fscal = rinv2*(vc - qqA*ewc*M_2_SQRTPI*exp(-ewcdr*ewcdr));
266 /* Use a fourth order series expansion for small ewcdr */
267 fscal = ewc*ewc*qqA*vr0*(2.0/3.0 - 0.4*ewcdr*ewcdr);
270 /* The force vector is obtained by multiplication with the
273 svmul(fscal, dx, df);
276 for (iv = 0; (iv < DIM); iv++)
278 for (jv = 0; (jv < DIM); jv++)
280 dxdf[iv][jv] += dx[iv]*df[jv];
292 /* Dipole correction on force */
293 if (dipole_coeff != 0)
295 for (j = 0; (j < DIM); j++)
297 f[i][j] -= dipcorrA[j]*chargeA[i];
302 else if (calc_excl_corr || dipole_coeff != 0)
304 for (i = start; (i < end); i++)
306 /* Initiate local variables (for this i-particle) to 0 */
307 qiA = chargeA[i]*one_4pi_eps;
308 qiB = chargeB[i]*one_4pi_eps;
313 i2 = excl->index[i+1];
315 /* Loop over excluded neighbours */
316 for (j = i1; (j < i2); j++)
321 qqA = qiA*chargeA[k];
322 qqB = qiB*chargeB[k];
323 if (qqA != 0.0 || qqB != 0.0)
325 qqL = L1*qqA + lambda*qqB;
326 rvec_sub(x[i], x[k], dx);
329 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
330 for (m = DIM-1; (m >= 0); m--)
332 if (dx[m] > 0.5*box[m][m])
334 rvec_dec(dx, box[m]);
336 else if (dx[m] < -0.5*box[m][m])
338 rvec_inc(dx, box[m]);
345 rinv = gmx_invsqrt(dr2);
348 v = gmx_erf(ewc*dr)*rinv;
351 fscal = rinv2*(vc-qqL*ewc*M_2_SQRTPI*exp(-ewc*ewc*dr2));
352 svmul(fscal, dx, df);
355 for (iv = 0; (iv < DIM); iv++)
357 for (jv = 0; (jv < DIM); jv++)
359 dxdf[iv][jv] += dx[iv]*df[jv];
362 dvdl_excl += (qqB - qqA)*v;
367 dvdl_excl += (qqB - qqA)*vr0;
373 /* Dipole correction on force */
374 if (dipole_coeff != 0)
376 for (j = 0; (j < DIM); j++)
378 f[i][j] -= L1*dipcorrA[j]*chargeA[i]
379 + lambda*dipcorrB[j]*chargeB[i];
384 for (iv = 0; (iv < DIM); iv++)
386 for (jv = 0; (jv < DIM); jv++)
388 vir[iv][jv] += 0.5*dxdf[iv][jv];
395 /* Global corrections only on master process */
396 if (MASTER(cr) && thread == 0)
398 for (q = 0; q < (bFreeEnergy ? 2 : 1); q++)
402 /* Self-energy correction */
403 Vself[q] = ewc*one_4pi_eps*fr->q2sum[q]*M_1_SQRTPI;
406 /* Apply surface dipole correction:
407 * correction = dipole_coeff * (dipole)^2
409 if (dipole_coeff != 0)
411 if (ewald_geometry == eewg3D)
413 Vdipole[q] = dipole_coeff*iprod(mutot[q], mutot[q]);
415 else if (ewald_geometry == eewg3DC)
417 Vdipole[q] = dipole_coeff*mutot[q][ZZ]*mutot[q][ZZ];
425 enercorr = Vdipole[0] - Vself[0] - Vexcl;
429 enercorr = L1*(Vdipole[0] - Vself[0])
430 + lambda*(Vdipole[1] - Vself[1])
432 *dvdlambda += Vdipole[1] - Vself[1]
433 - (Vdipole[0] - Vself[0]) - dvdl_excl;
438 fprintf(debug, "Long Range corrections for Ewald interactions:\n");
439 fprintf(debug, "start=%d,natoms=%d\n", start, end-start);
440 fprintf(debug, "q2sum = %g, Vself=%g\n",
441 L1*q2sumA+lambda*q2sumB, L1*Vself[0]+lambda*Vself[1]);
442 fprintf(debug, "Long Range correction: Vexcl=%g\n", Vexcl);
443 if (MASTER(cr) && thread == 0)
445 if (epsilon_surface > 0 || ewald_geometry == eewg3DC)
447 fprintf(debug, "Total dipole correction: Vdipole=%g\n",
448 L1*Vdipole[0]+lambda*Vdipole[1]);
453 /* Return the correction to the energy */
457 real ewald_charge_correction(t_commrec *cr, t_forcerec *fr, real lambda,
459 real *dvdlambda, tensor vir)
462 real vol, fac, qs2A, qs2B, vc, enercorr;
467 /* Apply charge correction */
468 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
470 fac = M_PI*ONE_4PI_EPS0/(fr->epsilon_r*2.0*vol*vol*sqr(fr->ewaldcoeff));
472 qs2A = fr->qsum[0]*fr->qsum[0];
473 qs2B = fr->qsum[1]*fr->qsum[1];
475 vc = (qs2A*(1 - lambda) + qs2B*lambda)*fac;
479 *dvdlambda += -vol*(qs2B - qs2A)*fac;
481 for (d = 0; d < DIM; d++)
488 fprintf(debug, "Total charge correction: Vcharge=%g\n", enercorr);