9353a24047aa105f840968e450cd0f975652f118
[alexxy/gromacs.git] / src / gromacs / mdtypes / interaction_const.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2020,2021, 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 "interaction_const.h"
38
39 #include <cstdio>
40
41 #include "gromacs/ewald/ewald_utils.h"
42 #include "gromacs/math/functions.h"
43 #include "gromacs/math/units.h"
44 #include "gromacs/mdlib/rf_util.h"
45 #include "gromacs/mdtypes/forcerec.h"
46 #include "gromacs/mdtypes/inputrec.h"
47 #include "gromacs/topology/topology.h"
48 #include "gromacs/utility/fatalerror.h"
49 #include "gromacs/utility/pleasecite.h"
50
51 interaction_const_t::SoftCoreParameters::SoftCoreParameters(const t_lambda& fepvals) :
52     alphaVdw(fepvals.sc_alpha),
53     alphaCoulomb(fepvals.bScCoul ? fepvals.sc_alpha : 0),
54     lambdaPower(fepvals.sc_power),
55     sigma6WithInvalidSigma(gmx::power6(fepvals.sc_sigma)),
56     sigma6Minimum(fepvals.bScCoul ? gmx::power6(fepvals.sc_sigma_min) : 0)
57 {
58     // This is checked during tpr reading, so we can assert here
59     GMX_RELEASE_ASSERT(fepvals.sc_r_power == 6.0, "We only support soft-core r-power 6");
60 }
61
62 /*! \brief Print Coulomb Ewald citations and set ewald coefficients */
63 static void initCoulombEwaldParameters(FILE*                fp,
64                                        const t_inputrec&    ir,
65                                        bool                 systemHasNetCharge,
66                                        interaction_const_t* ic)
67 {
68     if (!EEL_PME_EWALD(ir.coulombtype))
69     {
70         return;
71     }
72
73     if (fp)
74     {
75         fprintf(fp, "Will do PME sum in reciprocal space for electrostatic interactions.\n");
76
77         if (ir.coulombtype == CoulombInteractionType::P3mAD)
78         {
79             please_cite(fp, "Hockney1988");
80             please_cite(fp, "Ballenegger2012");
81         }
82         else
83         {
84             please_cite(fp, "Essmann95a");
85         }
86
87         if (ir.ewald_geometry == EwaldGeometry::ThreeDC)
88         {
89             if (fp)
90             {
91                 fprintf(fp,
92                         "Using the Ewald3DC correction for systems with a slab geometry%s.\n",
93                         systemHasNetCharge ? " and net charge" : "");
94             }
95             please_cite(fp, "In-Chul99a");
96             if (systemHasNetCharge)
97             {
98                 please_cite(fp, "Ballenegger2009");
99             }
100         }
101     }
102
103     ic->ewaldcoeff_q = calc_ewaldcoeff_q(ir.rcoulomb, ir.ewald_rtol);
104     if (fp)
105     {
106         fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for Ewald\n", 1 / ic->ewaldcoeff_q);
107     }
108
109     if (ic->coulomb_modifier == InteractionModifiers::PotShift)
110     {
111         GMX_RELEASE_ASSERT(ic->rcoulomb != 0, "Cutoff radius cannot be zero");
112         ic->sh_ewald = std::erfc(ic->ewaldcoeff_q * ic->rcoulomb) / ic->rcoulomb;
113     }
114     else
115     {
116         ic->sh_ewald = 0;
117     }
118 }
119
120 /*! \brief Print Van der Waals Ewald citations and set ewald coefficients */
121 static void initVdwEwaldParameters(FILE* fp, const t_inputrec& ir, interaction_const_t* ic)
122 {
123     if (!EVDW_PME(ir.vdwtype))
124     {
125         return;
126     }
127
128     if (fp)
129     {
130         fprintf(fp, "Will do PME sum in reciprocal space for LJ dispersion interactions.\n");
131         please_cite(fp, "Essmann95a");
132     }
133     ic->ewaldcoeff_lj = calc_ewaldcoeff_lj(ir.rvdw, ir.ewald_rtol_lj);
134     if (fp)
135     {
136         fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for LJ Ewald\n", 1 / ic->ewaldcoeff_lj);
137     }
138
139     if (ic->vdw_modifier == InteractionModifiers::PotShift)
140     {
141         real crc2       = gmx::square(ic->ewaldcoeff_lj * ic->rvdw);
142         ic->sh_lj_ewald = (std::exp(-crc2) * (1 + crc2 + 0.5 * crc2 * crc2) - 1) / gmx::power6(ic->rvdw);
143     }
144     else
145     {
146         ic->sh_lj_ewald = 0;
147     }
148 }
149
150 static real calcBuckinghamBMax(FILE* fplog, const gmx_mtop_t& mtop)
151 {
152     const t_atoms *at1, *at2;
153     int            i, j, tpi, tpj, ntypes;
154     real           b, bmin;
155
156     if (fplog)
157     {
158         fprintf(fplog, "Determining largest Buckingham b parameter for table\n");
159     }
160     ntypes = mtop.ffparams.atnr;
161
162     bmin            = -1;
163     real bham_b_max = 0;
164     for (size_t mt1 = 0; mt1 < mtop.moltype.size(); mt1++)
165     {
166         at1 = &mtop.moltype[mt1].atoms;
167         for (i = 0; (i < at1->nr); i++)
168         {
169             tpi = at1->atom[i].type;
170             if (tpi >= ntypes)
171             {
172                 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", i, tpi, ntypes);
173             }
174
175             for (size_t mt2 = mt1; mt2 < mtop.moltype.size(); mt2++)
176             {
177                 at2 = &mtop.moltype[mt2].atoms;
178                 for (j = 0; (j < at2->nr); j++)
179                 {
180                     tpj = at2->atom[j].type;
181                     if (tpj >= ntypes)
182                     {
183                         gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", j, tpj, ntypes);
184                     }
185                     b = mtop.ffparams.iparams[tpi * ntypes + tpj].bham.b;
186                     if (b > bham_b_max)
187                     {
188                         bham_b_max = b;
189                     }
190                     if ((b < bmin) || (bmin == -1))
191                     {
192                         bmin = b;
193                     }
194                 }
195             }
196         }
197     }
198     if (fplog)
199     {
200         fprintf(fplog, "Buckingham b parameters, min: %g, max: %g\n", bmin, bham_b_max);
201     }
202
203     return bham_b_max;
204 }
205
206 static void clear_force_switch_constants(shift_consts_t* sc)
207 {
208     sc->c2   = 0;
209     sc->c3   = 0;
210     sc->cpot = 0;
211 }
212
213 static void force_switch_constants(real p, real rsw, real rc, shift_consts_t* sc)
214 {
215     /* Here we determine the coefficient for shifting the force to zero
216      * between distance rsw and the cut-off rc.
217      * For a potential of r^-p, we have force p*r^-(p+1).
218      * But to save flops we absorb p in the coefficient.
219      * Thus we get:
220      * force/p   = r^-(p+1) + c2*r^2 + c3*r^3
221      * potential = r^-p + c2/3*r^3 + c3/4*r^4 + cpot
222      */
223     sc->c2   = ((p + 1) * rsw - (p + 4) * rc) / (pow(rc, p + 2) * gmx::square(rc - rsw));
224     sc->c3   = -((p + 1) * rsw - (p + 3) * rc) / (pow(rc, p + 2) * gmx::power3(rc - rsw));
225     sc->cpot = -pow(rc, -p) + p * sc->c2 / 3 * gmx::power3(rc - rsw)
226                + p * sc->c3 / 4 * gmx::power4(rc - rsw);
227 }
228
229 static void potential_switch_constants(real rsw, real rc, switch_consts_t* sc)
230 {
231     /* The switch function is 1 at rsw and 0 at rc.
232      * The derivative and second derivate are zero at both ends.
233      * rsw        = max(r - r_switch, 0)
234      * sw         = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
235      * dsw        = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
236      * force      = force*dsw - potential*sw
237      * potential *= sw
238      */
239     sc->c3 = -10 / gmx::power3(rc - rsw);
240     sc->c4 = 15 / gmx::power4(rc - rsw);
241     sc->c5 = -6 / gmx::power5(rc - rsw);
242 }
243
244
245 interaction_const_t init_interaction_const(FILE* fp, const t_inputrec& ir, const gmx_mtop_t& mtop, bool systemHasNetCharge)
246 {
247     interaction_const_t interactionConst;
248
249     interactionConst.coulombEwaldTables = std::make_unique<EwaldCorrectionTables>();
250     interactionConst.vdwEwaldTables     = std::make_unique<EwaldCorrectionTables>();
251
252     /* Lennard-Jones */
253     interactionConst.vdwtype         = ir.vdwtype;
254     interactionConst.vdw_modifier    = ir.vdw_modifier;
255     interactionConst.reppow          = mtop.ffparams.reppow;
256     interactionConst.rvdw            = cutoff_inf(ir.rvdw);
257     interactionConst.rvdw_switch     = ir.rvdw_switch;
258     interactionConst.ljpme_comb_rule = ir.ljpme_combination_rule;
259     interactionConst.useBuckingham   = (mtop.ffparams.functype[0] == F_BHAM);
260     if (interactionConst.useBuckingham)
261     {
262         interactionConst.buckinghamBMax = calcBuckinghamBMax(fp, mtop);
263     }
264
265     initVdwEwaldParameters(fp, ir, &interactionConst);
266
267     clear_force_switch_constants(&interactionConst.dispersion_shift);
268     clear_force_switch_constants(&interactionConst.repulsion_shift);
269
270     switch (interactionConst.vdw_modifier)
271     {
272         case InteractionModifiers::PotShift:
273             /* Only shift the potential, don't touch the force */
274             interactionConst.dispersion_shift.cpot = -1.0 / gmx::power6(interactionConst.rvdw);
275             interactionConst.repulsion_shift.cpot  = -1.0 / gmx::power12(interactionConst.rvdw);
276             break;
277         case InteractionModifiers::ForceSwitch:
278             /* Switch the force, switch and shift the potential */
279             force_switch_constants(
280                     6.0, interactionConst.rvdw_switch, interactionConst.rvdw, &interactionConst.dispersion_shift);
281             force_switch_constants(
282                     12.0, interactionConst.rvdw_switch, interactionConst.rvdw, &interactionConst.repulsion_shift);
283             break;
284         case InteractionModifiers::PotSwitch:
285             /* Switch the potential and force */
286             potential_switch_constants(
287                     interactionConst.rvdw_switch, interactionConst.rvdw, &interactionConst.vdw_switch);
288             break;
289         case InteractionModifiers::None:
290         case InteractionModifiers::ExactCutoff:
291             /* Nothing to do here */
292             break;
293         default: gmx_incons("unimplemented potential modifier");
294     }
295
296     /* Electrostatics */
297     interactionConst.eeltype          = ir.coulombtype;
298     interactionConst.coulomb_modifier = ir.coulomb_modifier;
299     interactionConst.rcoulomb         = cutoff_inf(ir.rcoulomb);
300     interactionConst.rcoulomb_switch  = ir.rcoulomb_switch;
301     interactionConst.epsilon_r        = ir.epsilon_r;
302
303     /* Set the Coulomb energy conversion factor */
304     if (interactionConst.epsilon_r != 0)
305     {
306         interactionConst.epsfac = gmx::c_one4PiEps0 / interactionConst.epsilon_r;
307     }
308     else
309     {
310         /* eps = 0 is infinite dieletric: no Coulomb interactions */
311         interactionConst.epsfac = 0;
312     }
313
314     /* Reaction-field */
315     if (EEL_RF(interactionConst.eeltype))
316     {
317         GMX_RELEASE_ASSERT(interactionConst.eeltype != CoulombInteractionType::GRFNotused,
318                            "GRF is no longer supported");
319         interactionConst.reactionFieldPermitivity = ir.epsilon_rf;
320         calc_rffac(fp,
321                    interactionConst.epsilon_r,
322                    interactionConst.reactionFieldPermitivity,
323                    interactionConst.rcoulomb,
324                    &interactionConst.reactionFieldCoefficient,
325                    &interactionConst.reactionFieldShift);
326     }
327     else
328     {
329         /* For plain cut-off we might use the reaction-field kernels */
330         interactionConst.reactionFieldPermitivity = interactionConst.epsilon_r;
331         interactionConst.reactionFieldCoefficient = 0;
332         if (ir.coulomb_modifier == InteractionModifiers::PotShift)
333         {
334             interactionConst.reactionFieldShift = 1 / interactionConst.rcoulomb;
335         }
336         else
337         {
338             interactionConst.reactionFieldShift = 0;
339         }
340     }
341
342     initCoulombEwaldParameters(fp, ir, systemHasNetCharge, &interactionConst);
343
344     if (fp != nullptr)
345     {
346         real dispersion_shift;
347
348         dispersion_shift = interactionConst.dispersion_shift.cpot;
349         if (EVDW_PME(interactionConst.vdwtype))
350         {
351             dispersion_shift -= interactionConst.sh_lj_ewald;
352         }
353         fprintf(fp,
354                 "Potential shift: LJ r^-12: %.3e r^-6: %.3e",
355                 interactionConst.repulsion_shift.cpot,
356                 dispersion_shift);
357
358         if (interactionConst.eeltype == CoulombInteractionType::Cut)
359         {
360             fprintf(fp, ", Coulomb %.e", -interactionConst.reactionFieldShift);
361         }
362         else if (EEL_PME(interactionConst.eeltype))
363         {
364             fprintf(fp, ", Ewald %.3e", -interactionConst.sh_ewald);
365         }
366         fprintf(fp, "\n");
367     }
368
369     if (ir.efep != FreeEnergyPerturbationType::No)
370     {
371         GMX_RELEASE_ASSERT(ir.fepvals, "ir.fepvals should be set wth free-energy");
372         interactionConst.softCoreParameters =
373                 std::make_unique<interaction_const_t::SoftCoreParameters>(*ir.fepvals);
374     }
375
376     return interactionConst;
377 }