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