d66280e64db08698502d1192bab4d2f55df24202
[alexxy/gromacs.git] / src / gromacs / mdlib / forcerec.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #include "gmxpre.h"
38
39 #include "forcerec.h"
40
41 #include "config.h"
42
43 #include <cassert>
44 #include <cmath>
45 #include <cstdlib>
46 #include <cstring>
47
48 #include <algorithm>
49 #include <memory>
50
51 #include "gromacs/commandline/filenm.h"
52 #include "gromacs/domdec/domdec.h"
53 #include "gromacs/domdec/domdec_struct.h"
54 #include "gromacs/ewald/ewald.h"
55 #include "gromacs/ewald/ewald_utils.h"
56 #include "gromacs/fileio/filetypes.h"
57 #include "gromacs/gmxlib/network.h"
58 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
59 #include "gromacs/gpu_utils/gpu_utils.h"
60 #include "gromacs/hardware/hw_info.h"
61 #include "gromacs/listed_forces/gpubonded.h"
62 #include "gromacs/listed_forces/manage_threading.h"
63 #include "gromacs/listed_forces/pairs.h"
64 #include "gromacs/math/functions.h"
65 #include "gromacs/math/units.h"
66 #include "gromacs/math/utilities.h"
67 #include "gromacs/math/vec.h"
68 #include "gromacs/mdlib/dispersioncorrection.h"
69 #include "gromacs/mdlib/force.h"
70 #include "gromacs/mdlib/forcerec_threading.h"
71 #include "gromacs/mdlib/gmx_omp_nthreads.h"
72 #include "gromacs/mdlib/md_support.h"
73 #include "gromacs/mdlib/qmmm.h"
74 #include "gromacs/mdlib/rf_util.h"
75 #include "gromacs/mdlib/wall.h"
76 #include "gromacs/mdtypes/commrec.h"
77 #include "gromacs/mdtypes/fcdata.h"
78 #include "gromacs/mdtypes/group.h"
79 #include "gromacs/mdtypes/iforceprovider.h"
80 #include "gromacs/mdtypes/inputrec.h"
81 #include "gromacs/mdtypes/md_enums.h"
82 #include "gromacs/nbnxm/gpu_data_mgmt.h"
83 #include "gromacs/nbnxm/nbnxm.h"
84 #include "gromacs/nbnxm/nbnxm_geometry.h"
85 #include "gromacs/pbcutil/ishift.h"
86 #include "gromacs/pbcutil/pbc.h"
87 #include "gromacs/tables/forcetable.h"
88 #include "gromacs/topology/mtop_util.h"
89 #include "gromacs/trajectory/trajectoryframe.h"
90 #include "gromacs/utility/cstringutil.h"
91 #include "gromacs/utility/exceptions.h"
92 #include "gromacs/utility/fatalerror.h"
93 #include "gromacs/utility/gmxassert.h"
94 #include "gromacs/utility/logger.h"
95 #include "gromacs/utility/physicalnodecommunicator.h"
96 #include "gromacs/utility/pleasecite.h"
97 #include "gromacs/utility/smalloc.h"
98 #include "gromacs/utility/strconvert.h"
99
100 static real *mk_nbfp(const gmx_ffparams_t *idef, gmx_bool bBHAM)
101 {
102     real *nbfp;
103     int   i, j, k, atnr;
104
105     atnr = idef->atnr;
106     if (bBHAM)
107     {
108         snew(nbfp, 3*atnr*atnr);
109         for (i = k = 0; (i < atnr); i++)
110         {
111             for (j = 0; (j < atnr); j++, k++)
112             {
113                 BHAMA(nbfp, atnr, i, j) = idef->iparams[k].bham.a;
114                 BHAMB(nbfp, atnr, i, j) = idef->iparams[k].bham.b;
115                 /* nbfp now includes the 6.0 derivative prefactor */
116                 BHAMC(nbfp, atnr, i, j) = idef->iparams[k].bham.c*6.0;
117             }
118         }
119     }
120     else
121     {
122         snew(nbfp, 2*atnr*atnr);
123         for (i = k = 0; (i < atnr); i++)
124         {
125             for (j = 0; (j < atnr); j++, k++)
126             {
127                 /* nbfp now includes the 6.0/12.0 derivative prefactors */
128                 C6(nbfp, atnr, i, j)   = idef->iparams[k].lj.c6*6.0;
129                 C12(nbfp, atnr, i, j)  = idef->iparams[k].lj.c12*12.0;
130             }
131         }
132     }
133
134     return nbfp;
135 }
136
137 static real *make_ljpme_c6grid(const gmx_ffparams_t *idef, t_forcerec *fr)
138 {
139     int        i, j, k, atnr;
140     real       c6, c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
141     real      *grid;
142
143     /* For LJ-PME simulations, we correct the energies with the reciprocal space
144      * inside of the cut-off. To do this the non-bonded kernels needs to have
145      * access to the C6-values used on the reciprocal grid in pme.c
146      */
147
148     atnr = idef->atnr;
149     snew(grid, 2*atnr*atnr);
150     for (i = k = 0; (i < atnr); i++)
151     {
152         for (j = 0; (j < atnr); j++, k++)
153         {
154             c6i  = idef->iparams[i*(atnr+1)].lj.c6;
155             c12i = idef->iparams[i*(atnr+1)].lj.c12;
156             c6j  = idef->iparams[j*(atnr+1)].lj.c6;
157             c12j = idef->iparams[j*(atnr+1)].lj.c12;
158             c6   = std::sqrt(c6i * c6j);
159             if (fr->ljpme_combination_rule == eljpmeLB
160                 && !gmx_numzero(c6) && !gmx_numzero(c12i) && !gmx_numzero(c12j))
161             {
162                 sigmai = gmx::sixthroot(c12i / c6i);
163                 sigmaj = gmx::sixthroot(c12j / c6j);
164                 epsi   = c6i * c6i / c12i;
165                 epsj   = c6j * c6j / c12j;
166                 c6     = std::sqrt(epsi * epsj) * gmx::power6(0.5*(sigmai+sigmaj));
167             }
168             /* Store the elements at the same relative positions as C6 in nbfp in order
169              * to simplify access in the kernels
170              */
171             grid[2*(atnr*i+j)] = c6*6.0;
172         }
173     }
174     return grid;
175 }
176
177 /* This routine sets fr->solvent_opt to the most common solvent in the
178  * system, e.g. esolSPC or esolTIP4P. It will also mark each charge group in
179  * the fr->solvent_type array with the correct type (or esolNO).
180  *
181  * Charge groups that fulfill the conditions but are not identical to the
182  * most common one will be marked as esolNO in the solvent_type array.
183  *
184  * TIP3p is identical to SPC for these purposes, so we call it
185  * SPC in the arrays (Apologies to Bill Jorgensen ;-)
186  *
187  * NOTE: QM particle should not
188  * become an optimized solvent. Not even if there is only one charge
189  * group in the Qm
190  */
191
192 typedef struct
193 {
194     int    model;
195     int    count;
196     int    vdwtype[4];
197     real   charge[4];
198 } solvent_parameters_t;
199
200 static void
201 check_solvent_cg(const gmx_moltype_t    *molt,
202                  int                     cg0,
203                  int                     nmol,
204                  const unsigned char    *qm_grpnr,
205                  const t_grps           *qm_grps,
206                  t_forcerec   *          fr,
207                  int                    *n_solvent_parameters,
208                  solvent_parameters_t  **solvent_parameters_p,
209                  int                     cginfo,
210                  int                    *cg_sp)
211 {
212     t_atom               *atom;
213     int                   j, k;
214     int                   j0, j1, nj;
215     gmx_bool              perturbed;
216     gmx_bool              has_vdw[4];
217     gmx_bool              match;
218     real                  tmp_charge[4]  = { 0.0 }; /* init to zero to make gcc 7 happy */
219     int                   tmp_vdwtype[4] = { 0 };   /* init to zero to make gcc 7 happy */
220     int                   tjA;
221     gmx_bool              qm;
222     solvent_parameters_t *solvent_parameters;
223
224     /* We use a list with parameters for each solvent type.
225      * Every time we discover a new molecule that fulfills the basic
226      * conditions for a solvent we compare with the previous entries
227      * in these lists. If the parameters are the same we just increment
228      * the counter for that type, and otherwise we create a new type
229      * based on the current molecule.
230      *
231      * Once we've finished going through all molecules we check which
232      * solvent is most common, and mark all those molecules while we
233      * clear the flag on all others.
234      */
235
236     solvent_parameters = *solvent_parameters_p;
237
238     /* Mark the cg first as non optimized */
239     *cg_sp = -1;
240
241     /* Check if this cg has no exclusions with atoms in other charge groups
242      * and all atoms inside the charge group excluded.
243      * We only have 3 or 4 atom solvent loops.
244      */
245     if (GET_CGINFO_EXCL_INTER(cginfo) ||
246         !GET_CGINFO_EXCL_INTRA(cginfo))
247     {
248         return;
249     }
250
251     /* Get the indices of the first atom in this charge group */
252     j0     = molt->cgs.index[cg0];
253     j1     = molt->cgs.index[cg0+1];
254
255     /* Number of atoms in our molecule */
256     nj     = j1 - j0;
257
258     if (debug)
259     {
260         fprintf(debug,
261                 "Moltype '%s': there are %d atoms in this charge group\n",
262                 *molt->name, nj);
263     }
264
265     /* Check if it could be an SPC (3 atoms) or TIP4p (4) water,
266      * otherwise skip it.
267      */
268     if (nj < 3 || nj > 4)
269     {
270         return;
271     }
272
273     /* Check if we are doing QM on this group */
274     qm = FALSE;
275     if (qm_grpnr != nullptr)
276     {
277         for (j = j0; j < j1 && !qm; j++)
278         {
279             qm = (qm_grpnr[j] < qm_grps->nr - 1);
280         }
281     }
282     /* Cannot use solvent optimization with QM */
283     if (qm)
284     {
285         return;
286     }
287
288     atom = molt->atoms.atom;
289
290     /* Still looks like a solvent, time to check parameters */
291
292     /* If it is perturbed (free energy) we can't use the solvent loops,
293      * so then we just skip to the next molecule.
294      */
295     perturbed = FALSE;
296
297     for (j = j0; j < j1 && !perturbed; j++)
298     {
299         perturbed = PERTURBED(atom[j]);
300     }
301
302     if (perturbed)
303     {
304         return;
305     }
306
307     /* Now it's only a question if the VdW and charge parameters
308      * are OK. Before doing the check we compare and see if they are
309      * identical to a possible previous solvent type.
310      * First we assign the current types and charges.
311      */
312     for (j = 0; j < nj; j++)
313     {
314         tmp_vdwtype[j] = atom[j0+j].type;
315         tmp_charge[j]  = atom[j0+j].q;
316     }
317
318     /* Does it match any previous solvent type? */
319     for (k = 0; k < *n_solvent_parameters; k++)
320     {
321         match = TRUE;
322
323
324         /* We can only match SPC with 3 atoms and TIP4p with 4 atoms */
325         if ( (solvent_parameters[k].model == esolSPC   && nj != 3)  ||
326              (solvent_parameters[k].model == esolTIP4P && nj != 4) )
327         {
328             match = FALSE;
329         }
330
331         /* Check that types & charges match for all atoms in molecule */
332         for (j = 0; j < nj && match; j++)
333         {
334             if (tmp_vdwtype[j] != solvent_parameters[k].vdwtype[j])
335             {
336                 match = FALSE;
337             }
338             if (tmp_charge[j] != solvent_parameters[k].charge[j])
339             {
340                 match = FALSE;
341             }
342         }
343         if (match)
344         {
345             /* Congratulations! We have a matched solvent.
346              * Flag it with this type for later processing.
347              */
348             *cg_sp = k;
349             solvent_parameters[k].count += nmol;
350
351             /* We are done with this charge group */
352             return;
353         }
354     }
355
356     /* If we get here, we have a tentative new solvent type.
357      * Before we add it we must check that it fulfills the requirements
358      * of the solvent optimized loops. First determine which atoms have
359      * VdW interactions.
360      */
361     for (j = 0; j < nj; j++)
362     {
363         has_vdw[j] = FALSE;
364         tjA        = tmp_vdwtype[j];
365
366         /* Go through all other tpes and see if any have non-zero
367          * VdW parameters when combined with this one.
368          */
369         for (k = 0; k < fr->ntype && (!has_vdw[j]); k++)
370         {
371             /* We already checked that the atoms weren't perturbed,
372              * so we only need to check state A now.
373              */
374             if (fr->bBHAM)
375             {
376                 has_vdw[j] = (has_vdw[j] ||
377                               (BHAMA(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
378                               (BHAMB(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
379                               (BHAMC(fr->nbfp, fr->ntype, tjA, k) != 0.0));
380             }
381             else
382             {
383                 /* Standard LJ */
384                 has_vdw[j] = (has_vdw[j] ||
385                               (C6(fr->nbfp, fr->ntype, tjA, k)  != 0.0) ||
386                               (C12(fr->nbfp, fr->ntype, tjA, k) != 0.0));
387             }
388         }
389     }
390
391     /* Now we know all we need to make the final check and assignment. */
392     if (nj == 3)
393     {
394         /* So, is it an SPC?
395          * For this we require thatn all atoms have charge,
396          * the charges on atom 2 & 3 should be the same, and only
397          * atom 1 might have VdW.
398          */
399         if (!has_vdw[1] &&
400             !has_vdw[2] &&
401             tmp_charge[0]  != 0 &&
402             tmp_charge[1]  != 0 &&
403             tmp_charge[2]  == tmp_charge[1])
404         {
405             srenew(solvent_parameters, *n_solvent_parameters+1);
406             solvent_parameters[*n_solvent_parameters].model = esolSPC;
407             solvent_parameters[*n_solvent_parameters].count = nmol;
408             for (k = 0; k < 3; k++)
409             {
410                 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
411                 solvent_parameters[*n_solvent_parameters].charge[k]  = tmp_charge[k];
412             }
413
414             *cg_sp = *n_solvent_parameters;
415             (*n_solvent_parameters)++;
416         }
417     }
418     else if (nj == 4)
419     {
420         /* Or could it be a TIP4P?
421          * For this we require thatn atoms 2,3,4 have charge, but not atom 1.
422          * Only atom 1 mght have VdW.
423          */
424         if (!has_vdw[1] &&
425             !has_vdw[2] &&
426             !has_vdw[3] &&
427             tmp_charge[0]  == 0 &&
428             tmp_charge[1]  != 0 &&
429             tmp_charge[2]  == tmp_charge[1] &&
430             tmp_charge[3]  != 0)
431         {
432             srenew(solvent_parameters, *n_solvent_parameters+1);
433             solvent_parameters[*n_solvent_parameters].model = esolTIP4P;
434             solvent_parameters[*n_solvent_parameters].count = nmol;
435             for (k = 0; k < 4; k++)
436             {
437                 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
438                 solvent_parameters[*n_solvent_parameters].charge[k]  = tmp_charge[k];
439             }
440
441             *cg_sp = *n_solvent_parameters;
442             (*n_solvent_parameters)++;
443         }
444     }
445
446     *solvent_parameters_p = solvent_parameters;
447 }
448
449 static void
450 check_solvent(FILE  *                fp,
451               const gmx_mtop_t  *    mtop,
452               t_forcerec  *          fr,
453               cginfo_mb_t           *cginfo_mb)
454 {
455     const t_block     *   cgs;
456     const gmx_moltype_t  *molt;
457     int                   mol, cg_mol, at_offset, am, cgm, i, nmol_ch, nmol;
458     int                   n_solvent_parameters;
459     solvent_parameters_t *solvent_parameters;
460     int                 **cg_sp;
461     int                   bestsp, bestsol;
462
463     if (debug)
464     {
465         fprintf(debug, "Going to determine what solvent types we have.\n");
466     }
467
468     n_solvent_parameters = 0;
469     solvent_parameters   = nullptr;
470     /* Allocate temporary array for solvent type */
471     snew(cg_sp, mtop->molblock.size());
472
473     at_offset = 0;
474     for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
475     {
476         molt = &mtop->moltype[mtop->molblock[mb].type];
477         cgs  = &molt->cgs;
478         /* Here we have to loop over all individual molecules
479          * because we need to check for QMMM particles.
480          */
481         snew(cg_sp[mb], cginfo_mb[mb].cg_mod);
482         nmol_ch = cginfo_mb[mb].cg_mod/cgs->nr;
483         nmol    = mtop->molblock[mb].nmol/nmol_ch;
484         for (mol = 0; mol < nmol_ch; mol++)
485         {
486             cgm = mol*cgs->nr;
487             am  = mol*cgs->index[cgs->nr];
488             for (cg_mol = 0; cg_mol < cgs->nr; cg_mol++)
489             {
490                 check_solvent_cg(molt, cg_mol, nmol,
491                                  mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics].empty() ?
492                                  nullptr : mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics].data()+at_offset+am,
493                                  &mtop->groups.groups[SimulationAtomGroupType::QuantumMechanics],
494                                  fr,
495                                  &n_solvent_parameters, &solvent_parameters,
496                                  cginfo_mb[mb].cginfo[cgm+cg_mol],
497                                  &cg_sp[mb][cgm+cg_mol]);
498             }
499         }
500         at_offset += cgs->index[cgs->nr];
501     }
502
503     /* Puh! We finished going through all charge groups.
504      * Now find the most common solvent model.
505      */
506
507     /* Most common solvent this far */
508     bestsp = -2;
509     for (i = 0; i < n_solvent_parameters; i++)
510     {
511         if (bestsp == -2 ||
512             solvent_parameters[i].count > solvent_parameters[bestsp].count)
513         {
514             bestsp = i;
515         }
516     }
517
518     if (bestsp >= 0)
519     {
520         bestsol = solvent_parameters[bestsp].model;
521     }
522     else
523     {
524         bestsol = esolNO;
525     }
526
527     fr->nWatMol = 0;
528     for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
529     {
530         cgs  = &mtop->moltype[mtop->molblock[mb].type].cgs;
531         nmol = (mtop->molblock[mb].nmol*cgs->nr)/cginfo_mb[mb].cg_mod;
532         for (i = 0; i < cginfo_mb[mb].cg_mod; i++)
533         {
534             if (cg_sp[mb][i] == bestsp)
535             {
536                 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i], bestsol);
537                 fr->nWatMol += nmol;
538             }
539             else
540             {
541                 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i], esolNO);
542             }
543         }
544         sfree(cg_sp[mb]);
545     }
546     sfree(cg_sp);
547
548     if (bestsol != esolNO && fp != nullptr)
549     {
550         fprintf(fp, "\nEnabling %s-like water optimization for %d molecules.\n\n",
551                 esol_names[bestsol],
552                 solvent_parameters[bestsp].count);
553     }
554
555     sfree(solvent_parameters);
556     fr->solvent_opt = bestsol;
557 }
558
559 enum {
560     acNONE = 0, acCONSTRAINT, acSETTLE
561 };
562
563 static cginfo_mb_t *init_cginfo_mb(FILE *fplog, const gmx_mtop_t *mtop,
564                                    t_forcerec *fr, gmx_bool bNoSolvOpt,
565                                    gmx_bool *bFEP_NonBonded,
566                                    gmx_bool *bExcl_IntraCGAll_InterCGNone)
567 {
568     const t_block        *cgs;
569     const t_blocka       *excl;
570     const gmx_moltype_t  *molt;
571     const gmx_molblock_t *molb;
572     cginfo_mb_t          *cginfo_mb;
573     gmx_bool             *type_VDW;
574     int                  *cginfo;
575     int                   cg_offset, a_offset;
576     int                   m, cg, a0, a1, gid, ai, j, aj, excl_nalloc;
577     int                  *a_con;
578     int                   ftype;
579     int                   ia;
580     gmx_bool              bId, *bExcl, bExclIntraAll, bExclInter, bHaveVDW, bHaveQ, bHavePerturbedAtoms;
581
582     snew(cginfo_mb, mtop->molblock.size());
583
584     snew(type_VDW, fr->ntype);
585     for (ai = 0; ai < fr->ntype; ai++)
586     {
587         type_VDW[ai] = FALSE;
588         for (j = 0; j < fr->ntype; j++)
589         {
590             type_VDW[ai] = type_VDW[ai] ||
591                 fr->bBHAM ||
592                 C6(fr->nbfp, fr->ntype, ai, j) != 0 ||
593                 C12(fr->nbfp, fr->ntype, ai, j) != 0;
594         }
595     }
596
597     *bFEP_NonBonded               = FALSE;
598     *bExcl_IntraCGAll_InterCGNone = TRUE;
599
600     excl_nalloc = 10;
601     snew(bExcl, excl_nalloc);
602     cg_offset = 0;
603     a_offset  = 0;
604     for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
605     {
606         molb = &mtop->molblock[mb];
607         molt = &mtop->moltype[molb->type];
608         cgs  = &molt->cgs;
609         excl = &molt->excls;
610
611         /* Check if the cginfo is identical for all molecules in this block.
612          * If so, we only need an array of the size of one molecule.
613          * Otherwise we make an array of #mol times #cgs per molecule.
614          */
615         bId = TRUE;
616         for (m = 0; m < molb->nmol; m++)
617         {
618             int am = m*cgs->index[cgs->nr];
619             for (cg = 0; cg < cgs->nr; cg++)
620             {
621                 a0 = cgs->index[cg];
622                 a1 = cgs->index[cg+1];
623                 if (getGroupType(mtop->groups, SimulationAtomGroupType::QuantumMechanics, a_offset+am+a0) !=
624                     getGroupType(mtop->groups, SimulationAtomGroupType::QuantumMechanics, a_offset   +a0))
625                 {
626                     bId = FALSE;
627                 }
628                 if (!mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics].empty())
629                 {
630                     for (ai = a0; ai < a1; ai++)
631                     {
632                         if (mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics][a_offset+am+ai] !=
633                             mtop->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics][a_offset   +ai])
634                         {
635                             bId = FALSE;
636                         }
637                     }
638                 }
639             }
640         }
641
642         cginfo_mb[mb].cg_start = cg_offset;
643         cginfo_mb[mb].cg_end   = cg_offset + molb->nmol*cgs->nr;
644         cginfo_mb[mb].cg_mod   = (bId ? 1 : molb->nmol)*cgs->nr;
645         snew(cginfo_mb[mb].cginfo, cginfo_mb[mb].cg_mod);
646         cginfo = cginfo_mb[mb].cginfo;
647
648         /* Set constraints flags for constrained atoms */
649         snew(a_con, molt->atoms.nr);
650         for (ftype = 0; ftype < F_NRE; ftype++)
651         {
652             if (interaction_function[ftype].flags & IF_CONSTRAINT)
653             {
654                 int nral;
655
656                 nral = NRAL(ftype);
657                 for (ia = 0; ia < molt->ilist[ftype].size(); ia += 1+nral)
658                 {
659                     int a;
660
661                     for (a = 0; a < nral; a++)
662                     {
663                         a_con[molt->ilist[ftype].iatoms[ia+1+a]] =
664                             (ftype == F_SETTLE ? acSETTLE : acCONSTRAINT);
665                     }
666                 }
667             }
668         }
669
670         for (m = 0; m < (bId ? 1 : molb->nmol); m++)
671         {
672             int cgm = m*cgs->nr;
673             int am  = m*cgs->index[cgs->nr];
674             for (cg = 0; cg < cgs->nr; cg++)
675             {
676                 a0 = cgs->index[cg];
677                 a1 = cgs->index[cg+1];
678
679                 /* Store the energy group in cginfo */
680                 gid = getGroupType(mtop->groups, SimulationAtomGroupType::EnergyOutput, a_offset+am+a0);
681                 SET_CGINFO_GID(cginfo[cgm+cg], gid);
682
683                 /* Check the intra/inter charge group exclusions */
684                 if (a1-a0 > excl_nalloc)
685                 {
686                     excl_nalloc = a1 - a0;
687                     srenew(bExcl, excl_nalloc);
688                 }
689                 /* bExclIntraAll: all intra cg interactions excluded
690                  * bExclInter:    any inter cg interactions excluded
691                  */
692                 bExclIntraAll       = TRUE;
693                 bExclInter          = FALSE;
694                 bHaveVDW            = FALSE;
695                 bHaveQ              = FALSE;
696                 bHavePerturbedAtoms = FALSE;
697                 for (ai = a0; ai < a1; ai++)
698                 {
699                     /* Check VDW and electrostatic interactions */
700                     bHaveVDW = bHaveVDW || (type_VDW[molt->atoms.atom[ai].type] ||
701                                             type_VDW[molt->atoms.atom[ai].typeB]);
702                     bHaveQ  = bHaveQ    || (molt->atoms.atom[ai].q != 0 ||
703                                             molt->atoms.atom[ai].qB != 0);
704
705                     bHavePerturbedAtoms = bHavePerturbedAtoms || (PERTURBED(molt->atoms.atom[ai]) != 0);
706
707                     /* Clear the exclusion list for atom ai */
708                     for (aj = a0; aj < a1; aj++)
709                     {
710                         bExcl[aj-a0] = FALSE;
711                     }
712                     /* Loop over all the exclusions of atom ai */
713                     for (j = excl->index[ai]; j < excl->index[ai+1]; j++)
714                     {
715                         aj = excl->a[j];
716                         if (aj < a0 || aj >= a1)
717                         {
718                             bExclInter = TRUE;
719                         }
720                         else
721                         {
722                             bExcl[aj-a0] = TRUE;
723                         }
724                     }
725                     /* Check if ai excludes a0 to a1 */
726                     for (aj = a0; aj < a1; aj++)
727                     {
728                         if (!bExcl[aj-a0])
729                         {
730                             bExclIntraAll = FALSE;
731                         }
732                     }
733
734                     switch (a_con[ai])
735                     {
736                         case acCONSTRAINT:
737                             SET_CGINFO_CONSTR(cginfo[cgm+cg]);
738                             break;
739                         case acSETTLE:
740                             SET_CGINFO_SETTLE(cginfo[cgm+cg]);
741                             break;
742                         default:
743                             break;
744                     }
745                 }
746                 if (bExclIntraAll)
747                 {
748                     SET_CGINFO_EXCL_INTRA(cginfo[cgm+cg]);
749                 }
750                 if (bExclInter)
751                 {
752                     SET_CGINFO_EXCL_INTER(cginfo[cgm+cg]);
753                 }
754                 if (a1 - a0 > MAX_CHARGEGROUP_SIZE)
755                 {
756                     /* The size in cginfo is currently only read with DD */
757                     gmx_fatal(FARGS, "A charge group has size %d which is larger than the limit of %d atoms", a1-a0, MAX_CHARGEGROUP_SIZE);
758                 }
759                 if (bHaveVDW)
760                 {
761                     SET_CGINFO_HAS_VDW(cginfo[cgm+cg]);
762                 }
763                 if (bHaveQ)
764                 {
765                     SET_CGINFO_HAS_Q(cginfo[cgm+cg]);
766                 }
767                 if (bHavePerturbedAtoms && fr->efep != efepNO)
768                 {
769                     SET_CGINFO_FEP(cginfo[cgm+cg]);
770                     *bFEP_NonBonded = TRUE;
771                 }
772                 /* Store the charge group size */
773                 SET_CGINFO_NATOMS(cginfo[cgm+cg], a1-a0);
774
775                 if (!bExclIntraAll || bExclInter)
776                 {
777                     *bExcl_IntraCGAll_InterCGNone = FALSE;
778                 }
779             }
780         }
781
782         sfree(a_con);
783
784         cg_offset += molb->nmol*cgs->nr;
785         a_offset  += molb->nmol*cgs->index[cgs->nr];
786     }
787     sfree(type_VDW);
788     sfree(bExcl);
789
790     /* the solvent optimizer is called after the QM is initialized,
791      * because we don't want to have the QM subsystemto become an
792      * optimized solvent
793      */
794
795     check_solvent(fplog, mtop, fr, cginfo_mb);
796
797     if (getenv("GMX_NO_SOLV_OPT"))
798     {
799         if (fplog)
800         {
801             fprintf(fplog, "Found environment variable GMX_NO_SOLV_OPT.\n"
802                     "Disabling all solvent optimization\n");
803         }
804         fr->solvent_opt = esolNO;
805     }
806     if (bNoSolvOpt)
807     {
808         fr->solvent_opt = esolNO;
809     }
810     if (!fr->solvent_opt)
811     {
812         for (size_t mb = 0; mb < mtop->molblock.size(); mb++)
813         {
814             for (cg = 0; cg < cginfo_mb[mb].cg_mod; cg++)
815             {
816                 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[cg], esolNO);
817             }
818         }
819     }
820
821     return cginfo_mb;
822 }
823
824 static std::vector<int> cginfo_expand(const int          nmb,
825                                       const cginfo_mb_t *cgi_mb)
826 {
827     const int        ncg = cgi_mb[nmb - 1].cg_end;
828
829     std::vector<int> cginfo(ncg);
830
831     int              mb = 0;
832     for (int cg = 0; cg < ncg; cg++)
833     {
834         while (cg >= cgi_mb[mb].cg_end)
835         {
836             mb++;
837         }
838         cginfo[cg] =
839             cgi_mb[mb].cginfo[(cg - cgi_mb[mb].cg_start) % cgi_mb[mb].cg_mod];
840     }
841
842     return cginfo;
843 }
844
845 static void done_cginfo_mb(cginfo_mb_t *cginfo_mb, int numMolBlocks)
846 {
847     if (cginfo_mb == nullptr)
848     {
849         return;
850     }
851     for (int mb = 0; mb < numMolBlocks; ++mb)
852     {
853         sfree(cginfo_mb[mb].cginfo);
854     }
855     sfree(cginfo_mb);
856 }
857
858 /* Sets the sum of charges (squared) and C6 in the system in fr.
859  * Returns whether the system has a net charge.
860  */
861 static bool set_chargesum(FILE *log, t_forcerec *fr, const gmx_mtop_t *mtop)
862 {
863     /*This now calculates sum for q and c6*/
864     double         qsum, q2sum, q, c6sum, c6;
865
866     qsum   = 0;
867     q2sum  = 0;
868     c6sum  = 0;
869     for (const gmx_molblock_t &molb : mtop->molblock)
870     {
871         int            nmol  = molb.nmol;
872         const t_atoms *atoms = &mtop->moltype[molb.type].atoms;
873         for (int i = 0; i < atoms->nr; i++)
874         {
875             q       = atoms->atom[i].q;
876             qsum   += nmol*q;
877             q2sum  += nmol*q*q;
878             c6      = mtop->ffparams.iparams[atoms->atom[i].type*(mtop->ffparams.atnr+1)].lj.c6;
879             c6sum  += nmol*c6;
880         }
881     }
882     fr->qsum[0]   = qsum;
883     fr->q2sum[0]  = q2sum;
884     fr->c6sum[0]  = c6sum;
885
886     if (fr->efep != efepNO)
887     {
888         qsum   = 0;
889         q2sum  = 0;
890         c6sum  = 0;
891         for (const gmx_molblock_t &molb : mtop->molblock)
892         {
893             int            nmol  = molb.nmol;
894             const t_atoms *atoms = &mtop->moltype[molb.type].atoms;
895             for (int i = 0; i < atoms->nr; i++)
896             {
897                 q       = atoms->atom[i].qB;
898                 qsum   += nmol*q;
899                 q2sum  += nmol*q*q;
900                 c6      = mtop->ffparams.iparams[atoms->atom[i].typeB*(mtop->ffparams.atnr+1)].lj.c6;
901                 c6sum  += nmol*c6;
902             }
903             fr->qsum[1]   = qsum;
904             fr->q2sum[1]  = q2sum;
905             fr->c6sum[1]  = c6sum;
906         }
907     }
908     else
909     {
910         fr->qsum[1]   = fr->qsum[0];
911         fr->q2sum[1]  = fr->q2sum[0];
912         fr->c6sum[1]  = fr->c6sum[0];
913     }
914     if (log)
915     {
916         if (fr->efep == efepNO)
917         {
918             fprintf(log, "System total charge: %.3f\n", fr->qsum[0]);
919         }
920         else
921         {
922             fprintf(log, "System total charge, top. A: %.3f top. B: %.3f\n",
923                     fr->qsum[0], fr->qsum[1]);
924         }
925     }
926
927     /* A cut-off of 1e-4 is used to catch rounding errors due to ascii input */
928     return (std::abs(fr->qsum[0]) > 1e-4 ||
929             std::abs(fr->qsum[1]) > 1e-4);
930 }
931
932 void update_forcerec(t_forcerec *fr, matrix box)
933 {
934     if (fr->ic->eeltype == eelGRF)
935     {
936         calc_rffac(nullptr, fr->ic->eeltype, fr->ic->epsilon_r, fr->ic->epsilon_rf,
937                    fr->ic->rcoulomb, fr->temp, fr->zsquare, box,
938                    &fr->ic->k_rf, &fr->ic->c_rf);
939     }
940 }
941
942 static real calcBuckinghamBMax(FILE *fplog, const gmx_mtop_t *mtop)
943 {
944     const t_atoms *at1, *at2;
945     int            i, j, tpi, tpj, ntypes;
946     real           b, bmin;
947
948     if (fplog)
949     {
950         fprintf(fplog, "Determining largest Buckingham b parameter for table\n");
951     }
952     ntypes = mtop->ffparams.atnr;
953
954     bmin            = -1;
955     real bham_b_max = 0;
956     for (size_t mt1 = 0; mt1 < mtop->moltype.size(); mt1++)
957     {
958         at1 = &mtop->moltype[mt1].atoms;
959         for (i = 0; (i < at1->nr); i++)
960         {
961             tpi = at1->atom[i].type;
962             if (tpi >= ntypes)
963             {
964                 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", i, tpi, ntypes);
965             }
966
967             for (size_t mt2 = mt1; mt2 < mtop->moltype.size(); mt2++)
968             {
969                 at2 = &mtop->moltype[mt2].atoms;
970                 for (j = 0; (j < at2->nr); j++)
971                 {
972                     tpj = at2->atom[j].type;
973                     if (tpj >= ntypes)
974                     {
975                         gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", j, tpj, ntypes);
976                     }
977                     b = mtop->ffparams.iparams[tpi*ntypes + tpj].bham.b;
978                     if (b > bham_b_max)
979                     {
980                         bham_b_max = b;
981                     }
982                     if ((b < bmin) || (bmin == -1))
983                     {
984                         bmin = b;
985                     }
986                 }
987             }
988         }
989     }
990     if (fplog)
991     {
992         fprintf(fplog, "Buckingham b parameters, min: %g, max: %g\n",
993                 bmin, bham_b_max);
994     }
995
996     return bham_b_max;
997 }
998
999 static void make_nbf_tables(FILE *fp,
1000                             const interaction_const_t *ic, real rtab,
1001                             const char *tabfn, char *eg1, char *eg2,
1002                             t_nblists *nbl)
1003 {
1004     char buf[STRLEN];
1005     int  i, j;
1006
1007     if (tabfn == nullptr)
1008     {
1009         if (debug)
1010         {
1011             fprintf(debug, "No table file name passed, can not read table, can not do non-bonded interactions\n");
1012         }
1013         return;
1014     }
1015
1016     sprintf(buf, "%s", tabfn);
1017     if (eg1 && eg2)
1018     {
1019         /* Append the two energy group names */
1020         sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "_%s_%s.%s",
1021                 eg1, eg2, ftp2ext(efXVG));
1022     }
1023     nbl->table_elec_vdw = make_tables(fp, ic, buf, rtab, 0);
1024     /* Copy the contents of the table to separate coulomb and LJ tables too,
1025      * to improve cache performance.
1026      */
1027     /* For performance reasons we want
1028      * the table data to be aligned to 16-byte.
1029      */
1030     nbl->table_elec                =
1031         new t_forcetable(GMX_TABLE_INTERACTION_ELEC,
1032                          nbl->table_elec_vdw->format);
1033     nbl->table_elec->r             = nbl->table_elec_vdw->r;
1034     nbl->table_elec->n             = nbl->table_elec_vdw->n;
1035     nbl->table_elec->scale         = nbl->table_elec_vdw->scale;
1036     nbl->table_elec->formatsize    = nbl->table_elec_vdw->formatsize;
1037     nbl->table_elec->ninteractions = 1;
1038     nbl->table_elec->stride        = nbl->table_elec->formatsize * nbl->table_elec->ninteractions;
1039     snew_aligned(nbl->table_elec->data, nbl->table_elec->stride*(nbl->table_elec->n+1), 32);
1040
1041     nbl->table_vdw                =
1042         new t_forcetable(GMX_TABLE_INTERACTION_VDWREP_VDWDISP,
1043                          nbl->table_elec_vdw->format);
1044     nbl->table_vdw->r             = nbl->table_elec_vdw->r;
1045     nbl->table_vdw->n             = nbl->table_elec_vdw->n;
1046     nbl->table_vdw->scale         = nbl->table_elec_vdw->scale;
1047     nbl->table_vdw->formatsize    = nbl->table_elec_vdw->formatsize;
1048     nbl->table_vdw->ninteractions = 2;
1049     nbl->table_vdw->stride        = nbl->table_vdw->formatsize * nbl->table_vdw->ninteractions;
1050     snew_aligned(nbl->table_vdw->data, nbl->table_vdw->stride*(nbl->table_vdw->n+1), 32);
1051
1052     /* NOTE: Using a single i-loop here leads to mix-up of data in table_vdw
1053      *       with (at least) gcc 6.2, 6.3 and 6.4 when compiled with -O3 and AVX
1054      */
1055     for (i = 0; i <= nbl->table_elec_vdw->n; i++)
1056     {
1057         for (j = 0; j < 4; j++)
1058         {
1059             nbl->table_elec->data[4*i+j] = nbl->table_elec_vdw->data[12*i+j];
1060         }
1061     }
1062     for (i = 0; i <= nbl->table_elec_vdw->n; i++)
1063     {
1064         for (j = 0; j < 8; j++)
1065         {
1066             nbl->table_vdw->data[8*i+j] = nbl->table_elec_vdw->data[12*i+4+j];
1067         }
1068     }
1069 }
1070
1071 /*!\brief If there's bonded interactions of type \c ftype1 or \c
1072  * ftype2 present in the topology, build an array of the number of
1073  * interactions present for each bonded interaction index found in the
1074  * topology.
1075  *
1076  * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
1077  * valid type with that parameter.
1078  *
1079  * \c count will be reallocated as necessary to fit the largest bonded
1080  * interaction index found, and its current size will be returned in
1081  * \c ncount. It will contain zero for every bonded interaction index
1082  * for which no interactions are present in the topology.
1083  */
1084 static void count_tables(int ftype1, int ftype2, const gmx_mtop_t *mtop,
1085                          int *ncount, int **count)
1086 {
1087     int                  ftype, i, j, tabnr;
1088
1089     // Loop over all moleculetypes
1090     for (const gmx_moltype_t &molt : mtop->moltype)
1091     {
1092         // Loop over all interaction types
1093         for (ftype = 0; ftype < F_NRE; ftype++)
1094         {
1095             // If the current interaction type is one of the types whose tables we're trying to count...
1096             if (ftype == ftype1 || ftype == ftype2)
1097             {
1098                 const InteractionList &il     = molt.ilist[ftype];
1099                 const int              stride = 1 + NRAL(ftype);
1100                 // ... and there are actually some interactions for this type
1101                 for (i = 0; i < il.size(); i += stride)
1102                 {
1103                     // Find out which table index the user wanted
1104                     tabnr = mtop->ffparams.iparams[il.iatoms[i]].tab.table;
1105                     if (tabnr < 0)
1106                     {
1107                         gmx_fatal(FARGS, "A bonded table number is smaller than 0: %d\n", tabnr);
1108                     }
1109                     // Make room for this index in the data structure
1110                     if (tabnr >= *ncount)
1111                     {
1112                         srenew(*count, tabnr+1);
1113                         for (j = *ncount; j < tabnr+1; j++)
1114                         {
1115                             (*count)[j] = 0;
1116                         }
1117                         *ncount = tabnr+1;
1118                     }
1119                     // Record that this table index is used and must have a valid file
1120                     (*count)[tabnr]++;
1121                 }
1122             }
1123         }
1124     }
1125 }
1126
1127 /*!\brief If there's bonded interactions of flavour \c tabext and type
1128  * \c ftype1 or \c ftype2 present in the topology, seek them in the
1129  * list of filenames passed to mdrun, and make bonded tables from
1130  * those files.
1131  *
1132  * \c ftype1 or \c ftype2 may be set to -1 to disable seeking for a
1133  * valid type with that parameter.
1134  *
1135  * A fatal error occurs if no matching filename is found.
1136  */
1137 static bondedtable_t *make_bonded_tables(FILE *fplog,
1138                                          int ftype1, int ftype2,
1139                                          const gmx_mtop_t *mtop,
1140                                          gmx::ArrayRef<const std::string> tabbfnm,
1141                                          const char *tabext)
1142 {
1143     int            ncount, *count;
1144     bondedtable_t *tab;
1145
1146     tab = nullptr;
1147
1148     ncount = 0;
1149     count  = nullptr;
1150     count_tables(ftype1, ftype2, mtop, &ncount, &count);
1151
1152     // Are there any relevant tabulated bond interactions?
1153     if (ncount > 0)
1154     {
1155         snew(tab, ncount);
1156         for (int i = 0; i < ncount; i++)
1157         {
1158             // Do any interactions exist that requires this table?
1159             if (count[i] > 0)
1160             {
1161                 // This pattern enforces the current requirement that
1162                 // table filenames end in a characteristic sequence
1163                 // before the file type extension, and avoids table 13
1164                 // being recognized and used for table 1.
1165                 std::string patternToFind = gmx::formatString("_%s%d.%s", tabext, i, ftp2ext(efXVG));
1166                 bool        madeTable     = false;
1167                 for (gmx::index j = 0; j < tabbfnm.ssize() && !madeTable; ++j)
1168                 {
1169                     if (gmx::endsWith(tabbfnm[j], patternToFind))
1170                     {
1171                         // Finally read the table from the file found
1172                         tab[i]    = make_bonded_table(fplog, tabbfnm[j].c_str(), NRAL(ftype1)-2);
1173                         madeTable = true;
1174                     }
1175                 }
1176                 if (!madeTable)
1177                 {
1178                     bool isPlural = (ftype2 != -1);
1179                     gmx_fatal(FARGS, "Tabulated interaction of type '%s%s%s' with index %d cannot be used because no table file whose name matched '%s' was passed via the gmx mdrun -tableb command-line option.",
1180                               interaction_function[ftype1].longname,
1181                               isPlural ? "' or '" : "",
1182                               isPlural ? interaction_function[ftype2].longname : "",
1183                               i,
1184                               patternToFind.c_str());
1185                 }
1186             }
1187         }
1188         sfree(count);
1189     }
1190
1191     return tab;
1192 }
1193
1194 void forcerec_set_ranges(t_forcerec *fr,
1195                          int ncg_home, int ncg_force,
1196                          int natoms_force,
1197                          int natoms_force_constr, int natoms_f_novirsum)
1198 {
1199     fr->cg0 = 0;
1200     fr->hcg = ncg_home;
1201
1202     /* fr->ncg_force is unused in the standard code,
1203      * but it can be useful for modified code dealing with charge groups.
1204      */
1205     fr->ncg_force           = ncg_force;
1206     fr->natoms_force        = natoms_force;
1207     fr->natoms_force_constr = natoms_force_constr;
1208
1209     if (fr->natoms_force_constr > fr->nalloc_force)
1210     {
1211         fr->nalloc_force = over_alloc_dd(fr->natoms_force_constr);
1212     }
1213
1214     if (fr->haveDirectVirialContributions)
1215     {
1216         fr->forceBufferForDirectVirialContributions->resize(natoms_f_novirsum);
1217     }
1218 }
1219
1220 static real cutoff_inf(real cutoff)
1221 {
1222     if (cutoff == 0)
1223     {
1224         cutoff = GMX_CUTOFF_INF;
1225     }
1226
1227     return cutoff;
1228 }
1229
1230 /*! \brief Print Coulomb Ewald citations and set ewald coefficients */
1231 static void initCoulombEwaldParameters(FILE *fp, const t_inputrec *ir,
1232                                        bool systemHasNetCharge,
1233                                        interaction_const_t *ic)
1234 {
1235     if (!EEL_PME_EWALD(ir->coulombtype))
1236     {
1237         return;
1238     }
1239
1240     if (fp)
1241     {
1242         fprintf(fp, "Will do PME sum in reciprocal space for electrostatic interactions.\n");
1243
1244         if (ir->coulombtype == eelP3M_AD)
1245         {
1246             please_cite(fp, "Hockney1988");
1247             please_cite(fp, "Ballenegger2012");
1248         }
1249         else
1250         {
1251             please_cite(fp, "Essmann95a");
1252         }
1253
1254         if (ir->ewald_geometry == eewg3DC)
1255         {
1256             if (fp)
1257             {
1258                 fprintf(fp, "Using the Ewald3DC correction for systems with a slab geometry%s.\n",
1259                         systemHasNetCharge ? " and net charge" : "");
1260             }
1261             please_cite(fp, "In-Chul99a");
1262             if (systemHasNetCharge)
1263             {
1264                 please_cite(fp, "Ballenegger2009");
1265             }
1266         }
1267     }
1268
1269     ic->ewaldcoeff_q = calc_ewaldcoeff_q(ir->rcoulomb, ir->ewald_rtol);
1270     if (fp)
1271     {
1272         fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for Ewald\n",
1273                 1/ic->ewaldcoeff_q);
1274     }
1275
1276     if (ic->coulomb_modifier == eintmodPOTSHIFT)
1277     {
1278         GMX_RELEASE_ASSERT(ic->rcoulomb != 0, "Cutoff radius cannot be zero");
1279         ic->sh_ewald = std::erfc(ic->ewaldcoeff_q*ic->rcoulomb) / ic->rcoulomb;
1280     }
1281     else
1282     {
1283         ic->sh_ewald = 0;
1284     }
1285 }
1286
1287 /*! \brief Print Van der Waals Ewald citations and set ewald coefficients */
1288 static void initVdwEwaldParameters(FILE *fp, const t_inputrec *ir,
1289                                    interaction_const_t *ic)
1290 {
1291     if (!EVDW_PME(ir->vdwtype))
1292     {
1293         return;
1294     }
1295
1296     if (fp)
1297     {
1298         fprintf(fp, "Will do PME sum in reciprocal space for LJ dispersion interactions.\n");
1299         please_cite(fp, "Essmann95a");
1300     }
1301     ic->ewaldcoeff_lj = calc_ewaldcoeff_lj(ir->rvdw, ir->ewald_rtol_lj);
1302     if (fp)
1303     {
1304         fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for LJ Ewald\n",
1305                 1/ic->ewaldcoeff_lj);
1306     }
1307
1308     if (ic->vdw_modifier == eintmodPOTSHIFT)
1309     {
1310         real crc2       = gmx::square(ic->ewaldcoeff_lj*ic->rvdw);
1311         ic->sh_lj_ewald = (std::exp(-crc2)*(1 + crc2 + 0.5*crc2*crc2) - 1)/gmx::power6(ic->rvdw);
1312     }
1313     else
1314     {
1315         ic->sh_lj_ewald = 0;
1316     }
1317 }
1318
1319 gmx_bool uses_simple_tables(int                       cutoff_scheme,
1320                             const nonbonded_verlet_t *nbv)
1321 {
1322     gmx_bool bUsesSimpleTables = TRUE;
1323
1324     switch (cutoff_scheme)
1325     {
1326         case ecutsGROUP:
1327             bUsesSimpleTables = TRUE;
1328             break;
1329         case ecutsVERLET:
1330             GMX_RELEASE_ASSERT(nullptr != nbv, "A non-bonded verlet object is required with the Verlet cutoff-scheme");
1331             bUsesSimpleTables = nbv->pairlistIsSimple();
1332             break;
1333         default:
1334             gmx_incons("unimplemented");
1335     }
1336     return bUsesSimpleTables;
1337 }
1338
1339 static void init_ewald_f_table(interaction_const_t *ic,
1340                                real                 rtab)
1341 {
1342     real maxr;
1343
1344     /* Get the Ewald table spacing based on Coulomb and/or LJ
1345      * Ewald coefficients and rtol.
1346      */
1347     ic->tabq_scale = ewald_spline3_table_scale(ic);
1348
1349     if (ic->cutoff_scheme == ecutsVERLET)
1350     {
1351         maxr = ic->rcoulomb;
1352     }
1353     else
1354     {
1355         maxr = std::max(ic->rcoulomb, rtab);
1356     }
1357     ic->tabq_size  = static_cast<int>(maxr*ic->tabq_scale) + 2;
1358
1359     sfree_aligned(ic->tabq_coul_FDV0);
1360     sfree_aligned(ic->tabq_coul_F);
1361     sfree_aligned(ic->tabq_coul_V);
1362
1363     sfree_aligned(ic->tabq_vdw_FDV0);
1364     sfree_aligned(ic->tabq_vdw_F);
1365     sfree_aligned(ic->tabq_vdw_V);
1366
1367     if (EEL_PME_EWALD(ic->eeltype))
1368     {
1369         /* Create the original table data in FDV0 */
1370         snew_aligned(ic->tabq_coul_FDV0, ic->tabq_size*4, 32);
1371         snew_aligned(ic->tabq_coul_F, ic->tabq_size, 32);
1372         snew_aligned(ic->tabq_coul_V, ic->tabq_size, 32);
1373         table_spline3_fill_ewald_lr(ic->tabq_coul_F, ic->tabq_coul_V, ic->tabq_coul_FDV0,
1374                                     ic->tabq_size, 1/ic->tabq_scale, ic->ewaldcoeff_q, v_q_ewald_lr);
1375     }
1376
1377     if (EVDW_PME(ic->vdwtype))
1378     {
1379         snew_aligned(ic->tabq_vdw_FDV0, ic->tabq_size*4, 32);
1380         snew_aligned(ic->tabq_vdw_F, ic->tabq_size, 32);
1381         snew_aligned(ic->tabq_vdw_V, ic->tabq_size, 32);
1382         table_spline3_fill_ewald_lr(ic->tabq_vdw_F, ic->tabq_vdw_V, ic->tabq_vdw_FDV0,
1383                                     ic->tabq_size, 1/ic->tabq_scale, ic->ewaldcoeff_lj, v_lj_ewald_lr);
1384     }
1385 }
1386
1387 void init_interaction_const_tables(FILE                *fp,
1388                                    interaction_const_t *ic,
1389                                    real                 rtab)
1390 {
1391     if (EEL_PME_EWALD(ic->eeltype) || EVDW_PME(ic->vdwtype))
1392     {
1393         init_ewald_f_table(ic, rtab);
1394
1395         if (fp != nullptr)
1396         {
1397             fprintf(fp, "Initialized non-bonded Ewald correction tables, spacing: %.2e size: %d\n\n",
1398                     1/ic->tabq_scale, ic->tabq_size);
1399         }
1400     }
1401 }
1402
1403 static void clear_force_switch_constants(shift_consts_t *sc)
1404 {
1405     sc->c2   = 0;
1406     sc->c3   = 0;
1407     sc->cpot = 0;
1408 }
1409
1410 static void force_switch_constants(real p,
1411                                    real rsw, real rc,
1412                                    shift_consts_t *sc)
1413 {
1414     /* Here we determine the coefficient for shifting the force to zero
1415      * between distance rsw and the cut-off rc.
1416      * For a potential of r^-p, we have force p*r^-(p+1).
1417      * But to save flops we absorb p in the coefficient.
1418      * Thus we get:
1419      * force/p   = r^-(p+1) + c2*r^2 + c3*r^3
1420      * potential = r^-p + c2/3*r^3 + c3/4*r^4 + cpot
1421      */
1422     sc->c2   =  ((p + 1)*rsw - (p + 4)*rc)/(pow(rc, p + 2)*gmx::square(rc - rsw));
1423     sc->c3   = -((p + 1)*rsw - (p + 3)*rc)/(pow(rc, p + 2)*gmx::power3(rc - rsw));
1424     sc->cpot = -pow(rc, -p) + p*sc->c2/3*gmx::power3(rc - rsw) + p*sc->c3/4*gmx::power4(rc - rsw);
1425 }
1426
1427 static void potential_switch_constants(real rsw, real rc,
1428                                        switch_consts_t *sc)
1429 {
1430     /* The switch function is 1 at rsw and 0 at rc.
1431      * The derivative and second derivate are zero at both ends.
1432      * rsw        = max(r - r_switch, 0)
1433      * sw         = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
1434      * dsw        = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
1435      * force      = force*dsw - potential*sw
1436      * potential *= sw
1437      */
1438     sc->c3 = -10/gmx::power3(rc - rsw);
1439     sc->c4 =  15/gmx::power4(rc - rsw);
1440     sc->c5 =  -6/gmx::power5(rc - rsw);
1441 }
1442
1443 /*! \brief Construct interaction constants
1444  *
1445  * This data is used (particularly) by search and force code for
1446  * short-range interactions. Many of these are constant for the whole
1447  * simulation; some are constant only after PME tuning completes.
1448  */
1449 static void
1450 init_interaction_const(FILE                       *fp,
1451                        interaction_const_t       **interaction_const,
1452                        const t_inputrec           *ir,
1453                        const gmx_mtop_t           *mtop,
1454                        bool                        systemHasNetCharge)
1455 {
1456     interaction_const_t *ic;
1457
1458     snew(ic, 1);
1459
1460     ic->cutoff_scheme   = ir->cutoff_scheme;
1461
1462     /* Just allocate something so we can free it */
1463     snew_aligned(ic->tabq_coul_FDV0, 16, 32);
1464     snew_aligned(ic->tabq_coul_F, 16, 32);
1465     snew_aligned(ic->tabq_coul_V, 16, 32);
1466
1467     /* Lennard-Jones */
1468     ic->vdwtype         = ir->vdwtype;
1469     ic->vdw_modifier    = ir->vdw_modifier;
1470     ic->reppow          = mtop->ffparams.reppow;
1471     ic->rvdw            = cutoff_inf(ir->rvdw);
1472     ic->rvdw_switch     = ir->rvdw_switch;
1473     ic->ljpme_comb_rule = ir->ljpme_combination_rule;
1474     ic->useBuckingham   = (mtop->ffparams.functype[0] == F_BHAM);
1475     if (ic->useBuckingham)
1476     {
1477         ic->buckinghamBMax = calcBuckinghamBMax(fp, mtop);
1478     }
1479
1480     initVdwEwaldParameters(fp, ir, ic);
1481
1482     clear_force_switch_constants(&ic->dispersion_shift);
1483     clear_force_switch_constants(&ic->repulsion_shift);
1484
1485     switch (ic->vdw_modifier)
1486     {
1487         case eintmodPOTSHIFT:
1488             /* Only shift the potential, don't touch the force */
1489             ic->dispersion_shift.cpot = -1.0/gmx::power6(ic->rvdw);
1490             ic->repulsion_shift.cpot  = -1.0/gmx::power12(ic->rvdw);
1491             break;
1492         case eintmodFORCESWITCH:
1493             /* Switch the force, switch and shift the potential */
1494             force_switch_constants(6.0, ic->rvdw_switch, ic->rvdw,
1495                                    &ic->dispersion_shift);
1496             force_switch_constants(12.0, ic->rvdw_switch, ic->rvdw,
1497                                    &ic->repulsion_shift);
1498             break;
1499         case eintmodPOTSWITCH:
1500             /* Switch the potential and force */
1501             potential_switch_constants(ic->rvdw_switch, ic->rvdw,
1502                                        &ic->vdw_switch);
1503             break;
1504         case eintmodNONE:
1505         case eintmodEXACTCUTOFF:
1506             /* Nothing to do here */
1507             break;
1508         default:
1509             gmx_incons("unimplemented potential modifier");
1510     }
1511
1512     ic->sh_invrc6 = -ic->dispersion_shift.cpot;
1513
1514     /* Electrostatics */
1515     ic->eeltype          = ir->coulombtype;
1516     ic->coulomb_modifier = ir->coulomb_modifier;
1517     ic->rcoulomb         = cutoff_inf(ir->rcoulomb);
1518     ic->rcoulomb_switch  = ir->rcoulomb_switch;
1519     ic->epsilon_r        = ir->epsilon_r;
1520
1521     /* Set the Coulomb energy conversion factor */
1522     if (ic->epsilon_r != 0)
1523     {
1524         ic->epsfac = ONE_4PI_EPS0/ic->epsilon_r;
1525     }
1526     else
1527     {
1528         /* eps = 0 is infinite dieletric: no Coulomb interactions */
1529         ic->epsfac = 0;
1530     }
1531
1532     /* Reaction-field */
1533     if (EEL_RF(ic->eeltype))
1534     {
1535         ic->epsilon_rf = ir->epsilon_rf;
1536         /* Generalized reaction field parameters are updated every step */
1537         if (ic->eeltype != eelGRF)
1538         {
1539             calc_rffac(fp, ic->eeltype, ic->epsilon_r, ic->epsilon_rf,
1540                        ic->rcoulomb, 0, 0, nullptr,
1541                        &ic->k_rf, &ic->c_rf);
1542         }
1543
1544         if (ir->cutoff_scheme == ecutsGROUP && ic->eeltype == eelRF_ZERO)
1545         {
1546             /* grompp should have done this, but this scheme is obsolete */
1547             ic->coulomb_modifier = eintmodEXACTCUTOFF;
1548         }
1549     }
1550     else
1551     {
1552         /* For plain cut-off we might use the reaction-field kernels */
1553         ic->epsilon_rf = ic->epsilon_r;
1554         ic->k_rf       = 0;
1555         if (ir->coulomb_modifier == eintmodPOTSHIFT)
1556         {
1557             ic->c_rf   = 1/ic->rcoulomb;
1558         }
1559         else
1560         {
1561             ic->c_rf   = 0;
1562         }
1563     }
1564
1565     initCoulombEwaldParameters(fp, ir, systemHasNetCharge, ic);
1566
1567     if (fp != nullptr)
1568     {
1569         real dispersion_shift;
1570
1571         dispersion_shift = ic->dispersion_shift.cpot;
1572         if (EVDW_PME(ic->vdwtype))
1573         {
1574             dispersion_shift -= ic->sh_lj_ewald;
1575         }
1576         fprintf(fp, "Potential shift: LJ r^-12: %.3e r^-6: %.3e",
1577                 ic->repulsion_shift.cpot, dispersion_shift);
1578
1579         if (ic->eeltype == eelCUT)
1580         {
1581             fprintf(fp, ", Coulomb %.e", -ic->c_rf);
1582         }
1583         else if (EEL_PME(ic->eeltype))
1584         {
1585             fprintf(fp, ", Ewald %.3e", -ic->sh_ewald);
1586         }
1587         fprintf(fp, "\n");
1588     }
1589
1590     *interaction_const = ic;
1591 }
1592
1593 static void
1594 done_interaction_const(interaction_const_t *interaction_const)
1595 {
1596     sfree_aligned(interaction_const->tabq_coul_FDV0);
1597     sfree_aligned(interaction_const->tabq_coul_F);
1598     sfree_aligned(interaction_const->tabq_coul_V);
1599     sfree(interaction_const);
1600 }
1601
1602 void init_forcerec(FILE                             *fp,
1603                    const gmx::MDLogger              &mdlog,
1604                    t_forcerec                       *fr,
1605                    t_fcdata                         *fcd,
1606                    const t_inputrec                 *ir,
1607                    const gmx_mtop_t                 *mtop,
1608                    const t_commrec                  *cr,
1609                    matrix                            box,
1610                    const char                       *tabfn,
1611                    const char                       *tabpfn,
1612                    gmx::ArrayRef<const std::string>  tabbfnm,
1613                    const gmx_hw_info_t              &hardwareInfo,
1614                    const gmx_device_info_t          *deviceInfo,
1615                    const bool                        useGpuForBonded,
1616                    gmx_bool                          bNoSolvOpt,
1617                    real                              print_force)
1618 {
1619     int            m, negp_pp, negptable, egi, egj;
1620     real           rtab;
1621     char          *env;
1622     double         dbl;
1623     const t_block *cgs;
1624     gmx_bool       bGenericKernelOnly;
1625     gmx_bool       needGroupSchemeTables, bSomeNormalNbListsAreInUse;
1626     gmx_bool       bFEP_NonBonded;
1627     int           *nm_ind, egp_flags;
1628
1629
1630     /* By default we turn SIMD kernels on, but it might be turned off further down... */
1631     fr->use_simd_kernels = TRUE;
1632
1633     if (check_box(ir->ePBC, box))
1634     {
1635         gmx_fatal(FARGS, "%s", check_box(ir->ePBC, box));
1636     }
1637
1638     /* Test particle insertion ? */
1639     if (EI_TPI(ir->eI))
1640     {
1641         /* Set to the size of the molecule to be inserted (the last one) */
1642         /* Because of old style topologies, we have to use the last cg
1643          * instead of the last molecule type.
1644          */
1645         cgs       = &mtop->moltype[mtop->molblock.back().type].cgs;
1646         fr->n_tpi = cgs->index[cgs->nr] - cgs->index[cgs->nr-1];
1647         gmx::RangePartitioning molecules = gmx_mtop_molecules(*mtop);
1648         if (fr->n_tpi != molecules.block(molecules.numBlocks() - 1).size())
1649         {
1650             gmx_fatal(FARGS, "The molecule to insert can not consist of multiple charge groups.\nMake it a single charge group.");
1651         }
1652     }
1653     else
1654     {
1655         fr->n_tpi = 0;
1656     }
1657
1658     if (ir->coulombtype == eelRF_NEC_UNSUPPORTED)
1659     {
1660         gmx_fatal(FARGS, "%s electrostatics is no longer supported",
1661                   eel_names[ir->coulombtype]);
1662     }
1663
1664     if (ir->bAdress)
1665     {
1666         gmx_fatal(FARGS, "AdResS simulations are no longer supported");
1667     }
1668     if (ir->useTwinRange)
1669     {
1670         gmx_fatal(FARGS, "Twin-range simulations are no longer supported");
1671     }
1672     /* Copy the user determined parameters */
1673     fr->userint1  = ir->userint1;
1674     fr->userint2  = ir->userint2;
1675     fr->userint3  = ir->userint3;
1676     fr->userint4  = ir->userint4;
1677     fr->userreal1 = ir->userreal1;
1678     fr->userreal2 = ir->userreal2;
1679     fr->userreal3 = ir->userreal3;
1680     fr->userreal4 = ir->userreal4;
1681
1682     /* Shell stuff */
1683     fr->fc_stepsize = ir->fc_stepsize;
1684
1685     /* Free energy */
1686     fr->efep        = ir->efep;
1687     fr->sc_alphavdw = ir->fepvals->sc_alpha;
1688     if (ir->fepvals->bScCoul)
1689     {
1690         fr->sc_alphacoul  = ir->fepvals->sc_alpha;
1691         fr->sc_sigma6_min = gmx::power6(ir->fepvals->sc_sigma_min);
1692     }
1693     else
1694     {
1695         fr->sc_alphacoul  = 0;
1696         fr->sc_sigma6_min = 0; /* only needed when bScCoul is on */
1697     }
1698     fr->sc_power      = ir->fepvals->sc_power;
1699     fr->sc_r_power    = ir->fepvals->sc_r_power;
1700     fr->sc_sigma6_def = gmx::power6(ir->fepvals->sc_sigma);
1701
1702     env = getenv("GMX_SCSIGMA_MIN");
1703     if (env != nullptr)
1704     {
1705         dbl = 0;
1706         sscanf(env, "%20lf", &dbl);
1707         fr->sc_sigma6_min = gmx::power6(dbl);
1708         if (fp)
1709         {
1710             fprintf(fp, "Setting the minimum soft core sigma to %g nm\n", dbl);
1711         }
1712     }
1713
1714     fr->bNonbonded = TRUE;
1715     if (getenv("GMX_NO_NONBONDED") != nullptr)
1716     {
1717         /* turn off non-bonded calculations */
1718         fr->bNonbonded = FALSE;
1719         GMX_LOG(mdlog.warning).asParagraph().appendText(
1720                 "Found environment variable GMX_NO_NONBONDED.\n"
1721                 "Disabling nonbonded calculations.");
1722     }
1723
1724     bGenericKernelOnly = FALSE;
1725
1726     /* We now check in the NS code whether a particular combination of interactions
1727      * can be used with water optimization, and disable it if that is not the case.
1728      */
1729
1730     if (getenv("GMX_NB_GENERIC") != nullptr)
1731     {
1732         if (fp != nullptr)
1733         {
1734             fprintf(fp,
1735                     "Found environment variable GMX_NB_GENERIC.\n"
1736                     "Disabling all interaction-specific nonbonded kernels, will only\n"
1737                     "use the slow generic ones in src/gmxlib/nonbonded/nb_generic.c\n\n");
1738         }
1739         bGenericKernelOnly = TRUE;
1740     }
1741
1742     if (bGenericKernelOnly)
1743     {
1744         bNoSolvOpt         = TRUE;
1745     }
1746
1747     if ( (getenv("GMX_DISABLE_SIMD_KERNELS") != nullptr) || (getenv("GMX_NOOPTIMIZEDKERNELS") != nullptr) )
1748     {
1749         fr->use_simd_kernels = FALSE;
1750         if (fp != nullptr)
1751         {
1752             fprintf(fp,
1753                     "\nFound environment variable GMX_DISABLE_SIMD_KERNELS.\n"
1754                     "Disabling the usage of any SIMD-specific non-bonded & bonded kernel routines\n"
1755                     "(e.g. SSE2/SSE4.1/AVX).\n\n");
1756         }
1757     }
1758
1759     fr->bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
1760
1761     /* Neighbour searching stuff */
1762     fr->cutoff_scheme = ir->cutoff_scheme;
1763     fr->bGrid         = (ir->ns_type == ensGRID);
1764     fr->ePBC          = ir->ePBC;
1765
1766     if (fr->cutoff_scheme == ecutsGROUP)
1767     {
1768         const char *note = "NOTE: This file uses the deprecated 'group' cutoff_scheme. This will be\n"
1769             "removed in a future release when 'verlet' supports all interaction forms.\n";
1770
1771         if (MASTER(cr))
1772         {
1773             fprintf(stderr, "\n%s\n", note);
1774         }
1775         if (fp != nullptr)
1776         {
1777             fprintf(fp, "\n%s\n", note);
1778         }
1779     }
1780
1781     /* Determine if we will do PBC for distances in bonded interactions */
1782     if (fr->ePBC == epbcNONE)
1783     {
1784         fr->bMolPBC = FALSE;
1785     }
1786     else
1787     {
1788         if (!DOMAINDECOMP(cr))
1789         {
1790             gmx_bool bSHAKE;
1791
1792             bSHAKE = (ir->eConstrAlg == econtSHAKE &&
1793                       (gmx_mtop_ftype_count(mtop, F_CONSTR) > 0 ||
1794                        gmx_mtop_ftype_count(mtop, F_CONSTRNC) > 0));
1795
1796             /* The group cut-off scheme and SHAKE assume charge groups
1797              * are whole, but not using molpbc is faster in most cases.
1798              * With intermolecular interactions we need PBC for calculating
1799              * distances between atoms in different molecules.
1800              */
1801             if ((fr->cutoff_scheme == ecutsGROUP || bSHAKE) &&
1802                 !mtop->bIntermolecularInteractions)
1803             {
1804                 fr->bMolPBC = ir->bPeriodicMols;
1805
1806                 if (bSHAKE && fr->bMolPBC)
1807                 {
1808                     gmx_fatal(FARGS, "SHAKE is not supported with periodic molecules");
1809                 }
1810             }
1811             else
1812             {
1813                 /* Not making molecules whole is faster in most cases,
1814                  * but With orientation restraints we need whole molecules.
1815                  */
1816                 fr->bMolPBC = (fcd->orires.nr == 0);
1817
1818                 if (getenv("GMX_USE_GRAPH") != nullptr)
1819                 {
1820                     fr->bMolPBC = FALSE;
1821                     if (fp)
1822                     {
1823                         GMX_LOG(mdlog.warning).asParagraph().appendText("GMX_USE_GRAPH is set, using the graph for bonded interactions");
1824                     }
1825
1826                     if (mtop->bIntermolecularInteractions)
1827                     {
1828                         GMX_LOG(mdlog.warning).asParagraph().appendText("WARNING: Molecules linked by intermolecular interactions have to reside in the same periodic image, otherwise artifacts will occur!");
1829                     }
1830                 }
1831
1832                 GMX_RELEASE_ASSERT(fr->bMolPBC || !mtop->bIntermolecularInteractions, "We need to use PBC within molecules with inter-molecular interactions");
1833
1834                 if (bSHAKE && fr->bMolPBC)
1835                 {
1836                     gmx_fatal(FARGS, "SHAKE is not properly supported with intermolecular interactions. For short simulations where linked molecules remain in the same periodic image, the environment variable GMX_USE_GRAPH can be used to override this check.\n");
1837                 }
1838             }
1839         }
1840         else
1841         {
1842             fr->bMolPBC = dd_bonded_molpbc(cr->dd, fr->ePBC);
1843         }
1844     }
1845
1846     fr->rc_scaling = ir->refcoord_scaling;
1847     copy_rvec(ir->posres_com, fr->posres_com);
1848     copy_rvec(ir->posres_comB, fr->posres_comB);
1849     fr->rlist                    = cutoff_inf(ir->rlist);
1850     fr->ljpme_combination_rule   = ir->ljpme_combination_rule;
1851
1852     /* This now calculates sum for q and c6*/
1853     bool systemHasNetCharge = set_chargesum(fp, fr, mtop);
1854
1855     /* fr->ic is used both by verlet and group kernels (to some extent) now */
1856     init_interaction_const(fp, &fr->ic, ir, mtop, systemHasNetCharge);
1857     init_interaction_const_tables(fp, fr->ic, ir->rlist + ir->tabext);
1858
1859     const interaction_const_t *ic = fr->ic;
1860
1861     /* TODO: Replace this Ewald table or move it into interaction_const_t */
1862     if (ir->coulombtype == eelEWALD)
1863     {
1864         init_ewald_tab(&(fr->ewald_table), ir, fp);
1865     }
1866
1867     /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
1868     switch (ic->eeltype)
1869     {
1870         case eelCUT:
1871             fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_COULOMB;
1872             break;
1873
1874         case eelRF:
1875         case eelGRF:
1876             fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
1877             break;
1878
1879         case eelRF_ZERO:
1880             fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
1881             GMX_RELEASE_ASSERT(ic->coulomb_modifier == eintmodEXACTCUTOFF, "With the group scheme RF-zero needs the exact cut-off modifier");
1882             break;
1883
1884         case eelSWITCH:
1885         case eelSHIFT:
1886         case eelUSER:
1887         case eelENCADSHIFT:
1888         case eelPMESWITCH:
1889         case eelPMEUSER:
1890         case eelPMEUSERSWITCH:
1891             fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
1892             break;
1893
1894         case eelPME:
1895         case eelP3M_AD:
1896         case eelEWALD:
1897             fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_EWALD;
1898             break;
1899
1900         default:
1901             gmx_fatal(FARGS, "Unsupported electrostatic interaction: %s", eel_names[ic->eeltype]);
1902     }
1903     fr->nbkernel_elec_modifier = ic->coulomb_modifier;
1904
1905     /* Vdw: Translate from mdp settings to kernel format */
1906     switch (ic->vdwtype)
1907     {
1908         case evdwCUT:
1909             if (fr->bBHAM)
1910             {
1911                 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_BUCKINGHAM;
1912             }
1913             else
1914             {
1915                 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LENNARDJONES;
1916             }
1917             break;
1918         case evdwPME:
1919             fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LJEWALD;
1920             break;
1921
1922         case evdwSWITCH:
1923         case evdwSHIFT:
1924         case evdwUSER:
1925         case evdwENCADSHIFT:
1926             fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
1927             break;
1928
1929         default:
1930             gmx_fatal(FARGS, "Unsupported vdw interaction: %s", evdw_names[ic->vdwtype]);
1931     }
1932     fr->nbkernel_vdw_modifier = ic->vdw_modifier;
1933
1934     if (ir->cutoff_scheme == ecutsGROUP)
1935     {
1936         fr->bvdwtab    = ((ic->vdwtype != evdwCUT || !gmx_within_tol(ic->reppow, 12.0, 10*GMX_DOUBLE_EPS))
1937                           && !EVDW_PME(ic->vdwtype));
1938         /* We have special kernels for standard Ewald and PME, but the pme-switch ones are tabulated above */
1939         fr->bcoultab   = !(ic->eeltype == eelCUT ||
1940                            ic->eeltype == eelEWALD ||
1941                            ic->eeltype == eelPME ||
1942                            ic->eeltype == eelP3M_AD ||
1943                            ic->eeltype == eelRF ||
1944                            ic->eeltype == eelRF_ZERO);
1945
1946         /* If the user absolutely wants different switch/shift settings for coul/vdw, it is likely
1947          * going to be faster to tabulate the interaction than calling the generic kernel.
1948          * However, if generic kernels have been requested we keep things analytically.
1949          */
1950         if (fr->nbkernel_elec_modifier == eintmodPOTSWITCH &&
1951             fr->nbkernel_vdw_modifier == eintmodPOTSWITCH &&
1952             !bGenericKernelOnly)
1953         {
1954             if ((ic->rcoulomb_switch != ic->rvdw_switch) || (ic->rcoulomb != ic->rvdw))
1955             {
1956                 fr->bcoultab = TRUE;
1957                 /* Once we tabulate electrostatics, we can use the switch function for LJ,
1958                  * which would otherwise need two tables.
1959                  */
1960             }
1961         }
1962         else if ((fr->nbkernel_elec_modifier == eintmodPOTSHIFT && fr->nbkernel_vdw_modifier == eintmodPOTSHIFT) ||
1963                  ((fr->nbkernel_elec_interaction == GMX_NBKERNEL_ELEC_REACTIONFIELD &&
1964                    fr->nbkernel_elec_modifier == eintmodEXACTCUTOFF &&
1965                    (fr->nbkernel_vdw_modifier == eintmodPOTSWITCH || fr->nbkernel_vdw_modifier == eintmodPOTSHIFT))))
1966         {
1967             if ((ic->rcoulomb != ic->rvdw) && (!bGenericKernelOnly))
1968             {
1969                 fr->bcoultab = TRUE;
1970             }
1971         }
1972
1973         if (fr->nbkernel_elec_modifier == eintmodFORCESWITCH)
1974         {
1975             fr->bcoultab = TRUE;
1976         }
1977         if (fr->nbkernel_vdw_modifier == eintmodFORCESWITCH)
1978         {
1979             fr->bvdwtab = TRUE;
1980         }
1981
1982         if (getenv("GMX_REQUIRE_TABLES"))
1983         {
1984             fr->bvdwtab  = TRUE;
1985             fr->bcoultab = TRUE;
1986         }
1987
1988         if (fp)
1989         {
1990             fprintf(fp, "Table routines are used for coulomb: %s\n",
1991                     gmx::boolToString(fr->bcoultab));
1992             fprintf(fp, "Table routines are used for vdw:     %s\n",
1993                     gmx::boolToString(fr->bvdwtab));
1994         }
1995
1996         if (fr->bvdwtab)
1997         {
1998             fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
1999             fr->nbkernel_vdw_modifier    = eintmodNONE;
2000         }
2001         if (fr->bcoultab)
2002         {
2003             fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
2004             fr->nbkernel_elec_modifier    = eintmodNONE;
2005         }
2006     }
2007
2008     if (ir->cutoff_scheme == ecutsVERLET)
2009     {
2010         if (!gmx_within_tol(ic->reppow, 12.0, 10*GMX_DOUBLE_EPS))
2011         {
2012             gmx_fatal(FARGS, "Cut-off scheme %s only supports LJ repulsion power 12", ecutscheme_names[ir->cutoff_scheme]);
2013         }
2014         /* Older tpr files can contain Coulomb user tables with the Verlet cutoff-scheme,
2015          * while mdrun does not (and never did) support this.
2016          */
2017         if (EEL_USER(fr->ic->eeltype))
2018         {
2019             gmx_fatal(FARGS, "Combination of %s and cutoff scheme %s is not supported",
2020                       eel_names[ir->coulombtype], ecutscheme_names[ir->cutoff_scheme]);
2021         }
2022
2023         fr->bvdwtab  = FALSE;
2024         fr->bcoultab = FALSE;
2025     }
2026
2027     /* 1-4 interaction electrostatics */
2028     fr->fudgeQQ = mtop->ffparams.fudgeQQ;
2029
2030     /* Parameters for generalized RF */
2031     fr->zsquare = 0.0;
2032     fr->temp    = 0.0;
2033
2034     if (ic->eeltype == eelGRF)
2035     {
2036         init_generalized_rf(fp, mtop, ir, fr);
2037     }
2038
2039     fr->haveDirectVirialContributions =
2040         (EEL_FULL(ic->eeltype) || EVDW_PME(ic->vdwtype) ||
2041          fr->forceProviders->hasForceProvider() ||
2042          gmx_mtop_ftype_count(mtop, F_POSRES) > 0 ||
2043          gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0 ||
2044          ir->nwall > 0 ||
2045          ir->bPull ||
2046          ir->bRot ||
2047          ir->bIMD);
2048
2049     if (fr->haveDirectVirialContributions)
2050     {
2051         fr->forceBufferForDirectVirialContributions = new std::vector<gmx::RVec>;
2052     }
2053
2054     if (fr->shift_vec == nullptr)
2055     {
2056         snew(fr->shift_vec, SHIFTS);
2057     }
2058
2059     if (fr->fshift == nullptr)
2060     {
2061         snew(fr->fshift, SHIFTS);
2062     }
2063
2064     if (fr->nbfp == nullptr)
2065     {
2066         fr->ntype = mtop->ffparams.atnr;
2067         fr->nbfp  = mk_nbfp(&mtop->ffparams, fr->bBHAM);
2068         if (EVDW_PME(ic->vdwtype))
2069         {
2070             fr->ljpme_c6grid  = make_ljpme_c6grid(&mtop->ffparams, fr);
2071         }
2072     }
2073
2074     /* Copy the energy group exclusions */
2075     fr->egp_flags = ir->opts.egp_flags;
2076
2077     /* Van der Waals stuff */
2078     if ((ic->vdwtype != evdwCUT) && (ic->vdwtype != evdwUSER) && !fr->bBHAM)
2079     {
2080         if (ic->rvdw_switch >= ic->rvdw)
2081         {
2082             gmx_fatal(FARGS, "rvdw_switch (%f) must be < rvdw (%f)",
2083                       ic->rvdw_switch, ic->rvdw);
2084         }
2085         if (fp)
2086         {
2087             fprintf(fp, "Using %s Lennard-Jones, switch between %g and %g nm\n",
2088                     (ic->eeltype == eelSWITCH) ? "switched" : "shifted",
2089                     ic->rvdw_switch, ic->rvdw);
2090         }
2091     }
2092
2093     if (fr->bBHAM && EVDW_PME(ic->vdwtype))
2094     {
2095         gmx_fatal(FARGS, "LJ PME not supported with Buckingham");
2096     }
2097
2098     if (fr->bBHAM && (ic->vdwtype == evdwSHIFT || ic->vdwtype == evdwSWITCH))
2099     {
2100         gmx_fatal(FARGS, "Switch/shift interaction not supported with Buckingham");
2101     }
2102
2103     if (fr->bBHAM && fr->cutoff_scheme == ecutsVERLET)
2104     {
2105         gmx_fatal(FARGS, "Verlet cutoff-scheme is not supported with Buckingham");
2106     }
2107
2108     if (fp && fr->cutoff_scheme == ecutsGROUP)
2109     {
2110         fprintf(fp, "Cut-off's:   NS: %g   Coulomb: %g   %s: %g\n",
2111                 fr->rlist, ic->rcoulomb, fr->bBHAM ? "BHAM" : "LJ", ic->rvdw);
2112     }
2113
2114     if (ir->implicit_solvent)
2115     {
2116         gmx_fatal(FARGS, "Implict solvation is no longer supported.");
2117     }
2118
2119     /* Construct tables for the group scheme. A little unnecessary to
2120      * make both vdw and coul tables sometimes, but what the
2121      * heck. Note that both cutoff schemes construct Ewald tables in
2122      * init_interaction_const_tables. */
2123     needGroupSchemeTables = (ir->cutoff_scheme == ecutsGROUP &&
2124                              (fr->bcoultab || fr->bvdwtab));
2125
2126     negp_pp   = ir->opts.ngener - ir->nwall;
2127     negptable = 0;
2128     if (!needGroupSchemeTables)
2129     {
2130         bSomeNormalNbListsAreInUse = TRUE;
2131     }
2132     else
2133     {
2134         bSomeNormalNbListsAreInUse = FALSE;
2135         for (egi = 0; egi < negp_pp; egi++)
2136         {
2137             for (egj = egi; egj < negp_pp; egj++)
2138             {
2139                 egp_flags = ir->opts.egp_flags[GID(egi, egj, ir->opts.ngener)];
2140                 if (!(egp_flags & EGP_EXCL))
2141                 {
2142                     if (egp_flags & EGP_TABLE)
2143                     {
2144                         negptable++;
2145                     }
2146                     else
2147                     {
2148                         bSomeNormalNbListsAreInUse = TRUE;
2149                     }
2150                 }
2151             }
2152         }
2153     }
2154
2155     /* This code automatically gives table length tabext without cut-off's,
2156      * in that case grompp should already have checked that we do not need
2157      * normal tables and we only generate tables for 1-4 interactions.
2158      */
2159     rtab = ir->rlist + ir->tabext;
2160
2161     if (needGroupSchemeTables)
2162     {
2163         /* make tables for ordinary interactions */
2164         if (bSomeNormalNbListsAreInUse)
2165         {
2166             make_nbf_tables(fp, ic, rtab, tabfn, nullptr, nullptr, nullptr);
2167             m = 1;
2168         }
2169         else
2170         {
2171             m = 0;
2172         }
2173         if (negptable > 0)
2174         {
2175             /* Read the special tables for certain energy group pairs */
2176             nm_ind = mtop->groups.groups[SimulationAtomGroupType::EnergyOutput].nm_ind;
2177             for (egi = 0; egi < negp_pp; egi++)
2178             {
2179                 for (egj = egi; egj < negp_pp; egj++)
2180                 {
2181                     egp_flags = ir->opts.egp_flags[GID(egi, egj, ir->opts.ngener)];
2182                     if ((egp_flags & EGP_TABLE) && !(egp_flags & EGP_EXCL))
2183                     {
2184                         /* Read the table file with the two energy groups names appended */
2185                         make_nbf_tables(fp, ic, rtab, tabfn,
2186                                         *mtop->groups.groupNames[nm_ind[egi]],
2187                                         *mtop->groups.groupNames[nm_ind[egj]],
2188                                         nullptr);
2189                         m++;
2190                     }
2191                 }
2192             }
2193         }
2194     }
2195
2196     /* We want to use unmodified tables for 1-4 coulombic
2197      * interactions, so we must in general have an extra set of
2198      * tables. */
2199     if (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 ||
2200         gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0 ||
2201         gmx_mtop_ftype_count(mtop, F_LJC_PAIRS_NB) > 0)
2202     {
2203         fr->pairsTable = make_tables(fp, ic, tabpfn, rtab,
2204                                      GMX_MAKETABLES_14ONLY);
2205     }
2206
2207     /* Wall stuff */
2208     fr->nwall = ir->nwall;
2209     if (ir->nwall && ir->wall_type == ewtTABLE)
2210     {
2211         make_wall_tables(fp, ir, tabfn, &mtop->groups, fr);
2212     }
2213
2214     if (fcd && !tabbfnm.empty())
2215     {
2216         // Need to catch std::bad_alloc
2217         // TODO Don't need to catch this here, when merging with master branch
2218         try
2219         {
2220             fcd->bondtab  = make_bonded_tables(fp,
2221                                                F_TABBONDS, F_TABBONDSNC,
2222                                                mtop, tabbfnm, "b");
2223             fcd->angletab = make_bonded_tables(fp,
2224                                                F_TABANGLES, -1,
2225                                                mtop, tabbfnm, "a");
2226             fcd->dihtab   = make_bonded_tables(fp,
2227                                                F_TABDIHS, -1,
2228                                                mtop, tabbfnm, "d");
2229         }
2230         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
2231     }
2232     else
2233     {
2234         if (debug)
2235         {
2236             fprintf(debug, "No fcdata or table file name passed, can not read table, can not do bonded interactions\n");
2237         }
2238     }
2239
2240     // QM/MM initialization if requested
2241     fr->bQMMM = ir->bQMMM;
2242     if (fr->bQMMM)
2243     {
2244         // Initialize QM/MM if supported
2245         if (GMX_QMMM)
2246         {
2247             GMX_LOG(mdlog.info).asParagraph().
2248                 appendText("Large parts of the QM/MM support is deprecated, and may be removed in a future "
2249                            "version. Please get in touch with the developers if you find the support useful, "
2250                            "as help is needed if the functionality is to continue to be available.");
2251             fr->qr = mk_QMMMrec();
2252             init_QMMMrec(cr, mtop, ir, fr);
2253         }
2254         else
2255         {
2256             gmx_incons("QM/MM was requested, but is only available when GROMACS "
2257                        "is configured with QM/MM support");
2258         }
2259     }
2260
2261     /* Set all the static charge group info */
2262     fr->cginfo_mb = init_cginfo_mb(fp, mtop, fr, bNoSolvOpt,
2263                                    &bFEP_NonBonded,
2264                                    &fr->bExcl_IntraCGAll_InterCGNone);
2265     if (!DOMAINDECOMP(cr))
2266     {
2267         fr->cginfo = cginfo_expand(mtop->molblock.size(), fr->cginfo_mb);
2268     }
2269
2270     if (!DOMAINDECOMP(cr))
2271     {
2272         forcerec_set_ranges(fr, ncg_mtop(mtop), ncg_mtop(mtop),
2273                             mtop->natoms, mtop->natoms, mtop->natoms);
2274     }
2275
2276     fr->print_force = print_force;
2277
2278     /* Initialize the thread working data for bonded interactions */
2279     fr->bondedThreading =
2280         init_bonded_threading(fp, mtop->groups.groups[SimulationAtomGroupType::EnergyOutput].nr);
2281
2282     fr->nthread_ewc = gmx_omp_nthreads_get(emntBonded);
2283     snew(fr->ewc_t, fr->nthread_ewc);
2284
2285     if (fr->cutoff_scheme == ecutsVERLET)
2286     {
2287         // We checked the cut-offs in grompp, but double-check here.
2288         // We have PME+LJcutoff kernels for rcoulomb>rvdw.
2289         if (EEL_PME_EWALD(ir->coulombtype) && ir->vdwtype == eelCUT)
2290         {
2291             GMX_RELEASE_ASSERT(ir->rcoulomb >= ir->rvdw, "With Verlet lists and PME we should have rcoulomb>=rvdw");
2292         }
2293         else
2294         {
2295             GMX_RELEASE_ASSERT(ir->rcoulomb == ir->rvdw, "With Verlet lists and no PME rcoulomb and rvdw should be identical");
2296         }
2297
2298         fr->nbv = Nbnxm::init_nb_verlet(mdlog, bFEP_NonBonded, ir, fr,
2299                                         cr, hardwareInfo, deviceInfo,
2300                                         mtop, box);
2301
2302         if (useGpuForBonded)
2303         {
2304             auto stream = DOMAINDECOMP(cr) ?
2305                 Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, Nbnxm::InteractionLocality::NonLocal) :
2306                 Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv, Nbnxm::InteractionLocality::Local);
2307             // TODO the heap allocation is only needed while
2308             // t_forcerec lacks a constructor.
2309             fr->gpuBonded = new gmx::GpuBonded(mtop->ffparams,
2310                                                stream);
2311         }
2312     }
2313
2314     if (ir->eDispCorr != edispcNO)
2315     {
2316         fr->dispersionCorrection =
2317             std::make_unique<DispersionCorrection>(*mtop, *ir, fr->bBHAM,
2318                                                    fr->ntype,
2319                                                    gmx::arrayRefFromArray(fr->nbfp, fr->ntype*fr->ntype*2),
2320                                                    *fr->ic, tabfn);
2321         fr->dispersionCorrection->print(mdlog);
2322     }
2323
2324     if (fp != nullptr)
2325     {
2326         /* Here we switch from using mdlog, which prints the newline before
2327          * the paragraph, to our old fprintf logging, which prints the newline
2328          * after the paragraph, so we should add a newline here.
2329          */
2330         fprintf(fp, "\n");
2331     }
2332 }
2333
2334 /* Frees GPU memory and sets a tMPI node barrier.
2335  *
2336  * Note that this function needs to be called even if GPUs are not used
2337  * in this run because the PME ranks have no knowledge of whether GPUs
2338  * are used or not, but all ranks need to enter the barrier below.
2339  * \todo Remove physical node barrier from this function after making sure
2340  * that it's not needed anymore (with a shared GPU run).
2341  */
2342 void free_gpu_resources(t_forcerec                          *fr,
2343                         const gmx::PhysicalNodeCommunicator &physicalNodeCommunicator)
2344 {
2345     bool isPPrankUsingGPU = (fr != nullptr) && (fr->nbv != nullptr) && fr->nbv->useGpu();
2346
2347     /* stop the GPU profiler (only CUDA) */
2348     stopGpuProfiler();
2349
2350     if (isPPrankUsingGPU)
2351     {
2352         /* Free data in GPU memory and pinned memory before destroying the GPU context */
2353         fr->nbv.reset();
2354
2355         delete fr->gpuBonded;
2356         fr->gpuBonded = nullptr;
2357     }
2358
2359     /* With tMPI we need to wait for all ranks to finish deallocation before
2360      * destroying the CUDA context in free_gpu() as some tMPI ranks may be sharing
2361      * GPU and context.
2362      *
2363      * This is not a concern in OpenCL where we use one context per rank which
2364      * is freed in nbnxn_gpu_free().
2365      *
2366      * Note: it is safe to not call the barrier on the ranks which do not use GPU,
2367      * but it is easier and more futureproof to call it on the whole node.
2368      */
2369     if (GMX_THREAD_MPI)
2370     {
2371         physicalNodeCommunicator.barrier();
2372     }
2373 }
2374
2375 void done_forcerec(t_forcerec *fr, int numMolBlocks)
2376 {
2377     if (fr == nullptr)
2378     {
2379         // PME-only ranks don't have a forcerec
2380         return;
2381     }
2382     done_cginfo_mb(fr->cginfo_mb, numMolBlocks);
2383     sfree(fr->nbfp);
2384     done_interaction_const(fr->ic);
2385     sfree(fr->shift_vec);
2386     sfree(fr->fshift);
2387     sfree(fr->ewc_t);
2388     tear_down_bonded_threading(fr->bondedThreading);
2389     GMX_RELEASE_ASSERT(fr->gpuBonded == nullptr, "Should have been deleted earlier, when used");
2390     fr->bondedThreading = nullptr;
2391     delete fr;
2392 }