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