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