Apply re-formatting to C++ in src/ tree.
[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/mdtypes/multipletimestepping.h"
86 #include "gromacs/nbnxm/nbnxm.h"
87 #include "gromacs/pbcutil/ishift.h"
88 #include "gromacs/pbcutil/pbc.h"
89 #include "gromacs/tables/forcetable.h"
90 #include "gromacs/topology/mtop_util.h"
91 #include "gromacs/trajectory/trajectoryframe.h"
92 #include "gromacs/utility/cstringutil.h"
93 #include "gromacs/utility/exceptions.h"
94 #include "gromacs/utility/fatalerror.h"
95 #include "gromacs/utility/gmxassert.h"
96 #include "gromacs/utility/logger.h"
97 #include "gromacs/utility/physicalnodecommunicator.h"
98 #include "gromacs/utility/pleasecite.h"
99 #include "gromacs/utility/smalloc.h"
100 #include "gromacs/utility/strconvert.h"
101
102 #include "gpuforcereduction.h"
103
104 ForceHelperBuffers::ForceHelperBuffers(bool haveDirectVirialContributions) :
105     haveDirectVirialContributions_(haveDirectVirialContributions)
106 {
107     shiftForces_.resize(SHIFTS);
108 }
109
110 void ForceHelperBuffers::resize(int numAtoms)
111 {
112     if (haveDirectVirialContributions_)
113     {
114         forceBufferForDirectVirialContributions_.resize(numAtoms);
115     }
116 }
117
118 static std::vector<real> mk_nbfp(const gmx_ffparams_t* idef, gmx_bool bBHAM)
119 {
120     std::vector<real> nbfp;
121     int               atnr;
122
123     atnr = idef->atnr;
124     if (bBHAM)
125     {
126         nbfp.resize(3 * atnr * atnr);
127         int k = 0;
128         for (int i = 0; (i < atnr); i++)
129         {
130             for (int j = 0; (j < atnr); j++, k++)
131             {
132                 BHAMA(nbfp, atnr, i, j) = idef->iparams[k].bham.a;
133                 BHAMB(nbfp, atnr, i, j) = idef->iparams[k].bham.b;
134                 /* nbfp now includes the 6.0 derivative prefactor */
135                 BHAMC(nbfp, atnr, i, j) = idef->iparams[k].bham.c * 6.0;
136             }
137         }
138     }
139     else
140     {
141         nbfp.resize(2 * atnr * atnr);
142         int k = 0;
143         for (int i = 0; (i < atnr); i++)
144         {
145             for (int j = 0; (j < atnr); j++, k++)
146             {
147                 /* nbfp now includes the 6.0/12.0 derivative prefactors */
148                 C6(nbfp, atnr, i, j)  = idef->iparams[k].lj.c6 * 6.0;
149                 C12(nbfp, atnr, i, j) = idef->iparams[k].lj.c12 * 12.0;
150             }
151         }
152     }
153
154     return nbfp;
155 }
156
157 static real* make_ljpme_c6grid(const gmx_ffparams_t* idef, t_forcerec* fr)
158 {
159     int   i, j, k, atnr;
160     real  c6, c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
161     real* grid;
162
163     /* For LJ-PME simulations, we correct the energies with the reciprocal space
164      * inside of the cut-off. To do this the non-bonded kernels needs to have
165      * access to the C6-values used on the reciprocal grid in pme.c
166      */
167
168     atnr = idef->atnr;
169     snew(grid, 2 * atnr * atnr);
170     for (i = k = 0; (i < atnr); i++)
171     {
172         for (j = 0; (j < atnr); j++, k++)
173         {
174             c6i  = idef->iparams[i * (atnr + 1)].lj.c6;
175             c12i = idef->iparams[i * (atnr + 1)].lj.c12;
176             c6j  = idef->iparams[j * (atnr + 1)].lj.c6;
177             c12j = idef->iparams[j * (atnr + 1)].lj.c12;
178             c6   = std::sqrt(c6i * c6j);
179             if (fr->ljpme_combination_rule == eljpmeLB && !gmx_numzero(c6) && !gmx_numzero(c12i)
180                 && !gmx_numzero(c12j))
181             {
182                 sigmai = gmx::sixthroot(c12i / c6i);
183                 sigmaj = gmx::sixthroot(c12j / c6j);
184                 epsi   = c6i * c6i / c12i;
185                 epsj   = c6j * c6j / c12j;
186                 c6     = std::sqrt(epsi * epsj) * gmx::power6(0.5 * (sigmai + sigmaj));
187             }
188             /* Store the elements at the same relative positions as C6 in nbfp in order
189              * to simplify access in the kernels
190              */
191             grid[2 * (atnr * i + j)] = c6 * 6.0;
192         }
193     }
194     return grid;
195 }
196
197 enum
198 {
199     acNONE = 0,
200     acCONSTRAINT,
201     acSETTLE
202 };
203
204 static std::vector<cginfo_mb_t> init_cginfo_mb(const gmx_mtop_t* mtop, const t_forcerec* fr)
205 {
206     gmx_bool* type_VDW;
207     int*      a_con;
208
209     snew(type_VDW, fr->ntype);
210     for (int ai = 0; ai < fr->ntype; ai++)
211     {
212         type_VDW[ai] = FALSE;
213         for (int j = 0; j < fr->ntype; j++)
214         {
215             type_VDW[ai] = type_VDW[ai] || fr->bBHAM || C6(fr->nbfp, fr->ntype, ai, j) != 0
216                            || C12(fr->nbfp, fr->ntype, ai, j) != 0;
217         }
218     }
219
220     std::vector<cginfo_mb_t> cginfoPerMolblock;
221     int                      a_offset = 0;
222     for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
223     {
224         const gmx_molblock_t& molb = mtop->molblock[mb];
225         const gmx_moltype_t&  molt = mtop->moltype[molb.type];
226         const auto&           excl = molt.excls;
227
228         /* Check if the cginfo is identical for all molecules in this block.
229          * If so, we only need an array of the size of one molecule.
230          * Otherwise we make an array of #mol times #cgs per molecule.
231          */
232         gmx_bool bId = TRUE;
233         for (int m = 0; m < molb.nmol; m++)
234         {
235             const int am = m * molt.atoms.nr;
236             for (int a = 0; a < molt.atoms.nr; a++)
237             {
238                 if (getGroupType(mtop->groups, SimulationAtomGroupType::QuantumMechanics, a_offset + am + a)
239                     != getGroupType(mtop->groups, SimulationAtomGroupType::QuantumMechanics, a_offset + a))
240                 {
241                     bId = FALSE;
242                 }
243                 if (!mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics].empty())
244                 {
245                     if (mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics][a_offset + am + a]
246                         != mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics][a_offset + a])
247                     {
248                         bId = FALSE;
249                     }
250                 }
251             }
252         }
253
254         cginfo_mb_t cginfo_mb;
255         cginfo_mb.cg_start = a_offset;
256         cginfo_mb.cg_end   = a_offset + molb.nmol * molt.atoms.nr;
257         cginfo_mb.cg_mod   = (bId ? 1 : molb.nmol) * molt.atoms.nr;
258         cginfo_mb.cginfo.resize(cginfo_mb.cg_mod);
259         gmx::ArrayRef<int> cginfo = cginfo_mb.cginfo;
260
261         /* Set constraints flags for constrained atoms */
262         snew(a_con, molt.atoms.nr);
263         for (int ftype = 0; ftype < F_NRE; ftype++)
264         {
265             if (interaction_function[ftype].flags & IF_CONSTRAINT)
266             {
267                 const int nral = NRAL(ftype);
268                 for (int ia = 0; ia < molt.ilist[ftype].size(); ia += 1 + nral)
269                 {
270                     int a;
271
272                     for (a = 0; a < nral; a++)
273                     {
274                         a_con[molt.ilist[ftype].iatoms[ia + 1 + a]] =
275                                 (ftype == F_SETTLE ? acSETTLE : acCONSTRAINT);
276                     }
277                 }
278             }
279         }
280
281         for (int m = 0; m < (bId ? 1 : molb.nmol); m++)
282         {
283             const int molculeOffsetInBlock = m * molt.atoms.nr;
284             for (int a = 0; a < molt.atoms.nr; a++)
285             {
286                 const t_atom& atom     = molt.atoms.atom[a];
287                 int&          atomInfo = cginfo[molculeOffsetInBlock + a];
288
289                 /* Store the energy group in cginfo */
290                 int gid = getGroupType(mtop->groups,
291                                        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,
602                               isPlural ? "' or '" : "",
603                               isPlural ? interaction_function[ftype2].longname : "",
604                               i,
605                               patternToFind.c_str());
606                 }
607             }
608         }
609         sfree(count);
610     }
611
612     return tab;
613 }
614
615 void forcerec_set_ranges(t_forcerec* fr, int natoms_force, int natoms_force_constr, int natoms_f_novirsum)
616 {
617     fr->natoms_force        = natoms_force;
618     fr->natoms_force_constr = natoms_force_constr;
619
620     for (auto& forceHelperBuffers : fr->forceHelperBuffers)
621     {
622         forceHelperBuffers.resize(natoms_f_novirsum);
623     }
624 }
625
626 static real cutoff_inf(real cutoff)
627 {
628     if (cutoff == 0)
629     {
630         cutoff = GMX_CUTOFF_INF;
631     }
632
633     return cutoff;
634 }
635
636 /*! \brief Print Coulomb Ewald citations and set ewald coefficients */
637 static void initCoulombEwaldParameters(FILE*                fp,
638                                        const t_inputrec*    ir,
639                                        bool                 systemHasNetCharge,
640                                        interaction_const_t* ic)
641 {
642     if (!EEL_PME_EWALD(ir->coulombtype))
643     {
644         return;
645     }
646
647     if (fp)
648     {
649         fprintf(fp, "Will do PME sum in reciprocal space for electrostatic interactions.\n");
650
651         if (ir->coulombtype == eelP3M_AD)
652         {
653             please_cite(fp, "Hockney1988");
654             please_cite(fp, "Ballenegger2012");
655         }
656         else
657         {
658             please_cite(fp, "Essmann95a");
659         }
660
661         if (ir->ewald_geometry == eewg3DC)
662         {
663             if (fp)
664             {
665                 fprintf(fp,
666                         "Using the Ewald3DC correction for systems with a slab geometry%s.\n",
667                         systemHasNetCharge ? " and net charge" : "");
668             }
669             please_cite(fp, "In-Chul99a");
670             if (systemHasNetCharge)
671             {
672                 please_cite(fp, "Ballenegger2009");
673             }
674         }
675     }
676
677     ic->ewaldcoeff_q = calc_ewaldcoeff_q(ir->rcoulomb, ir->ewald_rtol);
678     if (fp)
679     {
680         fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for Ewald\n", 1 / ic->ewaldcoeff_q);
681     }
682
683     if (ic->coulomb_modifier == eintmodPOTSHIFT)
684     {
685         GMX_RELEASE_ASSERT(ic->rcoulomb != 0, "Cutoff radius cannot be zero");
686         ic->sh_ewald = std::erfc(ic->ewaldcoeff_q * ic->rcoulomb) / ic->rcoulomb;
687     }
688     else
689     {
690         ic->sh_ewald = 0;
691     }
692 }
693
694 /*! \brief Print Van der Waals Ewald citations and set ewald coefficients */
695 static void initVdwEwaldParameters(FILE* fp, const t_inputrec* ir, interaction_const_t* ic)
696 {
697     if (!EVDW_PME(ir->vdwtype))
698     {
699         return;
700     }
701
702     if (fp)
703     {
704         fprintf(fp, "Will do PME sum in reciprocal space for LJ dispersion interactions.\n");
705         please_cite(fp, "Essmann95a");
706     }
707     ic->ewaldcoeff_lj = calc_ewaldcoeff_lj(ir->rvdw, ir->ewald_rtol_lj);
708     if (fp)
709     {
710         fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for LJ Ewald\n", 1 / ic->ewaldcoeff_lj);
711     }
712
713     if (ic->vdw_modifier == eintmodPOTSHIFT)
714     {
715         real crc2       = gmx::square(ic->ewaldcoeff_lj * ic->rvdw);
716         ic->sh_lj_ewald = (std::exp(-crc2) * (1 + crc2 + 0.5 * crc2 * crc2) - 1) / gmx::power6(ic->rvdw);
717     }
718     else
719     {
720         ic->sh_lj_ewald = 0;
721     }
722 }
723
724 /* Generate Coulomb and/or Van der Waals Ewald long-range correction tables
725  *
726  * Tables are generated for one or both, depending on if the pointers
727  * are non-null. The spacing for both table sets is the same and obeys
728  * both accuracy requirements, when relevant.
729  */
730 static void init_ewald_f_table(const interaction_const_t& ic,
731                                const real                 tableExtensionLength,
732                                EwaldCorrectionTables*     coulombTables,
733                                EwaldCorrectionTables*     vdwTables)
734 {
735     const bool useCoulombTable = (EEL_PME_EWALD(ic.eeltype) && coulombTables != nullptr);
736     const bool useVdwTable     = (EVDW_PME(ic.vdwtype) && vdwTables != nullptr);
737
738     /* Get the Ewald table spacing based on Coulomb and/or LJ
739      * Ewald coefficients and rtol.
740      */
741     const real tableScale = ewald_spline3_table_scale(ic, useCoulombTable, useVdwTable);
742
743     const bool havePerturbedNonbondeds = (ic.softCoreParameters != nullptr);
744
745     real tableLen = ic.rcoulomb;
746     if (useCoulombTable && havePerturbedNonbondeds && tableExtensionLength > 0.0)
747     {
748         /* TODO: Ideally this should also check if couple-intramol == no, but that isn't
749          * stored in ir. Grompp puts that info into an opts structure that doesn't make it into the tpr.
750          * The alternative is to look through all the exclusions and check if they come from
751          * couple-intramol == no. Meanwhile, always having larger tables should only affect
752          * memory consumption, not speed (barring cache issues).
753          */
754         tableLen = ic.rcoulomb + tableExtensionLength;
755     }
756     const int tableSize = static_cast<int>(tableLen * tableScale) + 2;
757
758     if (useCoulombTable)
759     {
760         *coulombTables =
761                 generateEwaldCorrectionTables(tableSize, tableScale, ic.ewaldcoeff_q, v_q_ewald_lr);
762     }
763
764     if (useVdwTable)
765     {
766         *vdwTables = generateEwaldCorrectionTables(tableSize, tableScale, ic.ewaldcoeff_lj, v_lj_ewald_lr);
767     }
768 }
769
770 void init_interaction_const_tables(FILE* fp, interaction_const_t* ic, const real tableExtensionLength)
771 {
772     if (EEL_PME_EWALD(ic->eeltype) || EVDW_PME(ic->vdwtype))
773     {
774         init_ewald_f_table(
775                 *ic, tableExtensionLength, ic->coulombEwaldTables.get(), ic->vdwEwaldTables.get());
776         if (fp != nullptr)
777         {
778             fprintf(fp,
779                     "Initialized non-bonded Ewald tables, spacing: %.2e size: %zu\n\n",
780                     1 / ic->coulombEwaldTables->scale,
781                     ic->coulombEwaldTables->tableF.size());
782         }
783     }
784 }
785
786 static void clear_force_switch_constants(shift_consts_t* sc)
787 {
788     sc->c2   = 0;
789     sc->c3   = 0;
790     sc->cpot = 0;
791 }
792
793 static void force_switch_constants(real p, real rsw, real rc, shift_consts_t* sc)
794 {
795     /* Here we determine the coefficient for shifting the force to zero
796      * between distance rsw and the cut-off rc.
797      * For a potential of r^-p, we have force p*r^-(p+1).
798      * But to save flops we absorb p in the coefficient.
799      * Thus we get:
800      * force/p   = r^-(p+1) + c2*r^2 + c3*r^3
801      * potential = r^-p + c2/3*r^3 + c3/4*r^4 + cpot
802      */
803     sc->c2   = ((p + 1) * rsw - (p + 4) * rc) / (pow(rc, p + 2) * gmx::square(rc - rsw));
804     sc->c3   = -((p + 1) * rsw - (p + 3) * rc) / (pow(rc, p + 2) * gmx::power3(rc - rsw));
805     sc->cpot = -pow(rc, -p) + p * sc->c2 / 3 * gmx::power3(rc - rsw)
806                + p * sc->c3 / 4 * gmx::power4(rc - rsw);
807 }
808
809 static void potential_switch_constants(real rsw, real rc, switch_consts_t* sc)
810 {
811     /* The switch function is 1 at rsw and 0 at rc.
812      * The derivative and second derivate are zero at both ends.
813      * rsw        = max(r - r_switch, 0)
814      * sw         = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
815      * dsw        = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
816      * force      = force*dsw - potential*sw
817      * potential *= sw
818      */
819     sc->c3 = -10 / gmx::power3(rc - rsw);
820     sc->c4 = 15 / gmx::power4(rc - rsw);
821     sc->c5 = -6 / gmx::power5(rc - rsw);
822 }
823
824 /*! \brief Construct interaction constants
825  *
826  * This data is used (particularly) by search and force code for
827  * short-range interactions. Many of these are constant for the whole
828  * simulation; some are constant only after PME tuning completes.
829  */
830 static void init_interaction_const(FILE*                 fp,
831                                    interaction_const_t** interaction_const,
832                                    const t_inputrec*     ir,
833                                    const gmx_mtop_t*     mtop,
834                                    bool                  systemHasNetCharge)
835 {
836     interaction_const_t* ic = new interaction_const_t;
837
838     ic->coulombEwaldTables = std::make_unique<EwaldCorrectionTables>();
839     ic->vdwEwaldTables     = std::make_unique<EwaldCorrectionTables>();
840
841     /* Lennard-Jones */
842     ic->vdwtype         = ir->vdwtype;
843     ic->vdw_modifier    = ir->vdw_modifier;
844     ic->reppow          = mtop->ffparams.reppow;
845     ic->rvdw            = cutoff_inf(ir->rvdw);
846     ic->rvdw_switch     = ir->rvdw_switch;
847     ic->ljpme_comb_rule = ir->ljpme_combination_rule;
848     ic->useBuckingham   = (mtop->ffparams.functype[0] == F_BHAM);
849     if (ic->useBuckingham)
850     {
851         ic->buckinghamBMax = calcBuckinghamBMax(fp, mtop);
852     }
853
854     initVdwEwaldParameters(fp, ir, ic);
855
856     clear_force_switch_constants(&ic->dispersion_shift);
857     clear_force_switch_constants(&ic->repulsion_shift);
858
859     switch (ic->vdw_modifier)
860     {
861         case eintmodPOTSHIFT:
862             /* Only shift the potential, don't touch the force */
863             ic->dispersion_shift.cpot = -1.0 / gmx::power6(ic->rvdw);
864             ic->repulsion_shift.cpot  = -1.0 / gmx::power12(ic->rvdw);
865             break;
866         case eintmodFORCESWITCH:
867             /* Switch the force, switch and shift the potential */
868             force_switch_constants(6.0, ic->rvdw_switch, ic->rvdw, &ic->dispersion_shift);
869             force_switch_constants(12.0, ic->rvdw_switch, ic->rvdw, &ic->repulsion_shift);
870             break;
871         case eintmodPOTSWITCH:
872             /* Switch the potential and force */
873             potential_switch_constants(ic->rvdw_switch, ic->rvdw, &ic->vdw_switch);
874             break;
875         case eintmodNONE:
876         case eintmodEXACTCUTOFF:
877             /* Nothing to do here */
878             break;
879         default: gmx_incons("unimplemented potential modifier");
880     }
881
882     /* Electrostatics */
883     ic->eeltype          = ir->coulombtype;
884     ic->coulomb_modifier = ir->coulomb_modifier;
885     ic->rcoulomb         = cutoff_inf(ir->rcoulomb);
886     ic->rcoulomb_switch  = ir->rcoulomb_switch;
887     ic->epsilon_r        = ir->epsilon_r;
888
889     /* Set the Coulomb energy conversion factor */
890     if (ic->epsilon_r != 0)
891     {
892         ic->epsfac = ONE_4PI_EPS0 / ic->epsilon_r;
893     }
894     else
895     {
896         /* eps = 0 is infinite dieletric: no Coulomb interactions */
897         ic->epsfac = 0;
898     }
899
900     /* Reaction-field */
901     if (EEL_RF(ic->eeltype))
902     {
903         GMX_RELEASE_ASSERT(ic->eeltype != eelGRF_NOTUSED, "GRF is no longer supported");
904         ic->epsilon_rf = ir->epsilon_rf;
905
906         calc_rffac(fp, ic->epsilon_r, ic->epsilon_rf, ic->rcoulomb, &ic->k_rf, &ic->c_rf);
907     }
908     else
909     {
910         /* For plain cut-off we might use the reaction-field kernels */
911         ic->epsilon_rf = ic->epsilon_r;
912         ic->k_rf       = 0;
913         if (ir->coulomb_modifier == eintmodPOTSHIFT)
914         {
915             ic->c_rf = 1 / ic->rcoulomb;
916         }
917         else
918         {
919             ic->c_rf = 0;
920         }
921     }
922
923     initCoulombEwaldParameters(fp, ir, systemHasNetCharge, ic);
924
925     if (fp != nullptr)
926     {
927         real dispersion_shift;
928
929         dispersion_shift = ic->dispersion_shift.cpot;
930         if (EVDW_PME(ic->vdwtype))
931         {
932             dispersion_shift -= ic->sh_lj_ewald;
933         }
934         fprintf(fp, "Potential shift: LJ r^-12: %.3e r^-6: %.3e", ic->repulsion_shift.cpot, dispersion_shift);
935
936         if (ic->eeltype == eelCUT)
937         {
938             fprintf(fp, ", Coulomb %.e", -ic->c_rf);
939         }
940         else if (EEL_PME(ic->eeltype))
941         {
942             fprintf(fp, ", Ewald %.3e", -ic->sh_ewald);
943         }
944         fprintf(fp, "\n");
945     }
946
947     if (ir->efep != efepNO)
948     {
949         GMX_RELEASE_ASSERT(ir->fepvals, "ir->fepvals should be set wth free-energy");
950         ic->softCoreParameters = std::make_unique<interaction_const_t::SoftCoreParameters>(*ir->fepvals);
951     }
952
953     *interaction_const = ic;
954 }
955
956 void init_forcerec(FILE*                            fp,
957                    const gmx::MDLogger&             mdlog,
958                    t_forcerec*                      fr,
959                    const t_inputrec*                ir,
960                    const gmx_mtop_t*                mtop,
961                    const t_commrec*                 cr,
962                    matrix                           box,
963                    const char*                      tabfn,
964                    const char*                      tabpfn,
965                    gmx::ArrayRef<const std::string> tabbfnm,
966                    real                             print_force)
967 {
968     /* The CMake default turns SIMD kernels on, but it might be turned off further down... */
969     fr->use_simd_kernels = GMX_USE_SIMD_KERNELS;
970
971     if (check_box(ir->pbcType, box))
972     {
973         gmx_fatal(FARGS, "%s", check_box(ir->pbcType, box));
974     }
975
976     /* Test particle insertion ? */
977     if (EI_TPI(ir->eI))
978     {
979         /* Set to the size of the molecule to be inserted (the last one) */
980         gmx::RangePartitioning molecules = gmx_mtop_molecules(*mtop);
981         fr->n_tpi                        = molecules.block(molecules.numBlocks() - 1).size();
982     }
983     else
984     {
985         fr->n_tpi = 0;
986     }
987
988     if (ir->coulombtype == eelRF_NEC_UNSUPPORTED || ir->coulombtype == eelGRF_NOTUSED)
989     {
990         gmx_fatal(FARGS, "%s electrostatics is no longer supported", eel_names[ir->coulombtype]);
991     }
992
993     if (ir->bAdress)
994     {
995         gmx_fatal(FARGS, "AdResS simulations are no longer supported");
996     }
997     if (ir->useTwinRange)
998     {
999         gmx_fatal(FARGS, "Twin-range simulations are no longer supported");
1000     }
1001     /* Copy the user determined parameters */
1002     fr->userint1  = ir->userint1;
1003     fr->userint2  = ir->userint2;
1004     fr->userint3  = ir->userint3;
1005     fr->userint4  = ir->userint4;
1006     fr->userreal1 = ir->userreal1;
1007     fr->userreal2 = ir->userreal2;
1008     fr->userreal3 = ir->userreal3;
1009     fr->userreal4 = ir->userreal4;
1010
1011     /* Shell stuff */
1012     fr->fc_stepsize = ir->fc_stepsize;
1013
1014     /* Free energy */
1015     fr->efep = ir->efep;
1016
1017     if ((getenv("GMX_DISABLE_SIMD_KERNELS") != nullptr) || (getenv("GMX_NOOPTIMIZEDKERNELS") != nullptr))
1018     {
1019         fr->use_simd_kernels = FALSE;
1020         if (fp != nullptr)
1021         {
1022             fprintf(fp,
1023                     "\nFound environment variable GMX_DISABLE_SIMD_KERNELS.\n"
1024                     "Disabling the usage of any SIMD-specific non-bonded & bonded kernel routines\n"
1025                     "(e.g. SSE2/SSE4.1/AVX).\n\n");
1026         }
1027     }
1028
1029     fr->bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
1030
1031     /* Neighbour searching stuff */
1032     fr->pbcType = ir->pbcType;
1033
1034     /* Determine if we will do PBC for distances in bonded interactions */
1035     if (fr->pbcType == PbcType::No)
1036     {
1037         fr->bMolPBC = FALSE;
1038     }
1039     else
1040     {
1041         const bool useEwaldSurfaceCorrection =
1042                 (EEL_PME_EWALD(ir->coulombtype) && ir->epsilon_surface != 0);
1043         const bool haveOrientationRestraints = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
1044         if (!DOMAINDECOMP(cr))
1045         {
1046             fr->bMolPBC = true;
1047
1048             if (useEwaldSurfaceCorrection || haveOrientationRestraints)
1049             {
1050                 fr->wholeMoleculeTransform =
1051                         std::make_unique<gmx::WholeMoleculeTransform>(*mtop, ir->pbcType);
1052             }
1053         }
1054         else
1055         {
1056             fr->bMolPBC = dd_bonded_molpbc(cr->dd, fr->pbcType);
1057
1058             /* With Ewald surface correction it is useful to support e.g. running water
1059              * in parallel with update groups.
1060              * With orientation restraints there is no sensible use case supported with DD.
1061              */
1062             if ((useEwaldSurfaceCorrection && !dd_moleculesAreAlwaysWhole(*cr->dd)) || haveOrientationRestraints)
1063             {
1064                 gmx_fatal(FARGS,
1065                           "You requested Ewald surface correction or orientation restraints, "
1066                           "but molecules are broken "
1067                           "over periodic boundary conditions by the domain decomposition. "
1068                           "Run without domain decomposition instead.");
1069             }
1070         }
1071
1072         if (useEwaldSurfaceCorrection)
1073         {
1074             GMX_RELEASE_ASSERT(!DOMAINDECOMP(cr) || dd_moleculesAreAlwaysWhole(*cr->dd),
1075                                "Molecules can not be broken by PBC with epsilon_surface > 0");
1076         }
1077     }
1078
1079     fr->rc_scaling = ir->refcoord_scaling;
1080     copy_rvec(ir->posres_com, fr->posres_com);
1081     copy_rvec(ir->posres_comB, fr->posres_comB);
1082     fr->rlist                  = cutoff_inf(ir->rlist);
1083     fr->ljpme_combination_rule = ir->ljpme_combination_rule;
1084
1085     /* This now calculates sum for q and c6*/
1086     bool systemHasNetCharge = set_chargesum(fp, fr, mtop);
1087
1088     /* fr->ic is used both by verlet and group kernels (to some extent) now */
1089     init_interaction_const(fp, &fr->ic, ir, mtop, systemHasNetCharge);
1090     init_interaction_const_tables(fp, fr->ic, ir->tabext);
1091
1092     const interaction_const_t* ic = fr->ic;
1093
1094     /* TODO: Replace this Ewald table or move it into interaction_const_t */
1095     if (ir->coulombtype == eelEWALD)
1096     {
1097         init_ewald_tab(&(fr->ewald_table), ir, fp);
1098     }
1099
1100     /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
1101     switch (ic->eeltype)
1102     {
1103         case eelCUT: fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_COULOMB; break;
1104
1105         case eelRF:
1106         case eelRF_ZERO: fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD; break;
1107
1108         case eelSWITCH:
1109         case eelSHIFT:
1110         case eelUSER:
1111         case eelPMESWITCH:
1112         case eelPMEUSER:
1113         case eelPMEUSERSWITCH:
1114             fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
1115             break;
1116
1117         case eelPME:
1118         case eelP3M_AD:
1119         case eelEWALD: fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_EWALD; break;
1120
1121         default:
1122             gmx_fatal(FARGS, "Unsupported electrostatic interaction: %s", eel_names[ic->eeltype]);
1123     }
1124     fr->nbkernel_elec_modifier = ic->coulomb_modifier;
1125
1126     /* Vdw: Translate from mdp settings to kernel format */
1127     switch (ic->vdwtype)
1128     {
1129         case evdwCUT:
1130             if (fr->bBHAM)
1131             {
1132                 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_BUCKINGHAM;
1133             }
1134             else
1135             {
1136                 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LENNARDJONES;
1137             }
1138             break;
1139         case evdwPME: fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LJEWALD; break;
1140
1141         case evdwSWITCH:
1142         case evdwSHIFT:
1143         case evdwUSER: fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE; break;
1144
1145         default: gmx_fatal(FARGS, "Unsupported vdw interaction: %s", evdw_names[ic->vdwtype]);
1146     }
1147     fr->nbkernel_vdw_modifier = ic->vdw_modifier;
1148
1149     if (!gmx_within_tol(ic->reppow, 12.0, 10 * GMX_DOUBLE_EPS))
1150     {
1151         gmx_fatal(FARGS, "Only LJ repulsion power 12 is supported");
1152     }
1153     /* Older tpr files can contain Coulomb user tables with the Verlet cutoff-scheme,
1154      * while mdrun does not (and never did) support this.
1155      */
1156     if (EEL_USER(fr->ic->eeltype))
1157     {
1158         gmx_fatal(FARGS, "Electrostatics type %s is currently not supported", eel_names[ir->coulombtype]);
1159     }
1160
1161     fr->bvdwtab  = FALSE;
1162     fr->bcoultab = FALSE;
1163
1164     /* 1-4 interaction electrostatics */
1165     fr->fudgeQQ = mtop->ffparams.fudgeQQ;
1166
1167     // Multiple time stepping
1168     fr->useMts = ir->useMts;
1169
1170     if (fr->useMts)
1171     {
1172         GMX_ASSERT(gmx::checkMtsRequirements(*ir).empty(),
1173                    "All MTS requirements should be met here");
1174     }
1175
1176     const bool haveDirectVirialContributionsFast =
1177             fr->forceProviders->hasForceProvider() || gmx_mtop_ftype_count(mtop, F_POSRES) > 0
1178             || gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0 || ir->nwall > 0 || ir->bPull || ir->bRot
1179             || ir->bIMD;
1180     const bool haveDirectVirialContributionsSlow = EEL_FULL(ic->eeltype) || EVDW_PME(ic->vdwtype);
1181     for (int i = 0; i < (fr->useMts ? 2 : 1); i++)
1182     {
1183         bool haveDirectVirialContributions =
1184                 (((!fr->useMts || i == 0) && haveDirectVirialContributionsFast)
1185                  || ((!fr->useMts || i == 1) && haveDirectVirialContributionsSlow));
1186         fr->forceHelperBuffers.emplace_back(haveDirectVirialContributions);
1187     }
1188
1189     if (fr->shift_vec == nullptr)
1190     {
1191         snew(fr->shift_vec, SHIFTS);
1192     }
1193
1194     if (fr->nbfp.empty())
1195     {
1196         fr->ntype = mtop->ffparams.atnr;
1197         fr->nbfp  = mk_nbfp(&mtop->ffparams, fr->bBHAM);
1198         if (EVDW_PME(ic->vdwtype))
1199         {
1200             fr->ljpme_c6grid = make_ljpme_c6grid(&mtop->ffparams, fr);
1201         }
1202     }
1203
1204     /* Copy the energy group exclusions */
1205     fr->egp_flags = ir->opts.egp_flags;
1206
1207     /* Van der Waals stuff */
1208     if ((ic->vdwtype != evdwCUT) && (ic->vdwtype != evdwUSER) && !fr->bBHAM)
1209     {
1210         if (ic->rvdw_switch >= ic->rvdw)
1211         {
1212             gmx_fatal(FARGS, "rvdw_switch (%f) must be < rvdw (%f)", ic->rvdw_switch, ic->rvdw);
1213         }
1214         if (fp)
1215         {
1216             fprintf(fp,
1217                     "Using %s Lennard-Jones, switch between %g and %g nm\n",
1218                     (ic->eeltype == eelSWITCH) ? "switched" : "shifted",
1219                     ic->rvdw_switch,
1220                     ic->rvdw);
1221         }
1222     }
1223
1224     if (fr->bBHAM && EVDW_PME(ic->vdwtype))
1225     {
1226         gmx_fatal(FARGS, "LJ PME not supported with Buckingham");
1227     }
1228
1229     if (fr->bBHAM && (ic->vdwtype == evdwSHIFT || ic->vdwtype == evdwSWITCH))
1230     {
1231         gmx_fatal(FARGS, "Switch/shift interaction not supported with Buckingham");
1232     }
1233
1234     if (fr->bBHAM)
1235     {
1236         gmx_fatal(FARGS, "The Verlet cutoff-scheme does not (yet) support Buckingham");
1237     }
1238
1239     if (ir->implicit_solvent)
1240     {
1241         gmx_fatal(FARGS, "Implict solvation is no longer supported.");
1242     }
1243
1244
1245     /* This code automatically gives table length tabext without cut-off's,
1246      * in that case grompp should already have checked that we do not need
1247      * normal tables and we only generate tables for 1-4 interactions.
1248      */
1249     real rtab = ir->rlist + ir->tabext;
1250
1251     /* We want to use unmodified tables for 1-4 coulombic
1252      * interactions, so we must in general have an extra set of
1253      * tables. */
1254     if (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 || gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0
1255         || gmx_mtop_ftype_count(mtop, F_LJC_PAIRS_NB) > 0)
1256     {
1257         fr->pairsTable = make_tables(fp, ic, tabpfn, rtab, GMX_MAKETABLES_14ONLY);
1258     }
1259
1260     /* Wall stuff */
1261     fr->nwall = ir->nwall;
1262     if (ir->nwall && ir->wall_type == ewtTABLE)
1263     {
1264         make_wall_tables(fp, ir, tabfn, &mtop->groups, fr);
1265     }
1266
1267     fr->fcdata = std::make_unique<t_fcdata>();
1268
1269     if (!tabbfnm.empty())
1270     {
1271         t_fcdata& fcdata = *fr->fcdata;
1272         // Need to catch std::bad_alloc
1273         // TODO Don't need to catch this here, when merging with master branch
1274         try
1275         {
1276             // TODO move these tables into a separate struct and store reference in ListedForces
1277             fcdata.bondtab  = make_bonded_tables(fp, F_TABBONDS, F_TABBONDSNC, mtop, tabbfnm, "b");
1278             fcdata.angletab = make_bonded_tables(fp, F_TABANGLES, -1, mtop, tabbfnm, "a");
1279             fcdata.dihtab   = make_bonded_tables(fp, F_TABDIHS, -1, mtop, tabbfnm, "d");
1280         }
1281         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1282     }
1283     else
1284     {
1285         if (debug)
1286         {
1287             fprintf(debug,
1288                     "No fcdata or table file name passed, can not read table, can not do bonded "
1289                     "interactions\n");
1290         }
1291     }
1292
1293     /* Initialize the thread working data for bonded interactions */
1294     if (fr->useMts)
1295     {
1296         // Add one ListedForces object for each MTS level
1297         bool isFirstLevel = true;
1298         for (const auto& mtsLevel : ir->mtsLevels)
1299         {
1300             ListedForces::InteractionSelection interactionSelection;
1301             const auto&                        forceGroups = mtsLevel.forceGroups;
1302             if (forceGroups[static_cast<int>(gmx::MtsForceGroups::Pair)])
1303             {
1304                 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Pairs));
1305             }
1306             if (forceGroups[static_cast<int>(gmx::MtsForceGroups::Dihedral)])
1307             {
1308                 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Dihedrals));
1309             }
1310             if (forceGroups[static_cast<int>(gmx::MtsForceGroups::Angle)])
1311             {
1312                 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Angles));
1313             }
1314             if (isFirstLevel)
1315             {
1316                 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Rest));
1317                 isFirstLevel = false;
1318             }
1319             fr->listedForces.emplace_back(mtop->ffparams,
1320                                           mtop->groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1321                                           gmx_omp_nthreads_get(emntBonded),
1322                                           interactionSelection,
1323                                           fp);
1324         }
1325     }
1326     else
1327     {
1328         // Add one ListedForces object with all listed interactions
1329         fr->listedForces.emplace_back(mtop->ffparams,
1330                                       mtop->groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1331                                       gmx_omp_nthreads_get(emntBonded),
1332                                       ListedForces::interactionSelectionAll(),
1333                                       fp);
1334     }
1335
1336     // QM/MM initialization if requested
1337     if (ir->bQMMM)
1338     {
1339         gmx_incons("QM/MM was requested, but is no longer available in GROMACS");
1340     }
1341
1342     /* Set all the static charge group info */
1343     fr->cginfo_mb = init_cginfo_mb(mtop, fr);
1344     if (!DOMAINDECOMP(cr))
1345     {
1346         fr->cginfo = cginfo_expand(mtop->molblock.size(), fr->cginfo_mb);
1347     }
1348
1349     if (!DOMAINDECOMP(cr))
1350     {
1351         forcerec_set_ranges(fr, mtop->natoms, mtop->natoms, mtop->natoms);
1352     }
1353
1354     fr->print_force = print_force;
1355
1356     fr->nthread_ewc = gmx_omp_nthreads_get(emntBonded);
1357     snew(fr->ewc_t, fr->nthread_ewc);
1358
1359     if (ir->eDispCorr != edispcNO)
1360     {
1361         fr->dispersionCorrection = std::make_unique<DispersionCorrection>(
1362                 *mtop, *ir, fr->bBHAM, fr->ntype, fr->nbfp, *fr->ic, tabfn);
1363         fr->dispersionCorrection->print(mdlog);
1364     }
1365
1366     if (fp != nullptr)
1367     {
1368         /* Here we switch from using mdlog, which prints the newline before
1369          * the paragraph, to our old fprintf logging, which prints the newline
1370          * after the paragraph, so we should add a newline here.
1371          */
1372         fprintf(fp, "\n");
1373     }
1374 }
1375
1376 t_forcerec::t_forcerec() = default;
1377
1378 t_forcerec::~t_forcerec()
1379 {
1380     /* Note: This code will disappear when types are converted to C++ */
1381     sfree(shift_vec);
1382     sfree(ewc_t);
1383 }