Fix cmake policy warning
[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                         const real *chargeA, const real *chargeB,
73                         const real *C6A, const real *C6B,
74                         const real *sigmaA, const real *sigmaB,
75                         const real *sigma3A, const 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     }
190     if (debug)
191     {
192         fprintf(debug, "dipcorr = %8.3f  %8.3f  %8.3f\n",
193                 dipcorrA[XX], dipcorrA[YY], dipcorrA[ZZ]);
194         fprintf(debug, "mutot   = %8.3f  %8.3f  %8.3f\n",
195                 mutot[0][XX], mutot[0][YY], mutot[0][ZZ]);
196     }
197     bNeedLongRangeCorrection = (calc_excl_corr || dipole_coeff != 0);
198     if (bNeedLongRangeCorrection && !bHaveChargeOrTypePerturbed)
199     {
200         for (i = start; (i < end); i++)
201         {
202             /* Initiate local variables (for this i-particle) to 0 */
203             qiA = chargeA[i]*one_4pi_eps;
204             if (vdwPme)
205             {
206                 c6Ai = C6A[i];
207                 if (bDoingLBRule)
208                 {
209                     c6Ai *= sigma3A[i];
210                 }
211             }
212             if (calc_excl_corr)
213             {
214                 i1  = excl->index[i];
215                 i2  = excl->index[i+1];
216
217                 /* Loop over excluded neighbours */
218                 for (j = i1; (j < i2); j++)
219                 {
220                     k = AA[j];
221                     /*
222                      * First we must test whether k <> i, and then,
223                      * because the exclusions are all listed twice i->k
224                      * and k->i we must select just one of the two.  As
225                      * a minor optimization we only compute forces when
226                      * the charges are non-zero.
227                      */
228                     if (k > i)
229                     {
230                         qqA = qiA*chargeA[k];
231                         if (vdwPme)
232                         {
233                             c6A  = c6Ai * C6A[k];
234                             if (bDoingLBRule)
235                             {
236                                 c6A *= gmx::power6(0.5*(sigmaA[i]+sigmaA[k]))*sigma3A[k];
237                             }
238                         }
239                         if (qqA != 0.0 || c6A != 0.0)
240                         {
241                             rvec_sub(x[i], x[k], dx);
242                             if (bMolPBC)
243                             {
244                                 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
245                                 for (m = DIM-1; (m >= 0); m--)
246                                 {
247                                     if (dx[m] > 0.5*box[m][m])
248                                     {
249                                         rvec_dec(dx, box[m]);
250                                     }
251                                     else if (dx[m] < -0.5*box[m][m])
252                                     {
253                                         rvec_inc(dx, box[m]);
254                                     }
255                                 }
256                             }
257                             dr2 = norm2(dx);
258                             /* Distance between two excluded particles
259                              * may be zero in the case of shells
260                              */
261                             if (dr2 != 0)
262                             {
263                                 rinv              = gmx::invsqrt(dr2);
264                                 rinv2             = rinv*rinv;
265                                 if (qqA != 0.0)
266                                 {
267                                     real dr, fscal;
268
269                                     dr       = 1.0/rinv;
270                                     ewcdr    = ewc_q*dr;
271                                     vc       = qqA*std::erf(ewcdr)*rinv;
272                                     Vexcl_q += vc;
273 #if GMX_DOUBLE
274                                     /* Relative accuracy at R_ERF_R_INACC of 3e-10 */
275 #define       R_ERF_R_INACC 0.006
276 #else
277                                     /* Relative accuracy at R_ERF_R_INACC of 2e-5 */
278 #define       R_ERF_R_INACC 0.1
279 #endif
280                                     /* fscal is the scalar force pre-multiplied by rinv,
281                                      * to normalise the relative position vector dx */
282                                     if (ewcdr > R_ERF_R_INACC)
283                                     {
284                                         fscal = rinv2*(vc - qqA*ewc_q*M_2_SQRTPI*std::exp(-ewcdr*ewcdr));
285                                     }
286                                     else
287                                     {
288                                         /* Use a fourth order series expansion for small ewcdr */
289                                         fscal = ewc_q*ewc_q*qqA*vr0_q*(2.0/3.0 - 0.4*ewcdr*ewcdr);
290                                     }
291
292                                     /* The force vector is obtained by multiplication with
293                                      * the relative position vector
294                                      */
295                                     svmul(fscal, dx, df);
296                                     rvec_inc(f[k], df);
297                                     rvec_dec(f[i], df);
298                                     for (iv = 0; (iv < DIM); iv++)
299                                     {
300                                         for (jv = 0; (jv < DIM); jv++)
301                                         {
302                                             dxdf_q[iv][jv] += dx[iv]*df[jv];
303                                         }
304                                     }
305                                 }
306
307                                 if (c6A != 0.0)
308                                 {
309                                     real fscal;
310
311                                     rinv6     = rinv2*rinv2*rinv2;
312                                     ewcdr2    = ewc_lj2*dr2;
313                                     ewcdr4    = ewcdr2*ewcdr2;
314                                     /* We get the excluded long-range contribution from -C6*(1-g(r))
315                                      * g(r) is also defined in the manual under LJ-PME
316                                      */
317                                     vc        = -c6A*rinv6*(1.0 - exp(-ewcdr2)*(1 + ewcdr2 + 0.5*ewcdr4));
318                                     Vexcl_lj += vc;
319                                     /* The force is the derivative of the potential vc.
320                                      * fscal is the scalar force pre-multiplied by rinv,
321                                      * to normalise the relative position vector dx */
322                                     fscal     = 6.0*vc*rinv2 + c6A*rinv6*exp(-ewcdr2)*ewc_lj2*ewcdr4;
323
324                                     /* The force vector is obtained by multiplication with
325                                      * the relative position vector
326                                      */
327                                     svmul(fscal, dx, df);
328                                     rvec_inc(f[k], df);
329                                     rvec_dec(f[i], df);
330                                     for (iv = 0; (iv < DIM); iv++)
331                                     {
332                                         for (jv = 0; (jv < DIM); jv++)
333                                         {
334                                             dxdf_lj[iv][jv] += dx[iv]*df[jv];
335                                         }
336                                     }
337                                 }
338                             }
339                             else
340                             {
341                                 Vexcl_q  += qqA*vr0_q;
342                                 Vexcl_lj += c6A*vr0_lj;
343                             }
344                         }
345                     }
346                 }
347             }
348             /* Dipole correction on force */
349             if (dipole_coeff != 0 && i < numAtomsLocal)
350             {
351                 for (j = 0; (j < DIM); j++)
352                 {
353                     f[i][j] -= dipcorrA[j]*chargeA[i];
354                 }
355                 if (chargecorr[0] != 0)
356                 {
357                     f[i][ZZ] += chargecorr[0]*chargeA[i]*x[i][ZZ];
358                 }
359             }
360         }
361     }
362     else if (bNeedLongRangeCorrection)
363     {
364         for (i = start; (i < end); i++)
365         {
366             /* Initiate local variables (for this i-particle) to 0 */
367             qiA = chargeA[i]*one_4pi_eps;
368             qiB = chargeB[i]*one_4pi_eps;
369             if (vdwPme)
370             {
371                 c6Ai = C6A[i];
372                 c6Bi = C6B[i];
373                 if (bDoingLBRule)
374                 {
375                     c6Ai *= sigma3A[i];
376                     c6Bi *= sigma3B[i];
377                 }
378             }
379             if (calc_excl_corr)
380             {
381                 i1  = excl->index[i];
382                 i2  = excl->index[i+1];
383
384                 /* Loop over excluded neighbours */
385                 for (j = i1; (j < i2); j++)
386                 {
387                     k = AA[j];
388                     if (k > i)
389                     {
390                         qqA = qiA*chargeA[k];
391                         qqB = qiB*chargeB[k];
392                         if (vdwPme)
393                         {
394                             c6A = c6Ai*C6A[k];
395                             c6B = c6Bi*C6B[k];
396                             if (bDoingLBRule)
397                             {
398                                 c6A *= gmx::power6(0.5*(sigmaA[i]+sigmaA[k]))*sigma3A[k];
399                                 c6B *= gmx::power6(0.5*(sigmaB[i]+sigmaB[k]))*sigma3B[k];
400                             }
401                         }
402                         if (qqA != 0.0 || qqB != 0.0 || c6A != 0.0 || c6B != 0.0)
403                         {
404                             real fscal;
405
406                             qqL   = L1_q*qqA + lambda_q*qqB;
407                             if (vdwPme)
408                             {
409                                 c6L = L1_lj*c6A + lambda_lj*c6B;
410                             }
411                             rvec_sub(x[i], x[k], dx);
412                             if (bMolPBC)
413                             {
414                                 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
415                                 for (m = DIM-1; (m >= 0); m--)
416                                 {
417                                     if (dx[m] > 0.5*box[m][m])
418                                     {
419                                         rvec_dec(dx, box[m]);
420                                     }
421                                     else if (dx[m] < -0.5*box[m][m])
422                                     {
423                                         rvec_inc(dx, box[m]);
424                                     }
425                                 }
426                             }
427                             dr2 = norm2(dx);
428                             if (dr2 != 0)
429                             {
430                                 rinv    = gmx::invsqrt(dr2);
431                                 rinv2   = rinv*rinv;
432                                 if (qqA != 0.0 || qqB != 0.0)
433                                 {
434                                     real dr;
435
436                                     dr           = 1.0/rinv;
437                                     v            = std::erf(ewc_q*dr)*rinv;
438                                     vc           = qqL*v;
439                                     Vexcl_q     += vc;
440                                     /* fscal is the scalar force pre-multiplied by rinv,
441                                      * to normalise the relative position vector dx */
442                                     fscal        = rinv2*(vc-qqL*ewc_q*M_2_SQRTPI*std::exp(-ewc_q*ewc_q*dr2));
443                                     dvdl_excl_q += (qqB - qqA)*v;
444
445                                     /* The force vector is obtained by multiplication with
446                                      * the relative position vector
447                                      */
448                                     svmul(fscal, dx, df);
449                                     rvec_inc(f[k], df);
450                                     rvec_dec(f[i], df);
451                                     for (iv = 0; (iv < DIM); iv++)
452                                     {
453                                         for (jv = 0; (jv < DIM); jv++)
454                                         {
455                                             dxdf_q[iv][jv] += dx[iv]*df[jv];
456                                         }
457                                     }
458                                 }
459
460                                 if ((c6A != 0.0 || c6B != 0.0) && vdwPme)
461                                 {
462                                     rinv6         = rinv2*rinv2*rinv2;
463                                     ewcdr2        = ewc_lj2*dr2;
464                                     ewcdr4        = ewcdr2*ewcdr2;
465                                     v             = -rinv6*(1.0 - exp(-ewcdr2)*(1 + ewcdr2 + 0.5*ewcdr4));
466                                     vc            = c6L*v;
467                                     Vexcl_lj     += vc;
468                                     /* fscal is the scalar force pre-multiplied by rinv,
469                                      * to normalise the relative position vector dx */
470                                     fscal         = 6.0*vc*rinv2 + c6L*rinv6*exp(-ewcdr2)*ewc_lj2*ewcdr4;
471                                     dvdl_excl_lj += (c6B - c6A)*v;
472
473                                     /* The force vector is obtained by multiplication with
474                                      * the relative position vector
475                                      */
476                                     svmul(fscal, dx, df);
477                                     rvec_inc(f[k], df);
478                                     rvec_dec(f[i], df);
479                                     for (iv = 0; (iv < DIM); iv++)
480                                     {
481                                         for (jv = 0; (jv < DIM); jv++)
482                                         {
483                                             dxdf_lj[iv][jv] += dx[iv]*df[jv];
484                                         }
485                                     }
486                                 }
487                             }
488                             else
489                             {
490                                 Vexcl_q      += qqL*vr0_q;
491                                 dvdl_excl_q  += (qqB - qqA)*vr0_q;
492                                 Vexcl_lj     += c6L*vr0_lj;
493                                 dvdl_excl_lj += (c6B - c6A)*vr0_lj;
494                             }
495                         }
496                     }
497                 }
498             }
499             /* Dipole correction on force */
500             if (dipole_coeff != 0 && i < numAtomsLocal)
501             {
502                 for (j = 0; (j < DIM); j++)
503                 {
504                     f[i][j] -= L1_q*dipcorrA[j]*chargeA[i]
505                         + lambda_q*dipcorrB[j]*chargeB[i];
506                 }
507                 if (chargecorr[0] != 0 || chargecorr[1] != 0)
508                 {
509                     f[i][ZZ] += (L1_q*chargecorr[0]*chargeA[i]
510                                  + lambda_q*chargecorr[1])*x[i][ZZ];
511                 }
512             }
513         }
514     }
515     for (iv = 0; (iv < DIM); iv++)
516     {
517         for (jv = 0; (jv < DIM); jv++)
518         {
519             vir_q[iv][jv]  += 0.5*dxdf_q[iv][jv];
520             vir_lj[iv][jv] += 0.5*dxdf_lj[iv][jv];
521         }
522     }
523
524     Vself_q[0]  = 0;
525     Vself_q[1]  = 0;
526     Vself_lj[0] = 0;
527     Vself_lj[1] = 0;
528
529     /* Global corrections only on master process */
530     if (MASTER(cr) && thread == 0)
531     {
532         for (q = 0; q < (bHaveChargeOrTypePerturbed ? 2 : 1); q++)
533         {
534             if (calc_excl_corr)
535             {
536                 /* Self-energy correction */
537                 Vself_q[q] = ewc_q*one_4pi_eps*fr->q2sum[q]*M_1_SQRTPI;
538                 if (vdwPme)
539                 {
540                     Vself_lj[q] =  fr->c6sum[q]*0.5*vr0_lj;
541                 }
542             }
543
544             /* Apply surface and charged surface dipole correction:
545              * correction = dipole_coeff * ( (dipole)^2
546              *              - qsum*sum_i q_i z_i^2 - qsum^2 * box_z^2 / 12 )
547              */
548             if (dipole_coeff != 0)
549             {
550                 if (ewald_geometry == eewg3D)
551                 {
552                     Vdipole[q] = dipole_coeff*iprod(mutot[q], mutot[q]);
553                 }
554                 else if (ewald_geometry == eewg3DC)
555                 {
556                     Vdipole[q] = dipole_coeff*mutot[q][ZZ]*mutot[q][ZZ];
557
558                     if (chargecorr[q] != 0)
559                     {
560                         /* Here we use a non thread-parallelized loop,
561                          * because this is the only loop over atoms for
562                          * energies and they need reduction (unlike forces).
563                          * We could implement a reduction over threads,
564                          * but this case is rarely used.
565                          */
566                         const real *qPtr   = (q == 0 ? chargeA : chargeB);
567                         real        sumQZ2 = 0;
568                         for (int i = 0; i < numAtomsLocal; i++)
569                         {
570                             sumQZ2 += qPtr[i]*x[i][ZZ]*x[i][ZZ];
571                         }
572                         Vdipole[q] -= dipole_coeff*fr->qsum[q]*(sumQZ2 + fr->qsum[q]*box[ZZ][ZZ]*box[ZZ][ZZ]/12);
573                     }
574                 }
575             }
576         }
577     }
578     if (!bHaveChargeOrTypePerturbed)
579     {
580         *Vcorr_q = Vdipole[0] - Vself_q[0] - Vexcl_q;
581         if (vdwPme)
582         {
583             *Vcorr_lj = -Vself_lj[0] - Vexcl_lj;
584         }
585     }
586     else
587     {
588         *Vcorr_q = L1_q*(Vdipole[0] - Vself_q[0])
589             + lambda_q*(Vdipole[1] - Vself_q[1])
590             - Vexcl_q;
591         *dvdlambda_q += Vdipole[1] - Vself_q[1]
592             - (Vdipole[0] - Vself_q[0]) - dvdl_excl_q;
593         if (vdwPme)
594         {
595             *Vcorr_lj      = -(L1_lj*Vself_lj[0] + lambda_lj*Vself_lj[1]) - Vexcl_lj;
596             *dvdlambda_lj += -Vself_lj[1] + Vself_lj[0] - dvdl_excl_lj;
597         }
598     }
599
600     if (debug)
601     {
602         fprintf(debug, "Long Range corrections for Ewald interactions:\n");
603         fprintf(debug, "q2sum = %g, Vself_q=%g c6sum = %g, Vself_lj=%g\n",
604                 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]);
605         fprintf(debug, "Electrostatic Long Range correction: Vexcl=%g\n", Vexcl_q);
606         fprintf(debug, "Lennard-Jones Long Range correction: Vexcl=%g\n", Vexcl_lj);
607         if (MASTER(cr) && thread == 0)
608         {
609             if (epsilon_surface > 0 || ewald_geometry == eewg3DC)
610             {
611                 fprintf(debug, "Total dipole correction: Vdipole=%g\n",
612                         L1_q*Vdipole[0]+lambda_q*Vdipole[1]);
613             }
614         }
615     }
616 }