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