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