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