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