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 double isp=0.564189583547756;
109 gmx_bool bFreeEnergy = (chargeB != NULL);
110 gmx_bool bMolPBC = fr->bMolPBC;
112 one_4pi_eps = ONE_4PI_EPS0/fr->epsilon_r;
113 vr0 = ewc*2/sqrt(M_PI);
124 /* Note that we have to transform back to gromacs units, since
125 * mu_tot contains the dipole in debye units (for output).
127 for(i=0; (i<DIM); i++) {
128 mutot[0][i] = mu_tot[0][i]*DEBYE2ENM;
129 mutot[1][i] = mu_tot[1][i]*DEBYE2ENM;
134 switch (ewald_geometry) {
136 if (epsilon_surface != 0) {
138 2*M_PI*ONE_4PI_EPS0/((2*epsilon_surface + fr->epsilon_r)*vol);
139 for(i=0; (i<DIM); i++) {
140 dipcorrA[i] = 2*dipole_coeff*mutot[0][i];
141 dipcorrB[i] = 2*dipole_coeff*mutot[1][i];
146 dipole_coeff = 2*M_PI*one_4pi_eps/vol;
147 dipcorrA[ZZ] = 2*dipole_coeff*mutot[0][ZZ];
148 dipcorrB[ZZ] = 2*dipole_coeff*mutot[1][ZZ];
151 gmx_incons("Unsupported Ewald geometry");
155 fprintf(debug,"dipcorr = %8.3f %8.3f %8.3f\n",
156 dipcorrA[XX],dipcorrA[YY],dipcorrA[ZZ]);
157 fprintf(debug,"mutot = %8.3f %8.3f %8.3f\n",
158 mutot[0][XX],mutot[0][YY],mutot[0][ZZ]);
162 if ((calc_excl_corr || dipole_coeff != 0) && !bFreeEnergy) {
163 for(i=start; (i<end); i++) {
164 /* Initiate local variables (for this i-particle) to 0 */
165 qiA = chargeA[i]*one_4pi_eps;
170 i2 = excl->index[i+1];
172 /* Loop over excluded neighbours */
173 for(j=i1; (j<i2); j++) {
176 * First we must test whether k <> i, and then, because the
177 * exclusions are all listed twice i->k and k->i we must select
178 * just one of the two.
179 * As a minor optimization we only compute forces when the charges
183 qqA = qiA*chargeA[k];
185 rvec_sub(x[i],x[k],dx);
187 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
188 for(m=DIM-1; (m>=0); m--) {
189 if (dx[m] > 0.5*box[m][m])
191 else if (dx[m] < -0.5*box[m][m])
196 /* Distance between two excluded particles may be zero in the
200 rinv = gmx_invsqrt(dr2);
213 Geps = eps*VFtab[nnn+2];
214 Heps2 = eps2*VFtab[nnn+3];
217 FF = Fp+Geps+2.0*Heps2;
222 fscal = vc*rinv2+fijC*tabscale*rinv;
223 /* End of tabulated interaction part */
226 /* This is the code you would want instead if not using
230 vc = qqA*gmx_erf(ewcdr)*rinv;
233 /* Relative accuracy at R_ERF_R_INACC of 3e-10 */
234 #define R_ERF_R_INACC 0.006
236 /* Relative accuracy at R_ERF_R_INACC of 2e-5 */
237 #define R_ERF_R_INACC 0.1
239 if (ewcdr > R_ERF_R_INACC) {
240 fscal = rinv2*(vc - 2.0*qqA*ewc*isp*exp(-ewcdr*ewcdr));
242 /* Use a fourth order series expansion for small ewcdr */
243 fscal = ewc*ewc*qqA*vr0*(2.0/3.0 - 0.4*ewcdr*ewcdr);
246 /* The force vector is obtained by multiplication with the
252 for(iv=0; (iv<DIM); iv++)
253 for(jv=0; (jv<DIM); jv++)
254 dxdf[iv][jv] += dx[iv]*df[jv];
262 /* Dipole correction on force */
263 if (dipole_coeff != 0) {
264 for(j=0; (j<DIM); j++)
265 f[i][j] -= dipcorrA[j]*chargeA[i];
268 } else if (calc_excl_corr || dipole_coeff != 0) {
269 for(i=start; (i<end); i++) {
270 /* Initiate local variables (for this i-particle) to 0 */
271 qiA = chargeA[i]*one_4pi_eps;
272 qiB = chargeB[i]*one_4pi_eps;
277 i2 = excl->index[i+1];
279 /* Loop over excluded neighbours */
280 for(j=i1; (j<i2); j++) {
283 qqA = qiA*chargeA[k];
284 qqB = qiB*chargeB[k];
285 if (qqA != 0.0 || qqB != 0.0) {
286 qqL = L1*qqA + lambda*qqB;
287 rvec_sub(x[i],x[k],dx);
289 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
290 for(m=DIM-1; (m>=0); m--) {
291 if (dx[m] > 0.5*box[m][m])
293 else if (dx[m] < -0.5*box[m][m])
299 rinv = gmx_invsqrt(dr2);
302 v = gmx_erf(ewc*dr)*rinv;
305 fscal = rinv2*(vc-2.0*qqL*ewc*isp*exp(-ewc*ewc*dr2));
309 for(iv=0; (iv<DIM); iv++)
310 for(jv=0; (jv<DIM); jv++)
311 dxdf[iv][jv] += dx[iv]*df[jv];
312 dvdl_excl += (qqB - qqA)*v;
315 dvdl_excl += (qqB - qqA)*vr0;
321 /* Dipole correction on force */
322 if (dipole_coeff != 0) {
323 for(j=0; (j<DIM); j++)
324 f[i][j] -= L1*dipcorrA[j]*chargeA[i]
325 + lambda*dipcorrB[j]*chargeB[i];
329 for(iv=0; (iv<DIM); iv++)
330 for(jv=0; (jv<DIM); jv++)
331 vir[iv][jv] += 0.5*dxdf[iv][jv];
336 /* Global corrections only on master process */
337 if (MASTER(cr) && thread == 0) {
338 for(q=0; q<(bFreeEnergy ? 2 : 1); q++) {
339 if (calc_excl_corr) {
340 /* Self-energy correction */
341 Vself[q] = ewc*one_4pi_eps*fr->q2sum[q]/sqrt(M_PI);
344 /* Apply surface dipole correction:
345 * correction = dipole_coeff * (dipole)^2
347 if (dipole_coeff != 0) {
348 if (ewald_geometry == eewg3D)
349 Vdipole[q] = dipole_coeff*iprod(mutot[q],mutot[q]);
350 else if (ewald_geometry == eewg3DC)
351 Vdipole[q] = dipole_coeff*mutot[q][ZZ]*mutot[q][ZZ];
357 enercorr = Vdipole[0] - Vself[0] - Vexcl;
359 enercorr = L1*(Vdipole[0] - Vself[0])
360 + lambda*(Vdipole[1] - Vself[1])
362 *dvdlambda += Vdipole[1] - Vself[1]
363 - (Vdipole[0] - Vself[0]) - dvdl_excl;
367 fprintf(debug,"Long Range corrections for Ewald interactions:\n");
368 fprintf(debug,"start=%d,natoms=%d\n",start,end-start);
369 fprintf(debug,"q2sum = %g, Vself=%g\n",
370 L1*q2sumA+lambda*q2sumB,L1*Vself[0]+lambda*Vself[1]);
371 fprintf(debug,"Long Range correction: Vexcl=%g\n",Vexcl);
372 if (MASTER(cr) && thread == 0) {
373 if (epsilon_surface > 0 || ewald_geometry == eewg3DC) {
374 fprintf(debug,"Total dipole correction: Vdipole=%g\n",
375 L1*Vdipole[0]+lambda*Vdipole[1]);
380 /* Return the correction to the energy */
384 real ewald_charge_correction(t_commrec *cr,t_forcerec *fr,real lambda,
386 real *dvdlambda,tensor vir)
389 real vol,fac,qs2A,qs2B,vc,enercorr;
394 /* Apply charge correction */
395 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
397 fac = M_PI*ONE_4PI_EPS0/(fr->epsilon_r*2.0*vol*vol*sqr(fr->ewaldcoeff));
399 qs2A = fr->qsum[0]*fr->qsum[0];
400 qs2B = fr->qsum[1]*fr->qsum[1];
402 vc = (qs2A*(1 - lambda) + qs2B*lambda)*fac;
406 *dvdlambda += -vol*(qs2B - qs2A)*fac;
415 fprintf(debug,"Total charge correction: Vcharge=%g\n",enercorr);