SYCL: Avoid using no_init read accessor in rocFFT
[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/gmx_omp_nthreads.h"
69 #include "gromacs/mdlib/md_support.h"
70 #include "gromacs/mdlib/wall.h"
71 #include "gromacs/mdlib/wholemoleculetransform.h"
72 #include "gromacs/mdtypes/commrec.h"
73 #include "gromacs/mdtypes/fcdata.h"
74 #include "gromacs/mdtypes/forcerec.h"
75 #include "gromacs/mdtypes/group.h"
76 #include "gromacs/mdtypes/iforceprovider.h"
77 #include "gromacs/mdtypes/inputrec.h"
78 #include "gromacs/mdtypes/interaction_const.h"
79 #include "gromacs/mdtypes/md_enums.h"
80 #include "gromacs/mdtypes/multipletimestepping.h"
81 #include "gromacs/mdtypes/nblist.h"
82 #include "gromacs/mdtypes/simulation_workload.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<gmx::AtomInfoWithinMoleculeBlock>
202 makeAtomInfoForEachMoleculeBlock(const gmx_mtop_t& mtop, 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<gmx::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         gmx::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                 int64_t& 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                 atomInfo = (atomInfo & ~gmx::sc_atomInfo_EnergyGroupIdMask) | 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:
315                         atomInfo |= gmx::sc_atomInfo_Constraint;
316                         break;
317                     case ConstraintTypeForAtom::Settle: atomInfo |= gmx::sc_atomInfo_Settle; break;
318                     default: break;
319                 }
320                 if (haveExclusions)
321                 {
322                     atomInfo |= gmx::sc_atomInfo_Exclusion;
323                 }
324                 if (bHaveVDW)
325                 {
326                     atomInfo |= gmx::sc_atomInfo_HasVdw;
327                 }
328                 if (bHaveQ)
329                 {
330                     atomInfo |= gmx::sc_atomInfo_HasCharge;
331                 }
332                 if (fr->efep != FreeEnergyPerturbationType::No)
333                 {
334                     if (PERTURBED(atom))
335                     {
336                         atomInfo |= gmx::sc_atomInfo_FreeEnergyPerturbation;
337                     }
338                     if (atomHasPerturbedChargeIn14Interaction(a, molt))
339                     {
340                         atomInfo |= gmx::sc_atomInfo_HasPerturbedChargeIn14Interaction;
341                     }
342                 }
343             }
344         }
345
346         atomInfoForEachMoleculeBlock.push_back(atomInfoOfMoleculeBlock);
347
348         indexOfFirstAtomInMoleculeBlock += molb.nmol * molt.atoms.nr;
349     }
350
351     return atomInfoForEachMoleculeBlock;
352 }
353
354 static std::vector<int64_t> expandAtomInfo(const int nmb,
355                                            gmx::ArrayRef<const gmx::AtomInfoWithinMoleculeBlock> atomInfoForEachMoleculeBlock)
356 {
357     const int numAtoms = atomInfoForEachMoleculeBlock[nmb - 1].indexOfLastAtomInMoleculeBlock;
358
359     std::vector<int64_t> atomInfo(numAtoms);
360
361     int mb = 0;
362     for (int a = 0; a < numAtoms; a++)
363     {
364         while (a >= atomInfoForEachMoleculeBlock[mb].indexOfLastAtomInMoleculeBlock)
365         {
366             mb++;
367         }
368         atomInfo[a] = atomInfoForEachMoleculeBlock[mb]
369                               .atomInfo[(a - atomInfoForEachMoleculeBlock[mb].indexOfFirstAtomInMoleculeBlock)
370                                         % atomInfoForEachMoleculeBlock[mb].atomInfo.size()];
371     }
372
373     return atomInfo;
374 }
375
376 /* Sets the sum of charges (squared) and C6 in the system in fr.
377  * Returns whether the system has a net charge.
378  */
379 static bool set_chargesum(FILE* log, t_forcerec* fr, const gmx_mtop_t& mtop)
380 {
381     /*This now calculates sum for q and c6*/
382     double qsum, q2sum, q, c6sum, c6;
383
384     qsum  = 0;
385     q2sum = 0;
386     c6sum = 0;
387     for (const gmx_molblock_t& molb : mtop.molblock)
388     {
389         int            nmol  = molb.nmol;
390         const t_atoms* atoms = &mtop.moltype[molb.type].atoms;
391         for (int i = 0; i < atoms->nr; i++)
392         {
393             q = atoms->atom[i].q;
394             qsum += nmol * q;
395             q2sum += nmol * q * q;
396             c6 = mtop.ffparams.iparams[atoms->atom[i].type * (mtop.ffparams.atnr + 1)].lj.c6;
397             c6sum += nmol * c6;
398         }
399     }
400     fr->qsum[0]  = qsum;
401     fr->q2sum[0] = q2sum;
402     fr->c6sum[0] = c6sum;
403
404     if (fr->efep != FreeEnergyPerturbationType::No)
405     {
406         qsum  = 0;
407         q2sum = 0;
408         c6sum = 0;
409         for (const gmx_molblock_t& molb : mtop.molblock)
410         {
411             int            nmol  = molb.nmol;
412             const t_atoms* atoms = &mtop.moltype[molb.type].atoms;
413             for (int i = 0; i < atoms->nr; i++)
414             {
415                 q = atoms->atom[i].qB;
416                 qsum += nmol * q;
417                 q2sum += nmol * q * q;
418                 c6 = mtop.ffparams.iparams[atoms->atom[i].typeB * (mtop.ffparams.atnr + 1)].lj.c6;
419                 c6sum += nmol * c6;
420             }
421             fr->qsum[1]  = qsum;
422             fr->q2sum[1] = q2sum;
423             fr->c6sum[1] = c6sum;
424         }
425     }
426     else
427     {
428         fr->qsum[1]  = fr->qsum[0];
429         fr->q2sum[1] = fr->q2sum[0];
430         fr->c6sum[1] = fr->c6sum[0];
431     }
432     if (log)
433     {
434         if (fr->efep == FreeEnergyPerturbationType::No)
435         {
436             fprintf(log, "System total charge: %.3f\n", fr->qsum[0]);
437         }
438         else
439         {
440             fprintf(log, "System total charge, top. A: %.3f top. B: %.3f\n", fr->qsum[0], fr->qsum[1]);
441         }
442     }
443
444     /* A cut-off of 1e-4 is used to catch rounding errors due to ascii input */
445     return (std::abs(fr->qsum[0]) > 1e-4 || std::abs(fr->qsum[1]) > 1e-4);
446 }
447
448 /*!\brief If there's bonded interactions of type \c ftype1 or \c
449  * ftype2 present in the topology, build an array of the number of
450  * interactions present for each bonded interaction index found in the
451  * topology.
452  *
453  * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
454  * valid type with that parameter.
455  *
456  * \c count will be reallocated as necessary to fit the largest bonded
457  * interaction index found, and its current size will be returned in
458  * \c ncount. It will contain zero for every bonded interaction index
459  * for which no interactions are present in the topology.
460  */
461 static void count_tables(int ftype1, int ftype2, const gmx_mtop_t& mtop, int* ncount, int** count)
462 {
463     int ftype, i, j, tabnr;
464
465     // Loop over all moleculetypes
466     for (const gmx_moltype_t& molt : mtop.moltype)
467     {
468         // Loop over all interaction types
469         for (ftype = 0; ftype < F_NRE; ftype++)
470         {
471             // If the current interaction type is one of the types whose tables we're trying to count...
472             if (ftype == ftype1 || ftype == ftype2)
473             {
474                 const InteractionList& il     = molt.ilist[ftype];
475                 const int              stride = 1 + NRAL(ftype);
476                 // ... and there are actually some interactions for this type
477                 for (i = 0; i < il.size(); i += stride)
478                 {
479                     // Find out which table index the user wanted
480                     tabnr = mtop.ffparams.iparams[il.iatoms[i]].tab.table;
481                     if (tabnr < 0)
482                     {
483                         gmx_fatal(FARGS, "A bonded table number is smaller than 0: %d\n", tabnr);
484                     }
485                     // Make room for this index in the data structure
486                     if (tabnr >= *ncount)
487                     {
488                         srenew(*count, tabnr + 1);
489                         for (j = *ncount; j < tabnr + 1; j++)
490                         {
491                             (*count)[j] = 0;
492                         }
493                         *ncount = tabnr + 1;
494                     }
495                     // Record that this table index is used and must have a valid file
496                     (*count)[tabnr]++;
497                 }
498             }
499         }
500     }
501 }
502
503 /*!\brief If there's bonded interactions of flavour \c tabext and type
504  * \c ftype1 or \c ftype2 present in the topology, seek them in the
505  * list of filenames passed to mdrun, and make bonded tables from
506  * those files.
507  *
508  * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
509  * valid type with that parameter.
510  *
511  * A fatal error occurs if no matching filename is found.
512  */
513 static std::vector<bondedtable_t> make_bonded_tables(FILE*                            fplog,
514                                                      int                              ftype1,
515                                                      int                              ftype2,
516                                                      const gmx_mtop_t&                mtop,
517                                                      gmx::ArrayRef<const std::string> tabbfnm,
518                                                      const char*                      tabext)
519 {
520     std::vector<bondedtable_t> tab;
521
522     int  ncount = 0;
523     int* count  = nullptr;
524     count_tables(ftype1, ftype2, mtop, &ncount, &count);
525
526     // Are there any relevant tabulated bond interactions?
527     if (ncount > 0)
528     {
529         tab.resize(ncount);
530         for (int i = 0; i < ncount; i++)
531         {
532             // Do any interactions exist that requires this table?
533             if (count[i] > 0)
534             {
535                 // This pattern enforces the current requirement that
536                 // table filenames end in a characteristic sequence
537                 // before the file type extension, and avoids table 13
538                 // being recognized and used for table 1.
539                 std::string patternToFind = gmx::formatString("_%s%d.%s", tabext, i, ftp2ext(efXVG));
540                 bool        madeTable     = false;
541                 for (gmx::index j = 0; j < tabbfnm.ssize() && !madeTable; ++j)
542                 {
543                     if (gmx::endsWith(tabbfnm[j], patternToFind))
544                     {
545                         // Finally read the table from the file found
546                         tab[i]    = make_bonded_table(fplog, tabbfnm[j].c_str(), NRAL(ftype1) - 2);
547                         madeTable = true;
548                     }
549                 }
550                 if (!madeTable)
551                 {
552                     bool isPlural = (ftype2 != -1);
553                     gmx_fatal(FARGS,
554                               "Tabulated interaction of type '%s%s%s' with index %d cannot be used "
555                               "because no table file whose name matched '%s' was passed via the "
556                               "gmx mdrun -tableb command-line option.",
557                               interaction_function[ftype1].longname,
558                               isPlural ? "' or '" : "",
559                               isPlural ? interaction_function[ftype2].longname : "",
560                               i,
561                               patternToFind.c_str());
562                 }
563             }
564         }
565         sfree(count);
566     }
567
568     return tab;
569 }
570
571 void forcerec_set_ranges(t_forcerec* fr, int natoms_force, int natoms_force_constr, int natoms_f_novirsum)
572 {
573     fr->natoms_force        = natoms_force;
574     fr->natoms_force_constr = natoms_force_constr;
575
576     for (auto& forceHelperBuffers : fr->forceHelperBuffers)
577     {
578         forceHelperBuffers.resize(natoms_f_novirsum);
579     }
580 }
581
582 /* Generate Coulomb and/or Van der Waals Ewald long-range correction tables
583  *
584  * Tables are generated for one or both, depending on if the pointers
585  * are non-null. The spacing for both table sets is the same and obeys
586  * both accuracy requirements, when relevant.
587  */
588 static void init_ewald_f_table(const interaction_const_t& ic,
589                                const real                 rlist,
590                                const real                 tabext,
591                                EwaldCorrectionTables*     coulombTables,
592                                EwaldCorrectionTables*     vdwTables)
593 {
594     const bool useCoulombTable = (EEL_PME_EWALD(ic.eeltype) && coulombTables != nullptr);
595     const bool useVdwTable     = (EVDW_PME(ic.vdwtype) && vdwTables != nullptr);
596
597     /* Get the Ewald table spacing based on Coulomb and/or LJ
598      * Ewald coefficients and rtol.
599      */
600     const real tableScale = ewald_spline3_table_scale(ic, useCoulombTable, useVdwTable);
601
602     const bool havePerturbedNonbondeds = (ic.softCoreParameters != nullptr);
603
604     real tableLen = ic.rcoulomb;
605     if ((useCoulombTable || useVdwTable) && havePerturbedNonbondeds && rlist + tabext > 0.0)
606     {
607         /* TODO: Ideally this should also check if couple-intramol == no, but that isn't
608          * stored in ir. Grompp puts that info into an opts structure that doesn't make it into the tpr.
609          * The alternative is to look through all the exclusions and check if they come from
610          * couple-intramol == no. Meanwhile, always having larger tables should only affect
611          * memory consumption, not speed (barring cache issues).
612          */
613         tableLen = rlist + tabext;
614     }
615     const int tableSize = static_cast<int>(tableLen * tableScale) + 2;
616
617     if (useCoulombTable)
618     {
619         *coulombTables =
620                 generateEwaldCorrectionTables(tableSize, tableScale, ic.ewaldcoeff_q, v_q_ewald_lr);
621     }
622
623     if (useVdwTable)
624     {
625         *vdwTables = generateEwaldCorrectionTables(tableSize, tableScale, ic.ewaldcoeff_lj, v_lj_ewald_lr);
626     }
627 }
628
629 void init_interaction_const_tables(FILE* fp, interaction_const_t* ic, const real rlist, const real tableExtensionLength)
630 {
631     if (EEL_PME_EWALD(ic->eeltype) || EVDW_PME(ic->vdwtype))
632     {
633         init_ewald_f_table(
634                 *ic, rlist, tableExtensionLength, ic->coulombEwaldTables.get(), ic->vdwEwaldTables.get());
635         if (fp != nullptr)
636         {
637             fprintf(fp,
638                     "Initialized non-bonded Ewald tables, spacing: %.2e size: %zu\n\n",
639                     1 / ic->coulombEwaldTables->scale,
640                     ic->coulombEwaldTables->tableF.size());
641         }
642     }
643 }
644
645 real cutoff_inf(real cutoff)
646 {
647     if (cutoff == 0)
648     {
649         cutoff = GMX_CUTOFF_INF;
650     }
651
652     return cutoff;
653 }
654
655 void init_forcerec(FILE*                            fplog,
656                    const gmx::MDLogger&             mdlog,
657                    const gmx::SimulationWorkload&   simulationWork,
658                    t_forcerec*                      forcerec,
659                    const t_inputrec&                inputrec,
660                    const gmx_mtop_t&                mtop,
661                    const t_commrec*                 commrec,
662                    matrix                           box,
663                    const char*                      tabfn,
664                    const char*                      tabpfn,
665                    gmx::ArrayRef<const std::string> tabbfnm,
666                    real                             print_force)
667 {
668     /* The CMake default turns SIMD kernels on, but it might be turned off further down... */
669     forcerec->use_simd_kernels = GMX_USE_SIMD_KERNELS;
670
671     if (check_box(inputrec.pbcType, box))
672     {
673         gmx_fatal(FARGS, "%s", check_box(inputrec.pbcType, box));
674     }
675
676     /* Test particle insertion ? */
677     if (EI_TPI(inputrec.eI))
678     {
679         /* Set to the size of the molecule to be inserted (the last one) */
680         gmx::RangePartitioning molecules = gmx_mtop_molecules(mtop);
681         forcerec->n_tpi                  = molecules.block(molecules.numBlocks() - 1).size();
682     }
683     else
684     {
685         forcerec->n_tpi = 0;
686     }
687
688     if (inputrec.coulombtype == CoulombInteractionType::RFNecUnsupported
689         || inputrec.coulombtype == CoulombInteractionType::GRFNotused)
690     {
691         gmx_fatal(FARGS, "%s electrostatics is no longer supported", enumValueToString(inputrec.coulombtype));
692     }
693
694     if (inputrec.bAdress)
695     {
696         gmx_fatal(FARGS, "AdResS simulations are no longer supported");
697     }
698     if (inputrec.useTwinRange)
699     {
700         gmx_fatal(FARGS, "Twin-range simulations are no longer supported");
701     }
702     /* Copy the user determined parameters */
703     forcerec->userint1  = inputrec.userint1;
704     forcerec->userint2  = inputrec.userint2;
705     forcerec->userint3  = inputrec.userint3;
706     forcerec->userint4  = inputrec.userint4;
707     forcerec->userreal1 = inputrec.userreal1;
708     forcerec->userreal2 = inputrec.userreal2;
709     forcerec->userreal3 = inputrec.userreal3;
710     forcerec->userreal4 = inputrec.userreal4;
711
712     /* Shell stuff */
713     forcerec->fc_stepsize = inputrec.fc_stepsize;
714
715     /* Free energy */
716     forcerec->efep = inputrec.efep;
717
718     if ((getenv("GMX_DISABLE_SIMD_KERNELS") != nullptr) || (getenv("GMX_NOOPTIMIZEDKERNELS") != nullptr))
719     {
720         forcerec->use_simd_kernels = FALSE;
721         if (fplog != nullptr)
722         {
723             fprintf(fplog,
724                     "\nFound environment variable GMX_DISABLE_SIMD_KERNELS.\n"
725                     "Disabling the usage of any SIMD-specific non-bonded & bonded kernel routines\n"
726                     "(e.g. SSE2/SSE4.1/AVX).\n\n");
727         }
728     }
729
730     forcerec->haveBuckingham = (mtop.ffparams.functype[0] == F_BHAM);
731
732     /* Neighbour searching stuff */
733     forcerec->pbcType = inputrec.pbcType;
734
735     /* Determine if we will do PBC for distances in bonded interactions */
736     if (forcerec->pbcType == PbcType::No)
737     {
738         forcerec->bMolPBC = FALSE;
739     }
740     else
741     {
742         forcerec->bMolPBC =
743                 (!haveDDAtomOrdering(*commrec) || dd_bonded_molpbc(*commrec->dd, forcerec->pbcType));
744
745         // Check and set up PBC for Ewald surface corrections or orientation restraints
746         const bool useEwaldSurfaceCorrection =
747                 (EEL_PME_EWALD(inputrec.coulombtype) && inputrec.epsilon_surface != 0);
748         const bool haveOrientationRestraints = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
749         const bool moleculesAreAlwaysWhole =
750                 (haveDDAtomOrdering(*commrec) && dd_moleculesAreAlwaysWhole(*commrec->dd));
751         // WholeMoleculeTransform is only supported with a single PP rank
752         if (!moleculesAreAlwaysWhole && !havePPDomainDecomposition(commrec)
753             && (useEwaldSurfaceCorrection || haveOrientationRestraints))
754         {
755             if (havePPDomainDecomposition(commrec))
756             {
757                 gmx_fatal(FARGS,
758                           "You requested Ewald surface correction or orientation restraints, "
759                           "but molecules are broken "
760                           "over periodic boundary conditions by the domain decomposition. "
761                           "Run without domain decomposition instead.");
762             }
763
764             forcerec->wholeMoleculeTransform = std::make_unique<gmx::WholeMoleculeTransform>(
765                     mtop, inputrec.pbcType, haveDDAtomOrdering(*commrec));
766         }
767
768         forcerec->bMolPBC =
769                 !haveDDAtomOrdering(*commrec) || dd_bonded_molpbc(*commrec->dd, forcerec->pbcType);
770
771         if (useEwaldSurfaceCorrection)
772         {
773             GMX_RELEASE_ASSERT(moleculesAreAlwaysWhole || forcerec->wholeMoleculeTransform,
774                                "Molecules can not be broken by PBC with epsilon_surface > 0");
775         }
776     }
777
778     forcerec->rc_scaling = inputrec.refcoord_scaling;
779     copy_rvec(inputrec.posres_com, forcerec->posres_com);
780     copy_rvec(inputrec.posres_comB, forcerec->posres_comB);
781     forcerec->rlist                  = cutoff_inf(inputrec.rlist);
782     forcerec->ljpme_combination_rule = inputrec.ljpme_combination_rule;
783
784     /* This now calculates sum for q and c6*/
785     bool systemHasNetCharge = set_chargesum(fplog, forcerec, mtop);
786
787     /* Make data structure used by kernels */
788     forcerec->ic = std::make_unique<interaction_const_t>(
789             init_interaction_const(fplog, inputrec, mtop, systemHasNetCharge));
790     init_interaction_const_tables(fplog, forcerec->ic.get(), forcerec->rlist, inputrec.tabext);
791
792     const interaction_const_t* interactionConst = forcerec->ic.get();
793
794     /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
795     switch (interactionConst->eeltype)
796     {
797         case CoulombInteractionType::Cut:
798             forcerec->nbkernel_elec_interaction = NbkernelElecType::Coulomb;
799             break;
800
801         case CoulombInteractionType::RF:
802         case CoulombInteractionType::RFZero:
803             forcerec->nbkernel_elec_interaction = NbkernelElecType::ReactionField;
804             break;
805
806         case CoulombInteractionType::Switch:
807         case CoulombInteractionType::Shift:
808         case CoulombInteractionType::User:
809         case CoulombInteractionType::PmeSwitch:
810         case CoulombInteractionType::PmeUser:
811         case CoulombInteractionType::PmeUserSwitch:
812             forcerec->nbkernel_elec_interaction = NbkernelElecType::CubicSplineTable;
813             break;
814
815         case CoulombInteractionType::Pme:
816         case CoulombInteractionType::P3mAD:
817         case CoulombInteractionType::Ewald:
818             forcerec->nbkernel_elec_interaction = NbkernelElecType::Ewald;
819             break;
820
821         default:
822             gmx_fatal(FARGS,
823                       "Unsupported electrostatic interaction: %s",
824                       enumValueToString(interactionConst->eeltype));
825     }
826     forcerec->nbkernel_elec_modifier = interactionConst->coulomb_modifier;
827
828     /* Vdw: Translate from mdp settings to kernel format */
829     switch (interactionConst->vdwtype)
830     {
831         case VanDerWaalsType::Cut:
832             if (forcerec->haveBuckingham)
833             {
834                 forcerec->nbkernel_vdw_interaction = NbkernelVdwType::Buckingham;
835             }
836             else
837             {
838                 forcerec->nbkernel_vdw_interaction = NbkernelVdwType::LennardJones;
839             }
840             break;
841         case VanDerWaalsType::Pme:
842             forcerec->nbkernel_vdw_interaction = NbkernelVdwType::LJEwald;
843             break;
844
845         case VanDerWaalsType::Switch:
846         case VanDerWaalsType::Shift:
847         case VanDerWaalsType::User:
848             forcerec->nbkernel_vdw_interaction = NbkernelVdwType::CubicSplineTable;
849             break;
850
851         default:
852             gmx_fatal(FARGS, "Unsupported vdw interaction: %s", enumValueToString(interactionConst->vdwtype));
853     }
854     forcerec->nbkernel_vdw_modifier = interactionConst->vdw_modifier;
855
856     if (!gmx_within_tol(interactionConst->reppow, 12.0, 10 * GMX_DOUBLE_EPS))
857     {
858         gmx_fatal(FARGS, "Only LJ repulsion power 12 is supported");
859     }
860     /* Older tpr files can contain Coulomb user tables with the Verlet cutoff-scheme,
861      * while mdrun does not (and never did) support this.
862      */
863     if (EEL_USER(forcerec->ic->eeltype))
864     {
865         gmx_fatal(FARGS,
866                   "Electrostatics type %s is currently not supported",
867                   enumValueToString(inputrec.coulombtype));
868     }
869
870     /* 1-4 interaction electrostatics */
871     forcerec->fudgeQQ = mtop.ffparams.fudgeQQ;
872
873     if (simulationWork.useMts)
874     {
875         GMX_ASSERT(gmx::checkMtsRequirements(inputrec).empty(),
876                    "All MTS requirements should be met here");
877     }
878
879     const bool haveDirectVirialContributionsFast =
880             forcerec->forceProviders->hasForceProvider() || gmx_mtop_ftype_count(mtop, F_POSRES) > 0
881             || gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0 || inputrec.nwall > 0 || inputrec.bPull
882             || inputrec.bRot || inputrec.bIMD;
883     const bool haveDirectVirialContributionsSlow =
884             EEL_FULL(interactionConst->eeltype) || EVDW_PME(interactionConst->vdwtype);
885     for (int i = 0; i < (simulationWork.useMts ? 2 : 1); i++)
886     {
887         bool haveDirectVirialContributions =
888                 (((!simulationWork.useMts || i == 0) && haveDirectVirialContributionsFast)
889                  || ((!simulationWork.useMts || i == 1) && haveDirectVirialContributionsSlow));
890         forcerec->forceHelperBuffers.emplace_back(haveDirectVirialContributions);
891     }
892
893     if (forcerec->shift_vec.empty())
894     {
895         forcerec->shift_vec.resize(gmx::c_numShiftVectors);
896     }
897
898     GMX_ASSERT(forcerec->nbfp.empty(), "The nonbonded force parameters should not be set up yet.");
899     forcerec->ntype = mtop.ffparams.atnr;
900     forcerec->nbfp  = makeNonBondedParameterLists(
901             mtop.ffparams.atnr, mtop.ffparams.iparams, forcerec->haveBuckingham);
902     if (EVDW_PME(interactionConst->vdwtype))
903     {
904         forcerec->ljpme_c6grid = makeLJPmeC6GridCorrectionParameters(
905                 mtop.ffparams.atnr, mtop.ffparams.iparams, forcerec->ljpme_combination_rule);
906     }
907
908     /* Copy the energy group exclusions */
909     forcerec->egp_flags = inputrec.opts.egp_flags;
910
911     /* Van der Waals stuff */
912     if ((interactionConst->vdwtype != VanDerWaalsType::Cut)
913         && (interactionConst->vdwtype != VanDerWaalsType::User) && !forcerec->haveBuckingham)
914     {
915         if (interactionConst->rvdw_switch >= interactionConst->rvdw)
916         {
917             gmx_fatal(FARGS,
918                       "rvdw_switch (%f) must be < rvdw (%f)",
919                       interactionConst->rvdw_switch,
920                       interactionConst->rvdw);
921         }
922         if (fplog)
923         {
924             fprintf(fplog,
925                     "Using %s Lennard-Jones, switch between %g and %g nm\n",
926                     (interactionConst->eeltype == CoulombInteractionType::Switch) ? "switched" : "shifted",
927                     interactionConst->rvdw_switch,
928                     interactionConst->rvdw);
929         }
930     }
931
932     if (forcerec->haveBuckingham && EVDW_PME(interactionConst->vdwtype))
933     {
934         gmx_fatal(FARGS, "LJ PME not supported with Buckingham");
935     }
936
937     if (forcerec->haveBuckingham
938         && (interactionConst->vdwtype == VanDerWaalsType::Shift
939             || interactionConst->vdwtype == VanDerWaalsType::Switch))
940     {
941         gmx_fatal(FARGS, "Switch/shift interaction not supported with Buckingham");
942     }
943
944     if (forcerec->haveBuckingham)
945     {
946         gmx_fatal(FARGS, "The Verlet cutoff-scheme does not (yet) support Buckingham");
947     }
948
949     if (inputrec.implicit_solvent)
950     {
951         gmx_fatal(FARGS, "Implict solvation is no longer supported.");
952     }
953
954
955     /* This code automatically gives table length tabext without cut-off's,
956      * in that case grompp should already have checked that we do not need
957      * normal tables and we only generate tables for 1-4 interactions.
958      */
959     real rtab = inputrec.rlist + inputrec.tabext;
960
961     /* We want to use unmodified tables for 1-4 coulombic
962      * interactions, so we must in general have an extra set of
963      * tables. */
964     if (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 || gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0
965         || gmx_mtop_ftype_count(mtop, F_LJC_PAIRS_NB) > 0)
966     {
967         forcerec->pairsTable = make_tables(fplog, interactionConst, tabpfn, rtab, GMX_MAKETABLES_14ONLY);
968     }
969
970     /* Wall stuff */
971     forcerec->nwall = inputrec.nwall;
972     if (inputrec.nwall && inputrec.wall_type == WallType::Table)
973     {
974         make_wall_tables(fplog, inputrec, tabfn, &mtop.groups, forcerec);
975     }
976
977     forcerec->fcdata = std::make_unique<t_fcdata>();
978
979     if (!tabbfnm.empty())
980     {
981         t_fcdata& fcdata = *forcerec->fcdata;
982         // Need to catch std::bad_alloc
983         // TODO Don't need to catch this here, when merging with master branch
984         try
985         {
986             // TODO move these tables into a separate struct and store reference in ListedForces
987             fcdata.bondtab = make_bonded_tables(fplog, F_TABBONDS, F_TABBONDSNC, mtop, tabbfnm, "b");
988             fcdata.angletab = make_bonded_tables(fplog, F_TABANGLES, -1, mtop, tabbfnm, "a");
989             fcdata.dihtab   = make_bonded_tables(fplog, F_TABDIHS, -1, mtop, tabbfnm, "d");
990         }
991         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
992     }
993     else
994     {
995         if (debug)
996         {
997             fprintf(debug,
998                     "No fcdata or table file name passed, can not read table, can not do bonded "
999                     "interactions\n");
1000         }
1001     }
1002
1003     /* Initialize the thread working data for bonded interactions */
1004     if (simulationWork.useMts)
1005     {
1006         // Add one ListedForces object for each MTS level
1007         bool isFirstLevel = true;
1008         for (const auto& mtsLevel : inputrec.mtsLevels)
1009         {
1010             ListedForces::InteractionSelection interactionSelection;
1011             const auto&                        forceGroups = mtsLevel.forceGroups;
1012             if (forceGroups[static_cast<int>(gmx::MtsForceGroups::Pair)])
1013             {
1014                 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Pairs));
1015             }
1016             if (forceGroups[static_cast<int>(gmx::MtsForceGroups::Dihedral)])
1017             {
1018                 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Dihedrals));
1019             }
1020             if (forceGroups[static_cast<int>(gmx::MtsForceGroups::Angle)])
1021             {
1022                 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Angles));
1023             }
1024             if (isFirstLevel)
1025             {
1026                 interactionSelection.set(static_cast<int>(ListedForces::InteractionGroup::Rest));
1027                 isFirstLevel = false;
1028             }
1029             forcerec->listedForces.emplace_back(
1030                     mtop.ffparams,
1031                     mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1032                     gmx_omp_nthreads_get(ModuleMultiThread::Bonded),
1033                     interactionSelection,
1034                     fplog);
1035         }
1036     }
1037     else
1038     {
1039         // Add one ListedForces object with all listed interactions
1040         forcerec->listedForces.emplace_back(
1041                 mtop.ffparams,
1042                 mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size(),
1043                 gmx_omp_nthreads_get(ModuleMultiThread::Bonded),
1044                 ListedForces::interactionSelectionAll(),
1045                 fplog);
1046     }
1047
1048     // QM/MM initialization if requested
1049     if (inputrec.bQMMM)
1050     {
1051         gmx_incons("QM/MM was requested, but is no longer available in GROMACS");
1052     }
1053
1054     /* Set all the static charge group info */
1055     forcerec->atomInfoForEachMoleculeBlock = makeAtomInfoForEachMoleculeBlock(mtop, forcerec);
1056     if (!haveDDAtomOrdering(*commrec))
1057     {
1058         forcerec->atomInfo = expandAtomInfo(mtop.molblock.size(), forcerec->atomInfoForEachMoleculeBlock);
1059     }
1060
1061     if (!haveDDAtomOrdering(*commrec))
1062     {
1063         forcerec_set_ranges(forcerec, mtop.natoms, mtop.natoms, mtop.natoms);
1064     }
1065
1066     forcerec->print_force = print_force;
1067
1068     if (inputrec.eDispCorr != DispersionCorrectionType::No)
1069     {
1070         forcerec->dispersionCorrection = std::make_unique<DispersionCorrection>(
1071                 mtop, inputrec, forcerec->haveBuckingham, forcerec->ntype, forcerec->nbfp, *forcerec->ic, tabfn);
1072         forcerec->dispersionCorrection->print(mdlog);
1073     }
1074
1075     if (fplog != nullptr)
1076     {
1077         /* Here we switch from using mdlog, which prints the newline before
1078          * the paragraph, to our old fprintf logging, which prints the newline
1079          * after the paragraph, so we should add a newline here.
1080          */
1081         fprintf(fplog, "\n");
1082     }
1083 }
1084
1085 t_forcerec::t_forcerec() = default;
1086
1087 t_forcerec::~t_forcerec() = default;