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