0b090cd656dd484e28c801246524c3d1d207c2cb
[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,2021, 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 <cmath>
44 #include <cstdlib>
45 #include <cstring>
46
47 #include <algorithm>
48 #include <memory>
49
50 #include "gromacs/commandline/filenm.h"
51 #include "gromacs/domdec/domdec.h"
52 #include "gromacs/domdec/domdec_struct.h"
53 #include "gromacs/ewald/ewald.h"
54 #include "gromacs/ewald/pme_pp_comm_gpu.h"
55 #include "gromacs/fileio/filetypes.h"
56 #include "gromacs/gmxlib/network.h"
57 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
58 #include "gromacs/gpu_utils/gpu_utils.h"
59 #include "gromacs/hardware/hw_info.h"
60 #include "gromacs/listed_forces/listed_forces_gpu.h"
61 #include "gromacs/listed_forces/listed_forces.h"
62 #include "gromacs/listed_forces/pairs.h"
63 #include "gromacs/math/functions.h"
64 #include "gromacs/math/utilities.h"
65 #include "gromacs/math/vec.h"
66 #include "gromacs/mdlib/dispersioncorrection.h"
67 #include "gromacs/mdlib/force.h"
68 #include "gromacs/mdlib/forcerec_threading.h"
69 #include "gromacs/mdlib/gmx_omp_nthreads.h"
70 #include "gromacs/mdlib/md_support.h"
71 #include "gromacs/mdlib/wall.h"
72 #include "gromacs/mdlib/wholemoleculetransform.h"
73 #include "gromacs/mdtypes/commrec.h"
74 #include "gromacs/mdtypes/fcdata.h"
75 #include "gromacs/mdtypes/forcerec.h"
76 #include "gromacs/mdtypes/group.h"
77 #include "gromacs/mdtypes/iforceprovider.h"
78 #include "gromacs/mdtypes/inputrec.h"
79 #include "gromacs/mdtypes/interaction_const.h"
80 #include "gromacs/mdtypes/md_enums.h"
81 #include "gromacs/mdtypes/multipletimestepping.h"
82 #include "gromacs/mdtypes/nblist.h"
83 #include "gromacs/nbnxm/nbnxm.h"
84 #include "gromacs/pbcutil/ishift.h"
85 #include "gromacs/pbcutil/pbc.h"
86 #include "gromacs/tables/forcetable.h"
87 #include "gromacs/topology/mtop_util.h"
88 #include "gromacs/topology/idef.h"
89 #include "gromacs/trajectory/trajectoryframe.h"
90 #include "gromacs/utility/cstringutil.h"
91 #include "gromacs/utility/exceptions.h"
92 #include "gromacs/utility/fatalerror.h"
93 #include "gromacs/utility/gmxassert.h"
94 #include "gromacs/utility/logger.h"
95 #include "gromacs/utility/physicalnodecommunicator.h"
96 #include "gromacs/utility/pleasecite.h"
97 #include "gromacs/utility/smalloc.h"
98 #include "gromacs/utility/strconvert.h"
99
100 #include "gpuforcereduction.h"
101
102 ForceHelperBuffers::ForceHelperBuffers(bool haveDirectVirialContributions) :
103     haveDirectVirialContributions_(haveDirectVirialContributions)
104 {
105     shiftForces_.resize(gmx::c_numShiftVectors);
106 }
107
108 void ForceHelperBuffers::resize(int numAtoms)
109 {
110     if (haveDirectVirialContributions_)
111     {
112         forceBufferForDirectVirialContributions_.resize(numAtoms);
113     }
114 }
115
116 std::vector<real> makeNonBondedParameterLists(const int                      numAtomTypes,
117                                               gmx::ArrayRef<const t_iparams> iparams,
118                                               bool                           useBuckinghamPotential)
119 {
120     std::vector<real> nbfp;
121
122     if (useBuckinghamPotential)
123     {
124         nbfp.resize(3 * numAtomTypes * numAtomTypes);
125         int k = 0;
126         for (int i = 0; (i < numAtomTypes); i++)
127         {
128             for (int j = 0; (j < numAtomTypes); j++, k++)
129             {
130                 BHAMA(nbfp, numAtomTypes, i, j) = iparams[k].bham.a;
131                 BHAMB(nbfp, numAtomTypes, i, j) = iparams[k].bham.b;
132                 /* nbfp now includes the 6.0 derivative prefactor */
133                 BHAMC(nbfp, numAtomTypes, i, j) = iparams[k].bham.c * 6.0;
134             }
135         }
136     }
137     else
138     {
139         nbfp.resize(2 * numAtomTypes * numAtomTypes);
140         int k = 0;
141         for (int i = 0; (i < numAtomTypes); i++)
142         {
143             for (int j = 0; (j < numAtomTypes); j++, k++)
144             {
145                 /* nbfp now includes the 6.0/12.0 derivative prefactors */
146                 C6(nbfp, numAtomTypes, i, j)  = iparams[k].lj.c6 * 6.0;
147                 C12(nbfp, numAtomTypes, i, j) = iparams[k].lj.c12 * 12.0;
148             }
149         }
150     }
151
152     return nbfp;
153 }
154
155 std::vector<real> makeLJPmeC6GridCorrectionParameters(const int                      numAtomTypes,
156                                                       gmx::ArrayRef<const t_iparams> iparams,
157                                                       LongRangeVdW ljpme_combination_rule)
158 {
159     /* For LJ-PME simulations, we correct the energies with the reciprocal space
160      * inside of the cut-off. To do this the non-bonded kernels needs to have
161      * access to the C6-values used on the reciprocal grid in pme.c
162      */
163
164     std::vector<real> grid(2 * numAtomTypes * numAtomTypes, 0.0);
165     int               k = 0;
166     for (int i = 0; (i < numAtomTypes); i++)
167     {
168         for (int j = 0; (j < numAtomTypes); j++, k++)
169         {
170             real c6i  = iparams[i * (numAtomTypes + 1)].lj.c6;
171             real c12i = iparams[i * (numAtomTypes + 1)].lj.c12;
172             real c6j  = iparams[j * (numAtomTypes + 1)].lj.c6;
173             real c12j = iparams[j * (numAtomTypes + 1)].lj.c12;
174             real c6   = std::sqrt(c6i * c6j);
175             if (ljpme_combination_rule == LongRangeVdW::LB && !gmx_numzero(c6) && !gmx_numzero(c12i)
176                 && !gmx_numzero(c12j))
177             {
178                 real sigmai = gmx::sixthroot(c12i / c6i);
179                 real sigmaj = gmx::sixthroot(c12j / c6j);
180                 real epsi   = c6i * c6i / c12i;
181                 real epsj   = c6j * c6j / c12j;
182                 c6          = std::sqrt(epsi * epsj) * gmx::power6(0.5 * (sigmai + sigmaj));
183             }
184             /* Store the elements at the same relative positions as C6 in nbfp in order
185              * to simplify access in the kernels
186              */
187             grid[2 * (numAtomTypes * i + j)] = c6 * 6.0;
188         }
189     }
190     return grid;
191 }
192
193 //! What kind of constraint affects an atom
194 enum class ConstraintTypeForAtom : int
195 {
196     None,       //!< No constraint active
197     Constraint, //!< F_CONSTR or F_CONSTRNC active
198     Settle,     //! F_SETTLE active
199 };
200
201 static std::vector<AtomInfoWithinMoleculeBlock> makeAtomInfoForEachMoleculeBlock(const gmx_mtop_t& mtop,
202                                                                                  const t_forcerec* fr)
203 {
204     std::vector<bool> atomUsesVdw(fr->ntype, false);
205     for (int ai = 0; ai < fr->ntype; ai++)
206     {
207         for (int j = 0; j < fr->ntype; j++)
208         {
209             atomUsesVdw[ai] = atomUsesVdw[ai] || fr->haveBuckingham || C6(fr->nbfp, fr->ntype, ai, j) != 0
210                               || C12(fr->nbfp, fr->ntype, ai, j) != 0;
211         }
212     }
213
214     std::vector<AtomInfoWithinMoleculeBlock> atomInfoForEachMoleculeBlock;
215     int                                      indexOfFirstAtomInMoleculeBlock = 0;
216     for (size_t mb = 0; mb < mtop.molblock.size(); mb++)
217     {
218         const gmx_molblock_t& molb = mtop.molblock[mb];
219         const gmx_moltype_t&  molt = mtop.moltype[molb.type];
220         const auto&           excl = molt.excls;
221
222         /* Check if all molecules in this block have identical
223          * atominfo. (That's true unless some kind of group- or
224          * distance-based algorithm is involved, e.g. QM/MM groups
225          * affecting multiple molecules within a block differently.)
226          * If so, we only need an array of the size of one molecule.
227          * Otherwise we make an array of #mol times #atoms per
228          * molecule.
229          */
230         bool allMoleculesWithinBlockAreIdentical = true;
231         for (int m = 0; m < molb.nmol; m++)
232         {
233             const int numAtomsInAllMolecules = m * molt.atoms.nr;
234             for (int a = 0; a < molt.atoms.nr; a++)
235             {
236                 if (getGroupType(mtop.groups,
237                                  SimulationAtomGroupType::QuantumMechanics,
238                                  indexOfFirstAtomInMoleculeBlock + numAtomsInAllMolecules + a)
239                     != getGroupType(mtop.groups,
240                                     SimulationAtomGroupType::QuantumMechanics,
241                                     indexOfFirstAtomInMoleculeBlock + a))
242                 {
243                     allMoleculesWithinBlockAreIdentical = false;
244                 }
245                 if (!mtop.groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics].empty())
246                 {
247                     if (mtop.groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics]
248                                                 [indexOfFirstAtomInMoleculeBlock + numAtomsInAllMolecules + a]
249                         != mtop.groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics]
250                                                    [indexOfFirstAtomInMoleculeBlock + a])
251                     {
252                         allMoleculesWithinBlockAreIdentical = false;
253                     }
254                 }
255             }
256         }
257
258         AtomInfoWithinMoleculeBlock atomInfoOfMoleculeBlock;
259         atomInfoOfMoleculeBlock.indexOfFirstAtomInMoleculeBlock = indexOfFirstAtomInMoleculeBlock;
260         atomInfoOfMoleculeBlock.indexOfLastAtomInMoleculeBlock =
261                 indexOfFirstAtomInMoleculeBlock + molb.nmol * molt.atoms.nr;
262         int atomInfoSize = (allMoleculesWithinBlockAreIdentical ? 1 : molb.nmol) * molt.atoms.nr;
263         atomInfoOfMoleculeBlock.atomInfo.resize(atomInfoSize);
264
265         /* Set constraints flags for constrained atoms */
266         std::vector<ConstraintTypeForAtom> constraintTypeOfAtom(molt.atoms.nr, ConstraintTypeForAtom::None);
267         for (int ftype = 0; ftype < F_NRE; ftype++)
268         {
269             if (interaction_function[ftype].flags & IF_CONSTRAINT)
270             {
271                 const int nral = NRAL(ftype);
272                 for (int ia = 0; ia < molt.ilist[ftype].size(); ia += 1 + nral)
273                 {
274                     for (int a = 0; a < nral; a++)
275                     {
276                         constraintTypeOfAtom[molt.ilist[ftype].iatoms[ia + 1 + a]] =
277                                 (ftype == F_SETTLE ? ConstraintTypeForAtom::Settle
278                                                    : ConstraintTypeForAtom::Constraint);
279                     }
280                 }
281             }
282         }
283
284         for (int m = 0; m < (allMoleculesWithinBlockAreIdentical ? 1 : molb.nmol); m++)
285         {
286             const int moleculeOffsetInBlock = m * molt.atoms.nr;
287             for (int a = 0; a < molt.atoms.nr; a++)
288             {
289                 const t_atom& atom = molt.atoms.atom[a];
290                 int& atomInfo      = atomInfoOfMoleculeBlock.atomInfo[moleculeOffsetInBlock + a];
291
292                 /* Store the energy group in atomInfo */
293                 int gid = getGroupType(mtop.groups,
294                                        SimulationAtomGroupType::EnergyOutput,
295                                        indexOfFirstAtomInMoleculeBlock + moleculeOffsetInBlock + a);
296                 SET_CGINFO_GID(atomInfo, gid);
297
298                 bool bHaveVDW = (atomUsesVdw[atom.type] || atomUsesVdw[atom.typeB]);
299                 bool bHaveQ   = (atom.q != 0 || atom.qB != 0);
300
301                 bool haveExclusions = false;
302                 /* Loop over all the exclusions of atom ai */
303                 for (const int j : excl[a])
304                 {
305                     if (j != a)
306                     {
307                         haveExclusions = true;
308                         break;
309                     }
310                 }
311
312                 switch (constraintTypeOfAtom[a])
313                 {
314                     case ConstraintTypeForAtom::Constraint: SET_CGINFO_CONSTR(atomInfo); break;
315                     case ConstraintTypeForAtom::Settle: SET_CGINFO_SETTLE(atomInfo); break;
316                     default: break;
317                 }
318                 if (haveExclusions)
319                 {
320                     SET_CGINFO_EXCL_INTER(atomInfo);
321                 }
322                 if (bHaveVDW)
323                 {
324                     SET_CGINFO_HAS_VDW(atomInfo);
325                 }
326                 if (bHaveQ)
327                 {
328                     SET_CGINFO_HAS_Q(atomInfo);
329                 }
330                 if (fr->efep != FreeEnergyPerturbationType::No && PERTURBED(atom))
331                 {
332                     SET_CGINFO_FEP(atomInfo);
333                 }
334             }
335         }
336
337         atomInfoForEachMoleculeBlock.push_back(atomInfoOfMoleculeBlock);
338
339         indexOfFirstAtomInMoleculeBlock += molb.nmol * molt.atoms.nr;
340     }
341
342     return atomInfoForEachMoleculeBlock;
343 }
344
345 static std::vector<int> expandAtomInfo(const int                                        nmb,
346                                        gmx::ArrayRef<const AtomInfoWithinMoleculeBlock> atomInfoForEachMoleculeBlock)
347 {
348     const int numAtoms = atomInfoForEachMoleculeBlock[nmb - 1].indexOfLastAtomInMoleculeBlock;
349
350     std::vector<int> atomInfo(numAtoms);
351
352     int mb = 0;
353     for (int a = 0; a < numAtoms; a++)
354     {
355         while (a >= atomInfoForEachMoleculeBlock[mb].indexOfLastAtomInMoleculeBlock)
356         {
357             mb++;
358         }
359         atomInfo[a] = atomInfoForEachMoleculeBlock[mb]
360                               .atomInfo[(a - atomInfoForEachMoleculeBlock[mb].indexOfFirstAtomInMoleculeBlock)
361                                         % atomInfoForEachMoleculeBlock[mb].atomInfo.size()];
362     }
363
364     return atomInfo;
365 }
366
367 /* Sets the sum of charges (squared) and C6 in the system in fr.
368  * Returns whether the system has a net charge.
369  */
370 static bool set_chargesum(FILE* log, t_forcerec* fr, const gmx_mtop_t& mtop)
371 {
372     /*This now calculates sum for q and c6*/
373     double qsum, q2sum, q, c6sum, c6;
374
375     qsum  = 0;
376     q2sum = 0;
377     c6sum = 0;
378     for (const gmx_molblock_t& molb : mtop.molblock)
379     {
380         int            nmol  = molb.nmol;
381         const t_atoms* atoms = &mtop.moltype[molb.type].atoms;
382         for (int i = 0; i < atoms->nr; i++)
383         {
384             q = atoms->atom[i].q;
385             qsum += nmol * q;
386             q2sum += nmol * q * q;
387             c6 = mtop.ffparams.iparams[atoms->atom[i].type * (mtop.ffparams.atnr + 1)].lj.c6;
388             c6sum += nmol * c6;
389         }
390     }
391     fr->qsum[0]  = qsum;
392     fr->q2sum[0] = q2sum;
393     fr->c6sum[0] = c6sum;
394
395     if (fr->efep != FreeEnergyPerturbationType::No)
396     {
397         qsum  = 0;
398         q2sum = 0;
399         c6sum = 0;
400         for (const gmx_molblock_t& molb : mtop.molblock)
401         {
402             int            nmol  = molb.nmol;
403             const t_atoms* atoms = &mtop.moltype[molb.type].atoms;
404             for (int i = 0; i < atoms->nr; i++)
405             {
406                 q = atoms->atom[i].qB;
407                 qsum += nmol * q;
408                 q2sum += nmol * q * q;
409                 c6 = mtop.ffparams.iparams[atoms->atom[i].typeB * (mtop.ffparams.atnr + 1)].lj.c6;
410                 c6sum += nmol * c6;
411             }
412             fr->qsum[1]  = qsum;
413             fr->q2sum[1] = q2sum;
414             fr->c6sum[1] = c6sum;
415         }
416     }
417     else
418     {
419         fr->qsum[1]  = fr->qsum[0];
420         fr->q2sum[1] = fr->q2sum[0];
421         fr->c6sum[1] = fr->c6sum[0];
422     }
423     if (log)
424     {
425         if (fr->efep == FreeEnergyPerturbationType::No)
426         {
427             fprintf(log, "System total charge: %.3f\n", fr->qsum[0]);
428         }
429         else
430         {
431             fprintf(log, "System total charge, top. A: %.3f top. B: %.3f\n", fr->qsum[0], fr->qsum[1]);
432         }
433     }
434
435     /* A cut-off of 1e-4 is used to catch rounding errors due to ascii input */
436     return (std::abs(fr->qsum[0]) > 1e-4 || std::abs(fr->qsum[1]) > 1e-4);
437 }
438
439 /*!\brief If there's bonded interactions of type \c ftype1 or \c
440  * ftype2 present in the topology, build an array of the number of
441  * interactions present for each bonded interaction index found in the
442  * topology.
443  *
444  * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
445  * valid type with that parameter.
446  *
447  * \c count will be reallocated as necessary to fit the largest bonded
448  * interaction index found, and its current size will be returned in
449  * \c ncount. It will contain zero for every bonded interaction index
450  * for which no interactions are present in the topology.
451  */
452 static void count_tables(int ftype1, int ftype2, const gmx_mtop_t& mtop, int* ncount, int** count)
453 {
454     int ftype, i, j, tabnr;
455
456     // Loop over all moleculetypes
457     for (const gmx_moltype_t& molt : mtop.moltype)
458     {
459         // Loop over all interaction types
460         for (ftype = 0; ftype < F_NRE; ftype++)
461         {
462             // If the current interaction type is one of the types whose tables we're trying to count...
463             if (ftype == ftype1 || ftype == ftype2)
464             {
465                 const InteractionList& il     = molt.ilist[ftype];
466                 const int              stride = 1 + NRAL(ftype);
467                 // ... and there are actually some interactions for this type
468                 for (i = 0; i < il.size(); i += stride)
469                 {
470                     // Find out which table index the user wanted
471                     tabnr = mtop.ffparams.iparams[il.iatoms[i]].tab.table;
472                     if (tabnr < 0)
473                     {
474                         gmx_fatal(FARGS, "A bonded table number is smaller than 0: %d\n", tabnr);
475                     }
476                     // Make room for this index in the data structure
477                     if (tabnr >= *ncount)
478                     {
479                         srenew(*count, tabnr + 1);
480                         for (j = *ncount; j < tabnr + 1; j++)
481                         {
482                             (*count)[j] = 0;
483                         }
484                         *ncount = tabnr + 1;
485                     }
486                     // Record that this table index is used and must have a valid file
487                     (*count)[tabnr]++;
488                 }
489             }
490         }
491     }
492 }
493
494 /*!\brief If there's bonded interactions of flavour \c tabext and type
495  * \c ftype1 or \c ftype2 present in the topology, seek them in the
496  * list of filenames passed to mdrun, and make bonded tables from
497  * those files.
498  *
499  * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
500  * valid type with that parameter.
501  *
502  * A fatal error occurs if no matching filename is found.
503  */
504 static std::vector<bondedtable_t> make_bonded_tables(FILE*                            fplog,
505                                                      int                              ftype1,
506                                                      int                              ftype2,
507                                                      const gmx_mtop_t&                mtop,
508                                                      gmx::ArrayRef<const std::string> tabbfnm,
509                                                      const char*                      tabext)
510 {
511     std::vector<bondedtable_t> tab;
512
513     int  ncount = 0;
514     int* count  = nullptr;
515     count_tables(ftype1, ftype2, mtop, &ncount, &count);
516
517     // Are there any relevant tabulated bond interactions?
518     if (ncount > 0)
519     {
520         tab.resize(ncount);
521         for (int i = 0; i < ncount; i++)
522         {
523             // Do any interactions exist that requires this table?
524             if (count[i] > 0)
525             {
526                 // This pattern enforces the current requirement that
527                 // table filenames end in a characteristic sequence
528                 // before the file type extension, and avoids table 13
529                 // being recognized and used for table 1.
530                 std::string patternToFind = gmx::formatString("_%s%d.%s", tabext, i, ftp2ext(efXVG));
531                 bool        madeTable     = false;
532                 for (gmx::index j = 0; j < tabbfnm.ssize() && !madeTable; ++j)
533                 {
534                     if (gmx::endsWith(tabbfnm[j], patternToFind))
535                     {
536                         // Finally read the table from the file found
537                         tab[i]    = make_bonded_table(fplog, tabbfnm[j].c_str(), NRAL(ftype1) - 2);
538                         madeTable = true;
539                     }
540                 }
541                 if (!madeTable)
542                 {
543                     bool isPlural = (ftype2 != -1);
544                     gmx_fatal(FARGS,
545                               "Tabulated interaction of type '%s%s%s' with index %d cannot be used "
546                               "because no table file whose name matched '%s' was passed via the "
547                               "gmx mdrun -tableb command-line option.",
548                               interaction_function[ftype1].longname,
549                               isPlural ? "' or '" : "",
550                               isPlural ? interaction_function[ftype2].longname : "",
551                               i,
552                               patternToFind.c_str());
553                 }
554             }
555         }
556         sfree(count);
557     }
558
559     return tab;
560 }
561
562 void forcerec_set_ranges(t_forcerec* fr, int natoms_force, int natoms_force_constr, int natoms_f_novirsum)
563 {
564     fr->natoms_force        = natoms_force;
565     fr->natoms_force_constr = natoms_force_constr;
566
567     for (auto& forceHelperBuffers : fr->forceHelperBuffers)
568     {
569         forceHelperBuffers.resize(natoms_f_novirsum);
570     }
571 }
572
573 /* Generate Coulomb and/or Van der Waals Ewald long-range correction tables
574  *
575  * Tables are generated for one or both, depending on if the pointers
576  * are non-null. The spacing for both table sets is the same and obeys
577  * both accuracy requirements, when relevant.
578  */
579 static void init_ewald_f_table(const interaction_const_t& ic,
580                                const real                 rlist,
581                                const real                 tabext,
582                                EwaldCorrectionTables*     coulombTables,
583                                EwaldCorrectionTables*     vdwTables)
584 {
585     const bool useCoulombTable = (EEL_PME_EWALD(ic.eeltype) && coulombTables != nullptr);
586     const bool useVdwTable     = (EVDW_PME(ic.vdwtype) && vdwTables != nullptr);
587
588     /* Get the Ewald table spacing based on Coulomb and/or LJ
589      * Ewald coefficients and rtol.
590      */
591     const real tableScale = ewald_spline3_table_scale(ic, useCoulombTable, useVdwTable);
592
593     const bool havePerturbedNonbondeds = (ic.softCoreParameters != nullptr);
594
595     real tableLen = ic.rcoulomb;
596     if ((useCoulombTable || useVdwTable) && havePerturbedNonbondeds && rlist + tabext > 0.0)
597     {
598         /* TODO: Ideally this should also check if couple-intramol == no, but that isn't
599          * stored in ir. Grompp puts that info into an opts structure that doesn't make it into the tpr.
600          * The alternative is to look through all the exclusions and check if they come from
601          * couple-intramol == no. Meanwhile, always having larger tables should only affect
602          * memory consumption, not speed (barring cache issues).
603          */
604         tableLen = rlist + tabext;
605     }
606     const int tableSize = static_cast<int>(tableLen * tableScale) + 2;
607
608     if (useCoulombTable)
609     {
610         *coulombTables =
611                 generateEwaldCorrectionTables(tableSize, tableScale, ic.ewaldcoeff_q, v_q_ewald_lr);
612     }
613
614     if (useVdwTable)
615     {
616         *vdwTables = generateEwaldCorrectionTables(tableSize, tableScale, ic.ewaldcoeff_lj, v_lj_ewald_lr);
617     }
618 }
619
620 void init_interaction_const_tables(FILE* fp, interaction_const_t* ic, const real rlist, const real tableExtensionLength)
621 {
622     if (EEL_PME_EWALD(ic->eeltype) || EVDW_PME(ic->vdwtype))
623     {
624         init_ewald_f_table(
625                 *ic, rlist, tableExtensionLength, ic->coulombEwaldTables.get(), ic->vdwEwaldTables.get());
626         if (fp != nullptr)
627         {
628             fprintf(fp,
629                     "Initialized non-bonded Ewald tables, spacing: %.2e size: %zu\n\n",
630                     1 / ic->coulombEwaldTables->scale,
631                     ic->coulombEwaldTables->tableF.size());
632         }
633     }
634 }
635
636 real cutoff_inf(real cutoff)
637 {
638     if (cutoff == 0)
639     {
640         cutoff = GMX_CUTOFF_INF;
641     }
642
643     return cutoff;
644 }
645
646 void init_forcerec(FILE*                            fplog,
647                    const gmx::MDLogger&             mdlog,
648                    t_forcerec*                      forcerec,
649                    const t_inputrec&                inputrec,
650                    const gmx_mtop_t&                mtop,
651                    const t_commrec*                 commrec,
652                    matrix                           box,
653                    const char*                      tabfn,
654                    const char*                      tabpfn,
655                    gmx::ArrayRef<const std::string> tabbfnm,
656                    real                             print_force)
657 {
658     /* The CMake default turns SIMD kernels on, but it might be turned off further down... */
659     forcerec->use_simd_kernels = GMX_USE_SIMD_KERNELS;
660
661     if (check_box(inputrec.pbcType, box))
662     {
663         gmx_fatal(FARGS, "%s", check_box(inputrec.pbcType, box));
664     }
665
666     /* Test particle insertion ? */
667     if (EI_TPI(inputrec.eI))
668     {
669         /* Set to the size of the molecule to be inserted (the last one) */
670         gmx::RangePartitioning molecules = gmx_mtop_molecules(mtop);
671         forcerec->n_tpi                  = molecules.block(molecules.numBlocks() - 1).size();
672     }
673     else
674     {
675         forcerec->n_tpi = 0;
676     }
677
678     if (inputrec.coulombtype == CoulombInteractionType::RFNecUnsupported
679         || inputrec.coulombtype == CoulombInteractionType::GRFNotused)
680     {
681         gmx_fatal(FARGS, "%s electrostatics is no longer supported", enumValueToString(inputrec.coulombtype));
682     }
683
684     if (inputrec.bAdress)
685     {
686         gmx_fatal(FARGS, "AdResS simulations are no longer supported");
687     }
688     if (inputrec.useTwinRange)
689     {
690         gmx_fatal(FARGS, "Twin-range simulations are no longer supported");
691     }
692     /* Copy the user determined parameters */
693     forcerec->userint1  = inputrec.userint1;
694     forcerec->userint2  = inputrec.userint2;
695     forcerec->userint3  = inputrec.userint3;
696     forcerec->userint4  = inputrec.userint4;
697     forcerec->userreal1 = inputrec.userreal1;
698     forcerec->userreal2 = inputrec.userreal2;
699     forcerec->userreal3 = inputrec.userreal3;
700     forcerec->userreal4 = inputrec.userreal4;
701
702     /* Shell stuff */
703     forcerec->fc_stepsize = inputrec.fc_stepsize;
704
705     /* Free energy */
706     forcerec->efep = inputrec.efep;
707
708     if ((getenv("GMX_DISABLE_SIMD_KERNELS") != nullptr) || (getenv("GMX_NOOPTIMIZEDKERNELS") != nullptr))
709     {
710         forcerec->use_simd_kernels = FALSE;
711         if (fplog != nullptr)
712         {
713             fprintf(fplog,
714                     "\nFound environment variable GMX_DISABLE_SIMD_KERNELS.\n"
715                     "Disabling the usage of any SIMD-specific non-bonded & bonded kernel routines\n"
716                     "(e.g. SSE2/SSE4.1/AVX).\n\n");
717         }
718     }
719
720     forcerec->haveBuckingham = (mtop.ffparams.functype[0] == F_BHAM);
721
722     /* Neighbour searching stuff */
723     forcerec->pbcType = inputrec.pbcType;
724
725     /* Determine if we will do PBC for distances in bonded interactions */
726     if (forcerec->pbcType == PbcType::No)
727     {
728         forcerec->bMolPBC = FALSE;
729     }
730     else
731     {
732         const bool useEwaldSurfaceCorrection =
733                 (EEL_PME_EWALD(inputrec.coulombtype) && inputrec.epsilon_surface != 0);
734         const bool haveOrientationRestraints = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
735         if (!DOMAINDECOMP(commrec))
736         {
737             forcerec->bMolPBC = true;
738
739             if (useEwaldSurfaceCorrection || haveOrientationRestraints)
740             {
741                 forcerec->wholeMoleculeTransform =
742                         std::make_unique<gmx::WholeMoleculeTransform>(mtop, inputrec.pbcType);
743             }
744         }
745         else
746         {
747             forcerec->bMolPBC = dd_bonded_molpbc(*commrec->dd, forcerec->pbcType);
748
749             /* With Ewald surface correction it is useful to support e.g. running water
750              * in parallel with update groups.
751              * With orientation restraints there is no sensible use case supported with DD.
752              */
753             if ((useEwaldSurfaceCorrection && !dd_moleculesAreAlwaysWhole(*commrec->dd))
754                 || haveOrientationRestraints)
755             {
756                 gmx_fatal(FARGS,
757                           "You requested Ewald surface correction or orientation restraints, "
758                           "but molecules are broken "
759                           "over periodic boundary conditions by the domain decomposition. "
760                           "Run without domain decomposition instead.");
761             }
762         }
763
764         if (useEwaldSurfaceCorrection)
765         {
766             GMX_RELEASE_ASSERT(!DOMAINDECOMP(commrec) || dd_moleculesAreAlwaysWhole(*commrec->dd),
767                                "Molecules can not be broken by PBC with epsilon_surface > 0");
768         }
769     }
770
771     forcerec->rc_scaling = inputrec.refcoord_scaling;
772     copy_rvec(inputrec.posres_com, forcerec->posres_com);
773     copy_rvec(inputrec.posres_comB, forcerec->posres_comB);
774     forcerec->rlist                  = cutoff_inf(inputrec.rlist);
775     forcerec->ljpme_combination_rule = inputrec.ljpme_combination_rule;
776
777     /* This now calculates sum for q and c6*/
778     bool systemHasNetCharge = set_chargesum(fplog, forcerec, mtop);
779
780     /* Make data structure used by kernels */
781     forcerec->ic = std::make_unique<interaction_const_t>(
782             init_interaction_const(fplog, inputrec, mtop, systemHasNetCharge));
783     init_interaction_const_tables(fplog, forcerec->ic.get(), forcerec->rlist, inputrec.tabext);
784
785     const interaction_const_t* interactionConst = forcerec->ic.get();
786
787     /* TODO: Replace this Ewald table or move it into interaction_const_t */
788     if (inputrec.coulombtype == CoulombInteractionType::Ewald)
789     {
790         forcerec->ewald_table = std::make_unique<gmx_ewald_tab_t>(inputrec, fplog);
791     }
792
793     /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
794     switch (interactionConst->eeltype)
795     {
796         case CoulombInteractionType::Cut:
797             forcerec->nbkernel_elec_interaction = NbkernelElecType::Coulomb;
798             break;
799
800         case CoulombInteractionType::RF:
801         case CoulombInteractionType::RFZero:
802             forcerec->nbkernel_elec_interaction = NbkernelElecType::ReactionField;
803             break;
804
805         case CoulombInteractionType::Switch:
806         case CoulombInteractionType::Shift:
807         case CoulombInteractionType::User:
808         case CoulombInteractionType::PmeSwitch:
809         case CoulombInteractionType::PmeUser:
810         case CoulombInteractionType::PmeUserSwitch:
811             forcerec->nbkernel_elec_interaction = NbkernelElecType::CubicSplineTable;
812             break;
813
814         case CoulombInteractionType::Pme:
815         case CoulombInteractionType::P3mAD:
816         case CoulombInteractionType::Ewald:
817             forcerec->nbkernel_elec_interaction = NbkernelElecType::Ewald;
818             break;
819
820         default:
821             gmx_fatal(FARGS,
822                       "Unsupported electrostatic interaction: %s",
823                       enumValueToString(interactionConst->eeltype));
824     }
825     forcerec->nbkernel_elec_modifier = interactionConst->coulomb_modifier;
826
827     /* Vdw: Translate from mdp settings to kernel format */
828     switch (interactionConst->vdwtype)
829     {
830         case VanDerWaalsType::Cut:
831             if (forcerec->haveBuckingham)
832             {
833                 forcerec->nbkernel_vdw_interaction = NbkernelVdwType::Buckingham;
834             }
835             else
836             {
837                 forcerec->nbkernel_vdw_interaction = NbkernelVdwType::LennardJones;
838             }
839             break;
840         case VanDerWaalsType::Pme:
841             forcerec->nbkernel_vdw_interaction = NbkernelVdwType::LJEwald;
842             break;
843
844         case VanDerWaalsType::Switch:
845         case VanDerWaalsType::Shift:
846         case VanDerWaalsType::User:
847             forcerec->nbkernel_vdw_interaction = NbkernelVdwType::CubicSplineTable;
848             break;
849
850         default:
851             gmx_fatal(FARGS, "Unsupported vdw interaction: %s", enumValueToString(interactionConst->vdwtype));
852     }
853     forcerec->nbkernel_vdw_modifier = interactionConst->vdw_modifier;
854
855     if (!gmx_within_tol(interactionConst->reppow, 12.0, 10 * GMX_DOUBLE_EPS))
856     {
857         gmx_fatal(FARGS, "Only LJ repulsion power 12 is supported");
858     }
859     /* Older tpr files can contain Coulomb user tables with the Verlet cutoff-scheme,
860      * while mdrun does not (and never did) support this.
861      */
862     if (EEL_USER(forcerec->ic->eeltype))
863     {
864         gmx_fatal(FARGS,
865                   "Electrostatics type %s is currently not supported",
866                   enumValueToString(inputrec.coulombtype));
867     }
868
869     /* 1-4 interaction electrostatics */
870     forcerec->fudgeQQ = mtop.ffparams.fudgeQQ;
871
872     // Multiple time stepping
873     forcerec->useMts = inputrec.useMts;
874
875     if (forcerec->useMts)
876     {
877         GMX_ASSERT(gmx::checkMtsRequirements(inputrec).empty(),
878                    "All MTS requirements should be met here");
879     }
880
881     const bool haveDirectVirialContributionsFast =
882             forcerec->forceProviders->hasForceProvider() || gmx_mtop_ftype_count(mtop, F_POSRES) > 0
883             || gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0 || inputrec.nwall > 0 || inputrec.bPull
884             || inputrec.bRot || inputrec.bIMD;
885     const bool haveDirectVirialContributionsSlow =
886             EEL_FULL(interactionConst->eeltype) || EVDW_PME(interactionConst->vdwtype);
887     for (int i = 0; i < (forcerec->useMts ? 2 : 1); i++)
888     {
889         bool haveDirectVirialContributions =
890                 (((!forcerec->useMts || i == 0) && haveDirectVirialContributionsFast)
891                  || ((!forcerec->useMts || i == 1) && haveDirectVirialContributionsSlow));
892         forcerec->forceHelperBuffers.emplace_back(haveDirectVirialContributions);
893     }
894
895     if (forcerec->shift_vec.empty())
896     {
897         forcerec->shift_vec.resize(gmx::c_numShiftVectors);
898     }
899
900     GMX_ASSERT(forcerec->nbfp.empty(), "The nonbonded force parameters should not be set up yet.");
901     forcerec->ntype = mtop.ffparams.atnr;
902     forcerec->nbfp  = makeNonBondedParameterLists(
903             mtop.ffparams.atnr, mtop.ffparams.iparams, forcerec->haveBuckingham);
904     if (EVDW_PME(interactionConst->vdwtype))
905     {
906         forcerec->ljpme_c6grid = makeLJPmeC6GridCorrectionParameters(
907                 mtop.ffparams.atnr, mtop.ffparams.iparams, forcerec->ljpme_combination_rule);
908     }
909
910     /* Copy the energy group exclusions */
911     forcerec->egp_flags = inputrec.opts.egp_flags;
912
913     /* Van der Waals stuff */
914     if ((interactionConst->vdwtype != VanDerWaalsType::Cut)
915         && (interactionConst->vdwtype != VanDerWaalsType::User) && !forcerec->haveBuckingham)
916     {
917         if (interactionConst->rvdw_switch >= interactionConst->rvdw)
918         {
919             gmx_fatal(FARGS,
920                       "rvdw_switch (%f) must be < rvdw (%f)",
921                       interactionConst->rvdw_switch,
922                       interactionConst->rvdw);
923         }
924         if (fplog)
925         {
926             fprintf(fplog,
927                     "Using %s Lennard-Jones, switch between %g and %g nm\n",
928                     (interactionConst->eeltype == CoulombInteractionType::Switch) ? "switched" : "shifted",
929                     interactionConst->rvdw_switch,
930                     interactionConst->rvdw);
931         }
932     }
933
934     if (forcerec->haveBuckingham && EVDW_PME(interactionConst->vdwtype))
935     {
936         gmx_fatal(FARGS, "LJ PME not supported with Buckingham");
937     }
938
939     if (forcerec->haveBuckingham
940         && (interactionConst->vdwtype == VanDerWaalsType::Shift
941             || interactionConst->vdwtype == VanDerWaalsType::Switch))
942     {
943         gmx_fatal(FARGS, "Switch/shift interaction not supported with Buckingham");
944     }
945
946     if (forcerec->haveBuckingham)
947     {
948         gmx_fatal(FARGS, "The Verlet cutoff-scheme does not (yet) support Buckingham");
949     }
950
951     if (inputrec.implicit_solvent)
952     {
953         gmx_fatal(FARGS, "Implict solvation is no longer supported.");
954     }
955
956
957     /* This code automatically gives table length tabext without cut-off's,
958      * in that case grompp should already have checked that we do not need
959      * normal tables and we only generate tables for 1-4 interactions.
960      */
961     real rtab = inputrec.rlist + inputrec.tabext;
962
963     /* We want to use unmodified tables for 1-4 coulombic
964      * interactions, so we must in general have an extra set of
965      * tables. */
966     if (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 || gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0
967         || gmx_mtop_ftype_count(mtop, F_LJC_PAIRS_NB) > 0)
968     {
969         forcerec->pairsTable = make_tables(fplog, interactionConst, tabpfn, rtab, GMX_MAKETABLES_14ONLY);
970     }
971
972     /* Wall stuff */
973     forcerec->nwall = inputrec.nwall;
974     if (inputrec.nwall && inputrec.wall_type == WallType::Table)
975     {
976         make_wall_tables(fplog, inputrec, tabfn, &mtop.groups, forcerec);
977     }
978
979     forcerec->fcdata = std::make_unique<t_fcdata>();
980
981     if (!tabbfnm.empty())
982     {
983         t_fcdata& fcdata = *forcerec->fcdata;
984         // Need to catch std::bad_alloc
985         // TODO Don't need to catch this here, when merging with master branch
986         try
987         {
988             // TODO move these tables into a separate struct and store reference in ListedForces
989             fcdata.bondtab = make_bonded_tables(fplog, F_TABBONDS, F_TABBONDSNC, mtop, tabbfnm, "b");
990             fcdata.angletab = make_bonded_tables(fplog, F_TABANGLES, -1, mtop, tabbfnm, "a");
991             fcdata.dihtab   = make_bonded_tables(fplog, F_TABDIHS, -1, mtop, tabbfnm, "d");
992         }
993         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
994     }
995     else
996     {
997         if (debug)
998         {
999             fprintf(debug,
1000                     "No fcdata or table file name passed, can not read table, can not do bonded "
1001                     "interactions\n");
1002         }
1003     }
1004
1005     /* Initialize the thread working data for bonded interactions */
1006     if (forcerec->useMts)
1007     {
1008         // Add one ListedForces object for each MTS level
1009         bool isFirstLevel = true;
1010         for (const auto& mtsLevel : inputrec.mtsLevels)
1011         {
1012             ListedForces::InteractionSelection interactionSelection;
1013             const auto&                        forceGroups = mtsLevel.forceGroups;
1014             if (forceGroups[static_cast<int>(gmx::MtsForceGroups::Pair)])
1015             {
1016                 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Pairs));
1017             }
1018             if (forceGroups[static_cast<int>(gmx::MtsForceGroups::Dihedral)])
1019             {
1020                 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Dihedrals));
1021             }
1022             if (forceGroups[static_cast<int>(gmx::MtsForceGroups::Angle)])
1023             {
1024                 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Angles));
1025             }
1026             if (isFirstLevel)
1027             {
1028                 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Rest));
1029                 isFirstLevel = false;
1030             }
1031             forcerec->listedForces.emplace_back(
1032                     mtop.ffparams,
1033                     mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1034                     gmx_omp_nthreads_get(ModuleMultiThread::Bonded),
1035                     interactionSelection,
1036                     fplog);
1037         }
1038     }
1039     else
1040     {
1041         // Add one ListedForces object with all listed interactions
1042         forcerec->listedForces.emplace_back(
1043                 mtop.ffparams,
1044                 mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1045                 gmx_omp_nthreads_get(ModuleMultiThread::Bonded),
1046                 ListedForces::interactionSelectionAll(),
1047                 fplog);
1048     }
1049
1050     // QM/MM initialization if requested
1051     if (inputrec.bQMMM)
1052     {
1053         gmx_incons("QM/MM was requested, but is no longer available in GROMACS");
1054     }
1055
1056     /* Set all the static charge group info */
1057     forcerec->atomInfoForEachMoleculeBlock = makeAtomInfoForEachMoleculeBlock(mtop, forcerec);
1058     if (!DOMAINDECOMP(commrec))
1059     {
1060         forcerec->atomInfo = expandAtomInfo(mtop.molblock.size(), forcerec->atomInfoForEachMoleculeBlock);
1061     }
1062
1063     if (!DOMAINDECOMP(commrec))
1064     {
1065         forcerec_set_ranges(forcerec, mtop.natoms, mtop.natoms, mtop.natoms);
1066     }
1067
1068     forcerec->print_force = print_force;
1069
1070     forcerec->nthread_ewc = gmx_omp_nthreads_get(ModuleMultiThread::Bonded);
1071     forcerec->ewc_t.resize(forcerec->nthread_ewc);
1072
1073     if (inputrec.eDispCorr != DispersionCorrectionType::No)
1074     {
1075         forcerec->dispersionCorrection = std::make_unique<DispersionCorrection>(
1076                 mtop, inputrec, forcerec->haveBuckingham, forcerec->ntype, forcerec->nbfp, *forcerec->ic, tabfn);
1077         forcerec->dispersionCorrection->print(mdlog);
1078     }
1079
1080     if (fplog != nullptr)
1081     {
1082         /* Here we switch from using mdlog, which prints the newline before
1083          * the paragraph, to our old fprintf logging, which prints the newline
1084          * after the paragraph, so we should add a newline here.
1085          */
1086         fprintf(fplog, "\n");
1087     }
1088 }
1089
1090 t_forcerec::t_forcerec() = default;
1091
1092 t_forcerec::~t_forcerec() = default;