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