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)
62 } while (gmx_erfc(x*rc) > dtol);
64 n=i+60; /* search tolerance is 2^-60 */
69 if (gmx_erfc(x*rc) > dtol)
79 real ewald_LRcorrection(FILE *fplog,
81 t_commrec *cr,t_forcerec *fr,
82 real *chargeA,real *chargeB,
83 t_blocka *excl,rvec x[],
84 matrix box,rvec mu_tot[],
85 int ewald_geometry,real epsilon_surface,
86 real lambda,real *dvdlambda,
87 real *vdip,real *vcharge)
89 int i,i1,i2,j,k,m,iv,jv,q;
91 double q2sumA,q2sumB,Vexcl,dvdl_excl; /* Necessary for precision */
93 real v,vc,qiA,qiB,dr,dr2,rinv,fscal,enercorr;
94 real VselfA,VselfB=0,Vcharge[2],Vdipole[2],rinv2,ewc=fr->ewaldcoeff,ewcdr;
95 rvec df,dx,mutot[2],dipcorrA,dipcorrB;
96 rvec *f=fr->f_novirsum;
98 real vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
99 real L1,dipole_coeff,qqA,qqB,qqL,vr0;
102 real tabscale=fr->tabscale;
103 real eps,eps2,VV,FF,F,Y,Geps,Heps2,Fp,fijC,r1t;
104 real *VFtab=fr->coulvdwtab;
107 double isp=0.564189583547756;
110 gmx_bool bFreeEnergy = (chargeB != NULL);
111 gmx_bool bMolPBC = fr->bMolPBC;
113 one_4pi_eps = ONE_4PI_EPS0/fr->epsilon_r;
114 vr0 = ewc*2/sqrt(M_PI);
127 /* Note that we have to transform back to gromacs units, since
128 * mu_tot contains the dipole in debye units (for output).
130 for(i=0; (i<DIM); i++) {
131 mutot[0][i] = mu_tot[0][i]*DEBYE2ENM;
132 mutot[1][i] = mu_tot[1][i]*DEBYE2ENM;
137 switch (ewald_geometry) {
139 if (epsilon_surface != 0) {
141 2*M_PI*ONE_4PI_EPS0/((2*epsilon_surface + fr->epsilon_r)*vol);
142 for(i=0; (i<DIM); i++) {
143 dipcorrA[i] = 2*dipole_coeff*mutot[0][i];
144 dipcorrB[i] = 2*dipole_coeff*mutot[1][i];
149 dipole_coeff = 2*M_PI*one_4pi_eps/vol;
150 dipcorrA[ZZ] = 2*dipole_coeff*mutot[0][ZZ];
151 dipcorrB[ZZ] = 2*dipole_coeff*mutot[1][ZZ];
154 gmx_incons("Unsupported Ewald geometry");
158 fprintf(debug,"dipcorr = %8.3f %8.3f %8.3f\n",
159 dipcorrA[XX],dipcorrA[YY],dipcorrA[ZZ]);
160 fprintf(debug,"mutot = %8.3f %8.3f %8.3f\n",
161 mutot[0][XX],mutot[0][YY],mutot[0][ZZ]);
164 if (DOMAINDECOMP(cr))
171 for(i=start; (i<niat); i++) {
172 /* Initiate local variables (for this i-particle) to 0 */
173 qiA = chargeA[i]*one_4pi_eps;
175 i2 = excl->index[i+1];
177 q2sumA += chargeA[i]*chargeA[i];
179 /* Loop over excluded neighbours */
180 for(j=i1; (j<i2); j++) {
183 * First we must test whether k <> i, and then, because the
184 * exclusions are all listed twice i->k and k->i we must select
185 * just one of the two.
186 * As a minor optimization we only compute forces when the charges
190 qqA = qiA*chargeA[k];
192 rvec_sub(x[i],x[k],dx);
194 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
195 for(m=DIM-1; (m>=0); m--) {
196 if (dx[m] > 0.5*box[m][m])
198 else if (dx[m] < -0.5*box[m][m])
203 /* Distance between two excluded particles may be zero in the
207 rinv = gmx_invsqrt(dr2);
220 Geps = eps*VFtab[nnn+2];
221 Heps2 = eps2*VFtab[nnn+3];
224 FF = Fp+Geps+2.0*Heps2;
229 fscal = vc*rinv2+fijC*tabscale*rinv;
230 /* End of tabulated interaction part */
233 /* This is the code you would want instead if not using
237 vc = qqA*gmx_erf(ewcdr)*rinv;
240 /* Relative accuracy at R_ERF_R_INACC of 3e-10 */
241 #define R_ERF_R_INACC 0.006
243 /* Relative accuracy at R_ERF_R_INACC of 2e-5 */
244 #define R_ERF_R_INACC 0.1
246 if (ewcdr > R_ERF_R_INACC) {
247 fscal = rinv2*(vc - 2.0*qqA*ewc*isp*exp(-ewcdr*ewcdr));
249 /* Use a fourth order series expansion for small ewcdr */
250 fscal = ewc*ewc*qqA*vr0*(2.0/3.0 - 0.4*ewcdr*ewcdr);
253 /* The force vector is obtained by multiplication with the
259 for(iv=0; (iv<DIM); iv++)
260 for(jv=0; (jv<DIM); jv++)
261 dxdf[iv][jv] += dx[iv]*df[jv];
268 /* Dipole correction on force */
269 if (dipole_coeff != 0) {
270 for(j=0; (j<DIM); j++)
271 f[i][j] -= dipcorrA[j]*chargeA[i];
275 for(i=start; (i<niat); i++) {
276 /* Initiate local variables (for this i-particle) to 0 */
277 qiA = chargeA[i]*one_4pi_eps;
278 qiB = chargeB[i]*one_4pi_eps;
280 i2 = excl->index[i+1];
282 q2sumA += chargeA[i]*chargeA[i];
283 q2sumB += chargeB[i]*chargeB[i];
286 /* Loop over excluded neighbours */
287 for(j=i1; (j<i2); j++) {
290 qqA = qiA*chargeA[k];
291 qqB = qiB*chargeB[k];
292 if (qqA != 0.0 || qqB != 0.0) {
293 qqL = L1*qqA + lambda*qqB;
294 rvec_sub(x[i],x[k],dx);
296 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
297 for(m=DIM-1; (m>=0); m--) {
298 if (dx[m] > 0.5*box[m][m])
300 else if (dx[m] < -0.5*box[m][m])
306 rinv = gmx_invsqrt(dr2);
309 v = gmx_erf(ewc*dr)*rinv;
312 fscal = rinv2*(vc-2.0*qqL*ewc*isp*exp(-ewc*ewc*dr2));
316 for(iv=0; (iv<DIM); iv++)
317 for(jv=0; (jv<DIM); jv++)
318 dxdf[iv][jv] += dx[iv]*df[jv];
319 dvdl_excl += (qqB - qqA)*v;
322 dvdl_excl += (qqB - qqA)*vr0;
327 /* Dipole correction on force */
328 if (dipole_coeff != 0) {
329 for(j=0; (j<DIM); j++)
330 f[i][j] -= L1*dipcorrA[j]*chargeA[i]
331 + lambda*dipcorrB[j]*chargeB[i];
335 for(iv=0; (iv<DIM); iv++)
336 for(jv=0; (jv<DIM); jv++)
337 fr->vir_el_recip[iv][jv] += 0.5*dxdf[iv][jv];
339 /* Global corrections only on master process */
341 for(q=0; q<(bFreeEnergy ? 2 : 1); q++) {
342 /* Apply charge correction */
343 /* use vc as a dummy variable */
344 vc = fr->qsum[q]*fr->qsum[q]*M_PI*one_4pi_eps/(2.0*vol*vol*ewc*ewc);
345 for(iv=0; (iv<DIM); iv++)
346 fr->vir_el_recip[iv][iv] +=
347 (bFreeEnergy ? (q==0 ? L1*vc : lambda*vc) : vc);
348 Vcharge[q] = -vol*vc;
350 /* Apply surface dipole correction:
351 * correction = dipole_coeff * (dipole)^2
353 if (dipole_coeff != 0) {
354 if (ewald_geometry == eewg3D)
355 Vdipole[q] = dipole_coeff*iprod(mutot[q],mutot[q]);
356 else if (ewald_geometry == eewg3DC)
357 Vdipole[q] = dipole_coeff*mutot[q][ZZ]*mutot[q][ZZ];
362 VselfA = ewc*one_4pi_eps*q2sumA/sqrt(M_PI);
365 *vcharge = Vcharge[0];
367 enercorr = *vcharge + *vdip - VselfA - Vexcl;
369 VselfB = ewc*one_4pi_eps*q2sumB/sqrt(M_PI);
370 *vcharge = L1*Vcharge[0] + lambda*Vcharge[1];
371 *vdip = L1*Vdipole[0] + lambda*Vdipole[1];
372 enercorr = *vcharge + *vdip - (L1*VselfA + lambda*VselfB) - Vexcl;
373 *dvdlambda += Vdipole[1] + Vcharge[1] - VselfB
374 - (Vdipole[0] + Vcharge[0] - VselfA) - dvdl_excl;
378 fprintf(debug,"Long Range corrections for Ewald interactions:\n");
379 fprintf(debug,"start=%d,natoms=%d\n",start,end-start);
380 fprintf(debug,"q2sum = %g, Vself=%g\n",
381 L1*q2sumA+lambda*q2sumB,L1*VselfA+lambda*VselfB);
382 fprintf(debug,"Long Range correction: Vexcl=%g\n",Vexcl);
384 fprintf(debug,"Total charge correction: Vcharge=%g\n",
385 L1*Vcharge[0]+lambda*Vcharge[1]);
386 if (epsilon_surface > 0 || ewald_geometry == eewg3DC) {
387 fprintf(debug,"Total dipole correction: Vdipole=%g\n",
388 L1*Vdipole[0]+lambda*Vdipole[1]);
393 /* Return the correction to the energy */