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