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