Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / mdlib / forcerec.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013-2019, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #include "gmxpre.h"
38
39 #include "forcerec.h"
40
41 #include "config.h"
42
43 #include <cassert>
44 #include <cmath>
45 #include <cstdlib>
46 #include <cstring>
47
48 #include <algorithm>
49 #include <memory>
50
51 #include "gromacs/commandline/filenm.h"
52 #include "gromacs/domdec/domdec.h"
53 #include "gromacs/domdec/domdec_struct.h"
54 #include "gromacs/ewald/ewald.h"
55 #include "gromacs/ewald/ewald_utils.h"
56 #include "gromacs/ewald/pme_pp_comm_gpu.h"
57 #include "gromacs/fileio/filetypes.h"
58 #include "gromacs/gmxlib/network.h"
59 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
60 #include "gromacs/gpu_utils/gpu_utils.h"
61 #include "gromacs/hardware/hw_info.h"
62 #include "gromacs/listed_forces/gpubonded.h"
63 #include "gromacs/listed_forces/manage_threading.h"
64 #include "gromacs/listed_forces/pairs.h"
65 #include "gromacs/math/functions.h"
66 #include "gromacs/math/units.h"
67 #include "gromacs/math/utilities.h"
68 #include "gromacs/math/vec.h"
69 #include "gromacs/mdlib/dispersioncorrection.h"
70 #include "gromacs/mdlib/force.h"
71 #include "gromacs/mdlib/forcerec_threading.h"
72 #include "gromacs/mdlib/gmx_omp_nthreads.h"
73 #include "gromacs/mdlib/md_support.h"
74 #include "gromacs/mdlib/qmmm.h"
75 #include "gromacs/mdlib/rf_util.h"
76 #include "gromacs/mdlib/wall.h"
77 #include "gromacs/mdtypes/commrec.h"
78 #include "gromacs/mdtypes/fcdata.h"
79 #include "gromacs/mdtypes/group.h"
80 #include "gromacs/mdtypes/iforceprovider.h"
81 #include "gromacs/mdtypes/inputrec.h"
82 #include "gromacs/mdtypes/md_enums.h"
83 #include "gromacs/nbnxm/gpu_data_mgmt.h"
84 #include "gromacs/nbnxm/nbnxm.h"
85 #include "gromacs/nbnxm/nbnxm_geometry.h"
86 #include "gromacs/pbcutil/ishift.h"
87 #include "gromacs/pbcutil/pbc.h"
88 #include "gromacs/tables/forcetable.h"
89 #include "gromacs/topology/mtop_util.h"
90 #include "gromacs/trajectory/trajectoryframe.h"
91 #include "gromacs/utility/cstringutil.h"
92 #include "gromacs/utility/exceptions.h"
93 #include "gromacs/utility/fatalerror.h"
94 #include "gromacs/utility/gmxassert.h"
95 #include "gromacs/utility/logger.h"
96 #include "gromacs/utility/physicalnodecommunicator.h"
97 #include "gromacs/utility/pleasecite.h"
98 #include "gromacs/utility/smalloc.h"
99 #include "gromacs/utility/strconvert.h"
100
101 /*! \brief environment variable to enable GPU P2P communication */
102 static const bool c_enableGpuPmePpComms =
103         (getenv("GMX_GPU_PME_PP_COMMS") != nullptr) && GMX_THREAD_MPI && (GMX_GPU == GMX_GPU_CUDA);
104
105 static real* mk_nbfp(const gmx_ffparams_t* idef, gmx_bool bBHAM)
106 {
107     real* nbfp;
108     int   i, j, k, atnr;
109
110     atnr = idef->atnr;
111     if (bBHAM)
112     {
113         snew(nbfp, 3 * atnr * atnr);
114         for (i = k = 0; (i < atnr); i++)
115         {
116             for (j = 0; (j < atnr); j++, k++)
117             {
118                 BHAMA(nbfp, atnr, i, j) = idef->iparams[k].bham.a;
119                 BHAMB(nbfp, atnr, i, j) = idef->iparams[k].bham.b;
120                 /* nbfp now includes the 6.0 derivative prefactor */
121                 BHAMC(nbfp, atnr, i, j) = idef->iparams[k].bham.c * 6.0;
122             }
123         }
124     }
125     else
126     {
127         snew(nbfp, 2 * atnr * atnr);
128         for (i = k = 0; (i < atnr); i++)
129         {
130             for (j = 0; (j < atnr); j++, k++)
131             {
132                 /* nbfp now includes the 6.0/12.0 derivative prefactors */
133                 C6(nbfp, atnr, i, j)  = idef->iparams[k].lj.c6 * 6.0;
134                 C12(nbfp, atnr, i, j) = idef->iparams[k].lj.c12 * 12.0;
135             }
136         }
137     }
138
139     return nbfp;
140 }
141
142 static real* make_ljpme_c6grid(const gmx_ffparams_t* idef, t_forcerec* fr)
143 {
144     int   i, j, k, atnr;
145     real  c6, c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
146     real* grid;
147
148     /* For LJ-PME simulations, we correct the energies with the reciprocal space
149      * inside of the cut-off. To do this the non-bonded kernels needs to have
150      * access to the C6-values used on the reciprocal grid in pme.c
151      */
152
153     atnr = idef->atnr;
154     snew(grid, 2 * atnr * atnr);
155     for (i = k = 0; (i < atnr); i++)
156     {
157         for (j = 0; (j < atnr); j++, k++)
158         {
159             c6i  = idef->iparams[i * (atnr + 1)].lj.c6;
160             c12i = idef->iparams[i * (atnr + 1)].lj.c12;
161             c6j  = idef->iparams[j * (atnr + 1)].lj.c6;
162             c12j = idef->iparams[j * (atnr + 1)].lj.c12;
163             c6   = std::sqrt(c6i * c6j);
164             if (fr->ljpme_combination_rule == eljpmeLB && !gmx_numzero(c6) && !gmx_numzero(c12i)
165                 && !gmx_numzero(c12j))
166             {
167                 sigmai = gmx::sixthroot(c12i / c6i);
168                 sigmaj = gmx::sixthroot(c12j / c6j);
169                 epsi   = c6i * c6i / c12i;
170                 epsj   = c6j * c6j / c12j;
171                 c6     = std::sqrt(epsi * epsj) * gmx::power6(0.5 * (sigmai + sigmaj));
172             }
173             /* Store the elements at the same relative positions as C6 in nbfp in order
174              * to simplify access in the kernels
175              */
176             grid[2 * (atnr * i + j)] = c6 * 6.0;
177         }
178     }
179     return grid;
180 }
181
182 enum
183 {
184     acNONE = 0,
185     acCONSTRAINT,
186     acSETTLE
187 };
188
189 static cginfo_mb_t* init_cginfo_mb(const gmx_mtop_t* mtop, const t_forcerec* fr, gmx_bool* bFEP_NonBonded)
190 {
191     cginfo_mb_t* cginfo_mb;
192     gmx_bool*    type_VDW;
193     int*         cginfo;
194     int*         a_con;
195
196     snew(cginfo_mb, mtop->molblock.size());
197
198     snew(type_VDW, fr->ntype);
199     for (int ai = 0; ai < fr->ntype; ai++)
200     {
201         type_VDW[ai] = FALSE;
202         for (int j = 0; j < fr->ntype; j++)
203         {
204             type_VDW[ai] = type_VDW[ai] || fr->bBHAM || C6(fr->nbfp, fr->ntype, ai, j) != 0
205                            || C12(fr->nbfp, fr->ntype, ai, j) != 0;
206         }
207     }
208
209     *bFEP_NonBonded = FALSE;
210
211     int a_offset = 0;
212     for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
213     {
214         const gmx_molblock_t& molb = mtop->molblock[mb];
215         const gmx_moltype_t&  molt = mtop->moltype[molb.type];
216         const t_blocka&       excl = molt.excls;
217
218         /* Check if the cginfo is identical for all molecules in this block.
219          * If so, we only need an array of the size of one molecule.
220          * Otherwise we make an array of #mol times #cgs per molecule.
221          */
222         gmx_bool bId = TRUE;
223         for (int m = 0; m < molb.nmol; m++)
224         {
225             const int am = m * molt.atoms.nr;
226             for (int a = 0; a < molt.atoms.nr; a++)
227             {
228                 if (getGroupType(mtop->groups, SimulationAtomGroupType::QuantumMechanics, a_offset + am + a)
229                     != getGroupType(mtop->groups, SimulationAtomGroupType::QuantumMechanics, a_offset + a))
230                 {
231                     bId = FALSE;
232                 }
233                 if (!mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics].empty())
234                 {
235                     if (mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics][a_offset + am + a]
236                         != mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics][a_offset + a])
237                     {
238                         bId = FALSE;
239                     }
240                 }
241             }
242         }
243
244         cginfo_mb[mb].cg_start = a_offset;
245         cginfo_mb[mb].cg_end   = a_offset + molb.nmol * molt.atoms.nr;
246         cginfo_mb[mb].cg_mod   = (bId ? 1 : molb.nmol) * molt.atoms.nr;
247         snew(cginfo_mb[mb].cginfo, cginfo_mb[mb].cg_mod);
248         cginfo = cginfo_mb[mb].cginfo;
249
250         /* Set constraints flags for constrained atoms */
251         snew(a_con, molt.atoms.nr);
252         for (int ftype = 0; ftype < F_NRE; ftype++)
253         {
254             if (interaction_function[ftype].flags & IF_CONSTRAINT)
255             {
256                 const int nral = NRAL(ftype);
257                 for (int ia = 0; ia < molt.ilist[ftype].size(); ia += 1 + nral)
258                 {
259                     int a;
260
261                     for (a = 0; a < nral; a++)
262                     {
263                         a_con[molt.ilist[ftype].iatoms[ia + 1 + a]] =
264                                 (ftype == F_SETTLE ? acSETTLE : acCONSTRAINT);
265                     }
266                 }
267             }
268         }
269
270         for (int m = 0; m < (bId ? 1 : molb.nmol); m++)
271         {
272             const int molculeOffsetInBlock = m * molt.atoms.nr;
273             for (int a = 0; a < molt.atoms.nr; a++)
274             {
275                 const t_atom& atom     = molt.atoms.atom[a];
276                 int&          atomInfo = cginfo[molculeOffsetInBlock + a];
277
278                 /* Store the energy group in cginfo */
279                 int gid = getGroupType(mtop->groups, SimulationAtomGroupType::EnergyOutput,
280                                        a_offset + molculeOffsetInBlock + a);
281                 SET_CGINFO_GID(atomInfo, gid);
282
283                 bool bHaveVDW = (type_VDW[atom.type] || type_VDW[atom.typeB]);
284                 bool bHaveQ   = (atom.q != 0 || atom.qB != 0);
285
286                 bool haveExclusions = false;
287                 /* Loop over all the exclusions of atom ai */
288                 for (int j = excl.index[a]; j < excl.index[a + 1]; j++)
289                 {
290                     if (excl.a[j] != a)
291                     {
292                         haveExclusions = true;
293                         break;
294                     }
295                 }
296
297                 switch (a_con[a])
298                 {
299                     case acCONSTRAINT: SET_CGINFO_CONSTR(atomInfo); break;
300                     case acSETTLE: SET_CGINFO_SETTLE(atomInfo); break;
301                     default: break;
302                 }
303                 if (haveExclusions)
304                 {
305                     SET_CGINFO_EXCL_INTER(atomInfo);
306                 }
307                 if (bHaveVDW)
308                 {
309                     SET_CGINFO_HAS_VDW(atomInfo);
310                 }
311                 if (bHaveQ)
312                 {
313                     SET_CGINFO_HAS_Q(atomInfo);
314                 }
315                 if (fr->efep != efepNO && PERTURBED(atom))
316                 {
317                     SET_CGINFO_FEP(atomInfo);
318                     *bFEP_NonBonded = TRUE;
319                 }
320             }
321         }
322
323         sfree(a_con);
324
325         a_offset += molb.nmol * molt.atoms.nr;
326     }
327     sfree(type_VDW);
328
329     return cginfo_mb;
330 }
331
332 static std::vector<int> cginfo_expand(const int nmb, const cginfo_mb_t* cgi_mb)
333 {
334     const int ncg = cgi_mb[nmb - 1].cg_end;
335
336     std::vector<int> cginfo(ncg);
337
338     int mb = 0;
339     for (int cg = 0; cg < ncg; cg++)
340     {
341         while (cg >= cgi_mb[mb].cg_end)
342         {
343             mb++;
344         }
345         cginfo[cg] = cgi_mb[mb].cginfo[(cg - cgi_mb[mb].cg_start) % cgi_mb[mb].cg_mod];
346     }
347
348     return cginfo;
349 }
350
351 static void done_cginfo_mb(cginfo_mb_t* cginfo_mb, int numMolBlocks)
352 {
353     if (cginfo_mb == nullptr)
354     {
355         return;
356     }
357     for (int mb = 0; mb < numMolBlocks; ++mb)
358     {
359         sfree(cginfo_mb[mb].cginfo);
360     }
361     sfree(cginfo_mb);
362 }
363
364 /* Sets the sum of charges (squared) and C6 in the system in fr.
365  * Returns whether the system has a net charge.
366  */
367 static bool set_chargesum(FILE* log, t_forcerec* fr, const gmx_mtop_t* mtop)
368 {
369     /*This now calculates sum for q and c6*/
370     double qsum, q2sum, q, c6sum, c6;
371
372     qsum  = 0;
373     q2sum = 0;
374     c6sum = 0;
375     for (const gmx_molblock_t& molb : mtop->molblock)
376     {
377         int            nmol  = molb.nmol;
378         const t_atoms* atoms = &mtop->moltype[molb.type].atoms;
379         for (int i = 0; i < atoms->nr; i++)
380         {
381             q = atoms->atom[i].q;
382             qsum += nmol * q;
383             q2sum += nmol * q * q;
384             c6 = mtop->ffparams.iparams[atoms->atom[i].type * (mtop->ffparams.atnr + 1)].lj.c6;
385             c6sum += nmol * c6;
386         }
387     }
388     fr->qsum[0]  = qsum;
389     fr->q2sum[0] = q2sum;
390     fr->c6sum[0] = c6sum;
391
392     if (fr->efep != efepNO)
393     {
394         qsum  = 0;
395         q2sum = 0;
396         c6sum = 0;
397         for (const gmx_molblock_t& molb : mtop->molblock)
398         {
399             int            nmol  = molb.nmol;
400             const t_atoms* atoms = &mtop->moltype[molb.type].atoms;
401             for (int i = 0; i < atoms->nr; i++)
402             {
403                 q = atoms->atom[i].qB;
404                 qsum += nmol * q;
405                 q2sum += nmol * q * q;
406                 c6 = mtop->ffparams.iparams[atoms->atom[i].typeB * (mtop->ffparams.atnr + 1)].lj.c6;
407                 c6sum += nmol * c6;
408             }
409             fr->qsum[1]  = qsum;
410             fr->q2sum[1] = q2sum;
411             fr->c6sum[1] = c6sum;
412         }
413     }
414     else
415     {
416         fr->qsum[1]  = fr->qsum[0];
417         fr->q2sum[1] = fr->q2sum[0];
418         fr->c6sum[1] = fr->c6sum[0];
419     }
420     if (log)
421     {
422         if (fr->efep == efepNO)
423         {
424             fprintf(log, "System total charge: %.3f\n", fr->qsum[0]);
425         }
426         else
427         {
428             fprintf(log, "System total charge, top. A: %.3f top. B: %.3f\n", fr->qsum[0], fr->qsum[1]);
429         }
430     }
431
432     /* A cut-off of 1e-4 is used to catch rounding errors due to ascii input */
433     return (std::abs(fr->qsum[0]) > 1e-4 || std::abs(fr->qsum[1]) > 1e-4);
434 }
435
436 static real calcBuckinghamBMax(FILE* fplog, const gmx_mtop_t* mtop)
437 {
438     const t_atoms *at1, *at2;
439     int            i, j, tpi, tpj, ntypes;
440     real           b, bmin;
441
442     if (fplog)
443     {
444         fprintf(fplog, "Determining largest Buckingham b parameter for table\n");
445     }
446     ntypes = mtop->ffparams.atnr;
447
448     bmin            = -1;
449     real bham_b_max = 0;
450     for (size_t mt1 = 0; mt1 < mtop->moltype.size(); mt1++)
451     {
452         at1 = &mtop->moltype[mt1].atoms;
453         for (i = 0; (i < at1->nr); i++)
454         {
455             tpi = at1->atom[i].type;
456             if (tpi >= ntypes)
457             {
458                 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", i, tpi, ntypes);
459             }
460
461             for (size_t mt2 = mt1; mt2 < mtop->moltype.size(); mt2++)
462             {
463                 at2 = &mtop->moltype[mt2].atoms;
464                 for (j = 0; (j < at2->nr); j++)
465                 {
466                     tpj = at2->atom[j].type;
467                     if (tpj >= ntypes)
468                     {
469                         gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", j, tpj, ntypes);
470                     }
471                     b = mtop->ffparams.iparams[tpi * ntypes + tpj].bham.b;
472                     if (b > bham_b_max)
473                     {
474                         bham_b_max = b;
475                     }
476                     if ((b < bmin) || (bmin == -1))
477                     {
478                         bmin = b;
479                     }
480                 }
481             }
482         }
483     }
484     if (fplog)
485     {
486         fprintf(fplog, "Buckingham b parameters, min: %g, max: %g\n", bmin, bham_b_max);
487     }
488
489     return bham_b_max;
490 }
491
492 /*!\brief If there's bonded interactions of type \c ftype1 or \c
493  * ftype2 present in the topology, build an array of the number of
494  * interactions present for each bonded interaction index found in the
495  * topology.
496  *
497  * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
498  * valid type with that parameter.
499  *
500  * \c count will be reallocated as necessary to fit the largest bonded
501  * interaction index found, and its current size will be returned in
502  * \c ncount. It will contain zero for every bonded interaction index
503  * for which no interactions are present in the topology.
504  */
505 static void count_tables(int ftype1, int ftype2, const gmx_mtop_t* mtop, int* ncount, int** count)
506 {
507     int ftype, i, j, tabnr;
508
509     // Loop over all moleculetypes
510     for (const gmx_moltype_t& molt : mtop->moltype)
511     {
512         // Loop over all interaction types
513         for (ftype = 0; ftype < F_NRE; ftype++)
514         {
515             // If the current interaction type is one of the types whose tables we're trying to count...
516             if (ftype == ftype1 || ftype == ftype2)
517             {
518                 const InteractionList& il     = molt.ilist[ftype];
519                 const int              stride = 1 + NRAL(ftype);
520                 // ... and there are actually some interactions for this type
521                 for (i = 0; i < il.size(); i += stride)
522                 {
523                     // Find out which table index the user wanted
524                     tabnr = mtop->ffparams.iparams[il.iatoms[i]].tab.table;
525                     if (tabnr < 0)
526                     {
527                         gmx_fatal(FARGS, "A bonded table number is smaller than 0: %d\n", tabnr);
528                     }
529                     // Make room for this index in the data structure
530                     if (tabnr >= *ncount)
531                     {
532                         srenew(*count, tabnr + 1);
533                         for (j = *ncount; j < tabnr + 1; j++)
534                         {
535                             (*count)[j] = 0;
536                         }
537                         *ncount = tabnr + 1;
538                     }
539                     // Record that this table index is used and must have a valid file
540                     (*count)[tabnr]++;
541                 }
542             }
543         }
544     }
545 }
546
547 /*!\brief If there's bonded interactions of flavour \c tabext and type
548  * \c ftype1 or \c ftype2 present in the topology, seek them in the
549  * list of filenames passed to mdrun, and make bonded tables from
550  * those files.
551  *
552  * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
553  * valid type with that parameter.
554  *
555  * A fatal error occurs if no matching filename is found.
556  */
557 static bondedtable_t* make_bonded_tables(FILE*                            fplog,
558                                          int                              ftype1,
559                                          int                              ftype2,
560                                          const gmx_mtop_t*                mtop,
561                                          gmx::ArrayRef<const std::string> tabbfnm,
562                                          const char*                      tabext)
563 {
564     int            ncount, *count;
565     bondedtable_t* tab;
566
567     tab = nullptr;
568
569     ncount = 0;
570     count  = nullptr;
571     count_tables(ftype1, ftype2, mtop, &ncount, &count);
572
573     // Are there any relevant tabulated bond interactions?
574     if (ncount > 0)
575     {
576         snew(tab, ncount);
577         for (int i = 0; i < ncount; i++)
578         {
579             // Do any interactions exist that requires this table?
580             if (count[i] > 0)
581             {
582                 // This pattern enforces the current requirement that
583                 // table filenames end in a characteristic sequence
584                 // before the file type extension, and avoids table 13
585                 // being recognized and used for table 1.
586                 std::string patternToFind = gmx::formatString("_%s%d.%s", tabext, i, ftp2ext(efXVG));
587                 bool        madeTable     = false;
588                 for (gmx::index j = 0; j < tabbfnm.ssize() && !madeTable; ++j)
589                 {
590                     if (gmx::endsWith(tabbfnm[j], patternToFind))
591                     {
592                         // Finally read the table from the file found
593                         tab[i]    = make_bonded_table(fplog, tabbfnm[j].c_str(), NRAL(ftype1) - 2);
594                         madeTable = true;
595                     }
596                 }
597                 if (!madeTable)
598                 {
599                     bool isPlural = (ftype2 != -1);
600                     gmx_fatal(FARGS,
601                               "Tabulated interaction of type '%s%s%s' with index %d cannot be used "
602                               "because no table file whose name matched '%s' was passed via the "
603                               "gmx mdrun -tableb command-line option.",
604                               interaction_function[ftype1].longname, isPlural ? "' or '" : "",
605                               isPlural ? interaction_function[ftype2].longname : "", i,
606                               patternToFind.c_str());
607                 }
608             }
609         }
610         sfree(count);
611     }
612
613     return tab;
614 }
615
616 void forcerec_set_ranges(t_forcerec* fr, int natoms_force, int natoms_force_constr, int natoms_f_novirsum)
617 {
618     fr->natoms_force        = natoms_force;
619     fr->natoms_force_constr = natoms_force_constr;
620
621     if (fr->natoms_force_constr > fr->nalloc_force)
622     {
623         fr->nalloc_force = over_alloc_dd(fr->natoms_force_constr);
624     }
625
626     if (fr->haveDirectVirialContributions)
627     {
628         fr->forceBufferForDirectVirialContributions.resize(natoms_f_novirsum);
629     }
630 }
631
632 static real cutoff_inf(real cutoff)
633 {
634     if (cutoff == 0)
635     {
636         cutoff = GMX_CUTOFF_INF;
637     }
638
639     return cutoff;
640 }
641
642 /*! \brief Print Coulomb Ewald citations and set ewald coefficients */
643 static void initCoulombEwaldParameters(FILE*                fp,
644                                        const t_inputrec*    ir,
645                                        bool                 systemHasNetCharge,
646                                        interaction_const_t* ic)
647 {
648     if (!EEL_PME_EWALD(ir->coulombtype))
649     {
650         return;
651     }
652
653     if (fp)
654     {
655         fprintf(fp, "Will do PME sum in reciprocal space for electrostatic interactions.\n");
656
657         if (ir->coulombtype == eelP3M_AD)
658         {
659             please_cite(fp, "Hockney1988");
660             please_cite(fp, "Ballenegger2012");
661         }
662         else
663         {
664             please_cite(fp, "Essmann95a");
665         }
666
667         if (ir->ewald_geometry == eewg3DC)
668         {
669             if (fp)
670             {
671                 fprintf(fp, "Using the Ewald3DC correction for systems with a slab geometry%s.\n",
672                         systemHasNetCharge ? " and net charge" : "");
673             }
674             please_cite(fp, "In-Chul99a");
675             if (systemHasNetCharge)
676             {
677                 please_cite(fp, "Ballenegger2009");
678             }
679         }
680     }
681
682     ic->ewaldcoeff_q = calc_ewaldcoeff_q(ir->rcoulomb, ir->ewald_rtol);
683     if (fp)
684     {
685         fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for Ewald\n", 1 / ic->ewaldcoeff_q);
686     }
687
688     if (ic->coulomb_modifier == eintmodPOTSHIFT)
689     {
690         GMX_RELEASE_ASSERT(ic->rcoulomb != 0, "Cutoff radius cannot be zero");
691         ic->sh_ewald = std::erfc(ic->ewaldcoeff_q * ic->rcoulomb) / ic->rcoulomb;
692     }
693     else
694     {
695         ic->sh_ewald = 0;
696     }
697 }
698
699 /*! \brief Print Van der Waals Ewald citations and set ewald coefficients */
700 static void initVdwEwaldParameters(FILE* fp, const t_inputrec* ir, interaction_const_t* ic)
701 {
702     if (!EVDW_PME(ir->vdwtype))
703     {
704         return;
705     }
706
707     if (fp)
708     {
709         fprintf(fp, "Will do PME sum in reciprocal space for LJ dispersion interactions.\n");
710         please_cite(fp, "Essmann95a");
711     }
712     ic->ewaldcoeff_lj = calc_ewaldcoeff_lj(ir->rvdw, ir->ewald_rtol_lj);
713     if (fp)
714     {
715         fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for LJ Ewald\n", 1 / ic->ewaldcoeff_lj);
716     }
717
718     if (ic->vdw_modifier == eintmodPOTSHIFT)
719     {
720         real crc2       = gmx::square(ic->ewaldcoeff_lj * ic->rvdw);
721         ic->sh_lj_ewald = (std::exp(-crc2) * (1 + crc2 + 0.5 * crc2 * crc2) - 1) / gmx::power6(ic->rvdw);
722     }
723     else
724     {
725         ic->sh_lj_ewald = 0;
726     }
727 }
728
729 /* Generate Coulomb and/or Van der Waals Ewald long-range correction tables
730  *
731  * Tables are generated for one or both, depending on if the pointers
732  * are non-null. The spacing for both table sets is the same and obeys
733  * both accuracy requirements, when relevant.
734  */
735 static void init_ewald_f_table(const interaction_const_t& ic,
736                                EwaldCorrectionTables*     coulombTables,
737                                EwaldCorrectionTables*     vdwTables)
738 {
739     const bool useCoulombTable = (EEL_PME_EWALD(ic.eeltype) && coulombTables != nullptr);
740     const bool useVdwTable     = (EVDW_PME(ic.vdwtype) && vdwTables != nullptr);
741
742     /* Get the Ewald table spacing based on Coulomb and/or LJ
743      * Ewald coefficients and rtol.
744      */
745     const real tableScale = ewald_spline3_table_scale(ic, useCoulombTable, useVdwTable);
746
747     const int tableSize = static_cast<int>(ic.rcoulomb * tableScale) + 2;
748
749     if (useCoulombTable)
750     {
751         *coulombTables =
752                 generateEwaldCorrectionTables(tableSize, tableScale, ic.ewaldcoeff_q, v_q_ewald_lr);
753     }
754
755     if (useVdwTable)
756     {
757         *vdwTables = generateEwaldCorrectionTables(tableSize, tableScale, ic.ewaldcoeff_lj, v_lj_ewald_lr);
758     }
759 }
760
761 void init_interaction_const_tables(FILE* fp, interaction_const_t* ic)
762 {
763     if (EEL_PME_EWALD(ic->eeltype))
764     {
765         init_ewald_f_table(*ic, ic->coulombEwaldTables.get(), nullptr);
766         if (fp != nullptr)
767         {
768             fprintf(fp, "Initialized non-bonded Coulomb Ewald tables, spacing: %.2e size: %zu\n\n",
769                     1 / ic->coulombEwaldTables->scale, ic->coulombEwaldTables->tableF.size());
770         }
771     }
772 }
773
774 static void clear_force_switch_constants(shift_consts_t* sc)
775 {
776     sc->c2   = 0;
777     sc->c3   = 0;
778     sc->cpot = 0;
779 }
780
781 static void force_switch_constants(real p, real rsw, real rc, shift_consts_t* sc)
782 {
783     /* Here we determine the coefficient for shifting the force to zero
784      * between distance rsw and the cut-off rc.
785      * For a potential of r^-p, we have force p*r^-(p+1).
786      * But to save flops we absorb p in the coefficient.
787      * Thus we get:
788      * force/p   = r^-(p+1) + c2*r^2 + c3*r^3
789      * potential = r^-p + c2/3*r^3 + c3/4*r^4 + cpot
790      */
791     sc->c2   = ((p + 1) * rsw - (p + 4) * rc) / (pow(rc, p + 2) * gmx::square(rc - rsw));
792     sc->c3   = -((p + 1) * rsw - (p + 3) * rc) / (pow(rc, p + 2) * gmx::power3(rc - rsw));
793     sc->cpot = -pow(rc, -p) + p * sc->c2 / 3 * gmx::power3(rc - rsw)
794                + p * sc->c3 / 4 * gmx::power4(rc - rsw);
795 }
796
797 static void potential_switch_constants(real rsw, real rc, switch_consts_t* sc)
798 {
799     /* The switch function is 1 at rsw and 0 at rc.
800      * The derivative and second derivate are zero at both ends.
801      * rsw        = max(r - r_switch, 0)
802      * sw         = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
803      * dsw        = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
804      * force      = force*dsw - potential*sw
805      * potential *= sw
806      */
807     sc->c3 = -10 / gmx::power3(rc - rsw);
808     sc->c4 = 15 / gmx::power4(rc - rsw);
809     sc->c5 = -6 / gmx::power5(rc - rsw);
810 }
811
812 /*! \brief Construct interaction constants
813  *
814  * This data is used (particularly) by search and force code for
815  * short-range interactions. Many of these are constant for the whole
816  * simulation; some are constant only after PME tuning completes.
817  */
818 static void init_interaction_const(FILE*                 fp,
819                                    interaction_const_t** interaction_const,
820                                    const t_inputrec*     ir,
821                                    const gmx_mtop_t*     mtop,
822                                    bool                  systemHasNetCharge)
823 {
824     interaction_const_t* ic = new interaction_const_t;
825
826     ic->cutoff_scheme = ir->cutoff_scheme;
827
828     ic->coulombEwaldTables = std::make_unique<EwaldCorrectionTables>();
829
830     /* Lennard-Jones */
831     ic->vdwtype         = ir->vdwtype;
832     ic->vdw_modifier    = ir->vdw_modifier;
833     ic->reppow          = mtop->ffparams.reppow;
834     ic->rvdw            = cutoff_inf(ir->rvdw);
835     ic->rvdw_switch     = ir->rvdw_switch;
836     ic->ljpme_comb_rule = ir->ljpme_combination_rule;
837     ic->useBuckingham   = (mtop->ffparams.functype[0] == F_BHAM);
838     if (ic->useBuckingham)
839     {
840         ic->buckinghamBMax = calcBuckinghamBMax(fp, mtop);
841     }
842
843     initVdwEwaldParameters(fp, ir, ic);
844
845     clear_force_switch_constants(&ic->dispersion_shift);
846     clear_force_switch_constants(&ic->repulsion_shift);
847
848     switch (ic->vdw_modifier)
849     {
850         case eintmodPOTSHIFT:
851             /* Only shift the potential, don't touch the force */
852             ic->dispersion_shift.cpot = -1.0 / gmx::power6(ic->rvdw);
853             ic->repulsion_shift.cpot  = -1.0 / gmx::power12(ic->rvdw);
854             break;
855         case eintmodFORCESWITCH:
856             /* Switch the force, switch and shift the potential */
857             force_switch_constants(6.0, ic->rvdw_switch, ic->rvdw, &ic->dispersion_shift);
858             force_switch_constants(12.0, ic->rvdw_switch, ic->rvdw, &ic->repulsion_shift);
859             break;
860         case eintmodPOTSWITCH:
861             /* Switch the potential and force */
862             potential_switch_constants(ic->rvdw_switch, ic->rvdw, &ic->vdw_switch);
863             break;
864         case eintmodNONE:
865         case eintmodEXACTCUTOFF:
866             /* Nothing to do here */
867             break;
868         default: gmx_incons("unimplemented potential modifier");
869     }
870
871     /* Electrostatics */
872     ic->eeltype          = ir->coulombtype;
873     ic->coulomb_modifier = ir->coulomb_modifier;
874     ic->rcoulomb         = cutoff_inf(ir->rcoulomb);
875     ic->rcoulomb_switch  = ir->rcoulomb_switch;
876     ic->epsilon_r        = ir->epsilon_r;
877
878     /* Set the Coulomb energy conversion factor */
879     if (ic->epsilon_r != 0)
880     {
881         ic->epsfac = ONE_4PI_EPS0 / ic->epsilon_r;
882     }
883     else
884     {
885         /* eps = 0 is infinite dieletric: no Coulomb interactions */
886         ic->epsfac = 0;
887     }
888
889     /* Reaction-field */
890     if (EEL_RF(ic->eeltype))
891     {
892         GMX_RELEASE_ASSERT(ic->eeltype != eelGRF_NOTUSED, "GRF is no longer supported");
893         ic->epsilon_rf = ir->epsilon_rf;
894
895         calc_rffac(fp, ic->epsilon_r, ic->epsilon_rf, ic->rcoulomb, &ic->k_rf, &ic->c_rf);
896     }
897     else
898     {
899         /* For plain cut-off we might use the reaction-field kernels */
900         ic->epsilon_rf = ic->epsilon_r;
901         ic->k_rf       = 0;
902         if (ir->coulomb_modifier == eintmodPOTSHIFT)
903         {
904             ic->c_rf = 1 / ic->rcoulomb;
905         }
906         else
907         {
908             ic->c_rf = 0;
909         }
910     }
911
912     initCoulombEwaldParameters(fp, ir, systemHasNetCharge, ic);
913
914     if (fp != nullptr)
915     {
916         real dispersion_shift;
917
918         dispersion_shift = ic->dispersion_shift.cpot;
919         if (EVDW_PME(ic->vdwtype))
920         {
921             dispersion_shift -= ic->sh_lj_ewald;
922         }
923         fprintf(fp, "Potential shift: LJ r^-12: %.3e r^-6: %.3e", ic->repulsion_shift.cpot, dispersion_shift);
924
925         if (ic->eeltype == eelCUT)
926         {
927             fprintf(fp, ", Coulomb %.e", -ic->c_rf);
928         }
929         else if (EEL_PME(ic->eeltype))
930         {
931             fprintf(fp, ", Ewald %.3e", -ic->sh_ewald);
932         }
933         fprintf(fp, "\n");
934     }
935
936     *interaction_const = ic;
937 }
938
939 void init_forcerec(FILE*                            fp,
940                    const gmx::MDLogger&             mdlog,
941                    t_forcerec*                      fr,
942                    t_fcdata*                        fcd,
943                    const t_inputrec*                ir,
944                    const gmx_mtop_t*                mtop,
945                    const t_commrec*                 cr,
946                    matrix                           box,
947                    const char*                      tabfn,
948                    const char*                      tabpfn,
949                    gmx::ArrayRef<const std::string> tabbfnm,
950                    const gmx_hw_info_t&             hardwareInfo,
951                    const gmx_device_info_t*         deviceInfo,
952                    const bool                       useGpuForBonded,
953                    const bool                       pmeOnlyRankUsesGpu,
954                    real                             print_force,
955                    gmx_wallcycle*                   wcycle)
956 {
957     real     rtab;
958     char*    env;
959     double   dbl;
960     gmx_bool bFEP_NonBonded;
961
962     /* By default we turn SIMD kernels on, but it might be turned off further down... */
963     fr->use_simd_kernels = TRUE;
964
965     if (check_box(ir->ePBC, box))
966     {
967         gmx_fatal(FARGS, "%s", check_box(ir->ePBC, box));
968     }
969
970     /* Test particle insertion ? */
971     if (EI_TPI(ir->eI))
972     {
973         /* Set to the size of the molecule to be inserted (the last one) */
974         gmx::RangePartitioning molecules = gmx_mtop_molecules(*mtop);
975         fr->n_tpi                        = molecules.block(molecules.numBlocks() - 1).size();
976     }
977     else
978     {
979         fr->n_tpi = 0;
980     }
981
982     if (ir->coulombtype == eelRF_NEC_UNSUPPORTED || ir->coulombtype == eelGRF_NOTUSED)
983     {
984         gmx_fatal(FARGS, "%s electrostatics is no longer supported", eel_names[ir->coulombtype]);
985     }
986
987     if (ir->bAdress)
988     {
989         gmx_fatal(FARGS, "AdResS simulations are no longer supported");
990     }
991     if (ir->useTwinRange)
992     {
993         gmx_fatal(FARGS, "Twin-range simulations are no longer supported");
994     }
995     /* Copy the user determined parameters */
996     fr->userint1  = ir->userint1;
997     fr->userint2  = ir->userint2;
998     fr->userint3  = ir->userint3;
999     fr->userint4  = ir->userint4;
1000     fr->userreal1 = ir->userreal1;
1001     fr->userreal2 = ir->userreal2;
1002     fr->userreal3 = ir->userreal3;
1003     fr->userreal4 = ir->userreal4;
1004
1005     /* Shell stuff */
1006     fr->fc_stepsize = ir->fc_stepsize;
1007
1008     /* Free energy */
1009     fr->efep        = ir->efep;
1010     fr->sc_alphavdw = ir->fepvals->sc_alpha;
1011     if (ir->fepvals->bScCoul)
1012     {
1013         fr->sc_alphacoul  = ir->fepvals->sc_alpha;
1014         fr->sc_sigma6_min = gmx::power6(ir->fepvals->sc_sigma_min);
1015     }
1016     else
1017     {
1018         fr->sc_alphacoul  = 0;
1019         fr->sc_sigma6_min = 0; /* only needed when bScCoul is on */
1020     }
1021     fr->sc_power      = ir->fepvals->sc_power;
1022     fr->sc_r_power    = ir->fepvals->sc_r_power;
1023     fr->sc_sigma6_def = gmx::power6(ir->fepvals->sc_sigma);
1024
1025     env = getenv("GMX_SCSIGMA_MIN");
1026     if (env != nullptr)
1027     {
1028         dbl = 0;
1029         sscanf(env, "%20lf", &dbl);
1030         fr->sc_sigma6_min = gmx::power6(dbl);
1031         if (fp)
1032         {
1033             fprintf(fp, "Setting the minimum soft core sigma to %g nm\n", dbl);
1034         }
1035     }
1036
1037     fr->bNonbonded = TRUE;
1038     if (getenv("GMX_NO_NONBONDED") != nullptr)
1039     {
1040         /* turn off non-bonded calculations */
1041         fr->bNonbonded = FALSE;
1042         GMX_LOG(mdlog.warning)
1043                 .asParagraph()
1044                 .appendText(
1045                         "Found environment variable GMX_NO_NONBONDED.\n"
1046                         "Disabling nonbonded calculations.");
1047     }
1048
1049     if ((getenv("GMX_DISABLE_SIMD_KERNELS") != nullptr) || (getenv("GMX_NOOPTIMIZEDKERNELS") != nullptr))
1050     {
1051         fr->use_simd_kernels = FALSE;
1052         if (fp != nullptr)
1053         {
1054             fprintf(fp,
1055                     "\nFound environment variable GMX_DISABLE_SIMD_KERNELS.\n"
1056                     "Disabling the usage of any SIMD-specific non-bonded & bonded kernel routines\n"
1057                     "(e.g. SSE2/SSE4.1/AVX).\n\n");
1058         }
1059     }
1060
1061     fr->bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
1062
1063     /* Neighbour searching stuff */
1064     fr->cutoff_scheme = ir->cutoff_scheme;
1065     fr->ePBC          = ir->ePBC;
1066
1067     /* Determine if we will do PBC for distances in bonded interactions */
1068     if (fr->ePBC == epbcNONE)
1069     {
1070         fr->bMolPBC = FALSE;
1071     }
1072     else
1073     {
1074         const bool useEwaldSurfaceCorrection =
1075                 (EEL_PME_EWALD(ir->coulombtype) && ir->epsilon_surface != 0);
1076         if (!DOMAINDECOMP(cr))
1077         {
1078             gmx_bool bSHAKE;
1079
1080             bSHAKE = (ir->eConstrAlg == econtSHAKE
1081                       && (gmx_mtop_ftype_count(mtop, F_CONSTR) > 0
1082                           || gmx_mtop_ftype_count(mtop, F_CONSTRNC) > 0));
1083
1084             /* The group cut-off scheme and SHAKE assume charge groups
1085              * are whole, but not using molpbc is faster in most cases.
1086              * With intermolecular interactions we need PBC for calculating
1087              * distances between atoms in different molecules.
1088              */
1089             if (bSHAKE && !mtop->bIntermolecularInteractions)
1090             {
1091                 fr->bMolPBC = ir->bPeriodicMols;
1092
1093                 if (bSHAKE && fr->bMolPBC)
1094                 {
1095                     gmx_fatal(FARGS, "SHAKE is not supported with periodic molecules");
1096                 }
1097             }
1098             else
1099             {
1100                 /* Not making molecules whole is faster in most cases,
1101                  * but with orientation restraints or non-tinfoil boundary
1102                  * conditions we need whole molecules.
1103                  */
1104                 fr->bMolPBC = (fcd->orires.nr == 0 && !useEwaldSurfaceCorrection);
1105
1106                 if (getenv("GMX_USE_GRAPH") != nullptr)
1107                 {
1108                     fr->bMolPBC = FALSE;
1109                     if (fp)
1110                     {
1111                         GMX_LOG(mdlog.warning)
1112                                 .asParagraph()
1113                                 .appendText(
1114                                         "GMX_USE_GRAPH is set, using the graph for bonded "
1115                                         "interactions");
1116                     }
1117
1118                     if (mtop->bIntermolecularInteractions)
1119                     {
1120                         GMX_LOG(mdlog.warning)
1121                                 .asParagraph()
1122                                 .appendText(
1123                                         "WARNING: Molecules linked by intermolecular interactions "
1124                                         "have to reside in the same periodic image, otherwise "
1125                                         "artifacts will occur!");
1126                     }
1127                 }
1128
1129                 GMX_RELEASE_ASSERT(
1130                         fr->bMolPBC || !mtop->bIntermolecularInteractions,
1131                         "We need to use PBC within molecules with inter-molecular interactions");
1132
1133                 if (bSHAKE && fr->bMolPBC)
1134                 {
1135                     gmx_fatal(FARGS,
1136                               "SHAKE is not properly supported with intermolecular interactions. "
1137                               "For short simulations where linked molecules remain in the same "
1138                               "periodic image, the environment variable GMX_USE_GRAPH can be used "
1139                               "to override this check.\n");
1140                 }
1141             }
1142         }
1143         else
1144         {
1145             fr->bMolPBC = dd_bonded_molpbc(cr->dd, fr->ePBC);
1146
1147             if (useEwaldSurfaceCorrection && !dd_moleculesAreAlwaysWhole(*cr->dd))
1148             {
1149                 gmx_fatal(FARGS,
1150                           "You requested dipole correction (epsilon_surface > 0), but molecules "
1151                           "are broken "
1152                           "over periodic boundary conditions by the domain decomposition. "
1153                           "Run without domain decomposition instead.");
1154             }
1155         }
1156
1157         if (useEwaldSurfaceCorrection)
1158         {
1159             GMX_RELEASE_ASSERT((!DOMAINDECOMP(cr) && !fr->bMolPBC)
1160                                        || (DOMAINDECOMP(cr) && dd_moleculesAreAlwaysWhole(*cr->dd)),
1161                                "Molecules can not be broken by PBC with epsilon_surface > 0");
1162         }
1163     }
1164
1165     fr->rc_scaling = ir->refcoord_scaling;
1166     copy_rvec(ir->posres_com, fr->posres_com);
1167     copy_rvec(ir->posres_comB, fr->posres_comB);
1168     fr->rlist                  = cutoff_inf(ir->rlist);
1169     fr->ljpme_combination_rule = ir->ljpme_combination_rule;
1170
1171     /* This now calculates sum for q and c6*/
1172     bool systemHasNetCharge = set_chargesum(fp, fr, mtop);
1173
1174     /* fr->ic is used both by verlet and group kernels (to some extent) now */
1175     init_interaction_const(fp, &fr->ic, ir, mtop, systemHasNetCharge);
1176     init_interaction_const_tables(fp, fr->ic);
1177
1178     const interaction_const_t* ic = fr->ic;
1179
1180     /* TODO: Replace this Ewald table or move it into interaction_const_t */
1181     if (ir->coulombtype == eelEWALD)
1182     {
1183         init_ewald_tab(&(fr->ewald_table), ir, fp);
1184     }
1185
1186     /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
1187     switch (ic->eeltype)
1188     {
1189         case eelCUT: fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_COULOMB; break;
1190
1191         case eelRF:
1192         case eelRF_ZERO: fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD; break;
1193
1194         case eelSWITCH:
1195         case eelSHIFT:
1196         case eelUSER:
1197         case eelENCADSHIFT:
1198         case eelPMESWITCH:
1199         case eelPMEUSER:
1200         case eelPMEUSERSWITCH:
1201             fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
1202             break;
1203
1204         case eelPME:
1205         case eelP3M_AD:
1206         case eelEWALD: fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_EWALD; break;
1207
1208         default:
1209             gmx_fatal(FARGS, "Unsupported electrostatic interaction: %s", eel_names[ic->eeltype]);
1210     }
1211     fr->nbkernel_elec_modifier = ic->coulomb_modifier;
1212
1213     /* Vdw: Translate from mdp settings to kernel format */
1214     switch (ic->vdwtype)
1215     {
1216         case evdwCUT:
1217             if (fr->bBHAM)
1218             {
1219                 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_BUCKINGHAM;
1220             }
1221             else
1222             {
1223                 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LENNARDJONES;
1224             }
1225             break;
1226         case evdwPME: fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LJEWALD; break;
1227
1228         case evdwSWITCH:
1229         case evdwSHIFT:
1230         case evdwUSER:
1231         case evdwENCADSHIFT:
1232             fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
1233             break;
1234
1235         default: gmx_fatal(FARGS, "Unsupported vdw interaction: %s", evdw_names[ic->vdwtype]);
1236     }
1237     fr->nbkernel_vdw_modifier = ic->vdw_modifier;
1238
1239     if (ir->cutoff_scheme == ecutsVERLET)
1240     {
1241         if (!gmx_within_tol(ic->reppow, 12.0, 10 * GMX_DOUBLE_EPS))
1242         {
1243             gmx_fatal(FARGS, "Cut-off scheme %s only supports LJ repulsion power 12",
1244                       ecutscheme_names[ir->cutoff_scheme]);
1245         }
1246         /* Older tpr files can contain Coulomb user tables with the Verlet cutoff-scheme,
1247          * while mdrun does not (and never did) support this.
1248          */
1249         if (EEL_USER(fr->ic->eeltype))
1250         {
1251             gmx_fatal(FARGS, "Combination of %s and cutoff scheme %s is not supported",
1252                       eel_names[ir->coulombtype], ecutscheme_names[ir->cutoff_scheme]);
1253         }
1254
1255         fr->bvdwtab  = FALSE;
1256         fr->bcoultab = FALSE;
1257     }
1258
1259     /* 1-4 interaction electrostatics */
1260     fr->fudgeQQ = mtop->ffparams.fudgeQQ;
1261
1262     fr->haveDirectVirialContributions =
1263             (EEL_FULL(ic->eeltype) || EVDW_PME(ic->vdwtype) || fr->forceProviders->hasForceProvider()
1264              || gmx_mtop_ftype_count(mtop, F_POSRES) > 0 || gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0
1265              || ir->nwall > 0 || ir->bPull || ir->bRot || ir->bIMD);
1266
1267     if (fr->shift_vec == nullptr)
1268     {
1269         snew(fr->shift_vec, SHIFTS);
1270     }
1271
1272     fr->shiftForces.resize(SHIFTS);
1273
1274     if (fr->nbfp == nullptr)
1275     {
1276         fr->ntype = mtop->ffparams.atnr;
1277         fr->nbfp  = mk_nbfp(&mtop->ffparams, fr->bBHAM);
1278         if (EVDW_PME(ic->vdwtype))
1279         {
1280             fr->ljpme_c6grid = make_ljpme_c6grid(&mtop->ffparams, fr);
1281         }
1282     }
1283
1284     /* Copy the energy group exclusions */
1285     fr->egp_flags = ir->opts.egp_flags;
1286
1287     /* Van der Waals stuff */
1288     if ((ic->vdwtype != evdwCUT) && (ic->vdwtype != evdwUSER) && !fr->bBHAM)
1289     {
1290         if (ic->rvdw_switch >= ic->rvdw)
1291         {
1292             gmx_fatal(FARGS, "rvdw_switch (%f) must be < rvdw (%f)", ic->rvdw_switch, ic->rvdw);
1293         }
1294         if (fp)
1295         {
1296             fprintf(fp, "Using %s Lennard-Jones, switch between %g and %g nm\n",
1297                     (ic->eeltype == eelSWITCH) ? "switched" : "shifted", ic->rvdw_switch, ic->rvdw);
1298         }
1299     }
1300
1301     if (fr->bBHAM && EVDW_PME(ic->vdwtype))
1302     {
1303         gmx_fatal(FARGS, "LJ PME not supported with Buckingham");
1304     }
1305
1306     if (fr->bBHAM && (ic->vdwtype == evdwSHIFT || ic->vdwtype == evdwSWITCH))
1307     {
1308         gmx_fatal(FARGS, "Switch/shift interaction not supported with Buckingham");
1309     }
1310
1311     if (fr->bBHAM && fr->cutoff_scheme == ecutsVERLET)
1312     {
1313         gmx_fatal(FARGS, "Verlet cutoff-scheme is not supported with Buckingham");
1314     }
1315
1316     if (fp && fr->cutoff_scheme == ecutsGROUP)
1317     {
1318         fprintf(fp, "Cut-off's:   NS: %g   Coulomb: %g   %s: %g\n", fr->rlist, ic->rcoulomb,
1319                 fr->bBHAM ? "BHAM" : "LJ", ic->rvdw);
1320     }
1321
1322     if (ir->implicit_solvent)
1323     {
1324         gmx_fatal(FARGS, "Implict solvation is no longer supported.");
1325     }
1326
1327
1328     /* This code automatically gives table length tabext without cut-off's,
1329      * in that case grompp should already have checked that we do not need
1330      * normal tables and we only generate tables for 1-4 interactions.
1331      */
1332     rtab = ir->rlist + ir->tabext;
1333
1334     /* We want to use unmodified tables for 1-4 coulombic
1335      * interactions, so we must in general have an extra set of
1336      * tables. */
1337     if (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 || gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0
1338         || gmx_mtop_ftype_count(mtop, F_LJC_PAIRS_NB) > 0)
1339     {
1340         fr->pairsTable = make_tables(fp, ic, tabpfn, rtab, GMX_MAKETABLES_14ONLY);
1341     }
1342
1343     /* Wall stuff */
1344     fr->nwall = ir->nwall;
1345     if (ir->nwall && ir->wall_type == ewtTABLE)
1346     {
1347         make_wall_tables(fp, ir, tabfn, &mtop->groups, fr);
1348     }
1349
1350     if (fcd && !tabbfnm.empty())
1351     {
1352         // Need to catch std::bad_alloc
1353         // TODO Don't need to catch this here, when merging with master branch
1354         try
1355         {
1356             fcd->bondtab  = make_bonded_tables(fp, F_TABBONDS, F_TABBONDSNC, mtop, tabbfnm, "b");
1357             fcd->angletab = make_bonded_tables(fp, F_TABANGLES, -1, mtop, tabbfnm, "a");
1358             fcd->dihtab   = make_bonded_tables(fp, F_TABDIHS, -1, mtop, tabbfnm, "d");
1359         }
1360         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1361     }
1362     else
1363     {
1364         if (debug)
1365         {
1366             fprintf(debug,
1367                     "No fcdata or table file name passed, can not read table, can not do bonded "
1368                     "interactions\n");
1369         }
1370     }
1371
1372     // QM/MM initialization if requested
1373     fr->bQMMM = ir->bQMMM;
1374     if (fr->bQMMM)
1375     {
1376         // Initialize QM/MM if supported
1377         if (GMX_QMMM)
1378         {
1379             GMX_LOG(mdlog.info)
1380                     .asParagraph()
1381                     .appendText(
1382                             "Large parts of the QM/MM support is deprecated, and may be removed in "
1383                             "a future "
1384                             "version. Please get in touch with the developers if you find the "
1385                             "support useful, "
1386                             "as help is needed if the functionality is to continue to be "
1387                             "available.");
1388             fr->qr = mk_QMMMrec();
1389             init_QMMMrec(cr, mtop, ir, fr);
1390         }
1391         else
1392         {
1393             gmx_incons(
1394                     "QM/MM was requested, but is only available when GROMACS "
1395                     "is configured with QM/MM support");
1396         }
1397     }
1398
1399     /* Set all the static charge group info */
1400     fr->cginfo_mb = init_cginfo_mb(mtop, fr, &bFEP_NonBonded);
1401     if (!DOMAINDECOMP(cr))
1402     {
1403         fr->cginfo = cginfo_expand(mtop->molblock.size(), fr->cginfo_mb);
1404     }
1405
1406     if (!DOMAINDECOMP(cr))
1407     {
1408         forcerec_set_ranges(fr, mtop->natoms, mtop->natoms, mtop->natoms);
1409     }
1410
1411     fr->print_force = print_force;
1412
1413     /* Initialize the thread working data for bonded interactions */
1414     fr->bondedThreading = init_bonded_threading(
1415             fp, mtop->groups.groups[SimulationAtomGroupType::EnergyOutput].size());
1416
1417     fr->nthread_ewc = gmx_omp_nthreads_get(emntBonded);
1418     snew(fr->ewc_t, fr->nthread_ewc);
1419
1420     if (fr->cutoff_scheme == ecutsVERLET)
1421     {
1422         // We checked the cut-offs in grompp, but double-check here.
1423         // We have PME+LJcutoff kernels for rcoulomb>rvdw.
1424         if (EEL_PME_EWALD(ir->coulombtype) && ir->vdwtype == eelCUT)
1425         {
1426             GMX_RELEASE_ASSERT(ir->rcoulomb >= ir->rvdw,
1427                                "With Verlet lists and PME we should have rcoulomb>=rvdw");
1428         }
1429         else
1430         {
1431             GMX_RELEASE_ASSERT(
1432                     ir->rcoulomb == ir->rvdw,
1433                     "With Verlet lists and no PME rcoulomb and rvdw should be identical");
1434         }
1435
1436         fr->nbv = Nbnxm::init_nb_verlet(mdlog, bFEP_NonBonded, ir, fr, cr, hardwareInfo, deviceInfo,
1437                                         mtop, box, wcycle);
1438
1439         if (useGpuForBonded)
1440         {
1441             auto stream = DOMAINDECOMP(cr)
1442                                   ? Nbnxm::gpu_get_command_stream(
1443                                             fr->nbv->gpu_nbv, gmx::InteractionLocality::NonLocal)
1444                                   : Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv,
1445                                                                   gmx::InteractionLocality::Local);
1446             // TODO the heap allocation is only needed while
1447             // t_forcerec lacks a constructor.
1448             fr->gpuBonded = new gmx::GpuBonded(mtop->ffparams, stream, wcycle);
1449         }
1450     }
1451
1452     if (ir->eDispCorr != edispcNO)
1453     {
1454         fr->dispersionCorrection = std::make_unique<DispersionCorrection>(
1455                 *mtop, *ir, fr->bBHAM, fr->ntype,
1456                 gmx::arrayRefFromArray(fr->nbfp, fr->ntype * fr->ntype * 2), *fr->ic, tabfn);
1457         fr->dispersionCorrection->print(mdlog);
1458     }
1459
1460     if (fp != nullptr)
1461     {
1462         /* Here we switch from using mdlog, which prints the newline before
1463          * the paragraph, to our old fprintf logging, which prints the newline
1464          * after the paragraph, so we should add a newline here.
1465          */
1466         fprintf(fp, "\n");
1467     }
1468
1469     if (pmeOnlyRankUsesGpu && c_enableGpuPmePpComms)
1470     {
1471         fr->pmePpCommGpu = std::make_unique<gmx::PmePpCommGpu>(cr->mpi_comm_mysim, cr->dd->pme_nodeid);
1472     }
1473 }
1474
1475 t_forcerec::t_forcerec() = default;
1476
1477 t_forcerec::~t_forcerec() = default;
1478
1479 /* Frees GPU memory and sets a tMPI node barrier.
1480  *
1481  * Note that this function needs to be called even if GPUs are not used
1482  * in this run because the PME ranks have no knowledge of whether GPUs
1483  * are used or not, but all ranks need to enter the barrier below.
1484  * \todo Remove physical node barrier from this function after making sure
1485  * that it's not needed anymore (with a shared GPU run).
1486  */
1487 void free_gpu_resources(t_forcerec*                          fr,
1488                         const gmx::PhysicalNodeCommunicator& physicalNodeCommunicator,
1489                         const gmx_gpu_info_t&                gpu_info)
1490 {
1491     bool isPPrankUsingGPU = (fr != nullptr) && (fr->nbv != nullptr) && fr->nbv->useGpu();
1492
1493     /* stop the GPU profiler (only CUDA) */
1494     if (gpu_info.n_dev > 0)
1495     {
1496         stopGpuProfiler();
1497     }
1498
1499     if (isPPrankUsingGPU)
1500     {
1501         /* Free data in GPU memory and pinned memory before destroying the GPU context */
1502         fr->nbv.reset();
1503
1504         delete fr->gpuBonded;
1505         fr->gpuBonded = nullptr;
1506     }
1507
1508     /* With tMPI we need to wait for all ranks to finish deallocation before
1509      * destroying the CUDA context in free_gpu() as some tMPI ranks may be sharing
1510      * GPU and context.
1511      *
1512      * This is not a concern in OpenCL where we use one context per rank which
1513      * is freed in nbnxn_gpu_free().
1514      *
1515      * Note: it is safe to not call the barrier on the ranks which do not use GPU,
1516      * but it is easier and more futureproof to call it on the whole node.
1517      */
1518     if (GMX_THREAD_MPI)
1519     {
1520         physicalNodeCommunicator.barrier();
1521     }
1522 }
1523
1524 void done_forcerec(t_forcerec* fr, int numMolBlocks)
1525 {
1526     if (fr == nullptr)
1527     {
1528         // PME-only ranks don't have a forcerec
1529         return;
1530     }
1531     done_cginfo_mb(fr->cginfo_mb, numMolBlocks);
1532     sfree(fr->nbfp);
1533     delete fr->ic;
1534     sfree(fr->shift_vec);
1535     sfree(fr->ewc_t);
1536     tear_down_bonded_threading(fr->bondedThreading);
1537     GMX_RELEASE_ASSERT(fr->gpuBonded == nullptr, "Should have been deleted earlier, when used");
1538     fr->bondedThreading = nullptr;
1539     delete fr;
1540 }