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