Version bumps after new release
[alexxy/gromacs.git] / src / gromacs / gmxlib / ewald_util.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <stdio.h>
42 #include <math.h>
43 #include "gromacs/math/utilities.h"
44 #include "typedefs.h"
45 #include "types/commrec.h"
46 #include "vec.h"
47 #include "coulomb.h"
48 #include "gromacs/utility/smalloc.h"
49 #include "physics.h"
50 #include "txtdump.h"
51 #include "gromacs/fileio/futil.h"
52 #include "names.h"
53 #include "macros.h"
54
55 real calc_ewaldcoeff_q(real rc, real dtol)
56 {
57     real x = 5, low, high;
58     int  n, i = 0;
59
60
61     do
62     {
63         i++;
64         x *= 2;
65     }
66     while (gmx_erfc(x*rc) > dtol);
67
68     n    = i+60; /* search tolerance is 2^-60 */
69     low  = 0;
70     high = x;
71     for (i = 0; i < n; i++)
72     {
73         x = (low+high)/2;
74         if (gmx_erfc(x*rc) > dtol)
75         {
76             low = x;
77         }
78         else
79         {
80             high = x;
81         }
82     }
83     return x;
84 }
85
86 static real ewald_function_lj(real x, real rc)
87 {
88     real xrc, xrc2, xrc4, factor;
89     xrc  = x*rc;
90     xrc2 = xrc*xrc;
91     xrc4 = xrc2*xrc2;
92 #ifdef GMX_DOUBLE
93     factor = exp(-xrc2)*(1 + xrc2 + xrc4/2.0);
94 #else
95     factor = expf(-xrc2)*(1 + xrc2 + xrc4/2.0);
96 #endif
97
98     return factor;
99 }
100
101 real calc_ewaldcoeff_lj(real rc, real dtol)
102 {
103     real x = 5, low, high;
104     int  n, i = 0;
105
106     do
107     {
108         i++;
109         x *= 2.0;
110     }
111     while (ewald_function_lj(x, rc) > dtol);
112
113     n    = i + 60; /* search tolerance is 2^-60 */
114     low  = 0;
115     high = x;
116     for (i = 0; i < n; ++i)
117     {
118         x = (low + high) / 2.0;
119         if (ewald_function_lj(x, rc) > dtol)
120         {
121             low = x;
122         }
123         else
124         {
125             high = x;
126         }
127     }
128     return x;
129 }
130
131 /* There's nothing special to do here if just masses are perturbed,
132  * but if either charge or type is perturbed then the implementation
133  * requires that B states are defined for both charge and type, and
134  * does not optimize for the cases where only one changes.
135  *
136  * The parameter vectors for B states are left undefined in atoms2md()
137  * when either FEP is inactive, or when there are no mass/charge/type
138  * perturbations. The parameter vectors for LJ-PME are likewise
139  * undefined when LJ-PME is not active. This works because
140  * bHaveChargeOrTypePerturbed handles the control flow. */
141 void ewald_LRcorrection(int start, int end,
142                         t_commrec *cr, int thread, t_forcerec *fr,
143                         real *chargeA, real *chargeB,
144                         real *C6A, real *C6B,
145                         real *sigmaA, real *sigmaB,
146                         real *sigma3A, real *sigma3B,
147                         gmx_bool bHaveChargeOrTypePerturbed,
148                         gmx_bool calc_excl_corr,
149                         t_blocka *excl, rvec x[],
150                         matrix box, rvec mu_tot[],
151                         int ewald_geometry, real epsilon_surface,
152                         rvec *f, tensor vir_q, tensor vir_lj,
153                         real *Vcorr_q, real *Vcorr_lj,
154                         real lambda_q, real lambda_lj,
155                         real *dvdlambda_q, real *dvdlambda_lj)
156 {
157     int         i, i1, i2, j, k, m, iv, jv, q;
158     atom_id    *AA;
159     double      Vexcl_q, dvdl_excl_q, dvdl_excl_lj; /* Necessary for precision */
160     double      Vexcl_lj;
161     real        one_4pi_eps;
162     real        v, vc, qiA, qiB, dr2, rinv, enercorr;
163     real        Vself_q[2], Vself_lj[2], Vdipole[2], rinv2, ewc_q = fr->ewaldcoeff_q, ewcdr;
164     real        ewc_lj = fr->ewaldcoeff_lj, ewc_lj2 = ewc_lj * ewc_lj;
165     real        c6Ai   = 0, c6Bi = 0, c6A = 0, c6B = 0, ewcdr2, ewcdr4, c6L = 0, rinv6;
166     rvec        df, dx, mutot[2], dipcorrA, dipcorrB;
167     tensor      dxdf_q, dxdf_lj;
168     real        vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
169     real        L1_q, L1_lj, dipole_coeff, qqA, qqB, qqL, vr0_q, vr0_lj = 0;
170     gmx_bool    bMolPBC      = fr->bMolPBC;
171     gmx_bool    bDoingLBRule = (fr->ljpme_combination_rule == eljpmeLB);
172
173     /* This routine can be made faster by using tables instead of analytical interactions
174      * However, that requires a thorough verification that they are correct in all cases.
175      */
176
177     one_4pi_eps   = ONE_4PI_EPS0/fr->epsilon_r;
178     vr0_q         = ewc_q*M_2_SQRTPI;
179     if (EVDW_PME(fr->vdwtype))
180     {
181         vr0_lj    = -pow(ewc_lj, 6)/6.0;
182     }
183
184     AA           = excl->a;
185     Vexcl_q      = 0;
186     Vexcl_lj     = 0;
187     dvdl_excl_q  = 0;
188     dvdl_excl_lj = 0;
189     Vdipole[0]   = 0;
190     Vdipole[1]   = 0;
191     L1_q         = 1.0-lambda_q;
192     L1_lj        = 1.0-lambda_lj;
193     /* Note that we have to transform back to gromacs units, since
194      * mu_tot contains the dipole in debye units (for output).
195      */
196     for (i = 0; (i < DIM); i++)
197     {
198         mutot[0][i] = mu_tot[0][i]*DEBYE2ENM;
199         mutot[1][i] = mu_tot[1][i]*DEBYE2ENM;
200         dipcorrA[i] = 0;
201         dipcorrB[i] = 0;
202     }
203     dipole_coeff = 0;
204     switch (ewald_geometry)
205     {
206         case eewg3D:
207             if (epsilon_surface != 0)
208             {
209                 dipole_coeff =
210                     2*M_PI*ONE_4PI_EPS0/((2*epsilon_surface + fr->epsilon_r)*vol);
211                 for (i = 0; (i < DIM); i++)
212                 {
213                     dipcorrA[i] = 2*dipole_coeff*mutot[0][i];
214                     dipcorrB[i] = 2*dipole_coeff*mutot[1][i];
215                 }
216             }
217             break;
218         case eewg3DC:
219             dipole_coeff = 2*M_PI*one_4pi_eps/vol;
220             dipcorrA[ZZ] = 2*dipole_coeff*mutot[0][ZZ];
221             dipcorrB[ZZ] = 2*dipole_coeff*mutot[1][ZZ];
222             break;
223         default:
224             gmx_incons("Unsupported Ewald geometry");
225             break;
226     }
227     if (debug)
228     {
229         fprintf(debug, "dipcorr = %8.3f  %8.3f  %8.3f\n",
230                 dipcorrA[XX], dipcorrA[YY], dipcorrA[ZZ]);
231         fprintf(debug, "mutot   = %8.3f  %8.3f  %8.3f\n",
232                 mutot[0][XX], mutot[0][YY], mutot[0][ZZ]);
233     }
234     clear_mat(dxdf_q);
235     if (EVDW_PME(fr->vdwtype))
236     {
237         clear_mat(dxdf_lj);
238     }
239     if ((calc_excl_corr || dipole_coeff != 0) && !bHaveChargeOrTypePerturbed)
240     {
241         for (i = start; (i < end); i++)
242         {
243             /* Initiate local variables (for this i-particle) to 0 */
244             qiA = chargeA[i]*one_4pi_eps;
245             if (EVDW_PME(fr->vdwtype))
246             {
247                 c6Ai = C6A[i];
248                 if (bDoingLBRule)
249                 {
250                     c6Ai *= sigma3A[i];
251                 }
252             }
253             if (calc_excl_corr)
254             {
255                 i1  = excl->index[i];
256                 i2  = excl->index[i+1];
257
258                 /* Loop over excluded neighbours */
259                 for (j = i1; (j < i2); j++)
260                 {
261                     k = AA[j];
262                     /*
263                      * First we must test whether k <> i, and then,
264                      * because the exclusions are all listed twice i->k
265                      * and k->i we must select just one of the two.  As
266                      * a minor optimization we only compute forces when
267                      * the charges are non-zero.
268                      */
269                     if (k > i)
270                     {
271                         qqA = qiA*chargeA[k];
272                         if (EVDW_PME(fr->vdwtype))
273                         {
274                             c6A  = c6Ai * C6A[k];
275                             if (bDoingLBRule)
276                             {
277                                 c6A *= pow(0.5*(sigmaA[i]+sigmaA[k]), 6)*sigma3A[k];
278                             }
279                         }
280                         if (qqA != 0.0 || c6A != 0.0)
281                         {
282                             real fscal;
283
284                             fscal = 0;
285                             rvec_sub(x[i], x[k], dx);
286                             if (bMolPBC)
287                             {
288                                 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
289                                 for (m = DIM-1; (m >= 0); m--)
290                                 {
291                                     if (dx[m] > 0.5*box[m][m])
292                                     {
293                                         rvec_dec(dx, box[m]);
294                                     }
295                                     else if (dx[m] < -0.5*box[m][m])
296                                     {
297                                         rvec_inc(dx, box[m]);
298                                     }
299                                 }
300                             }
301                             dr2 = norm2(dx);
302                             /* Distance between two excluded particles
303                              * may be zero in the case of shells
304                              */
305                             if (dr2 != 0)
306                             {
307                                 rinv              = gmx_invsqrt(dr2);
308                                 rinv2             = rinv*rinv;
309                                 if (qqA != 0.0)
310                                 {
311                                     real dr;
312
313                                     dr       = 1.0/rinv;
314                                     ewcdr    = ewc_q*dr;
315                                     vc       = qqA*gmx_erf(ewcdr)*rinv;
316                                     Vexcl_q += vc;
317 #ifdef GMX_DOUBLE
318                                     /* Relative accuracy at R_ERF_R_INACC of 3e-10 */
319 #define       R_ERF_R_INACC 0.006
320 #else
321                                     /* Relative accuracy at R_ERF_R_INACC of 2e-5 */
322 #define       R_ERF_R_INACC 0.1
323 #endif
324                                     /* fscal is the scalar force pre-multiplied by rinv,
325                                      * to normalise the relative position vector dx */
326                                     if (ewcdr > R_ERF_R_INACC)
327                                     {
328                                         fscal = rinv2*(vc - qqA*ewc_q*M_2_SQRTPI*exp(-ewcdr*ewcdr));
329                                     }
330                                     else
331                                     {
332                                         /* Use a fourth order series expansion for small ewcdr */
333                                         fscal = ewc_q*ewc_q*qqA*vr0_q*(2.0/3.0 - 0.4*ewcdr*ewcdr);
334                                     }
335
336                                     /* The force vector is obtained by multiplication with
337                                      * the relative position vector
338                                      */
339                                     svmul(fscal, dx, df);
340                                     rvec_inc(f[k], df);
341                                     rvec_dec(f[i], df);
342                                     for (iv = 0; (iv < DIM); iv++)
343                                     {
344                                         for (jv = 0; (jv < DIM); jv++)
345                                         {
346                                             dxdf_q[iv][jv] += dx[iv]*df[jv];
347                                         }
348                                     }
349                                 }
350
351                                 if (c6A != 0.0)
352                                 {
353                                     rinv6     = rinv2*rinv2*rinv2;
354                                     ewcdr2    = ewc_lj2*dr2;
355                                     ewcdr4    = ewcdr2*ewcdr2;
356                                     /* We get the excluded long-range contribution from -C6*(1-g(r))
357                                      * g(r) is also defined in the manual under LJ-PME
358                                      */
359                                     vc        = -c6A*rinv6*(1.0 - exp(-ewcdr2)*(1 + ewcdr2 + 0.5*ewcdr4));
360                                     Vexcl_lj += vc;
361                                     /* The force is the derivative of the potential vc.
362                                      * fscal is the scalar force pre-multiplied by rinv,
363                                      * to normalise the relative position vector dx */
364                                     fscal     = 6.0*vc*rinv2 + c6A*rinv6*exp(-ewcdr2)*ewc_lj2*ewcdr4;
365
366                                     /* The force vector is obtained by multiplication with
367                                      * the relative position vector
368                                      */
369                                     svmul(fscal, dx, df);
370                                     rvec_inc(f[k], df);
371                                     rvec_dec(f[i], df);
372                                     for (iv = 0; (iv < DIM); iv++)
373                                     {
374                                         for (jv = 0; (jv < DIM); jv++)
375                                         {
376                                             dxdf_lj[iv][jv] += dx[iv]*df[jv];
377                                         }
378                                     }
379                                 }
380                             }
381                             else
382                             {
383                                 Vexcl_q  += qqA*vr0_q;
384                                 Vexcl_lj += c6A*vr0_lj;
385                             }
386                         }
387                     }
388                 }
389             }
390             /* Dipole correction on force */
391             if (dipole_coeff != 0)
392             {
393                 for (j = 0; (j < DIM); j++)
394                 {
395                     f[i][j] -= dipcorrA[j]*chargeA[i];
396                 }
397             }
398         }
399     }
400     else if (calc_excl_corr || dipole_coeff != 0)
401     {
402         for (i = start; (i < end); i++)
403         {
404             /* Initiate local variables (for this i-particle) to 0 */
405             qiA = chargeA[i]*one_4pi_eps;
406             qiB = chargeB[i]*one_4pi_eps;
407             if (EVDW_PME(fr->vdwtype))
408             {
409                 c6Ai = C6A[i];
410                 c6Bi = C6B[i];
411                 if (bDoingLBRule)
412                 {
413                     c6Ai *= sigma3A[i];
414                     c6Bi *= sigma3B[i];
415                 }
416             }
417             if (calc_excl_corr)
418             {
419                 i1  = excl->index[i];
420                 i2  = excl->index[i+1];
421
422                 /* Loop over excluded neighbours */
423                 for (j = i1; (j < i2); j++)
424                 {
425                     k = AA[j];
426                     if (k > i)
427                     {
428                         qqA = qiA*chargeA[k];
429                         qqB = qiB*chargeB[k];
430                         if (EVDW_PME(fr->vdwtype))
431                         {
432                             c6A = c6Ai*C6A[k];
433                             c6B = c6Bi*C6B[k];
434                             if (bDoingLBRule)
435                             {
436                                 c6A *= pow(0.5*(sigmaA[i]+sigmaA[k]), 6)*sigma3A[k];
437                                 c6B *= pow(0.5*(sigmaB[i]+sigmaB[k]), 6)*sigma3B[k];
438                             }
439                         }
440                         if (qqA != 0.0 || qqB != 0.0 || c6A != 0.0 || c6B != 0.0)
441                         {
442                             real fscal;
443
444                             fscal = 0;
445                             qqL   = L1_q*qqA + lambda_q*qqB;
446                             if (EVDW_PME(fr->vdwtype))
447                             {
448                                 c6L = L1_lj*c6A + lambda_lj*c6B;
449                             }
450                             rvec_sub(x[i], x[k], dx);
451                             if (bMolPBC)
452                             {
453                                 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
454                                 for (m = DIM-1; (m >= 0); m--)
455                                 {
456                                     if (dx[m] > 0.5*box[m][m])
457                                     {
458                                         rvec_dec(dx, box[m]);
459                                     }
460                                     else if (dx[m] < -0.5*box[m][m])
461                                     {
462                                         rvec_inc(dx, box[m]);
463                                     }
464                                 }
465                             }
466                             dr2 = norm2(dx);
467                             if (dr2 != 0)
468                             {
469                                 rinv    = gmx_invsqrt(dr2);
470                                 rinv2   = rinv*rinv;
471                                 if (qqA != 0.0 || qqB != 0.0)
472                                 {
473                                     real dr;
474
475                                     dr           = 1.0/rinv;
476                                     v            = gmx_erf(ewc_q*dr)*rinv;
477                                     vc           = qqL*v;
478                                     Vexcl_q     += vc;
479                                     /* fscal is the scalar force pre-multiplied by rinv,
480                                      * to normalise the relative position vector dx */
481                                     fscal        = rinv2*(vc-qqL*ewc_q*M_2_SQRTPI*exp(-ewc_q*ewc_q*dr2));
482                                     dvdl_excl_q += (qqB - qqA)*v;
483
484                                     /* The force vector is obtained by multiplication with
485                                      * the relative position vector
486                                      */
487                                     svmul(fscal, dx, df);
488                                     rvec_inc(f[k], df);
489                                     rvec_dec(f[i], df);
490                                     for (iv = 0; (iv < DIM); iv++)
491                                     {
492                                         for (jv = 0; (jv < DIM); jv++)
493                                         {
494                                             dxdf_q[iv][jv] += dx[iv]*df[jv];
495                                         }
496                                     }
497                                 }
498
499                                 if ((c6A != 0.0 || c6B != 0.0) && EVDW_PME(fr->vdwtype))
500                                 {
501                                     rinv6         = rinv2*rinv2*rinv2;
502                                     ewcdr2        = ewc_lj2*dr2;
503                                     ewcdr4        = ewcdr2*ewcdr2;
504                                     v             = -rinv6*(1.0 - exp(-ewcdr2)*(1 + ewcdr2 + 0.5*ewcdr4));
505                                     vc            = c6L*v;
506                                     Vexcl_lj     += vc;
507                                     /* fscal is the scalar force pre-multiplied by rinv,
508                                      * to normalise the relative position vector dx */
509                                     fscal         = 6.0*vc*rinv2 + c6L*rinv6*exp(-ewcdr2)*ewc_lj2*ewcdr4;
510                                     dvdl_excl_lj += (c6B - c6A)*v;
511
512                                     /* The force vector is obtained by multiplication with
513                                      * the relative position vector
514                                      */
515                                     svmul(fscal, dx, df);
516                                     rvec_inc(f[k], df);
517                                     rvec_dec(f[i], df);
518                                     for (iv = 0; (iv < DIM); iv++)
519                                     {
520                                         for (jv = 0; (jv < DIM); jv++)
521                                         {
522                                             dxdf_lj[iv][jv] += dx[iv]*df[jv];
523                                         }
524                                     }
525                                 }
526                             }
527                             else
528                             {
529                                 Vexcl_q      += qqL*vr0_q;
530                                 dvdl_excl_q  += (qqB - qqA)*vr0_q;
531                                 Vexcl_lj     += c6L*vr0_lj;
532                                 dvdl_excl_lj += (c6B - c6A)*vr0_lj;
533                             }
534                         }
535                     }
536                 }
537             }
538             /* Dipole correction on force */
539             if (dipole_coeff != 0)
540             {
541                 for (j = 0; (j < DIM); j++)
542                 {
543                     f[i][j] -= L1_q*dipcorrA[j]*chargeA[i]
544                         + lambda_q*dipcorrB[j]*chargeB[i];
545                 }
546             }
547         }
548     }
549     for (iv = 0; (iv < DIM); iv++)
550     {
551         for (jv = 0; (jv < DIM); jv++)
552         {
553             vir_q[iv][jv]  += 0.5*dxdf_q[iv][jv];
554             vir_lj[iv][jv] += 0.5*dxdf_lj[iv][jv];
555         }
556     }
557
558     Vself_q[0]  = 0;
559     Vself_q[1]  = 0;
560     Vself_lj[0] = 0;
561     Vself_lj[1] = 0;
562
563     /* Global corrections only on master process */
564     if (MASTER(cr) && thread == 0)
565     {
566         for (q = 0; q < (bHaveChargeOrTypePerturbed ? 2 : 1); q++)
567         {
568             if (calc_excl_corr)
569             {
570                 /* Self-energy correction */
571                 Vself_q[q] = ewc_q*one_4pi_eps*fr->q2sum[q]*M_1_SQRTPI;
572                 if (EVDW_PME(fr->vdwtype))
573                 {
574                     Vself_lj[q] =  fr->c6sum[q]*0.5*vr0_lj;
575                 }
576             }
577
578             /* Apply surface dipole correction:
579              * correction = dipole_coeff * (dipole)^2
580              */
581             if (dipole_coeff != 0)
582             {
583                 if (ewald_geometry == eewg3D)
584                 {
585                     Vdipole[q] = dipole_coeff*iprod(mutot[q], mutot[q]);
586                 }
587                 else if (ewald_geometry == eewg3DC)
588                 {
589                     Vdipole[q] = dipole_coeff*mutot[q][ZZ]*mutot[q][ZZ];
590                 }
591             }
592         }
593     }
594     if (!bHaveChargeOrTypePerturbed)
595     {
596         *Vcorr_q = Vdipole[0] - Vself_q[0] - Vexcl_q;
597         if (EVDW_PME(fr->vdwtype))
598         {
599             *Vcorr_lj = -Vself_lj[0] - Vexcl_lj;
600         }
601     }
602     else
603     {
604         *Vcorr_q = L1_q*(Vdipole[0] - Vself_q[0])
605             + lambda_q*(Vdipole[1] - Vself_q[1])
606             - Vexcl_q;
607         *dvdlambda_q += Vdipole[1] - Vself_q[1]
608             - (Vdipole[0] - Vself_q[0]) - dvdl_excl_q;
609         if (EVDW_PME(fr->vdwtype))
610         {
611             *Vcorr_lj      = -(L1_lj*Vself_lj[0] + lambda_lj*Vself_lj[1]) - Vexcl_lj;
612             *dvdlambda_lj += -Vself_lj[1] + Vself_lj[0] - dvdl_excl_lj;
613         }
614     }
615
616     if (debug)
617     {
618         fprintf(debug, "Long Range corrections for Ewald interactions:\n");
619         fprintf(debug, "start=%d,natoms=%d\n", start, end-start);
620         fprintf(debug, "q2sum = %g, Vself_q=%g c6sum = %g, Vself_lj=%g\n",
621                 L1_q*fr->q2sum[0]+lambda_q*fr->q2sum[1], L1_q*Vself_q[0]+lambda_q*Vself_q[1], L1_lj*fr->c6sum[0]+lambda_lj*fr->c6sum[1], L1_lj*Vself_lj[0]+lambda_lj*Vself_lj[1]);
622         fprintf(debug, "Electrostatic Long Range correction: Vexcl=%g\n", Vexcl_q);
623         fprintf(debug, "Lennard-Jones Long Range correction: Vexcl=%g\n", Vexcl_lj);
624         if (MASTER(cr) && thread == 0)
625         {
626             if (epsilon_surface > 0 || ewald_geometry == eewg3DC)
627             {
628                 fprintf(debug, "Total dipole correction: Vdipole=%g\n",
629                         L1_q*Vdipole[0]+lambda_q*Vdipole[1]);
630             }
631         }
632     }
633 }
634
635 real ewald_charge_correction(t_commrec *cr, t_forcerec *fr, real lambda,
636                              matrix box,
637                              real *dvdlambda, tensor vir)
638
639 {
640     real vol, fac, qs2A, qs2B, vc, enercorr;
641     int  d;
642
643     if (MASTER(cr))
644     {
645         /* Apply charge correction */
646         vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
647
648         fac = M_PI*ONE_4PI_EPS0/(fr->epsilon_r*2.0*vol*vol*sqr(fr->ewaldcoeff_q));
649
650         qs2A = fr->qsum[0]*fr->qsum[0];
651         qs2B = fr->qsum[1]*fr->qsum[1];
652
653         vc = (qs2A*(1 - lambda) + qs2B*lambda)*fac;
654
655         enercorr = -vol*vc;
656
657         *dvdlambda += -vol*(qs2B - qs2A)*fac;
658
659         for (d = 0; d < DIM; d++)
660         {
661             vir[d][d] += vc;
662         }
663
664         if (debug)
665         {
666             fprintf(debug, "Total charge correction: Vcharge=%g\n", enercorr);
667         }
668     }
669     else
670     {
671         enercorr = 0;
672     }
673
674     return enercorr;
675 }