Use ListOfLists in gmx_mtop_t and gmx_localtop_t
[alexxy/gromacs.git] / src / gromacs / mdlib / dispersioncorrection.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2019, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 #include "gmxpre.h"
36
37 #include "dispersioncorrection.h"
38
39 #include <cstdio>
40
41 #include "gromacs/math/units.h"
42 #include "gromacs/math/utilities.h"
43 #include "gromacs/math/vec.h"
44 #include "gromacs/mdtypes/forcerec.h"
45 #include "gromacs/mdtypes/inputrec.h"
46 #include "gromacs/mdtypes/md_enums.h"
47 #include "gromacs/mdtypes/nblist.h"
48 #include "gromacs/tables/forcetable.h"
49 #include "gromacs/topology/mtop_util.h"
50 #include "gromacs/topology/topology.h"
51 #include "gromacs/utility/arrayref.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/logger.h"
54
55 /* Implementation here to avoid other files needing to include the file that defines t_nblists */
56 DispersionCorrection::InteractionParams::~InteractionParams() = default;
57
58 /* Returns a matrix, as flat list, of combination rule combined LJ parameters */
59 static std::vector<real> mk_nbfp_combination_rule(const gmx_ffparams_t& ffparams, const int comb_rule)
60 {
61     const int atnr = ffparams.atnr;
62
63     std::vector<real> nbfp(atnr * atnr * 2);
64
65     for (int i = 0; i < atnr; ++i)
66     {
67         for (int j = 0; j < atnr; ++j)
68         {
69             const real c6i  = ffparams.iparams[i * (atnr + 1)].lj.c6;
70             const real c12i = ffparams.iparams[i * (atnr + 1)].lj.c12;
71             const real c6j  = ffparams.iparams[j * (atnr + 1)].lj.c6;
72             const real c12j = ffparams.iparams[j * (atnr + 1)].lj.c12;
73             real       c6   = std::sqrt(c6i * c6j);
74             real       c12  = std::sqrt(c12i * c12j);
75             if (comb_rule == eCOMB_ARITHMETIC && !gmx_numzero(c6) && !gmx_numzero(c12))
76             {
77                 const real sigmai = gmx::sixthroot(c12i / c6i);
78                 const real sigmaj = gmx::sixthroot(c12j / c6j);
79                 const real epsi   = c6i * c6i / c12i;
80                 const real epsj   = c6j * c6j / c12j;
81                 const real sigma  = 0.5 * (sigmai + sigmaj);
82                 const real eps    = std::sqrt(epsi * epsj);
83                 c6                = eps * gmx::power6(sigma);
84                 c12               = eps * gmx::power12(sigma);
85             }
86             C6(nbfp, atnr, i, j)  = c6 * 6.0;
87             C12(nbfp, atnr, i, j) = c12 * 12.0;
88         }
89     }
90
91     return nbfp;
92 }
93
94 /* Returns the A-topology atom type when aOrB=0, the B-topology atom type when aOrB=1 */
95 static int atomtypeAOrB(const t_atom& atom, int aOrB)
96 {
97     if (aOrB == 0)
98     {
99         return atom.type;
100     }
101     else
102     {
103         return atom.typeB;
104     }
105 }
106
107 DispersionCorrection::TopologyParams::TopologyParams(const gmx_mtop_t&         mtop,
108                                                      const t_inputrec&         inputrec,
109                                                      bool                      useBuckingham,
110                                                      int                       numAtomTypes,
111                                                      gmx::ArrayRef<const real> nonbondedForceParameters)
112 {
113     const int      ntp   = numAtomTypes;
114     const gmx_bool bBHAM = useBuckingham;
115
116     gmx::ArrayRef<const real> nbfp = nonbondedForceParameters;
117     std::vector<real>         nbfp_comb;
118     /* For LJ-PME, we want to correct for the difference between the
119      * actual C6 values and the C6 values used by the LJ-PME based on
120      * combination rules. */
121     if (EVDW_PME(inputrec.vdwtype))
122     {
123         nbfp_comb = mk_nbfp_combination_rule(
124                 mtop.ffparams,
125                 (inputrec.ljpme_combination_rule == eljpmeLB) ? eCOMB_ARITHMETIC : eCOMB_GEOMETRIC);
126         for (int tpi = 0; tpi < ntp; ++tpi)
127         {
128             for (int tpj = 0; tpj < ntp; ++tpj)
129             {
130                 C6(nbfp_comb, ntp, tpi, tpj) = C6(nbfp, ntp, tpi, tpj) - C6(nbfp_comb, ntp, tpi, tpj);
131                 C12(nbfp_comb, ntp, tpi, tpj) = C12(nbfp, ntp, tpi, tpj);
132             }
133         }
134         nbfp = nbfp_comb;
135     }
136
137     for (int q = 0; q < (inputrec.efep == efepNO ? 1 : 2); q++)
138     {
139         double csix    = 0;
140         double ctwelve = 0;
141         int    npair   = 0;
142         int    nexcl   = 0;
143         if (!EI_TPI(inputrec.eI))
144         {
145             numAtomsForDensity_ = mtop.natoms;
146             numCorrections_     = 0.5 * mtop.natoms;
147
148             /* Count the types so we avoid natoms^2 operations */
149             std::vector<int> typecount(ntp);
150             gmx_mtop_count_atomtypes(&mtop, q, typecount.data());
151
152             for (int tpi = 0; tpi < ntp; tpi++)
153             {
154                 for (int tpj = tpi; tpj < ntp; tpj++)
155                 {
156                     const int iCount = typecount[tpi];
157                     const int jCount = typecount[tpj];
158                     int       npair_ij;
159                     if (tpi != tpj)
160                     {
161                         npair_ij = iCount * jCount;
162                     }
163                     else
164                     {
165                         npair_ij = iCount * (iCount - 1) / 2;
166                     }
167                     if (bBHAM)
168                     {
169                         /* nbfp now includes the 6.0 derivative prefactor */
170                         csix += npair_ij * BHAMC(nbfp, ntp, tpi, tpj) / 6.0;
171                     }
172                     else
173                     {
174                         /* nbfp now includes the 6.0/12.0 derivative prefactors */
175                         csix += npair_ij * C6(nbfp, ntp, tpi, tpj) / 6.0;
176                         ctwelve += npair_ij * C12(nbfp, ntp, tpi, tpj) / 12.0;
177                     }
178                     npair += npair_ij;
179                 }
180             }
181             /* Subtract the excluded pairs.
182              * The main reason for substracting exclusions is that in some cases
183              * some combinations might never occur and the parameters could have
184              * any value. These unused values should not influence the dispersion
185              * correction.
186              */
187             for (const gmx_molblock_t& molb : mtop.molblock)
188             {
189                 const int      nmol  = molb.nmol;
190                 const t_atoms* atoms = &mtop.moltype[molb.type].atoms;
191                 const auto&    excl  = mtop.moltype[molb.type].excls;
192                 for (int i = 0; (i < atoms->nr); i++)
193                 {
194                     const int tpi = atomtypeAOrB(atoms->atom[i], q);
195                     for (const int k : excl[i])
196                     {
197                         if (k > i)
198                         {
199                             const int tpj = atomtypeAOrB(atoms->atom[k], q);
200                             if (bBHAM)
201                             {
202                                 /* nbfp now includes the 6.0 derivative prefactor */
203                                 csix -= nmol * BHAMC(nbfp, ntp, tpi, tpj) / 6.0;
204                             }
205                             else
206                             {
207                                 /* nbfp now includes the 6.0/12.0 derivative prefactors */
208                                 csix -= nmol * C6(nbfp, ntp, tpi, tpj) / 6.0;
209                                 ctwelve -= nmol * C12(nbfp, ntp, tpi, tpj) / 12.0;
210                             }
211                             nexcl += molb.nmol;
212                         }
213                     }
214                 }
215             }
216         }
217         else
218         {
219             const t_atoms& atoms_tpi = mtop.moltype[mtop.molblock.back().type].atoms;
220
221             /* Only correct for the interaction of the test particle
222              * with the rest of the system.
223              */
224             numAtomsForDensity_ = mtop.natoms - atoms_tpi.nr;
225             numCorrections_     = atoms_tpi.nr;
226
227             npair = 0;
228             for (size_t mb = 0; mb < mtop.molblock.size(); mb++)
229             {
230                 const gmx_molblock_t& molb  = mtop.molblock[mb];
231                 const t_atoms&        atoms = mtop.moltype[molb.type].atoms;
232                 for (int j = 0; j < atoms.nr; j++)
233                 {
234                     int nmolc = molb.nmol;
235                     /* Remove the interaction of the test charge group
236                      * with itself.
237                      */
238                     if (mb == mtop.molblock.size() - 1)
239                     {
240                         nmolc--;
241
242                         if (mb == 0 && molb.nmol == 1)
243                         {
244                             gmx_fatal(FARGS,
245                                       "Old format tpr with TPI, please generate a new tpr file");
246                         }
247                     }
248                     const int tpj = atomtypeAOrB(atoms.atom[j], q);
249                     for (int i = 0; i < atoms_tpi.nr; i++)
250                     {
251                         const int tpi = atomtypeAOrB(atoms_tpi.atom[i], q);
252                         if (bBHAM)
253                         {
254                             /* nbfp now includes the 6.0 derivative prefactor */
255                             csix += nmolc * BHAMC(nbfp, ntp, tpi, tpj) / 6.0;
256                         }
257                         else
258                         {
259                             /* nbfp now includes the 6.0/12.0 derivative prefactors */
260                             csix += nmolc * C6(nbfp, ntp, tpi, tpj) / 6.0;
261                             ctwelve += nmolc * C12(nbfp, ntp, tpi, tpj) / 12.0;
262                         }
263                         npair += nmolc;
264                     }
265                 }
266             }
267         }
268         if (npair - nexcl <= 0)
269         {
270             csix    = 0;
271             ctwelve = 0;
272         }
273         else
274         {
275             csix /= npair - nexcl;
276             ctwelve /= npair - nexcl;
277         }
278         if (debug)
279         {
280             fprintf(debug, "Counted %d exclusions\n", nexcl);
281             fprintf(debug, "Average C6 parameter is: %10g\n", csix);
282             fprintf(debug, "Average C12 parameter is: %10g\n", ctwelve);
283         }
284         avcsix_[q]    = csix;
285         avctwelve_[q] = ctwelve;
286     }
287 }
288
289 static void integrate_table(const real vdwtab[],
290                             const real scale,
291                             const int  offstart,
292                             const int  rstart,
293                             const int  rend,
294                             double*    enerout,
295                             double*    virout)
296 {
297     const double invscale  = 1.0 / scale;
298     const double invscale2 = invscale * invscale;
299     const double invscale3 = invscale * invscale2;
300
301     /* Following summation derived from cubic spline definition,
302      * Numerical Recipies in C, second edition, p. 113-116.  Exact for
303      * the cubic spline.  We first calculate the negative of the
304      * energy from rvdw to rvdw_switch, assuming that g(r)=1, and then
305      * add the more standard, abrupt cutoff correction to that result,
306      * yielding the long-range correction for a switched function.  We
307      * perform both the pressure and energy loops at the same time for
308      * simplicity, as the computational cost is low. */
309
310     /* Since the dispersion table has been scaled down a factor
311      * 6.0 and the repulsion a factor 12.0 to compensate for the
312      * c6/c12 parameters inside nbfp[] being scaled up (to save
313      * flops in kernels), we need to correct for this.
314      */
315     const double tabfactor = (offstart == 0 ? 6.0 : 12.0);
316
317     double enersum = 0.0;
318     double virsum  = 0.0;
319     for (int ri = rstart; ri < rend; ++ri)
320     {
321         const double r  = ri * invscale;
322         const double ea = invscale3;
323         const double eb = 2.0 * invscale2 * r;
324         const double ec = invscale * r * r;
325
326         const double pa = invscale3;
327         const double pb = 3.0 * invscale2 * r;
328         const double pc = 3.0 * invscale * r * r;
329         const double pd = r * r * r;
330
331         /* this "8" is from the packing in the vdwtab array - perhaps
332            should be defined? */
333
334         const int    offset = 8 * ri + offstart;
335         const double y0     = vdwtab[offset];
336         const double f      = vdwtab[offset + 1];
337         const double g      = vdwtab[offset + 2];
338         const double h      = vdwtab[offset + 3];
339
340         enersum += y0 * (ea / 3 + eb / 2 + ec) + f * (ea / 4 + eb / 3 + ec / 2)
341                    + g * (ea / 5 + eb / 4 + ec / 3) + h * (ea / 6 + eb / 5 + ec / 4);
342         virsum += f * (pa / 4 + pb / 3 + pc / 2 + pd) + 2 * g * (pa / 5 + pb / 4 + pc / 3 + pd / 2)
343                   + 3 * h * (pa / 6 + pb / 5 + pc / 4 + pd / 3);
344     }
345     *enerout = 4.0 * M_PI * enersum * tabfactor;
346     *virout  = 4.0 * M_PI * virsum * tabfactor;
347 }
348
349 /* Struct for storing and passing energy or virial corrections */
350 struct InteractionCorrection
351 {
352     real dispersion = 0;
353     real repulsion  = 0;
354 };
355
356 /* Adds the energy and virial corrections beyond the cut-off */
357 static void addCorrectionBeyondCutoff(InteractionCorrection* energy,
358                                       InteractionCorrection* virial,
359                                       const double           cutoffDistance)
360 {
361     const double rc3 = cutoffDistance * cutoffDistance * cutoffDistance;
362     const double rc9 = rc3 * rc3 * rc3;
363
364     energy->dispersion += -4.0 * M_PI / (3.0 * rc3);
365     energy->repulsion += 4.0 * M_PI / (9.0 * rc9);
366     virial->dispersion += 8.0 * M_PI / rc3;
367     virial->repulsion += -16.0 * M_PI / (3.0 * rc9);
368 }
369
370 void DispersionCorrection::setInteractionParameters(InteractionParams*         iParams,
371                                                     const interaction_const_t& ic,
372                                                     const char*                tableFileName)
373 {
374     /* We only need to set the tables at first call, i.e. tableFileName!=nullptr
375      * or when we changed the cut-off with LJ-PME tuning.
376      */
377     if (tableFileName || EVDW_PME(ic.vdwtype))
378     {
379         iParams->dispersionCorrectionTable_ =
380                 makeDispersionCorrectionTable(nullptr, &ic, ic.rvdw, tableFileName);
381     }
382
383     InteractionCorrection energy;
384     InteractionCorrection virial;
385
386     if ((ic.vdw_modifier == eintmodPOTSHIFT) || (ic.vdw_modifier == eintmodPOTSWITCH)
387         || (ic.vdw_modifier == eintmodFORCESWITCH) || (ic.vdwtype == evdwSHIFT)
388         || (ic.vdwtype == evdwSWITCH))
389     {
390         if (((ic.vdw_modifier == eintmodPOTSWITCH) || (ic.vdw_modifier == eintmodFORCESWITCH)
391              || (ic.vdwtype == evdwSWITCH))
392             && ic.rvdw_switch == 0)
393         {
394             gmx_fatal(FARGS,
395                       "With dispersion correction rvdw-switch can not be zero "
396                       "for vdw-type = %s",
397                       evdw_names[ic.vdwtype]);
398         }
399
400         GMX_ASSERT(iParams->dispersionCorrectionTable_, "We need an initialized table");
401
402         /* TODO This code depends on the logic in tables.c that
403            constructs the table layout, which should be made
404            explicit in future cleanup. */
405         GMX_ASSERT(iParams->dispersionCorrectionTable_->interaction == GMX_TABLE_INTERACTION_VDWREP_VDWDISP,
406                    "Dispersion-correction code needs a table with both repulsion and dispersion "
407                    "terms");
408         const real  scale  = iParams->dispersionCorrectionTable_->scale;
409         const real* vdwtab = iParams->dispersionCorrectionTable_->data;
410
411         /* Round the cut-offs to exact table values for precision */
412         int ri0 = static_cast<int>(std::floor(ic.rvdw_switch * scale));
413         int ri1 = static_cast<int>(std::ceil(ic.rvdw * scale));
414
415         /* The code below has some support for handling force-switching, i.e.
416          * when the force (instead of potential) is switched over a limited
417          * region. This leads to a constant shift in the potential inside the
418          * switching region, which we can handle by adding a constant energy
419          * term in the force-switch case just like when we do potential-shift.
420          *
421          * For now this is not enabled, but to keep the functionality in the
422          * code we check separately for switch and shift. When we do force-switch
423          * the shifting point is rvdw_switch, while it is the cutoff when we
424          * have a classical potential-shift.
425          *
426          * For a pure potential-shift the potential has a constant shift
427          * all the way out to the cutoff, and that is it. For other forms
428          * we need to calculate the constant shift up to the point where we
429          * start modifying the potential.
430          */
431         ri0 = (ic.vdw_modifier == eintmodPOTSHIFT) ? ri1 : ri0;
432
433         const double r0  = ri0 / scale;
434         const double rc3 = r0 * r0 * r0;
435         const double rc9 = rc3 * rc3 * rc3;
436
437         if ((ic.vdw_modifier == eintmodFORCESWITCH) || (ic.vdwtype == evdwSHIFT))
438         {
439             /* Determine the constant energy shift below rvdw_switch.
440              * Table has a scale factor since we have scaled it down to compensate
441              * for scaling-up c6/c12 with the derivative factors to save flops in analytical kernels.
442              */
443             iParams->enershiftsix_ = static_cast<real>(-1.0 / (rc3 * rc3)) - 6.0 * vdwtab[8 * ri0];
444             iParams->enershifttwelve_ = static_cast<real>(1.0 / (rc9 * rc3)) - 12.0 * vdwtab[8 * ri0 + 4];
445         }
446         else if (ic.vdw_modifier == eintmodPOTSHIFT)
447         {
448             iParams->enershiftsix_    = static_cast<real>(-1.0 / (rc3 * rc3));
449             iParams->enershifttwelve_ = static_cast<real>(1.0 / (rc9 * rc3));
450         }
451
452         /* Add the constant part from 0 to rvdw_switch.
453          * This integration from 0 to rvdw_switch overcounts the number
454          * of interactions by 1, as it also counts the self interaction.
455          * We will correct for this later.
456          */
457         energy.dispersion += 4.0 * M_PI * iParams->enershiftsix_ * rc3 / 3.0;
458         energy.repulsion += 4.0 * M_PI * iParams->enershifttwelve_ * rc3 / 3.0;
459
460         /* Calculate the contribution in the range [r0,r1] where we
461          * modify the potential. For a pure potential-shift modifier we will
462          * have ri0==ri1, and there will not be any contribution here.
463          */
464         double enersum = 0;
465         double virsum  = 0;
466         integrate_table(vdwtab, scale, 0, ri0, ri1, &enersum, &virsum);
467         energy.dispersion -= enersum;
468         virial.dispersion -= virsum;
469         integrate_table(vdwtab, scale, 4, ri0, ri1, &enersum, &virsum);
470         energy.repulsion -= enersum;
471         virial.repulsion -= virsum;
472
473
474         /* Alright: Above we compensated by REMOVING the parts outside r0
475          * corresponding to the ideal VdW 1/r6 and /r12 potentials.
476          *
477          * Regardless of whether r0 is the point where we start switching,
478          * or the cutoff where we calculated the constant shift, we include
479          * all the parts we are missing out to infinity from r0 by
480          * calculating the analytical dispersion correction.
481          */
482         addCorrectionBeyondCutoff(&energy, &virial, r0);
483     }
484     else if (ic.vdwtype == evdwCUT || EVDW_PME(ic.vdwtype) || ic.vdwtype == evdwUSER)
485     {
486         /* Note that with LJ-PME, the dispersion correction is multiplied
487          * by the difference between the actual C6 and the value of C6
488          * that would produce the combination rule.
489          * This means the normal energy and virial difference formulas
490          * can be used here.
491          */
492
493         const double rc3 = ic.rvdw * ic.rvdw * ic.rvdw;
494         const double rc9 = rc3 * rc3 * rc3;
495         if (ic.vdw_modifier == eintmodPOTSHIFT)
496         {
497             /* Contribution within the cut-off */
498             energy.dispersion += -4.0 * M_PI / (3.0 * rc3);
499             energy.repulsion += 4.0 * M_PI / (3.0 * rc9);
500         }
501         /* Contribution beyond the cut-off */
502         addCorrectionBeyondCutoff(&energy, &virial, ic.rvdw);
503     }
504     else
505     {
506         gmx_fatal(FARGS, "Dispersion correction is not implemented for vdw-type = %s",
507                   evdw_names[ic.vdwtype]);
508     }
509
510     iParams->enerdiffsix_    = energy.dispersion;
511     iParams->enerdifftwelve_ = energy.repulsion;
512     /* The 0.5 is due to the Gromacs definition of the virial */
513     iParams->virdiffsix_    = 0.5 * virial.dispersion;
514     iParams->virdifftwelve_ = 0.5 * virial.repulsion;
515 }
516
517 DispersionCorrection::DispersionCorrection(const gmx_mtop_t&          mtop,
518                                            const t_inputrec&          inputrec,
519                                            bool                       useBuckingham,
520                                            int                        numAtomTypes,
521                                            gmx::ArrayRef<const real>  nonbondedForceParameters,
522                                            const interaction_const_t& ic,
523                                            const char*                tableFileName) :
524     eDispCorr_(inputrec.eDispCorr),
525     vdwType_(inputrec.vdwtype),
526     eFep_(inputrec.efep),
527     topParams_(mtop, inputrec, useBuckingham, numAtomTypes, nonbondedForceParameters)
528 {
529     if (eDispCorr_ != edispcNO)
530     {
531         GMX_RELEASE_ASSERT(tableFileName, "Need a table file name");
532
533         setInteractionParameters(&iParams_, ic, tableFileName);
534     }
535 }
536
537 bool DispersionCorrection::correctFullInteraction() const
538 {
539     return (eDispCorr_ == edispcAllEner || eDispCorr_ == edispcAllEnerPres);
540 }
541
542 void DispersionCorrection::print(const gmx::MDLogger& mdlog) const
543 {
544     if (topParams_.avcsix_[0] == 0 && topParams_.avctwelve_[0] == 0)
545     {
546         GMX_LOG(mdlog.warning)
547                 .asParagraph()
548                 .appendText("WARNING: There are no atom pairs for dispersion correction");
549     }
550     else if (vdwType_ == evdwUSER)
551     {
552         GMX_LOG(mdlog.warning)
553                 .asParagraph()
554                 .appendText("WARNING: using dispersion correction with user tables\n");
555     }
556
557     std::string text = gmx::formatString("Long Range LJ corr.: <C6> %10.4e", topParams_.avcsix_[0]);
558     if (correctFullInteraction())
559     {
560         text += gmx::formatString(" <C12> %10.4e", topParams_.avctwelve_[0]);
561     }
562     GMX_LOG(mdlog.info).appendText(text);
563 }
564
565 void DispersionCorrection::setParameters(const interaction_const_t& ic)
566 {
567     if (eDispCorr_ != edispcNO)
568     {
569         setInteractionParameters(&iParams_, ic, nullptr);
570     }
571 }
572
573 DispersionCorrection::Correction DispersionCorrection::calculate(const matrix box, const real lambda) const
574
575 {
576     Correction corr;
577
578     if (eDispCorr_ == edispcNO)
579     {
580         return corr;
581     }
582
583     const bool bCorrAll  = correctFullInteraction();
584     const bool bCorrPres = (eDispCorr_ == edispcEnerPres || eDispCorr_ == edispcAllEnerPres);
585
586     const real invvol  = 1 / det(box);
587     const real density = topParams_.numAtomsForDensity_ * invvol;
588     const real numCorr = topParams_.numCorrections_;
589
590     real avcsix;
591     real avctwelve;
592     if (eFep_ == efepNO)
593     {
594         avcsix    = topParams_.avcsix_[0];
595         avctwelve = topParams_.avctwelve_[0];
596     }
597     else
598     {
599         avcsix    = (1 - lambda) * topParams_.avcsix_[0] + lambda * topParams_.avcsix_[1];
600         avctwelve = (1 - lambda) * topParams_.avctwelve_[0] + lambda * topParams_.avctwelve_[1];
601     }
602
603     const real enerdiff = numCorr * (density * iParams_.enerdiffsix_ - iParams_.enershiftsix_);
604     corr.energy += avcsix * enerdiff;
605     real dvdlambda = 0;
606     if (eFep_ != efepNO)
607     {
608         dvdlambda += (topParams_.avcsix_[1] - topParams_.avcsix_[0]) * enerdiff;
609     }
610     if (bCorrAll)
611     {
612         const real enerdiff = numCorr * (density * iParams_.enerdifftwelve_ - iParams_.enershifttwelve_);
613         corr.energy += avctwelve * enerdiff;
614         if (eFep_ != efepNO)
615         {
616             dvdlambda += (topParams_.avctwelve_[1] - topParams_.avctwelve_[0]) * enerdiff;
617         }
618     }
619
620     if (bCorrPres)
621     {
622         corr.virial = numCorr * density * avcsix * iParams_.virdiffsix_ / 3.0;
623         if (eDispCorr_ == edispcAllEnerPres)
624         {
625             corr.virial += numCorr * density * avctwelve * iParams_.virdifftwelve_ / 3.0;
626         }
627         /* The factor 2 is because of the Gromacs virial definition */
628         corr.pressure = -2.0 * invvol * corr.virial * PRESFAC;
629     }
630
631     if (eFep_ != efepNO)
632     {
633         corr.dvdl += dvdlambda;
634     }
635
636     return corr;
637 }