Bug-fix in long-range correction for LJ-PME
[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 "vec.h"
46 #include "coulomb.h"
47 #include "smalloc.h"
48 #include "physics.h"
49 #include "txtdump.h"
50 #include "gromacs/fileio/futil.h"
51 #include "names.h"
52 #include "macros.h"
53
54 real calc_ewaldcoeff_q(real rc, real dtol)
55 {
56     real x = 5, low, high;
57     int  n, i = 0;
58
59
60     do
61     {
62         i++;
63         x *= 2;
64     }
65     while (gmx_erfc(x*rc) > dtol);
66
67     n    = i+60; /* search tolerance is 2^-60 */
68     low  = 0;
69     high = x;
70     for (i = 0; i < n; i++)
71     {
72         x = (low+high)/2;
73         if (gmx_erfc(x*rc) > dtol)
74         {
75             low = x;
76         }
77         else
78         {
79             high = x;
80         }
81     }
82     return x;
83 }
84
85 static real ewald_function_lj(real x, real rc)
86 {
87     real xrc, xrc2, xrc4, factor;
88     xrc  = x*rc;
89     xrc2 = xrc*xrc;
90     xrc4 = xrc2*xrc2;
91 #ifdef GMX_DOUBLE
92     factor = exp(-xrc2)*(1 + xrc2 + xrc4/2.0);
93 #else
94     factor = expf(-xrc2)*(1 + xrc2 + xrc4/2.0);
95 #endif
96
97     return factor;
98 }
99
100 real calc_ewaldcoeff_lj(real rc, real dtol)
101 {
102     real x = 5, low, high;
103     int  n, i = 0;
104
105     do
106     {
107         i++;
108         x *= 2.0;
109     }
110     while (ewald_function_lj(x, rc) > dtol);
111
112     n    = i + 60; /* search tolerance is 2^-60 */
113     low  = 0;
114     high = x;
115     for (i = 0; i < n; ++i)
116     {
117         x = (low + high) / 2.0;
118         if (ewald_function_lj(x, rc) > dtol)
119         {
120             low = x;
121         }
122         else
123         {
124             high = x;
125         }
126     }
127     return x;
128 }
129
130 void ewald_LRcorrection(int start, int end,
131                         t_commrec *cr, int thread, t_forcerec *fr,
132                         real *chargeA, real *chargeB,
133                         real *C6A, real *C6B,
134                         real *sigmaA, real *sigmaB,
135                         real *sigma3A, real *sigma3B,
136                         gmx_bool calc_excl_corr,
137                         t_blocka *excl, rvec x[],
138                         matrix box, rvec mu_tot[],
139                         int ewald_geometry, real epsilon_surface,
140                         rvec *f, tensor vir_q, tensor vir_lj,
141                         real *Vcorr_q, real *Vcorr_lj,
142                         real lambda_q, real lambda_lj,
143                         real *dvdlambda_q, real *dvdlambda_lj)
144 {
145     int         i, i1, i2, j, k, m, iv, jv, q;
146     atom_id    *AA;
147     double      Vexcl_q, dvdl_excl_q, dvdl_excl_lj; /* Necessary for precision */
148     double      Vexcl_lj;
149     real        one_4pi_eps;
150     real        v, vc, qiA, qiB, dr, dr2, rinv, fscal, enercorr;
151     real        Vself_q[2], Vself_lj[2], Vdipole[2], rinv2, ewc_q = fr->ewaldcoeff_q, ewcdr;
152     real        ewc_lj = fr->ewaldcoeff_lj;
153     real        c6Ai   = 0, c6Bi = 0, c6A = 0, c6B = 0, ewcdr2, ewcdr4, ewcdr5, c6L = 0, rinv6;
154     rvec        df, dx, mutot[2], dipcorrA, dipcorrB;
155     tensor      dxdf_q, dxdf_lj;
156     real        vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
157     real        L1_q, L1_lj, dipole_coeff, qqA, qqB, qqL, vr0_q, vr0_lj = 0;
158     gmx_bool    bFreeEnergy  = (chargeB != NULL);
159     gmx_bool    bMolPBC      = fr->bMolPBC;
160     gmx_bool    bDoingLBRule = (fr->ljpme_combination_rule == eljpmeLB);
161
162     /* This routine can be made faster by using tables instead of analytical interactions
163      * However, that requires a thorough verification that they are correct in all cases.
164      */
165
166     one_4pi_eps   = ONE_4PI_EPS0/fr->epsilon_r;
167     vr0_q         = ewc_q*M_2_SQRTPI;
168     if (EVDW_PME(fr->vdwtype))
169     {
170         vr0_lj        = -pow(ewc_lj, 6)/6.0;
171     }
172
173     AA           = excl->a;
174     Vexcl_q      = 0;
175     Vexcl_lj     = 0;
176     dvdl_excl_q  = 0;
177     dvdl_excl_lj = 0;
178     Vdipole[0]   = 0;
179     Vdipole[1]   = 0;
180     L1_q         = 1.0-lambda_q;
181     L1_lj        = 1.0-lambda_lj;
182     /* Note that we have to transform back to gromacs units, since
183      * mu_tot contains the dipole in debye units (for output).
184      */
185     for (i = 0; (i < DIM); i++)
186     {
187         mutot[0][i] = mu_tot[0][i]*DEBYE2ENM;
188         mutot[1][i] = mu_tot[1][i]*DEBYE2ENM;
189         dipcorrA[i] = 0;
190         dipcorrB[i] = 0;
191     }
192     dipole_coeff = 0;
193     switch (ewald_geometry)
194     {
195         case eewg3D:
196             if (epsilon_surface != 0)
197             {
198                 dipole_coeff =
199                     2*M_PI*ONE_4PI_EPS0/((2*epsilon_surface + fr->epsilon_r)*vol);
200                 for (i = 0; (i < DIM); i++)
201                 {
202                     dipcorrA[i] = 2*dipole_coeff*mutot[0][i];
203                     dipcorrB[i] = 2*dipole_coeff*mutot[1][i];
204                 }
205             }
206             break;
207         case eewg3DC:
208             dipole_coeff = 2*M_PI*one_4pi_eps/vol;
209             dipcorrA[ZZ] = 2*dipole_coeff*mutot[0][ZZ];
210             dipcorrB[ZZ] = 2*dipole_coeff*mutot[1][ZZ];
211             break;
212         default:
213             gmx_incons("Unsupported Ewald geometry");
214             break;
215     }
216     if (debug)
217     {
218         fprintf(debug, "dipcorr = %8.3f  %8.3f  %8.3f\n",
219                 dipcorrA[XX], dipcorrA[YY], dipcorrA[ZZ]);
220         fprintf(debug, "mutot   = %8.3f  %8.3f  %8.3f\n",
221                 mutot[0][XX], mutot[0][YY], mutot[0][ZZ]);
222     }
223     clear_mat(dxdf_q);
224     if (EVDW_PME(fr->vdwtype))
225     {
226         clear_mat(dxdf_lj);
227     }
228     if ((calc_excl_corr || dipole_coeff != 0 || EVDW_PME(fr->vdwtype)) && !bFreeEnergy)
229     {
230         for (i = start; (i < end); i++)
231         {
232             /* Initiate local variables (for this i-particle) to 0 */
233             qiA = chargeA[i]*one_4pi_eps;
234             if (EVDW_PME(fr->vdwtype))
235             {
236                 c6Ai = C6A[i];
237                 if (bDoingLBRule)
238                 {
239                     c6Ai *= sigma3A[i];
240                 }
241             }
242             if (calc_excl_corr)
243             {
244                 i1  = excl->index[i];
245                 i2  = excl->index[i+1];
246
247                 /* Loop over excluded neighbours */
248                 for (j = i1; (j < i2); j++)
249                 {
250                     k = AA[j];
251                     /*
252                      * First we must test whether k <> i, and then,
253                      * because the exclusions are all listed twice i->k
254                      * and k->i we must select just one of the two.  As
255                      * a minor optimization we only compute forces when
256                      * the charges are non-zero.
257                      */
258                     if (k > i)
259                     {
260                         qqA = qiA*chargeA[k];
261                         if (EVDW_PME(fr->vdwtype))
262                         {
263                             c6A  = c6Ai * C6A[k];
264                             if (bDoingLBRule)
265                             {
266                                 c6A *= pow(0.5*(sigmaA[i]+sigmaA[k]), 6)*sigma3A[k];
267                             }
268                         }
269                         if (qqA != 0.0 || c6A != 0.0)
270                         {
271                             rvec_sub(x[i], x[k], dx);
272                             fscal = 0;
273                             if (bMolPBC)
274                             {
275                                 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
276                                 for (m = DIM-1; (m >= 0); m--)
277                                 {
278                                     if (dx[m] > 0.5*box[m][m])
279                                     {
280                                         rvec_dec(dx, box[m]);
281                                     }
282                                     else if (dx[m] < -0.5*box[m][m])
283                                     {
284                                         rvec_inc(dx, box[m]);
285                                     }
286                                 }
287                             }
288                             dr2 = norm2(dx);
289                             /* Distance between two excluded particles
290                              * may be zero in the case of shells
291                              */
292                             if (dr2 != 0)
293                             {
294                                 rinv              = gmx_invsqrt(dr2);
295                                 rinv2             = rinv*rinv;
296                                 dr                = 1.0/rinv;
297                                 if (qqA != 0.0)
298                                 {
299                                     ewcdr    = ewc_q*dr;
300                                     vc       = qqA*gmx_erf(ewcdr)*rinv;
301                                     Vexcl_q += vc;
302 #ifdef GMX_DOUBLE
303                                     /* Relative accuracy at R_ERF_R_INACC of 3e-10 */
304 #define       R_ERF_R_INACC 0.006
305 #else
306                                     /* Relative accuracy at R_ERF_R_INACC of 2e-5 */
307 #define       R_ERF_R_INACC 0.1
308 #endif
309                                     if (ewcdr > R_ERF_R_INACC)
310                                     {
311                                         fscal = rinv2*(vc - qqA*ewc_q*M_2_SQRTPI*exp(-ewcdr*ewcdr));
312                                     }
313                                     else
314                                     {
315                                         /* Use a fourth order series expansion for small ewcdr */
316                                         fscal = ewc_q*ewc_q*qqA*vr0_q*(2.0/3.0 - 0.4*ewcdr*ewcdr);
317                                     }
318
319                                     /* The force vector is obtained by multiplication with
320                                      * the distance vector
321                                      */
322                                     svmul(fscal, dx, df);
323                                     rvec_inc(f[k], df);
324                                     rvec_dec(f[i], df);
325                                     for (iv = 0; (iv < DIM); iv++)
326                                     {
327                                         for (jv = 0; (jv < DIM); jv++)
328                                         {
329                                             dxdf_q[iv][jv] += dx[iv]*df[jv];
330                                         }
331                                     }
332                                 }
333
334                                 if (c6A != 0.0)
335                                 {
336                                     rinv6     = rinv2*rinv2*rinv2;
337                                     ewcdr     = ewc_lj*dr;
338                                     ewcdr2    = ewcdr*ewcdr;
339                                     ewcdr4    = ewcdr2*ewcdr2;
340                                     ewcdr5    = ewcdr4*ewcdr;
341                                     /* We get the excluded long-range contribution from -C6*(1-g(r))
342                                      * g(r) is also defined in the manual under LJ-PME
343                                      */
344                                     vc        = -c6A*rinv6*(1.0 - exp(-ewcdr2)*(1 + ewcdr2 + 0.5*ewcdr4));
345                                     Vexcl_lj += vc;
346                                     /* The force is the derivative of the potential vc */
347                                     fscal     = 6.0*vc*rinv + c6A*rinv6*exp(-ewcdr2)*ewc_lj*ewcdr5;
348
349                                     /* The force vector is obtained by multiplication with
350                                      * the distance vector. Here we also divide by r, to account
351                                      * for the length of dx
352                                      */
353                                     svmul(fscal*rinv, dx, df);
354                                     rvec_inc(f[k], df);
355                                     rvec_dec(f[i], df);
356                                     for (iv = 0; (iv < DIM); iv++)
357                                     {
358                                         for (jv = 0; (jv < DIM); jv++)
359                                         {
360                                             dxdf_lj[iv][jv] += dx[iv]*df[jv];
361                                         }
362                                     }
363                                 }
364                             }
365                             else
366                             {
367                                 Vexcl_q  += qqA*vr0_q;
368                                 Vexcl_lj += c6A*vr0_lj;
369                             }
370                         }
371                     }
372                 }
373             }
374             /* Dipole correction on force */
375             if (dipole_coeff != 0)
376             {
377                 for (j = 0; (j < DIM); j++)
378                 {
379                     f[i][j] -= dipcorrA[j]*chargeA[i];
380                 }
381             }
382         }
383     }
384     else if (calc_excl_corr || dipole_coeff != 0 || EVDW_PME(fr->vdwtype))
385     {
386         for (i = start; (i < end); i++)
387         {
388             /* Initiate local variables (for this i-particle) to 0 */
389             qiA = chargeA[i]*one_4pi_eps;
390             qiB = chargeB[i]*one_4pi_eps;
391             if (EVDW_PME(fr->vdwtype))
392             {
393                 c6Ai = C6A[i];
394                 c6Bi = C6B[i];
395                 if (bDoingLBRule)
396                 {
397                     c6Ai *= sigma3A[i];
398                     c6Bi *= sigma3B[i];
399                 }
400             }
401             if (calc_excl_corr)
402             {
403                 i1  = excl->index[i];
404                 i2  = excl->index[i+1];
405
406                 /* Loop over excluded neighbours */
407                 for (j = i1; (j < i2); j++)
408                 {
409                     k = AA[j];
410                     if (k > i)
411                     {
412                         qqA = qiA*chargeA[k];
413                         qqB = qiB*chargeB[k];
414                         if (EVDW_PME(fr->vdwtype))
415                         {
416                             c6A = c6Ai*C6A[k];
417                             c6B = c6Bi*C6B[k];
418                             if (bDoingLBRule)
419                             {
420                                 c6A *= pow(0.5*(sigmaA[i]+sigmaA[k]), 6)*sigma3A[k];
421                                 c6B *= pow(0.5*(sigmaB[i]+sigmaB[k]), 6)*sigma3B[k];
422                             }
423                         }
424                         if (qqA != 0.0 || qqB != 0.0 || c6A != 0.0 || c6B != 0.0)
425                         {
426                             fscal = 0;
427                             qqL   = L1_q*qqA + lambda_q*qqB;
428                             if (EVDW_PME(fr->vdwtype))
429                             {
430                                 c6L = L1_lj*c6A + lambda_lj*c6B;
431                             }
432                             rvec_sub(x[i], x[k], dx);
433                             if (bMolPBC)
434                             {
435                                 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
436                                 for (m = DIM-1; (m >= 0); m--)
437                                 {
438                                     if (dx[m] > 0.5*box[m][m])
439                                     {
440                                         rvec_dec(dx, box[m]);
441                                     }
442                                     else if (dx[m] < -0.5*box[m][m])
443                                     {
444                                         rvec_inc(dx, box[m]);
445                                     }
446                                 }
447                             }
448                             dr2 = norm2(dx);
449                             if (dr2 != 0)
450                             {
451                                 rinv    = gmx_invsqrt(dr2);
452                                 rinv2   = rinv*rinv;
453                                 dr      = 1.0/rinv;
454                                 if (qqA != 0.0 || qqB != 0.0)
455                                 {
456                                     v            = gmx_erf(ewc_q*dr)*rinv;
457                                     vc           = qqL*v;
458                                     Vexcl_q     += vc;
459                                     fscal        = rinv2*(vc-qqL*ewc_q*M_2_SQRTPI*exp(-ewc_q*ewc_q*dr2));
460                                     dvdl_excl_q += (qqB - qqA)*v;
461
462                                     svmul(fscal, dx, df);
463                                     rvec_inc(f[k], df);
464                                     rvec_dec(f[i], df);
465                                     for (iv = 0; (iv < DIM); iv++)
466                                     {
467                                         for (jv = 0; (jv < DIM); jv++)
468                                         {
469                                             dxdf_q[iv][jv] += dx[iv]*df[jv];
470                                         }
471                                     }
472                                 }
473
474                                 if ((c6A != 0.0 || c6B != 0.0) && EVDW_PME(fr->vdwtype))
475                                 {
476                                     rinv6         = rinv2*rinv2*rinv2;
477                                     ewcdr         = ewc_lj*dr;
478                                     ewcdr2        = ewcdr*ewcdr;
479                                     ewcdr4        = ewcdr2*ewcdr2;
480                                     ewcdr5        = ewcdr4*ewcdr;
481                                     v             = -rinv6*(1.0 - exp(-ewcdr2)*(1 + ewcdr2 + 0.5*ewcdr4));
482                                     vc            = c6L*v;
483                                     Vexcl_lj     += vc;
484                                     fscal         = 6.0*vc*rinv + c6L*rinv6*exp(-ewcdr2)*ewc_lj*ewcdr5;
485                                     dvdl_excl_lj += (c6B - c6A)*v;
486
487                                     /* The force vector is obtained by multiplication with the
488                                      * distance vector. Here we also divide by r, to account
489                                      * for the length of dx
490                                      */
491                                     svmul(fscal*rinv, dx, df);
492                                     rvec_inc(f[k], df);
493                                     rvec_dec(f[i], df);
494                                     for (iv = 0; (iv < DIM); iv++)
495                                     {
496                                         for (jv = 0; (jv < DIM); jv++)
497                                         {
498                                             dxdf_lj[iv][jv] += dx[iv]*df[jv];
499                                         }
500                                     }
501                                 }
502                             }
503                             else
504                             {
505                                 Vexcl_q      += qqL*vr0_q;
506                                 dvdl_excl_q  += (qqB - qqA)*vr0_q;
507                                 Vexcl_lj     += c6L*vr0_lj;
508                                 dvdl_excl_lj += (c6B - c6A)*vr0_lj;
509                             }
510                         }
511                     }
512                 }
513             }
514             /* Dipole correction on force */
515             if (dipole_coeff != 0)
516             {
517                 for (j = 0; (j < DIM); j++)
518                 {
519                     f[i][j] -= L1_q*dipcorrA[j]*chargeA[i]
520                         + lambda_q*dipcorrB[j]*chargeB[i];
521                 }
522             }
523         }
524     }
525     for (iv = 0; (iv < DIM); iv++)
526     {
527         for (jv = 0; (jv < DIM); jv++)
528         {
529             vir_q[iv][jv]  += 0.5*dxdf_q[iv][jv];
530             vir_lj[iv][jv] += 0.5*dxdf_lj[iv][jv];
531         }
532     }
533
534     Vself_q[0]  = 0;
535     Vself_q[1]  = 0;
536     Vself_lj[0] = 0;
537     Vself_lj[1] = 0;
538
539     /* Global corrections only on master process */
540     if (MASTER(cr) && thread == 0)
541     {
542         for (q = 0; q < (bFreeEnergy ? 2 : 1); q++)
543         {
544             if (calc_excl_corr)
545             {
546                 /* Self-energy correction */
547                 Vself_q[q] = ewc_q*one_4pi_eps*fr->q2sum[q]*M_1_SQRTPI;
548                 if (EVDW_PME(fr->vdwtype))
549                 {
550                     Vself_lj[q] = -pow(ewc_lj, 6)*fr->c6sum[q];
551                 }
552             }
553
554             /* Apply surface dipole correction:
555              * correction = dipole_coeff * (dipole)^2
556              */
557             if (dipole_coeff != 0)
558             {
559                 if (ewald_geometry == eewg3D)
560                 {
561                     Vdipole[q] = dipole_coeff*iprod(mutot[q], mutot[q]);
562                 }
563                 else if (ewald_geometry == eewg3DC)
564                 {
565                     Vdipole[q] = dipole_coeff*mutot[q][ZZ]*mutot[q][ZZ];
566                 }
567             }
568         }
569     }
570     if (!bFreeEnergy)
571     {
572         *Vcorr_q = Vdipole[0] - Vself_q[0] - Vexcl_q;
573         if (EVDW_PME(fr->vdwtype))
574         {
575             *Vcorr_lj = -Vself_lj[0] - Vexcl_lj;
576         }
577     }
578     else
579     {
580         *Vcorr_q = L1_q*(Vdipole[0] - Vself_q[0])
581             + lambda_q*(Vdipole[1] - Vself_q[1])
582             - Vexcl_q;
583         *dvdlambda_q += Vdipole[1] - Vself_q[1]
584             - (Vdipole[0] - Vself_q[0]) - dvdl_excl_q;
585         if (EVDW_PME(fr->vdwtype))
586         {
587             *Vcorr_lj      = -(L1_lj*Vself_lj[0] + lambda_lj*Vself_lj[1]) - Vexcl_lj;
588             *dvdlambda_lj += -Vself_lj[1] + Vself_lj[0] - dvdl_excl_lj;
589         }
590     }
591
592     if (debug)
593     {
594         fprintf(debug, "Long Range corrections for Ewald interactions:\n");
595         fprintf(debug, "start=%d,natoms=%d\n", start, end-start);
596         fprintf(debug, "q2sum = %g, Vself_q=%g c6sum = %g, Vself_lj=%g\n",
597                 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]);
598         fprintf(debug, "Electrostatic Long Range correction: Vexcl=%g\n", Vexcl_q);
599         fprintf(debug, "Lennard-Jones Long Range correction: Vexcl=%g\n", Vexcl_lj);
600         if (MASTER(cr) && thread == 0)
601         {
602             if (epsilon_surface > 0 || ewald_geometry == eewg3DC)
603             {
604                 fprintf(debug, "Total dipole correction: Vdipole=%g\n",
605                         L1_q*Vdipole[0]+lambda_q*Vdipole[1]);
606             }
607         }
608     }
609 }
610
611 real ewald_charge_correction(t_commrec *cr, t_forcerec *fr, real lambda,
612                              matrix box,
613                              real *dvdlambda, tensor vir)
614
615 {
616     real vol, fac, qs2A, qs2B, vc, enercorr;
617     int  d;
618
619     if (MASTER(cr))
620     {
621         /* Apply charge correction */
622         vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
623
624         fac = M_PI*ONE_4PI_EPS0/(fr->epsilon_r*2.0*vol*vol*sqr(fr->ewaldcoeff_q));
625
626         qs2A = fr->qsum[0]*fr->qsum[0];
627         qs2B = fr->qsum[1]*fr->qsum[1];
628
629         vc = (qs2A*(1 - lambda) + qs2B*lambda)*fac;
630
631         enercorr = -vol*vc;
632
633         *dvdlambda += -vol*(qs2B - qs2A)*fac;
634
635         for (d = 0; d < DIM; d++)
636         {
637             vir[d][d] += vc;
638         }
639
640         if (debug)
641         {
642             fprintf(debug, "Total charge correction: Vcharge=%g\n", enercorr);
643         }
644     }
645     else
646     {
647         enercorr = 0;
648     }
649
650     return enercorr;
651 }