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,int thread,t_forcerec *fr,
82 real *chargeA,real *chargeB,
83 gmx_bool calc_excl_corr,
84 t_blocka *excl,rvec x[],
85 matrix box,rvec mu_tot[],
86 int ewald_geometry,real epsilon_surface,
88 real lambda,real *dvdlambda)
90 int i,i1,i2,j,k,m,iv,jv,q;
92 double q2sumA,q2sumB,Vexcl,dvdl_excl; /* Necessary for precision */
94 real v,vc,qiA,qiB,dr,dr2,rinv,fscal,enercorr;
95 real Vself[2],Vdipole[2],rinv2,ewc=fr->ewaldcoeff,ewcdr;
96 rvec df,dx,mutot[2],dipcorrA,dipcorrB;
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 gmx_bool bFreeEnergy = (chargeB != NULL);
108 gmx_bool bMolPBC = fr->bMolPBC;
110 one_4pi_eps = ONE_4PI_EPS0/fr->epsilon_r;
111 vr0 = ewc*M_2_SQRTPI;
122 /* Note that we have to transform back to gromacs units, since
123 * mu_tot contains the dipole in debye units (for output).
125 for(i=0; (i<DIM); i++) {
126 mutot[0][i] = mu_tot[0][i]*DEBYE2ENM;
127 mutot[1][i] = mu_tot[1][i]*DEBYE2ENM;
132 switch (ewald_geometry) {
134 if (epsilon_surface != 0) {
136 2*M_PI*ONE_4PI_EPS0/((2*epsilon_surface + fr->epsilon_r)*vol);
137 for(i=0; (i<DIM); i++) {
138 dipcorrA[i] = 2*dipole_coeff*mutot[0][i];
139 dipcorrB[i] = 2*dipole_coeff*mutot[1][i];
144 dipole_coeff = 2*M_PI*one_4pi_eps/vol;
145 dipcorrA[ZZ] = 2*dipole_coeff*mutot[0][ZZ];
146 dipcorrB[ZZ] = 2*dipole_coeff*mutot[1][ZZ];
149 gmx_incons("Unsupported Ewald geometry");
153 fprintf(debug,"dipcorr = %8.3f %8.3f %8.3f\n",
154 dipcorrA[XX],dipcorrA[YY],dipcorrA[ZZ]);
155 fprintf(debug,"mutot = %8.3f %8.3f %8.3f\n",
156 mutot[0][XX],mutot[0][YY],mutot[0][ZZ]);
160 if ((calc_excl_corr || dipole_coeff != 0) && !bFreeEnergy) {
161 for(i=start; (i<end); i++) {
162 /* Initiate local variables (for this i-particle) to 0 */
163 qiA = chargeA[i]*one_4pi_eps;
168 i2 = excl->index[i+1];
170 /* Loop over excluded neighbours */
171 for(j=i1; (j<i2); j++) {
174 * First we must test whether k <> i, and then, because the
175 * exclusions are all listed twice i->k and k->i we must select
176 * just one of the two.
177 * As a minor optimization we only compute forces when the charges
181 qqA = qiA*chargeA[k];
183 rvec_sub(x[i],x[k],dx);
185 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
186 for(m=DIM-1; (m>=0); m--) {
187 if (dx[m] > 0.5*box[m][m])
189 else if (dx[m] < -0.5*box[m][m])
194 /* Distance between two excluded particles may be zero in the
198 rinv = gmx_invsqrt(dr2);
211 Geps = eps*VFtab[nnn+2];
212 Heps2 = eps2*VFtab[nnn+3];
215 FF = Fp+Geps+2.0*Heps2;
220 fscal = vc*rinv2+fijC*tabscale*rinv;
221 /* End of tabulated interaction part */
224 /* This is the code you would want instead if not using
228 vc = qqA*gmx_erf(ewcdr)*rinv;
231 /* Relative accuracy at R_ERF_R_INACC of 3e-10 */
232 #define R_ERF_R_INACC 0.006
234 /* Relative accuracy at R_ERF_R_INACC of 2e-5 */
235 #define R_ERF_R_INACC 0.1
237 if (ewcdr > R_ERF_R_INACC) {
238 fscal = rinv2*(vc - qqA*ewc*M_2_SQRTPI*exp(-ewcdr*ewcdr));
240 /* Use a fourth order series expansion for small ewcdr */
241 fscal = ewc*ewc*qqA*vr0*(2.0/3.0 - 0.4*ewcdr*ewcdr);
244 /* The force vector is obtained by multiplication with the
250 for(iv=0; (iv<DIM); iv++)
251 for(jv=0; (jv<DIM); jv++)
252 dxdf[iv][jv] += dx[iv]*df[jv];
260 /* Dipole correction on force */
261 if (dipole_coeff != 0) {
262 for(j=0; (j<DIM); j++)
263 f[i][j] -= dipcorrA[j]*chargeA[i];
266 } else if (calc_excl_corr || dipole_coeff != 0) {
267 for(i=start; (i<end); i++) {
268 /* Initiate local variables (for this i-particle) to 0 */
269 qiA = chargeA[i]*one_4pi_eps;
270 qiB = chargeB[i]*one_4pi_eps;
275 i2 = excl->index[i+1];
277 /* Loop over excluded neighbours */
278 for(j=i1; (j<i2); j++) {
281 qqA = qiA*chargeA[k];
282 qqB = qiB*chargeB[k];
283 if (qqA != 0.0 || qqB != 0.0) {
284 qqL = L1*qqA + lambda*qqB;
285 rvec_sub(x[i],x[k],dx);
287 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
288 for(m=DIM-1; (m>=0); m--) {
289 if (dx[m] > 0.5*box[m][m])
291 else if (dx[m] < -0.5*box[m][m])
297 rinv = gmx_invsqrt(dr2);
300 v = gmx_erf(ewc*dr)*rinv;
303 fscal = rinv2*(vc-qqL*ewc*M_2_SQRTPI*exp(-ewc*ewc*dr2));
307 for(iv=0; (iv<DIM); iv++)
308 for(jv=0; (jv<DIM); jv++)
309 dxdf[iv][jv] += dx[iv]*df[jv];
310 dvdl_excl += (qqB - qqA)*v;
313 dvdl_excl += (qqB - qqA)*vr0;
319 /* Dipole correction on force */
320 if (dipole_coeff != 0) {
321 for(j=0; (j<DIM); j++)
322 f[i][j] -= L1*dipcorrA[j]*chargeA[i]
323 + lambda*dipcorrB[j]*chargeB[i];
327 for(iv=0; (iv<DIM); iv++)
328 for(jv=0; (jv<DIM); jv++)
329 vir[iv][jv] += 0.5*dxdf[iv][jv];
334 /* Global corrections only on master process */
335 if (MASTER(cr) && thread == 0) {
336 for(q=0; q<(bFreeEnergy ? 2 : 1); q++) {
337 if (calc_excl_corr) {
338 /* Self-energy correction */
339 Vself[q] = ewc*one_4pi_eps*fr->q2sum[q]*M_1_SQRTPI;
342 /* Apply surface dipole correction:
343 * correction = dipole_coeff * (dipole)^2
345 if (dipole_coeff != 0) {
346 if (ewald_geometry == eewg3D)
347 Vdipole[q] = dipole_coeff*iprod(mutot[q],mutot[q]);
348 else if (ewald_geometry == eewg3DC)
349 Vdipole[q] = dipole_coeff*mutot[q][ZZ]*mutot[q][ZZ];
355 enercorr = Vdipole[0] - Vself[0] - Vexcl;
357 enercorr = L1*(Vdipole[0] - Vself[0])
358 + lambda*(Vdipole[1] - Vself[1])
360 *dvdlambda += Vdipole[1] - Vself[1]
361 - (Vdipole[0] - Vself[0]) - dvdl_excl;
365 fprintf(debug,"Long Range corrections for Ewald interactions:\n");
366 fprintf(debug,"start=%d,natoms=%d\n",start,end-start);
367 fprintf(debug,"q2sum = %g, Vself=%g\n",
368 L1*q2sumA+lambda*q2sumB,L1*Vself[0]+lambda*Vself[1]);
369 fprintf(debug,"Long Range correction: Vexcl=%g\n",Vexcl);
370 if (MASTER(cr) && thread == 0) {
371 if (epsilon_surface > 0 || ewald_geometry == eewg3DC) {
372 fprintf(debug,"Total dipole correction: Vdipole=%g\n",
373 L1*Vdipole[0]+lambda*Vdipole[1]);
378 /* Return the correction to the energy */
382 real ewald_charge_correction(t_commrec *cr,t_forcerec *fr,real lambda,
384 real *dvdlambda,tensor vir)
387 real vol,fac,qs2A,qs2B,vc,enercorr;
392 /* Apply charge correction */
393 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
395 fac = M_PI*ONE_4PI_EPS0/(fr->epsilon_r*2.0*vol*vol*sqr(fr->ewaldcoeff));
397 qs2A = fr->qsum[0]*fr->qsum[0];
398 qs2B = fr->qsum[1]*fr->qsum[1];
400 vc = (qs2A*(1 - lambda) + qs2B*lambda)*fac;
404 *dvdlambda += -vol*(qs2B - qs2A)*fac;
413 fprintf(debug,"Total charge correction: Vcharge=%g\n",enercorr);