Merge branch release-2020 into master
[alexxy/gromacs.git] / src / gromacs / mdlib / forcerec.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-2019,2020, 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 "forcerec.h"
40
41 #include <cassert>
42 #include <cmath>
43 #include <cstdlib>
44 #include <cstring>
45
46 #include <algorithm>
47 #include <memory>
48
49 #include "gromacs/commandline/filenm.h"
50 #include "gromacs/domdec/domdec.h"
51 #include "gromacs/domdec/domdec_struct.h"
52 #include "gromacs/ewald/ewald.h"
53 #include "gromacs/ewald/ewald_utils.h"
54 #include "gromacs/ewald/pme_pp_comm_gpu.h"
55 #include "gromacs/fileio/filetypes.h"
56 #include "gromacs/gmxlib/network.h"
57 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
58 #include "gromacs/gpu_utils/gpu_utils.h"
59 #include "gromacs/hardware/hw_info.h"
60 #include "gromacs/listed_forces/gpubonded.h"
61 #include "gromacs/listed_forces/manage_threading.h"
62 #include "gromacs/listed_forces/pairs.h"
63 #include "gromacs/math/functions.h"
64 #include "gromacs/math/units.h"
65 #include "gromacs/math/utilities.h"
66 #include "gromacs/math/vec.h"
67 #include "gromacs/mdlib/dispersioncorrection.h"
68 #include "gromacs/mdlib/force.h"
69 #include "gromacs/mdlib/forcerec_threading.h"
70 #include "gromacs/mdlib/gmx_omp_nthreads.h"
71 #include "gromacs/mdlib/md_support.h"
72 #include "gromacs/mdlib/qmmm.h"
73 #include "gromacs/mdlib/rf_util.h"
74 #include "gromacs/mdlib/wall.h"
75 #include "gromacs/mdtypes/commrec.h"
76 #include "gromacs/mdtypes/fcdata.h"
77 #include "gromacs/mdtypes/forcerec.h"
78 #include "gromacs/mdtypes/group.h"
79 #include "gromacs/mdtypes/iforceprovider.h"
80 #include "gromacs/mdtypes/inputrec.h"
81 #include "gromacs/mdtypes/interaction_const.h"
82 #include "gromacs/mdtypes/md_enums.h"
83 #include "gromacs/nbnxm/gpu_data_mgmt.h"
84 #include "gromacs/nbnxm/nbnxm.h"
85 #include "gromacs/nbnxm/nbnxm_geometry.h"
86 #include "gromacs/pbcutil/ishift.h"
87 #include "gromacs/pbcutil/pbc.h"
88 #include "gromacs/tables/forcetable.h"
89 #include "gromacs/topology/mtop_util.h"
90 #include "gromacs/trajectory/trajectoryframe.h"
91 #include "gromacs/utility/cstringutil.h"
92 #include "gromacs/utility/exceptions.h"
93 #include "gromacs/utility/fatalerror.h"
94 #include "gromacs/utility/gmxassert.h"
95 #include "gromacs/utility/logger.h"
96 #include "gromacs/utility/physicalnodecommunicator.h"
97 #include "gromacs/utility/pleasecite.h"
98 #include "gromacs/utility/smalloc.h"
99 #include "gromacs/utility/strconvert.h"
100
101 static std::vector<real> mk_nbfp(const gmx_ffparams_t* idef, gmx_bool bBHAM)
102 {
103     std::vector<real> nbfp;
104     int               atnr;
105
106     atnr = idef->atnr;
107     if (bBHAM)
108     {
109         nbfp.resize(3 * atnr * atnr);
110         int k = 0;
111         for (int i = 0; (i < atnr); i++)
112         {
113             for (int j = 0; (j < atnr); j++, k++)
114             {
115                 BHAMA(nbfp, atnr, i, j) = idef->iparams[k].bham.a;
116                 BHAMB(nbfp, atnr, i, j) = idef->iparams[k].bham.b;
117                 /* nbfp now includes the 6.0 derivative prefactor */
118                 BHAMC(nbfp, atnr, i, j) = idef->iparams[k].bham.c * 6.0;
119             }
120         }
121     }
122     else
123     {
124         nbfp.resize(2 * atnr * atnr);
125         int k = 0;
126         for (int i = 0; (i < atnr); i++)
127         {
128             for (int j = 0; (j < atnr); j++, k++)
129             {
130                 /* nbfp now includes the 6.0/12.0 derivative prefactors */
131                 C6(nbfp, atnr, i, j)  = idef->iparams[k].lj.c6 * 6.0;
132                 C12(nbfp, atnr, i, j) = idef->iparams[k].lj.c12 * 12.0;
133             }
134         }
135     }
136
137     return nbfp;
138 }
139
140 static real* make_ljpme_c6grid(const gmx_ffparams_t* idef, t_forcerec* fr)
141 {
142     int   i, j, k, atnr;
143     real  c6, c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
144     real* grid;
145
146     /* For LJ-PME simulations, we correct the energies with the reciprocal space
147      * inside of the cut-off. To do this the non-bonded kernels needs to have
148      * access to the C6-values used on the reciprocal grid in pme.c
149      */
150
151     atnr = idef->atnr;
152     snew(grid, 2 * atnr * atnr);
153     for (i = k = 0; (i < atnr); i++)
154     {
155         for (j = 0; (j < atnr); j++, k++)
156         {
157             c6i  = idef->iparams[i * (atnr + 1)].lj.c6;
158             c12i = idef->iparams[i * (atnr + 1)].lj.c12;
159             c6j  = idef->iparams[j * (atnr + 1)].lj.c6;
160             c12j = idef->iparams[j * (atnr + 1)].lj.c12;
161             c6   = std::sqrt(c6i * c6j);
162             if (fr->ljpme_combination_rule == eljpmeLB && !gmx_numzero(c6) && !gmx_numzero(c12i)
163                 && !gmx_numzero(c12j))
164             {
165                 sigmai = gmx::sixthroot(c12i / c6i);
166                 sigmaj = gmx::sixthroot(c12j / c6j);
167                 epsi   = c6i * c6i / c12i;
168                 epsj   = c6j * c6j / c12j;
169                 c6     = std::sqrt(epsi * epsj) * gmx::power6(0.5 * (sigmai + sigmaj));
170             }
171             /* Store the elements at the same relative positions as C6 in nbfp in order
172              * to simplify access in the kernels
173              */
174             grid[2 * (atnr * i + j)] = c6 * 6.0;
175         }
176     }
177     return grid;
178 }
179
180 enum
181 {
182     acNONE = 0,
183     acCONSTRAINT,
184     acSETTLE
185 };
186
187 static std::vector<cginfo_mb_t> init_cginfo_mb(const gmx_mtop_t* mtop, const t_forcerec* fr)
188 {
189     gmx_bool* type_VDW;
190     int*      a_con;
191
192     snew(type_VDW, fr->ntype);
193     for (int ai = 0; ai < fr->ntype; ai++)
194     {
195         type_VDW[ai] = FALSE;
196         for (int j = 0; j < fr->ntype; j++)
197         {
198             type_VDW[ai] = type_VDW[ai] || fr->bBHAM || C6(fr->nbfp, fr->ntype, ai, j) != 0
199                            || C12(fr->nbfp, fr->ntype, ai, j) != 0;
200         }
201     }
202
203     std::vector<cginfo_mb_t> cginfoPerMolblock;
204     int                      a_offset = 0;
205     for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
206     {
207         const gmx_molblock_t& molb = mtop->molblock[mb];
208         const gmx_moltype_t&  molt = mtop->moltype[molb.type];
209         const auto&           excl = molt.excls;
210
211         /* Check if the cginfo is identical for all molecules in this block.
212          * If so, we only need an array of the size of one molecule.
213          * Otherwise we make an array of #mol times #cgs per molecule.
214          */
215         gmx_bool bId = TRUE;
216         for (int m = 0; m < molb.nmol; m++)
217         {
218             const int am = m * molt.atoms.nr;
219             for (int a = 0; a < molt.atoms.nr; a++)
220             {
221                 if (getGroupType(mtop->groups, SimulationAtomGroupType::QuantumMechanics, a_offset + am + a)
222                     != getGroupType(mtop->groups, SimulationAtomGroupType::QuantumMechanics, a_offset + a))
223                 {
224                     bId = FALSE;
225                 }
226                 if (!mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics].empty())
227                 {
228                     if (mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics][a_offset + am + a]
229                         != mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics][a_offset + a])
230                     {
231                         bId = FALSE;
232                     }
233                 }
234             }
235         }
236
237         cginfo_mb_t cginfo_mb;
238         cginfo_mb.cg_start = a_offset;
239         cginfo_mb.cg_end   = a_offset + molb.nmol * molt.atoms.nr;
240         cginfo_mb.cg_mod   = (bId ? 1 : molb.nmol) * molt.atoms.nr;
241         cginfo_mb.cginfo.resize(cginfo_mb.cg_mod);
242         gmx::ArrayRef<int> cginfo = cginfo_mb.cginfo;
243
244         /* Set constraints flags for constrained atoms */
245         snew(a_con, molt.atoms.nr);
246         for (int ftype = 0; ftype < F_NRE; ftype++)
247         {
248             if (interaction_function[ftype].flags & IF_CONSTRAINT)
249             {
250                 const int nral = NRAL(ftype);
251                 for (int ia = 0; ia < molt.ilist[ftype].size(); ia += 1 + nral)
252                 {
253                     int a;
254
255                     for (a = 0; a < nral; a++)
256                     {
257                         a_con[molt.ilist[ftype].iatoms[ia + 1 + a]] =
258                                 (ftype == F_SETTLE ? acSETTLE : acCONSTRAINT);
259                     }
260                 }
261             }
262         }
263
264         for (int m = 0; m < (bId ? 1 : molb.nmol); m++)
265         {
266             const int molculeOffsetInBlock = m * molt.atoms.nr;
267             for (int a = 0; a < molt.atoms.nr; a++)
268             {
269                 const t_atom& atom     = molt.atoms.atom[a];
270                 int&          atomInfo = cginfo[molculeOffsetInBlock + a];
271
272                 /* Store the energy group in cginfo */
273                 int gid = getGroupType(mtop->groups, SimulationAtomGroupType::EnergyOutput,
274                                        a_offset + molculeOffsetInBlock + a);
275                 SET_CGINFO_GID(atomInfo, gid);
276
277                 bool bHaveVDW = (type_VDW[atom.type] || type_VDW[atom.typeB]);
278                 bool bHaveQ   = (atom.q != 0 || atom.qB != 0);
279
280                 bool haveExclusions = false;
281                 /* Loop over all the exclusions of atom ai */
282                 for (const int j : excl[a])
283                 {
284                     if (j != a)
285                     {
286                         haveExclusions = true;
287                         break;
288                     }
289                 }
290
291                 switch (a_con[a])
292                 {
293                     case acCONSTRAINT: SET_CGINFO_CONSTR(atomInfo); break;
294                     case acSETTLE: SET_CGINFO_SETTLE(atomInfo); break;
295                     default: break;
296                 }
297                 if (haveExclusions)
298                 {
299                     SET_CGINFO_EXCL_INTER(atomInfo);
300                 }
301                 if (bHaveVDW)
302                 {
303                     SET_CGINFO_HAS_VDW(atomInfo);
304                 }
305                 if (bHaveQ)
306                 {
307                     SET_CGINFO_HAS_Q(atomInfo);
308                 }
309                 if (fr->efep != efepNO && PERTURBED(atom))
310                 {
311                     SET_CGINFO_FEP(atomInfo);
312                 }
313             }
314         }
315
316         sfree(a_con);
317
318         cginfoPerMolblock.push_back(cginfo_mb);
319
320         a_offset += molb.nmol * molt.atoms.nr;
321     }
322     sfree(type_VDW);
323
324     return cginfoPerMolblock;
325 }
326
327 static std::vector<int> cginfo_expand(const int nmb, gmx::ArrayRef<const cginfo_mb_t> cgi_mb)
328 {
329     const int ncg = cgi_mb[nmb - 1].cg_end;
330
331     std::vector<int> cginfo(ncg);
332
333     int mb = 0;
334     for (int cg = 0; cg < ncg; cg++)
335     {
336         while (cg >= cgi_mb[mb].cg_end)
337         {
338             mb++;
339         }
340         cginfo[cg] = cgi_mb[mb].cginfo[(cg - cgi_mb[mb].cg_start) % cgi_mb[mb].cg_mod];
341     }
342
343     return cginfo;
344 }
345
346 /* Sets the sum of charges (squared) and C6 in the system in fr.
347  * Returns whether the system has a net charge.
348  */
349 static bool set_chargesum(FILE* log, t_forcerec* fr, const gmx_mtop_t* mtop)
350 {
351     /*This now calculates sum for q and c6*/
352     double qsum, q2sum, q, c6sum, c6;
353
354     qsum  = 0;
355     q2sum = 0;
356     c6sum = 0;
357     for (const gmx_molblock_t& molb : mtop->molblock)
358     {
359         int            nmol  = molb.nmol;
360         const t_atoms* atoms = &mtop->moltype[molb.type].atoms;
361         for (int i = 0; i < atoms->nr; i++)
362         {
363             q = atoms->atom[i].q;
364             qsum += nmol * q;
365             q2sum += nmol * q * q;
366             c6 = mtop->ffparams.iparams[atoms->atom[i].type * (mtop->ffparams.atnr + 1)].lj.c6;
367             c6sum += nmol * c6;
368         }
369     }
370     fr->qsum[0]  = qsum;
371     fr->q2sum[0] = q2sum;
372     fr->c6sum[0] = c6sum;
373
374     if (fr->efep != efepNO)
375     {
376         qsum  = 0;
377         q2sum = 0;
378         c6sum = 0;
379         for (const gmx_molblock_t& molb : mtop->molblock)
380         {
381             int            nmol  = molb.nmol;
382             const t_atoms* atoms = &mtop->moltype[molb.type].atoms;
383             for (int i = 0; i < atoms->nr; i++)
384             {
385                 q = atoms->atom[i].qB;
386                 qsum += nmol * q;
387                 q2sum += nmol * q * q;
388                 c6 = mtop->ffparams.iparams[atoms->atom[i].typeB * (mtop->ffparams.atnr + 1)].lj.c6;
389                 c6sum += nmol * c6;
390             }
391             fr->qsum[1]  = qsum;
392             fr->q2sum[1] = q2sum;
393             fr->c6sum[1] = c6sum;
394         }
395     }
396     else
397     {
398         fr->qsum[1]  = fr->qsum[0];
399         fr->q2sum[1] = fr->q2sum[0];
400         fr->c6sum[1] = fr->c6sum[0];
401     }
402     if (log)
403     {
404         if (fr->efep == efepNO)
405         {
406             fprintf(log, "System total charge: %.3f\n", fr->qsum[0]);
407         }
408         else
409         {
410             fprintf(log, "System total charge, top. A: %.3f top. B: %.3f\n", fr->qsum[0], fr->qsum[1]);
411         }
412     }
413
414     /* A cut-off of 1e-4 is used to catch rounding errors due to ascii input */
415     return (std::abs(fr->qsum[0]) > 1e-4 || std::abs(fr->qsum[1]) > 1e-4);
416 }
417
418 static real calcBuckinghamBMax(FILE* fplog, const gmx_mtop_t* mtop)
419 {
420     const t_atoms *at1, *at2;
421     int            i, j, tpi, tpj, ntypes;
422     real           b, bmin;
423
424     if (fplog)
425     {
426         fprintf(fplog, "Determining largest Buckingham b parameter for table\n");
427     }
428     ntypes = mtop->ffparams.atnr;
429
430     bmin            = -1;
431     real bham_b_max = 0;
432     for (size_t mt1 = 0; mt1 < mtop->moltype.size(); mt1++)
433     {
434         at1 = &mtop->moltype[mt1].atoms;
435         for (i = 0; (i < at1->nr); i++)
436         {
437             tpi = at1->atom[i].type;
438             if (tpi >= ntypes)
439             {
440                 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", i, tpi, ntypes);
441             }
442
443             for (size_t mt2 = mt1; mt2 < mtop->moltype.size(); mt2++)
444             {
445                 at2 = &mtop->moltype[mt2].atoms;
446                 for (j = 0; (j < at2->nr); j++)
447                 {
448                     tpj = at2->atom[j].type;
449                     if (tpj >= ntypes)
450                     {
451                         gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", j, tpj, ntypes);
452                     }
453                     b = mtop->ffparams.iparams[tpi * ntypes + tpj].bham.b;
454                     if (b > bham_b_max)
455                     {
456                         bham_b_max = b;
457                     }
458                     if ((b < bmin) || (bmin == -1))
459                     {
460                         bmin = b;
461                     }
462                 }
463             }
464         }
465     }
466     if (fplog)
467     {
468         fprintf(fplog, "Buckingham b parameters, min: %g, max: %g\n", bmin, bham_b_max);
469     }
470
471     return bham_b_max;
472 }
473
474 /*!\brief If there's bonded interactions of type \c ftype1 or \c
475  * ftype2 present in the topology, build an array of the number of
476  * interactions present for each bonded interaction index found in the
477  * topology.
478  *
479  * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
480  * valid type with that parameter.
481  *
482  * \c count will be reallocated as necessary to fit the largest bonded
483  * interaction index found, and its current size will be returned in
484  * \c ncount. It will contain zero for every bonded interaction index
485  * for which no interactions are present in the topology.
486  */
487 static void count_tables(int ftype1, int ftype2, const gmx_mtop_t* mtop, int* ncount, int** count)
488 {
489     int ftype, i, j, tabnr;
490
491     // Loop over all moleculetypes
492     for (const gmx_moltype_t& molt : mtop->moltype)
493     {
494         // Loop over all interaction types
495         for (ftype = 0; ftype < F_NRE; ftype++)
496         {
497             // If the current interaction type is one of the types whose tables we're trying to count...
498             if (ftype == ftype1 || ftype == ftype2)
499             {
500                 const InteractionList& il     = molt.ilist[ftype];
501                 const int              stride = 1 + NRAL(ftype);
502                 // ... and there are actually some interactions for this type
503                 for (i = 0; i < il.size(); i += stride)
504                 {
505                     // Find out which table index the user wanted
506                     tabnr = mtop->ffparams.iparams[il.iatoms[i]].tab.table;
507                     if (tabnr < 0)
508                     {
509                         gmx_fatal(FARGS, "A bonded table number is smaller than 0: %d\n", tabnr);
510                     }
511                     // Make room for this index in the data structure
512                     if (tabnr >= *ncount)
513                     {
514                         srenew(*count, tabnr + 1);
515                         for (j = *ncount; j < tabnr + 1; j++)
516                         {
517                             (*count)[j] = 0;
518                         }
519                         *ncount = tabnr + 1;
520                     }
521                     // Record that this table index is used and must have a valid file
522                     (*count)[tabnr]++;
523                 }
524             }
525         }
526     }
527 }
528
529 /*!\brief If there's bonded interactions of flavour \c tabext and type
530  * \c ftype1 or \c ftype2 present in the topology, seek them in the
531  * list of filenames passed to mdrun, and make bonded tables from
532  * those files.
533  *
534  * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
535  * valid type with that parameter.
536  *
537  * A fatal error occurs if no matching filename is found.
538  */
539 static bondedtable_t* make_bonded_tables(FILE*                            fplog,
540                                          int                              ftype1,
541                                          int                              ftype2,
542                                          const gmx_mtop_t*                mtop,
543                                          gmx::ArrayRef<const std::string> tabbfnm,
544                                          const char*                      tabext)
545 {
546     int            ncount, *count;
547     bondedtable_t* tab;
548
549     tab = nullptr;
550
551     ncount = 0;
552     count  = nullptr;
553     count_tables(ftype1, ftype2, mtop, &ncount, &count);
554
555     // Are there any relevant tabulated bond interactions?
556     if (ncount > 0)
557     {
558         snew(tab, ncount);
559         for (int i = 0; i < ncount; i++)
560         {
561             // Do any interactions exist that requires this table?
562             if (count[i] > 0)
563             {
564                 // This pattern enforces the current requirement that
565                 // table filenames end in a characteristic sequence
566                 // before the file type extension, and avoids table 13
567                 // being recognized and used for table 1.
568                 std::string patternToFind = gmx::formatString("_%s%d.%s", tabext, i, ftp2ext(efXVG));
569                 bool        madeTable     = false;
570                 for (gmx::index j = 0; j < tabbfnm.ssize() && !madeTable; ++j)
571                 {
572                     if (gmx::endsWith(tabbfnm[j], patternToFind))
573                     {
574                         // Finally read the table from the file found
575                         tab[i]    = make_bonded_table(fplog, tabbfnm[j].c_str(), NRAL(ftype1) - 2);
576                         madeTable = true;
577                     }
578                 }
579                 if (!madeTable)
580                 {
581                     bool isPlural = (ftype2 != -1);
582                     gmx_fatal(FARGS,
583                               "Tabulated interaction of type '%s%s%s' with index %d cannot be used "
584                               "because no table file whose name matched '%s' was passed via the "
585                               "gmx mdrun -tableb command-line option.",
586                               interaction_function[ftype1].longname, isPlural ? "' or '" : "",
587                               isPlural ? interaction_function[ftype2].longname : "", i,
588                               patternToFind.c_str());
589                 }
590             }
591         }
592         sfree(count);
593     }
594
595     return tab;
596 }
597
598 void forcerec_set_ranges(t_forcerec* fr, int natoms_force, int natoms_force_constr, int natoms_f_novirsum)
599 {
600     fr->natoms_force        = natoms_force;
601     fr->natoms_force_constr = natoms_force_constr;
602
603     if (fr->natoms_force_constr > fr->nalloc_force)
604     {
605         fr->nalloc_force = over_alloc_dd(fr->natoms_force_constr);
606     }
607
608     if (fr->haveDirectVirialContributions)
609     {
610         fr->forceBufferForDirectVirialContributions.resize(natoms_f_novirsum);
611     }
612 }
613
614 static real cutoff_inf(real cutoff)
615 {
616     if (cutoff == 0)
617     {
618         cutoff = GMX_CUTOFF_INF;
619     }
620
621     return cutoff;
622 }
623
624 /*! \brief Print Coulomb Ewald citations and set ewald coefficients */
625 static void initCoulombEwaldParameters(FILE*                fp,
626                                        const t_inputrec*    ir,
627                                        bool                 systemHasNetCharge,
628                                        interaction_const_t* ic)
629 {
630     if (!EEL_PME_EWALD(ir->coulombtype))
631     {
632         return;
633     }
634
635     if (fp)
636     {
637         fprintf(fp, "Will do PME sum in reciprocal space for electrostatic interactions.\n");
638
639         if (ir->coulombtype == eelP3M_AD)
640         {
641             please_cite(fp, "Hockney1988");
642             please_cite(fp, "Ballenegger2012");
643         }
644         else
645         {
646             please_cite(fp, "Essmann95a");
647         }
648
649         if (ir->ewald_geometry == eewg3DC)
650         {
651             if (fp)
652             {
653                 fprintf(fp, "Using the Ewald3DC correction for systems with a slab geometry%s.\n",
654                         systemHasNetCharge ? " and net charge" : "");
655             }
656             please_cite(fp, "In-Chul99a");
657             if (systemHasNetCharge)
658             {
659                 please_cite(fp, "Ballenegger2009");
660             }
661         }
662     }
663
664     ic->ewaldcoeff_q = calc_ewaldcoeff_q(ir->rcoulomb, ir->ewald_rtol);
665     if (fp)
666     {
667         fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for Ewald\n", 1 / ic->ewaldcoeff_q);
668     }
669
670     if (ic->coulomb_modifier == eintmodPOTSHIFT)
671     {
672         GMX_RELEASE_ASSERT(ic->rcoulomb != 0, "Cutoff radius cannot be zero");
673         ic->sh_ewald = std::erfc(ic->ewaldcoeff_q * ic->rcoulomb) / ic->rcoulomb;
674     }
675     else
676     {
677         ic->sh_ewald = 0;
678     }
679 }
680
681 /*! \brief Print Van der Waals Ewald citations and set ewald coefficients */
682 static void initVdwEwaldParameters(FILE* fp, const t_inputrec* ir, interaction_const_t* ic)
683 {
684     if (!EVDW_PME(ir->vdwtype))
685     {
686         return;
687     }
688
689     if (fp)
690     {
691         fprintf(fp, "Will do PME sum in reciprocal space for LJ dispersion interactions.\n");
692         please_cite(fp, "Essmann95a");
693     }
694     ic->ewaldcoeff_lj = calc_ewaldcoeff_lj(ir->rvdw, ir->ewald_rtol_lj);
695     if (fp)
696     {
697         fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for LJ Ewald\n", 1 / ic->ewaldcoeff_lj);
698     }
699
700     if (ic->vdw_modifier == eintmodPOTSHIFT)
701     {
702         real crc2       = gmx::square(ic->ewaldcoeff_lj * ic->rvdw);
703         ic->sh_lj_ewald = (std::exp(-crc2) * (1 + crc2 + 0.5 * crc2 * crc2) - 1) / gmx::power6(ic->rvdw);
704     }
705     else
706     {
707         ic->sh_lj_ewald = 0;
708     }
709 }
710
711 /* Generate Coulomb and/or Van der Waals Ewald long-range correction tables
712  *
713  * Tables are generated for one or both, depending on if the pointers
714  * are non-null. The spacing for both table sets is the same and obeys
715  * both accuracy requirements, when relevant.
716  */
717 static void init_ewald_f_table(const interaction_const_t& ic,
718                                EwaldCorrectionTables*     coulombTables,
719                                EwaldCorrectionTables*     vdwTables)
720 {
721     const bool useCoulombTable = (EEL_PME_EWALD(ic.eeltype) && coulombTables != nullptr);
722     const bool useVdwTable     = (EVDW_PME(ic.vdwtype) && vdwTables != nullptr);
723
724     /* Get the Ewald table spacing based on Coulomb and/or LJ
725      * Ewald coefficients and rtol.
726      */
727     const real tableScale = ewald_spline3_table_scale(ic, useCoulombTable, useVdwTable);
728
729     const int tableSize = static_cast<int>(ic.rcoulomb * tableScale) + 2;
730
731     if (useCoulombTable)
732     {
733         *coulombTables =
734                 generateEwaldCorrectionTables(tableSize, tableScale, ic.ewaldcoeff_q, v_q_ewald_lr);
735     }
736
737     if (useVdwTable)
738     {
739         *vdwTables = generateEwaldCorrectionTables(tableSize, tableScale, ic.ewaldcoeff_lj, v_lj_ewald_lr);
740     }
741 }
742
743 void init_interaction_const_tables(FILE* fp, interaction_const_t* ic)
744 {
745     if (EEL_PME_EWALD(ic->eeltype))
746     {
747         init_ewald_f_table(*ic, ic->coulombEwaldTables.get(), nullptr);
748         if (fp != nullptr)
749         {
750             fprintf(fp, "Initialized non-bonded Coulomb Ewald tables, spacing: %.2e size: %zu\n\n",
751                     1 / ic->coulombEwaldTables->scale, ic->coulombEwaldTables->tableF.size());
752         }
753     }
754 }
755
756 static void clear_force_switch_constants(shift_consts_t* sc)
757 {
758     sc->c2   = 0;
759     sc->c3   = 0;
760     sc->cpot = 0;
761 }
762
763 static void force_switch_constants(real p, real rsw, real rc, shift_consts_t* sc)
764 {
765     /* Here we determine the coefficient for shifting the force to zero
766      * between distance rsw and the cut-off rc.
767      * For a potential of r^-p, we have force p*r^-(p+1).
768      * But to save flops we absorb p in the coefficient.
769      * Thus we get:
770      * force/p   = r^-(p+1) + c2*r^2 + c3*r^3
771      * potential = r^-p + c2/3*r^3 + c3/4*r^4 + cpot
772      */
773     sc->c2   = ((p + 1) * rsw - (p + 4) * rc) / (pow(rc, p + 2) * gmx::square(rc - rsw));
774     sc->c3   = -((p + 1) * rsw - (p + 3) * rc) / (pow(rc, p + 2) * gmx::power3(rc - rsw));
775     sc->cpot = -pow(rc, -p) + p * sc->c2 / 3 * gmx::power3(rc - rsw)
776                + p * sc->c3 / 4 * gmx::power4(rc - rsw);
777 }
778
779 static void potential_switch_constants(real rsw, real rc, switch_consts_t* sc)
780 {
781     /* The switch function is 1 at rsw and 0 at rc.
782      * The derivative and second derivate are zero at both ends.
783      * rsw        = max(r - r_switch, 0)
784      * sw         = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
785      * dsw        = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
786      * force      = force*dsw - potential*sw
787      * potential *= sw
788      */
789     sc->c3 = -10 / gmx::power3(rc - rsw);
790     sc->c4 = 15 / gmx::power4(rc - rsw);
791     sc->c5 = -6 / gmx::power5(rc - rsw);
792 }
793
794 /*! \brief Construct interaction constants
795  *
796  * This data is used (particularly) by search and force code for
797  * short-range interactions. Many of these are constant for the whole
798  * simulation; some are constant only after PME tuning completes.
799  */
800 static void init_interaction_const(FILE*                 fp,
801                                    interaction_const_t** interaction_const,
802                                    const t_inputrec*     ir,
803                                    const gmx_mtop_t*     mtop,
804                                    bool                  systemHasNetCharge)
805 {
806     interaction_const_t* ic = new interaction_const_t;
807
808     ic->cutoff_scheme = ir->cutoff_scheme;
809
810     ic->coulombEwaldTables = std::make_unique<EwaldCorrectionTables>();
811
812     /* Lennard-Jones */
813     ic->vdwtype         = ir->vdwtype;
814     ic->vdw_modifier    = ir->vdw_modifier;
815     ic->reppow          = mtop->ffparams.reppow;
816     ic->rvdw            = cutoff_inf(ir->rvdw);
817     ic->rvdw_switch     = ir->rvdw_switch;
818     ic->ljpme_comb_rule = ir->ljpme_combination_rule;
819     ic->useBuckingham   = (mtop->ffparams.functype[0] == F_BHAM);
820     if (ic->useBuckingham)
821     {
822         ic->buckinghamBMax = calcBuckinghamBMax(fp, mtop);
823     }
824
825     initVdwEwaldParameters(fp, ir, ic);
826
827     clear_force_switch_constants(&ic->dispersion_shift);
828     clear_force_switch_constants(&ic->repulsion_shift);
829
830     switch (ic->vdw_modifier)
831     {
832         case eintmodPOTSHIFT:
833             /* Only shift the potential, don't touch the force */
834             ic->dispersion_shift.cpot = -1.0 / gmx::power6(ic->rvdw);
835             ic->repulsion_shift.cpot  = -1.0 / gmx::power12(ic->rvdw);
836             break;
837         case eintmodFORCESWITCH:
838             /* Switch the force, switch and shift the potential */
839             force_switch_constants(6.0, ic->rvdw_switch, ic->rvdw, &ic->dispersion_shift);
840             force_switch_constants(12.0, ic->rvdw_switch, ic->rvdw, &ic->repulsion_shift);
841             break;
842         case eintmodPOTSWITCH:
843             /* Switch the potential and force */
844             potential_switch_constants(ic->rvdw_switch, ic->rvdw, &ic->vdw_switch);
845             break;
846         case eintmodNONE:
847         case eintmodEXACTCUTOFF:
848             /* Nothing to do here */
849             break;
850         default: gmx_incons("unimplemented potential modifier");
851     }
852
853     /* Electrostatics */
854     ic->eeltype          = ir->coulombtype;
855     ic->coulomb_modifier = ir->coulomb_modifier;
856     ic->rcoulomb         = cutoff_inf(ir->rcoulomb);
857     ic->rcoulomb_switch  = ir->rcoulomb_switch;
858     ic->epsilon_r        = ir->epsilon_r;
859
860     /* Set the Coulomb energy conversion factor */
861     if (ic->epsilon_r != 0)
862     {
863         ic->epsfac = ONE_4PI_EPS0 / ic->epsilon_r;
864     }
865     else
866     {
867         /* eps = 0 is infinite dieletric: no Coulomb interactions */
868         ic->epsfac = 0;
869     }
870
871     /* Reaction-field */
872     if (EEL_RF(ic->eeltype))
873     {
874         GMX_RELEASE_ASSERT(ic->eeltype != eelGRF_NOTUSED, "GRF is no longer supported");
875         ic->epsilon_rf = ir->epsilon_rf;
876
877         calc_rffac(fp, ic->epsilon_r, ic->epsilon_rf, ic->rcoulomb, &ic->k_rf, &ic->c_rf);
878     }
879     else
880     {
881         /* For plain cut-off we might use the reaction-field kernels */
882         ic->epsilon_rf = ic->epsilon_r;
883         ic->k_rf       = 0;
884         if (ir->coulomb_modifier == eintmodPOTSHIFT)
885         {
886             ic->c_rf = 1 / ic->rcoulomb;
887         }
888         else
889         {
890             ic->c_rf = 0;
891         }
892     }
893
894     initCoulombEwaldParameters(fp, ir, systemHasNetCharge, ic);
895
896     if (fp != nullptr)
897     {
898         real dispersion_shift;
899
900         dispersion_shift = ic->dispersion_shift.cpot;
901         if (EVDW_PME(ic->vdwtype))
902         {
903             dispersion_shift -= ic->sh_lj_ewald;
904         }
905         fprintf(fp, "Potential shift: LJ r^-12: %.3e r^-6: %.3e", ic->repulsion_shift.cpot, dispersion_shift);
906
907         if (ic->eeltype == eelCUT)
908         {
909             fprintf(fp, ", Coulomb %.e", -ic->c_rf);
910         }
911         else if (EEL_PME(ic->eeltype))
912         {
913             fprintf(fp, ", Ewald %.3e", -ic->sh_ewald);
914         }
915         fprintf(fp, "\n");
916     }
917
918     *interaction_const = ic;
919 }
920
921 bool areMoleculesDistributedOverPbc(const t_inputrec& ir, const gmx_mtop_t& mtop, const gmx::MDLogger& mdlog)
922 {
923     bool       areMoleculesDistributedOverPbc = false;
924     const bool useEwaldSurfaceCorrection = (EEL_PME_EWALD(ir.coulombtype) && ir.epsilon_surface != 0);
925
926     const bool bSHAKE =
927             (ir.eConstrAlg == econtSHAKE
928              && (gmx_mtop_ftype_count(mtop, F_CONSTR) > 0 || gmx_mtop_ftype_count(mtop, F_CONSTRNC) > 0));
929
930     /* The group cut-off scheme and SHAKE assume charge groups
931      * are whole, but not using molpbc is faster in most cases.
932      * With intermolecular interactions we need PBC for calculating
933      * distances between atoms in different molecules.
934      */
935     if (bSHAKE && !mtop.bIntermolecularInteractions)
936     {
937         areMoleculesDistributedOverPbc = ir.bPeriodicMols;
938
939         if (areMoleculesDistributedOverPbc)
940         {
941             gmx_fatal(FARGS, "SHAKE is not supported with periodic molecules");
942         }
943     }
944     else
945     {
946         /* Not making molecules whole is faster in most cases,
947          * but with orientation restraints or non-tinfoil boundary
948          * conditions we need whole molecules.
949          */
950         areMoleculesDistributedOverPbc =
951                 (gmx_mtop_ftype_count(mtop, F_ORIRES) == 0 && !useEwaldSurfaceCorrection);
952
953         if (getenv("GMX_USE_GRAPH") != nullptr)
954         {
955             areMoleculesDistributedOverPbc = false;
956
957             GMX_LOG(mdlog.warning)
958                     .asParagraph()
959                     .appendText(
960                             "GMX_USE_GRAPH is set, using the graph for bonded "
961                             "interactions");
962
963             if (mtop.bIntermolecularInteractions)
964             {
965                 GMX_LOG(mdlog.warning)
966                         .asParagraph()
967                         .appendText(
968                                 "WARNING: Molecules linked by intermolecular interactions "
969                                 "have to reside in the same periodic image, otherwise "
970                                 "artifacts will occur!");
971             }
972         }
973
974         GMX_RELEASE_ASSERT(areMoleculesDistributedOverPbc || !mtop.bIntermolecularInteractions,
975                            "We need to use PBC within molecules with inter-molecular interactions");
976
977         if (bSHAKE && areMoleculesDistributedOverPbc)
978         {
979             gmx_fatal(FARGS,
980                       "SHAKE is not properly supported with intermolecular interactions. "
981                       "For short simulations where linked molecules remain in the same "
982                       "periodic image, the environment variable GMX_USE_GRAPH can be used "
983                       "to override this check.\n");
984         }
985     }
986
987     return areMoleculesDistributedOverPbc;
988 }
989
990 void init_forcerec(FILE*                            fp,
991                    const gmx::MDLogger&             mdlog,
992                    t_forcerec*                      fr,
993                    t_fcdata*                        fcd,
994                    const t_inputrec*                ir,
995                    const gmx_mtop_t*                mtop,
996                    const t_commrec*                 cr,
997                    matrix                           box,
998                    const char*                      tabfn,
999                    const char*                      tabpfn,
1000                    gmx::ArrayRef<const std::string> tabbfnm,
1001                    real                             print_force)
1002 {
1003     /* By default we turn SIMD kernels on, but it might be turned off further down... */
1004     fr->use_simd_kernels = TRUE;
1005
1006     if (check_box(ir->pbcType, box))
1007     {
1008         gmx_fatal(FARGS, "%s", check_box(ir->pbcType, box));
1009     }
1010
1011     /* Test particle insertion ? */
1012     if (EI_TPI(ir->eI))
1013     {
1014         /* Set to the size of the molecule to be inserted (the last one) */
1015         gmx::RangePartitioning molecules = gmx_mtop_molecules(*mtop);
1016         fr->n_tpi                        = molecules.block(molecules.numBlocks() - 1).size();
1017     }
1018     else
1019     {
1020         fr->n_tpi = 0;
1021     }
1022
1023     if (ir->coulombtype == eelRF_NEC_UNSUPPORTED || ir->coulombtype == eelGRF_NOTUSED)
1024     {
1025         gmx_fatal(FARGS, "%s electrostatics is no longer supported", eel_names[ir->coulombtype]);
1026     }
1027
1028     if (ir->bAdress)
1029     {
1030         gmx_fatal(FARGS, "AdResS simulations are no longer supported");
1031     }
1032     if (ir->useTwinRange)
1033     {
1034         gmx_fatal(FARGS, "Twin-range simulations are no longer supported");
1035     }
1036     /* Copy the user determined parameters */
1037     fr->userint1  = ir->userint1;
1038     fr->userint2  = ir->userint2;
1039     fr->userint3  = ir->userint3;
1040     fr->userint4  = ir->userint4;
1041     fr->userreal1 = ir->userreal1;
1042     fr->userreal2 = ir->userreal2;
1043     fr->userreal3 = ir->userreal3;
1044     fr->userreal4 = ir->userreal4;
1045
1046     /* Shell stuff */
1047     fr->fc_stepsize = ir->fc_stepsize;
1048
1049     /* Free energy */
1050     fr->efep        = ir->efep;
1051     fr->sc_alphavdw = ir->fepvals->sc_alpha;
1052     if (ir->fepvals->bScCoul)
1053     {
1054         fr->sc_alphacoul  = ir->fepvals->sc_alpha;
1055         fr->sc_sigma6_min = gmx::power6(ir->fepvals->sc_sigma_min);
1056     }
1057     else
1058     {
1059         fr->sc_alphacoul  = 0;
1060         fr->sc_sigma6_min = 0; /* only needed when bScCoul is on */
1061     }
1062     fr->sc_power      = ir->fepvals->sc_power;
1063     fr->sc_r_power    = ir->fepvals->sc_r_power;
1064     fr->sc_sigma6_def = gmx::power6(ir->fepvals->sc_sigma);
1065
1066     char* env = getenv("GMX_SCSIGMA_MIN");
1067     if (env != nullptr)
1068     {
1069         double dbl = 0;
1070         sscanf(env, "%20lf", &dbl);
1071         fr->sc_sigma6_min = gmx::power6(dbl);
1072         if (fp)
1073         {
1074             fprintf(fp, "Setting the minimum soft core sigma to %g nm\n", dbl);
1075         }
1076     }
1077
1078     fr->bNonbonded = TRUE;
1079     if (getenv("GMX_NO_NONBONDED") != nullptr)
1080     {
1081         /* turn off non-bonded calculations */
1082         fr->bNonbonded = FALSE;
1083         GMX_LOG(mdlog.warning)
1084                 .asParagraph()
1085                 .appendText(
1086                         "Found environment variable GMX_NO_NONBONDED.\n"
1087                         "Disabling nonbonded calculations.");
1088     }
1089
1090     if ((getenv("GMX_DISABLE_SIMD_KERNELS") != nullptr) || (getenv("GMX_NOOPTIMIZEDKERNELS") != nullptr))
1091     {
1092         fr->use_simd_kernels = FALSE;
1093         if (fp != nullptr)
1094         {
1095             fprintf(fp,
1096                     "\nFound environment variable GMX_DISABLE_SIMD_KERNELS.\n"
1097                     "Disabling the usage of any SIMD-specific non-bonded & bonded kernel routines\n"
1098                     "(e.g. SSE2/SSE4.1/AVX).\n\n");
1099         }
1100     }
1101
1102     fr->bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
1103
1104     /* Neighbour searching stuff */
1105     fr->cutoff_scheme = ir->cutoff_scheme;
1106     fr->pbcType       = ir->pbcType;
1107
1108     /* Determine if we will do PBC for distances in bonded interactions */
1109     if (fr->pbcType == PbcType::No)
1110     {
1111         fr->bMolPBC = FALSE;
1112     }
1113     else
1114     {
1115         const bool useEwaldSurfaceCorrection =
1116                 (EEL_PME_EWALD(ir->coulombtype) && ir->epsilon_surface != 0);
1117         if (!DOMAINDECOMP(cr))
1118         {
1119             fr->bMolPBC = areMoleculesDistributedOverPbc(*ir, *mtop, mdlog);
1120             // The assert below is equivalent to fcd->orires.nr != gmx_mtop_ftype_count(mtop, F_ORIRES)
1121             GMX_RELEASE_ASSERT(!fr->bMolPBC || fcd->orires.nr == 0,
1122                                "Molecules broken over PBC exist in a simulation including "
1123                                "orientation restraints. "
1124                                "This likely means that the global topology and the force constant "
1125                                "data have gotten out of sync.");
1126         }
1127         else
1128         {
1129             fr->bMolPBC = dd_bonded_molpbc(cr->dd, fr->pbcType);
1130
1131             if (useEwaldSurfaceCorrection && !dd_moleculesAreAlwaysWhole(*cr->dd))
1132             {
1133                 gmx_fatal(FARGS,
1134                           "You requested dipole correction (epsilon_surface > 0), but molecules "
1135                           "are broken "
1136                           "over periodic boundary conditions by the domain decomposition. "
1137                           "Run without domain decomposition instead.");
1138             }
1139         }
1140
1141         if (useEwaldSurfaceCorrection)
1142         {
1143             GMX_RELEASE_ASSERT((!DOMAINDECOMP(cr) && !fr->bMolPBC)
1144                                        || (DOMAINDECOMP(cr) && dd_moleculesAreAlwaysWhole(*cr->dd)),
1145                                "Molecules can not be broken by PBC with epsilon_surface > 0");
1146         }
1147     }
1148
1149     fr->rc_scaling = ir->refcoord_scaling;
1150     copy_rvec(ir->posres_com, fr->posres_com);
1151     copy_rvec(ir->posres_comB, fr->posres_comB);
1152     fr->rlist                  = cutoff_inf(ir->rlist);
1153     fr->ljpme_combination_rule = ir->ljpme_combination_rule;
1154
1155     /* This now calculates sum for q and c6*/
1156     bool systemHasNetCharge = set_chargesum(fp, fr, mtop);
1157
1158     /* fr->ic is used both by verlet and group kernels (to some extent) now */
1159     init_interaction_const(fp, &fr->ic, ir, mtop, systemHasNetCharge);
1160     init_interaction_const_tables(fp, fr->ic);
1161
1162     const interaction_const_t* ic = fr->ic;
1163
1164     /* TODO: Replace this Ewald table or move it into interaction_const_t */
1165     if (ir->coulombtype == eelEWALD)
1166     {
1167         init_ewald_tab(&(fr->ewald_table), ir, fp);
1168     }
1169
1170     /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
1171     switch (ic->eeltype)
1172     {
1173         case eelCUT: fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_COULOMB; break;
1174
1175         case eelRF:
1176         case eelRF_ZERO: fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD; break;
1177
1178         case eelSWITCH:
1179         case eelSHIFT:
1180         case eelUSER:
1181         case eelENCADSHIFT:
1182         case eelPMESWITCH:
1183         case eelPMEUSER:
1184         case eelPMEUSERSWITCH:
1185             fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
1186             break;
1187
1188         case eelPME:
1189         case eelP3M_AD:
1190         case eelEWALD: fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_EWALD; break;
1191
1192         default:
1193             gmx_fatal(FARGS, "Unsupported electrostatic interaction: %s", eel_names[ic->eeltype]);
1194     }
1195     fr->nbkernel_elec_modifier = ic->coulomb_modifier;
1196
1197     /* Vdw: Translate from mdp settings to kernel format */
1198     switch (ic->vdwtype)
1199     {
1200         case evdwCUT:
1201             if (fr->bBHAM)
1202             {
1203                 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_BUCKINGHAM;
1204             }
1205             else
1206             {
1207                 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LENNARDJONES;
1208             }
1209             break;
1210         case evdwPME: fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LJEWALD; break;
1211
1212         case evdwSWITCH:
1213         case evdwSHIFT:
1214         case evdwUSER:
1215         case evdwENCADSHIFT:
1216             fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
1217             break;
1218
1219         default: gmx_fatal(FARGS, "Unsupported vdw interaction: %s", evdw_names[ic->vdwtype]);
1220     }
1221     fr->nbkernel_vdw_modifier = ic->vdw_modifier;
1222
1223     if (ir->cutoff_scheme == ecutsVERLET)
1224     {
1225         if (!gmx_within_tol(ic->reppow, 12.0, 10 * GMX_DOUBLE_EPS))
1226         {
1227             gmx_fatal(FARGS, "Cut-off scheme %s only supports LJ repulsion power 12",
1228                       ecutscheme_names[ir->cutoff_scheme]);
1229         }
1230         /* Older tpr files can contain Coulomb user tables with the Verlet cutoff-scheme,
1231          * while mdrun does not (and never did) support this.
1232          */
1233         if (EEL_USER(fr->ic->eeltype))
1234         {
1235             gmx_fatal(FARGS, "Combination of %s and cutoff scheme %s is not supported",
1236                       eel_names[ir->coulombtype], ecutscheme_names[ir->cutoff_scheme]);
1237         }
1238
1239         fr->bvdwtab  = FALSE;
1240         fr->bcoultab = FALSE;
1241     }
1242
1243     /* 1-4 interaction electrostatics */
1244     fr->fudgeQQ = mtop->ffparams.fudgeQQ;
1245
1246     fr->haveDirectVirialContributions =
1247             (EEL_FULL(ic->eeltype) || EVDW_PME(ic->vdwtype) || fr->forceProviders->hasForceProvider()
1248              || gmx_mtop_ftype_count(mtop, F_POSRES) > 0 || gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0
1249              || ir->nwall > 0 || ir->bPull || ir->bRot || ir->bIMD);
1250
1251     if (fr->shift_vec == nullptr)
1252     {
1253         snew(fr->shift_vec, SHIFTS);
1254     }
1255
1256     fr->shiftForces.resize(SHIFTS);
1257
1258     if (fr->nbfp.empty())
1259     {
1260         fr->ntype = mtop->ffparams.atnr;
1261         fr->nbfp  = mk_nbfp(&mtop->ffparams, fr->bBHAM);
1262         if (EVDW_PME(ic->vdwtype))
1263         {
1264             fr->ljpme_c6grid = make_ljpme_c6grid(&mtop->ffparams, fr);
1265         }
1266     }
1267
1268     /* Copy the energy group exclusions */
1269     fr->egp_flags = ir->opts.egp_flags;
1270
1271     /* Van der Waals stuff */
1272     if ((ic->vdwtype != evdwCUT) && (ic->vdwtype != evdwUSER) && !fr->bBHAM)
1273     {
1274         if (ic->rvdw_switch >= ic->rvdw)
1275         {
1276             gmx_fatal(FARGS, "rvdw_switch (%f) must be < rvdw (%f)", ic->rvdw_switch, ic->rvdw);
1277         }
1278         if (fp)
1279         {
1280             fprintf(fp, "Using %s Lennard-Jones, switch between %g and %g nm\n",
1281                     (ic->eeltype == eelSWITCH) ? "switched" : "shifted", ic->rvdw_switch, ic->rvdw);
1282         }
1283     }
1284
1285     if (fr->bBHAM && EVDW_PME(ic->vdwtype))
1286     {
1287         gmx_fatal(FARGS, "LJ PME not supported with Buckingham");
1288     }
1289
1290     if (fr->bBHAM && (ic->vdwtype == evdwSHIFT || ic->vdwtype == evdwSWITCH))
1291     {
1292         gmx_fatal(FARGS, "Switch/shift interaction not supported with Buckingham");
1293     }
1294
1295     if (fr->bBHAM && fr->cutoff_scheme == ecutsVERLET)
1296     {
1297         gmx_fatal(FARGS, "Verlet cutoff-scheme is not supported with Buckingham");
1298     }
1299
1300     if (fp && fr->cutoff_scheme == ecutsGROUP)
1301     {
1302         fprintf(fp, "Cut-off's:   NS: %g   Coulomb: %g   %s: %g\n", fr->rlist, ic->rcoulomb,
1303                 fr->bBHAM ? "BHAM" : "LJ", ic->rvdw);
1304     }
1305
1306     if (ir->implicit_solvent)
1307     {
1308         gmx_fatal(FARGS, "Implict solvation is no longer supported.");
1309     }
1310
1311
1312     /* This code automatically gives table length tabext without cut-off's,
1313      * in that case grompp should already have checked that we do not need
1314      * normal tables and we only generate tables for 1-4 interactions.
1315      */
1316     real rtab = ir->rlist + ir->tabext;
1317
1318     /* We want to use unmodified tables for 1-4 coulombic
1319      * interactions, so we must in general have an extra set of
1320      * tables. */
1321     if (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 || gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0
1322         || gmx_mtop_ftype_count(mtop, F_LJC_PAIRS_NB) > 0)
1323     {
1324         fr->pairsTable = make_tables(fp, ic, tabpfn, rtab, GMX_MAKETABLES_14ONLY);
1325     }
1326
1327     /* Wall stuff */
1328     fr->nwall = ir->nwall;
1329     if (ir->nwall && ir->wall_type == ewtTABLE)
1330     {
1331         make_wall_tables(fp, ir, tabfn, &mtop->groups, fr);
1332     }
1333
1334     if (fcd && !tabbfnm.empty())
1335     {
1336         // Need to catch std::bad_alloc
1337         // TODO Don't need to catch this here, when merging with master branch
1338         try
1339         {
1340             fcd->bondtab  = make_bonded_tables(fp, F_TABBONDS, F_TABBONDSNC, mtop, tabbfnm, "b");
1341             fcd->angletab = make_bonded_tables(fp, F_TABANGLES, -1, mtop, tabbfnm, "a");
1342             fcd->dihtab   = make_bonded_tables(fp, F_TABDIHS, -1, mtop, tabbfnm, "d");
1343         }
1344         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1345     }
1346     else
1347     {
1348         if (debug)
1349         {
1350             fprintf(debug,
1351                     "No fcdata or table file name passed, can not read table, can not do bonded "
1352                     "interactions\n");
1353         }
1354     }
1355
1356     // QM/MM initialization if requested
1357     fr->bQMMM = ir->bQMMM;
1358     if (fr->bQMMM)
1359     {
1360         // Initialize QM/MM if supported
1361         if (GMX_QMMM)
1362         {
1363             GMX_LOG(mdlog.info)
1364                     .asParagraph()
1365                     .appendText(
1366                             "Large parts of the QM/MM support is deprecated, and may be removed in "
1367                             "a future "
1368                             "version. Please get in touch with the developers if you find the "
1369                             "support useful, "
1370                             "as help is needed if the functionality is to continue to be "
1371                             "available.");
1372             fr->qr = mk_QMMMrec();
1373             init_QMMMrec(cr, mtop, ir, fr);
1374         }
1375         else
1376         {
1377             gmx_incons(
1378                     "QM/MM was requested, but is only available when GROMACS "
1379                     "is configured with QM/MM support");
1380         }
1381     }
1382
1383     /* Set all the static charge group info */
1384     fr->cginfo_mb = init_cginfo_mb(mtop, fr);
1385     if (!DOMAINDECOMP(cr))
1386     {
1387         fr->cginfo = cginfo_expand(mtop->molblock.size(), fr->cginfo_mb);
1388     }
1389
1390     if (!DOMAINDECOMP(cr))
1391     {
1392         forcerec_set_ranges(fr, mtop->natoms, mtop->natoms, mtop->natoms);
1393     }
1394
1395     fr->print_force = print_force;
1396
1397     /* Initialize the thread working data for bonded interactions */
1398     fr->bondedThreading = init_bonded_threading(
1399             fp, mtop->groups.groups[SimulationAtomGroupType::EnergyOutput].size());
1400
1401     fr->nthread_ewc = gmx_omp_nthreads_get(emntBonded);
1402     snew(fr->ewc_t, fr->nthread_ewc);
1403
1404     if (ir->eDispCorr != edispcNO)
1405     {
1406         fr->dispersionCorrection = std::make_unique<DispersionCorrection>(
1407                 *mtop, *ir, fr->bBHAM, fr->ntype, fr->nbfp, *fr->ic, tabfn);
1408         fr->dispersionCorrection->print(mdlog);
1409     }
1410
1411     if (fp != nullptr)
1412     {
1413         /* Here we switch from using mdlog, which prints the newline before
1414          * the paragraph, to our old fprintf logging, which prints the newline
1415          * after the paragraph, so we should add a newline here.
1416          */
1417         fprintf(fp, "\n");
1418     }
1419 }
1420
1421 t_forcerec::t_forcerec() = default;
1422
1423 t_forcerec::~t_forcerec()
1424 {
1425     /* Note: This code will disappear when types are converted to C++ */
1426     sfree(shift_vec);
1427     sfree(ewc_t);
1428     tear_down_bonded_threading(bondedThreading);
1429 }