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