738fce07e49d01c2f3b6596726e8285dabc98d6b
[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, 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 "config.h"
40
41 #include <assert.h>
42 #include <math.h>
43 #include <stdlib.h>
44 #include <string.h>
45
46 #include <algorithm>
47
48 #include "gromacs/domdec/domdec.h"
49 #include "gromacs/ewald/ewald.h"
50 #include "gromacs/gmxlib/gpu_utils/gpu_utils.h"
51 #include "gromacs/legacyheaders/copyrite.h"
52 #include "gromacs/legacyheaders/force.h"
53 #include "gromacs/legacyheaders/gmx_detect_hardware.h"
54 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
55 #include "gromacs/legacyheaders/inputrec.h"
56 #include "gromacs/legacyheaders/macros.h"
57 #include "gromacs/legacyheaders/md_logging.h"
58 #include "gromacs/legacyheaders/md_support.h"
59 #include "gromacs/legacyheaders/names.h"
60 #include "gromacs/legacyheaders/network.h"
61 #include "gromacs/legacyheaders/nonbonded.h"
62 #include "gromacs/legacyheaders/ns.h"
63 #include "gromacs/legacyheaders/qmmm.h"
64 #include "gromacs/legacyheaders/tables.h"
65 #include "gromacs/legacyheaders/txtdump.h"
66 #include "gromacs/legacyheaders/typedefs.h"
67 #include "gromacs/legacyheaders/types/commrec.h"
68 #include "gromacs/listed-forces/manage-threading.h"
69 #include "gromacs/math/calculate-ewald-splitting-coefficient.h"
70 #include "gromacs/math/units.h"
71 #include "gromacs/math/utilities.h"
72 #include "gromacs/math/vec.h"
73 #include "gromacs/mdlib/forcerec-threading.h"
74 #include "gromacs/mdlib/nb_verlet.h"
75 #include "gromacs/mdlib/nbnxn_atomdata.h"
76 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
77 #include "gromacs/mdlib/nbnxn_search.h"
78 #include "gromacs/mdlib/nbnxn_simd.h"
79 #include "gromacs/pbcutil/ishift.h"
80 #include "gromacs/pbcutil/pbc.h"
81 #include "gromacs/simd/simd.h"
82 #include "gromacs/topology/mtop_util.h"
83 #include "gromacs/utility/fatalerror.h"
84 #include "gromacs/utility/smalloc.h"
85
86 #include "nbnxn_gpu_jit_support.h"
87
88 t_forcerec *mk_forcerec(void)
89 {
90     t_forcerec *fr;
91
92     snew(fr, 1);
93
94     return fr;
95 }
96
97 #ifdef DEBUG
98 static void pr_nbfp(FILE *fp, real *nbfp, gmx_bool bBHAM, int atnr)
99 {
100     int i, j;
101
102     for (i = 0; (i < atnr); i++)
103     {
104         for (j = 0; (j < atnr); j++)
105         {
106             fprintf(fp, "%2d - %2d", i, j);
107             if (bBHAM)
108             {
109                 fprintf(fp, "  a=%10g, b=%10g, c=%10g\n", BHAMA(nbfp, atnr, i, j),
110                         BHAMB(nbfp, atnr, i, j), BHAMC(nbfp, atnr, i, j)/6.0);
111             }
112             else
113             {
114                 fprintf(fp, "  c6=%10g, c12=%10g\n", C6(nbfp, atnr, i, j)/6.0,
115                         C12(nbfp, atnr, i, j)/12.0);
116             }
117         }
118     }
119 }
120 #endif
121
122 static real *mk_nbfp(const gmx_ffparams_t *idef, gmx_bool bBHAM)
123 {
124     real *nbfp;
125     int   i, j, k, atnr;
126
127     atnr = idef->atnr;
128     if (bBHAM)
129     {
130         snew(nbfp, 3*atnr*atnr);
131         for (i = k = 0; (i < atnr); i++)
132         {
133             for (j = 0; (j < atnr); j++, k++)
134             {
135                 BHAMA(nbfp, atnr, i, j) = idef->iparams[k].bham.a;
136                 BHAMB(nbfp, atnr, i, j) = idef->iparams[k].bham.b;
137                 /* nbfp now includes the 6.0 derivative prefactor */
138                 BHAMC(nbfp, atnr, i, j) = idef->iparams[k].bham.c*6.0;
139             }
140         }
141     }
142     else
143     {
144         snew(nbfp, 2*atnr*atnr);
145         for (i = k = 0; (i < atnr); i++)
146         {
147             for (j = 0; (j < atnr); j++, k++)
148             {
149                 /* nbfp now includes the 6.0/12.0 derivative prefactors */
150                 C6(nbfp, atnr, i, j)   = idef->iparams[k].lj.c6*6.0;
151                 C12(nbfp, atnr, i, j)  = idef->iparams[k].lj.c12*12.0;
152             }
153         }
154     }
155
156     return nbfp;
157 }
158
159 static real *make_ljpme_c6grid(const gmx_ffparams_t *idef, t_forcerec *fr)
160 {
161     int        i, j, k, atnr;
162     real       c6, c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
163     real      *grid;
164     const real oneOverSix = 1.0 / 6.0;
165
166     /* For LJ-PME simulations, we correct the energies with the reciprocal space
167      * inside of the cut-off. To do this the non-bonded kernels needs to have
168      * access to the C6-values used on the reciprocal grid in pme.c
169      */
170
171     atnr = idef->atnr;
172     snew(grid, 2*atnr*atnr);
173     for (i = k = 0; (i < atnr); i++)
174     {
175         for (j = 0; (j < atnr); j++, k++)
176         {
177             c6i  = idef->iparams[i*(atnr+1)].lj.c6;
178             c12i = idef->iparams[i*(atnr+1)].lj.c12;
179             c6j  = idef->iparams[j*(atnr+1)].lj.c6;
180             c12j = idef->iparams[j*(atnr+1)].lj.c12;
181             c6   = sqrt(c6i * c6j);
182             if (fr->ljpme_combination_rule == eljpmeLB
183                 && !gmx_numzero(c6) && !gmx_numzero(c12i) && !gmx_numzero(c12j))
184             {
185                 sigmai = pow(c12i / c6i, oneOverSix);
186                 sigmaj = pow(c12j / c6j, oneOverSix);
187                 epsi   = c6i * c6i / c12i;
188                 epsj   = c6j * c6j / c12j;
189                 c6     = sqrt(epsi * epsj) * pow(0.5*(sigmai+sigmaj), 6);
190             }
191             /* Store the elements at the same relative positions as C6 in nbfp in order
192              * to simplify access in the kernels
193              */
194             grid[2*(atnr*i+j)] = c6*6.0;
195         }
196     }
197     return grid;
198 }
199
200 static real *mk_nbfp_combination_rule(const gmx_ffparams_t *idef, int comb_rule)
201 {
202     real      *nbfp;
203     int        i, j, atnr;
204     real       c6i, c6j, c12i, c12j, epsi, epsj, sigmai, sigmaj;
205     real       c6, c12;
206     const real oneOverSix = 1.0 / 6.0;
207
208     atnr = idef->atnr;
209     snew(nbfp, 2*atnr*atnr);
210     for (i = 0; i < atnr; ++i)
211     {
212         for (j = 0; j < atnr; ++j)
213         {
214             c6i  = idef->iparams[i*(atnr+1)].lj.c6;
215             c12i = idef->iparams[i*(atnr+1)].lj.c12;
216             c6j  = idef->iparams[j*(atnr+1)].lj.c6;
217             c12j = idef->iparams[j*(atnr+1)].lj.c12;
218             c6   = sqrt(c6i  * c6j);
219             c12  = sqrt(c12i * c12j);
220             if (comb_rule == eCOMB_ARITHMETIC
221                 && !gmx_numzero(c6) && !gmx_numzero(c12))
222             {
223                 sigmai = pow(c12i / c6i, oneOverSix);
224                 sigmaj = pow(c12j / c6j, oneOverSix);
225                 epsi   = c6i * c6i / c12i;
226                 epsj   = c6j * c6j / c12j;
227                 c6     = sqrt(epsi * epsj) * pow(0.5*(sigmai+sigmaj), 6);
228                 c12    = sqrt(epsi * epsj) * pow(0.5*(sigmai+sigmaj), 12);
229             }
230             C6(nbfp, atnr, i, j)   = c6*6.0;
231             C12(nbfp, atnr, i, j)  = c12*12.0;
232         }
233     }
234     return nbfp;
235 }
236
237 /* This routine sets fr->solvent_opt to the most common solvent in the
238  * system, e.g. esolSPC or esolTIP4P. It will also mark each charge group in
239  * the fr->solvent_type array with the correct type (or esolNO).
240  *
241  * Charge groups that fulfill the conditions but are not identical to the
242  * most common one will be marked as esolNO in the solvent_type array.
243  *
244  * TIP3p is identical to SPC for these purposes, so we call it
245  * SPC in the arrays (Apologies to Bill Jorgensen ;-)
246  *
247  * NOTE: QM particle should not
248  * become an optimized solvent. Not even if there is only one charge
249  * group in the Qm
250  */
251
252 typedef struct
253 {
254     int    model;
255     int    count;
256     int    vdwtype[4];
257     real   charge[4];
258 } solvent_parameters_t;
259
260 static void
261 check_solvent_cg(const gmx_moltype_t    *molt,
262                  int                     cg0,
263                  int                     nmol,
264                  const unsigned char    *qm_grpnr,
265                  const t_grps           *qm_grps,
266                  t_forcerec   *          fr,
267                  int                    *n_solvent_parameters,
268                  solvent_parameters_t  **solvent_parameters_p,
269                  int                     cginfo,
270                  int                    *cg_sp)
271 {
272     t_atom               *atom;
273     int                   j, k;
274     int                   j0, j1, nj;
275     gmx_bool              perturbed;
276     gmx_bool              has_vdw[4];
277     gmx_bool              match;
278     real                  tmp_charge[4]  = { 0.0 }; /* init to zero to make gcc4.8 happy */
279     int                   tmp_vdwtype[4] = { 0 };   /* init to zero to make gcc4.8 happy */
280     int                   tjA;
281     gmx_bool              qm;
282     solvent_parameters_t *solvent_parameters;
283
284     /* We use a list with parameters for each solvent type.
285      * Every time we discover a new molecule that fulfills the basic
286      * conditions for a solvent we compare with the previous entries
287      * in these lists. If the parameters are the same we just increment
288      * the counter for that type, and otherwise we create a new type
289      * based on the current molecule.
290      *
291      * Once we've finished going through all molecules we check which
292      * solvent is most common, and mark all those molecules while we
293      * clear the flag on all others.
294      */
295
296     solvent_parameters = *solvent_parameters_p;
297
298     /* Mark the cg first as non optimized */
299     *cg_sp = -1;
300
301     /* Check if this cg has no exclusions with atoms in other charge groups
302      * and all atoms inside the charge group excluded.
303      * We only have 3 or 4 atom solvent loops.
304      */
305     if (GET_CGINFO_EXCL_INTER(cginfo) ||
306         !GET_CGINFO_EXCL_INTRA(cginfo))
307     {
308         return;
309     }
310
311     /* Get the indices of the first atom in this charge group */
312     j0     = molt->cgs.index[cg0];
313     j1     = molt->cgs.index[cg0+1];
314
315     /* Number of atoms in our molecule */
316     nj     = j1 - j0;
317
318     if (debug)
319     {
320         fprintf(debug,
321                 "Moltype '%s': there are %d atoms in this charge group\n",
322                 *molt->name, nj);
323     }
324
325     /* Check if it could be an SPC (3 atoms) or TIP4p (4) water,
326      * otherwise skip it.
327      */
328     if (nj < 3 || nj > 4)
329     {
330         return;
331     }
332
333     /* Check if we are doing QM on this group */
334     qm = FALSE;
335     if (qm_grpnr != NULL)
336     {
337         for (j = j0; j < j1 && !qm; j++)
338         {
339             qm = (qm_grpnr[j] < qm_grps->nr - 1);
340         }
341     }
342     /* Cannot use solvent optimization with QM */
343     if (qm)
344     {
345         return;
346     }
347
348     atom = molt->atoms.atom;
349
350     /* Still looks like a solvent, time to check parameters */
351
352     /* If it is perturbed (free energy) we can't use the solvent loops,
353      * so then we just skip to the next molecule.
354      */
355     perturbed = FALSE;
356
357     for (j = j0; j < j1 && !perturbed; j++)
358     {
359         perturbed = PERTURBED(atom[j]);
360     }
361
362     if (perturbed)
363     {
364         return;
365     }
366
367     /* Now it's only a question if the VdW and charge parameters
368      * are OK. Before doing the check we compare and see if they are
369      * identical to a possible previous solvent type.
370      * First we assign the current types and charges.
371      */
372     for (j = 0; j < nj; j++)
373     {
374         tmp_vdwtype[j] = atom[j0+j].type;
375         tmp_charge[j]  = atom[j0+j].q;
376     }
377
378     /* Does it match any previous solvent type? */
379     for (k = 0; k < *n_solvent_parameters; k++)
380     {
381         match = TRUE;
382
383
384         /* We can only match SPC with 3 atoms and TIP4p with 4 atoms */
385         if ( (solvent_parameters[k].model == esolSPC   && nj != 3)  ||
386              (solvent_parameters[k].model == esolTIP4P && nj != 4) )
387         {
388             match = FALSE;
389         }
390
391         /* Check that types & charges match for all atoms in molecule */
392         for (j = 0; j < nj && match == TRUE; j++)
393         {
394             if (tmp_vdwtype[j] != solvent_parameters[k].vdwtype[j])
395             {
396                 match = FALSE;
397             }
398             if (tmp_charge[j] != solvent_parameters[k].charge[j])
399             {
400                 match = FALSE;
401             }
402         }
403         if (match == TRUE)
404         {
405             /* Congratulations! We have a matched solvent.
406              * Flag it with this type for later processing.
407              */
408             *cg_sp = k;
409             solvent_parameters[k].count += nmol;
410
411             /* We are done with this charge group */
412             return;
413         }
414     }
415
416     /* If we get here, we have a tentative new solvent type.
417      * Before we add it we must check that it fulfills the requirements
418      * of the solvent optimized loops. First determine which atoms have
419      * VdW interactions.
420      */
421     for (j = 0; j < nj; j++)
422     {
423         has_vdw[j] = FALSE;
424         tjA        = tmp_vdwtype[j];
425
426         /* Go through all other tpes and see if any have non-zero
427          * VdW parameters when combined with this one.
428          */
429         for (k = 0; k < fr->ntype && (has_vdw[j] == FALSE); k++)
430         {
431             /* We already checked that the atoms weren't perturbed,
432              * so we only need to check state A now.
433              */
434             if (fr->bBHAM)
435             {
436                 has_vdw[j] = (has_vdw[j] ||
437                               (BHAMA(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
438                               (BHAMB(fr->nbfp, fr->ntype, tjA, k) != 0.0) ||
439                               (BHAMC(fr->nbfp, fr->ntype, tjA, k) != 0.0));
440             }
441             else
442             {
443                 /* Standard LJ */
444                 has_vdw[j] = (has_vdw[j] ||
445                               (C6(fr->nbfp, fr->ntype, tjA, k)  != 0.0) ||
446                               (C12(fr->nbfp, fr->ntype, tjA, k) != 0.0));
447             }
448         }
449     }
450
451     /* Now we know all we need to make the final check and assignment. */
452     if (nj == 3)
453     {
454         /* So, is it an SPC?
455          * For this we require thatn all atoms have charge,
456          * the charges on atom 2 & 3 should be the same, and only
457          * atom 1 might have VdW.
458          */
459         if (has_vdw[1] == FALSE &&
460             has_vdw[2] == FALSE &&
461             tmp_charge[0]  != 0 &&
462             tmp_charge[1]  != 0 &&
463             tmp_charge[2]  == tmp_charge[1])
464         {
465             srenew(solvent_parameters, *n_solvent_parameters+1);
466             solvent_parameters[*n_solvent_parameters].model = esolSPC;
467             solvent_parameters[*n_solvent_parameters].count = nmol;
468             for (k = 0; k < 3; k++)
469             {
470                 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
471                 solvent_parameters[*n_solvent_parameters].charge[k]  = tmp_charge[k];
472             }
473
474             *cg_sp = *n_solvent_parameters;
475             (*n_solvent_parameters)++;
476         }
477     }
478     else if (nj == 4)
479     {
480         /* Or could it be a TIP4P?
481          * For this we require thatn atoms 2,3,4 have charge, but not atom 1.
482          * Only atom 1 mght have VdW.
483          */
484         if (has_vdw[1] == FALSE &&
485             has_vdw[2] == FALSE &&
486             has_vdw[3] == FALSE &&
487             tmp_charge[0]  == 0 &&
488             tmp_charge[1]  != 0 &&
489             tmp_charge[2]  == tmp_charge[1] &&
490             tmp_charge[3]  != 0)
491         {
492             srenew(solvent_parameters, *n_solvent_parameters+1);
493             solvent_parameters[*n_solvent_parameters].model = esolTIP4P;
494             solvent_parameters[*n_solvent_parameters].count = nmol;
495             for (k = 0; k < 4; k++)
496             {
497                 solvent_parameters[*n_solvent_parameters].vdwtype[k] = tmp_vdwtype[k];
498                 solvent_parameters[*n_solvent_parameters].charge[k]  = tmp_charge[k];
499             }
500
501             *cg_sp = *n_solvent_parameters;
502             (*n_solvent_parameters)++;
503         }
504     }
505
506     *solvent_parameters_p = solvent_parameters;
507 }
508
509 static void
510 check_solvent(FILE  *                fp,
511               const gmx_mtop_t  *    mtop,
512               t_forcerec  *          fr,
513               cginfo_mb_t           *cginfo_mb)
514 {
515     const t_block     *   cgs;
516     const gmx_moltype_t  *molt;
517     int                   mb, mol, cg_mol, at_offset, am, cgm, i, nmol_ch, nmol;
518     int                   n_solvent_parameters;
519     solvent_parameters_t *solvent_parameters;
520     int                 **cg_sp;
521     int                   bestsp, bestsol;
522
523     if (debug)
524     {
525         fprintf(debug, "Going to determine what solvent types we have.\n");
526     }
527
528     n_solvent_parameters = 0;
529     solvent_parameters   = NULL;
530     /* Allocate temporary array for solvent type */
531     snew(cg_sp, mtop->nmolblock);
532
533     at_offset = 0;
534     for (mb = 0; mb < mtop->nmolblock; mb++)
535     {
536         molt = &mtop->moltype[mtop->molblock[mb].type];
537         cgs  = &molt->cgs;
538         /* Here we have to loop over all individual molecules
539          * because we need to check for QMMM particles.
540          */
541         snew(cg_sp[mb], cginfo_mb[mb].cg_mod);
542         nmol_ch = cginfo_mb[mb].cg_mod/cgs->nr;
543         nmol    = mtop->molblock[mb].nmol/nmol_ch;
544         for (mol = 0; mol < nmol_ch; mol++)
545         {
546             cgm = mol*cgs->nr;
547             am  = mol*cgs->index[cgs->nr];
548             for (cg_mol = 0; cg_mol < cgs->nr; cg_mol++)
549             {
550                 check_solvent_cg(molt, cg_mol, nmol,
551                                  mtop->groups.grpnr[egcQMMM] ?
552                                  mtop->groups.grpnr[egcQMMM]+at_offset+am : 0,
553                                  &mtop->groups.grps[egcQMMM],
554                                  fr,
555                                  &n_solvent_parameters, &solvent_parameters,
556                                  cginfo_mb[mb].cginfo[cgm+cg_mol],
557                                  &cg_sp[mb][cgm+cg_mol]);
558             }
559         }
560         at_offset += cgs->index[cgs->nr];
561     }
562
563     /* Puh! We finished going through all charge groups.
564      * Now find the most common solvent model.
565      */
566
567     /* Most common solvent this far */
568     bestsp = -2;
569     for (i = 0; i < n_solvent_parameters; i++)
570     {
571         if (bestsp == -2 ||
572             solvent_parameters[i].count > solvent_parameters[bestsp].count)
573         {
574             bestsp = i;
575         }
576     }
577
578     if (bestsp >= 0)
579     {
580         bestsol = solvent_parameters[bestsp].model;
581     }
582     else
583     {
584         bestsol = esolNO;
585     }
586
587     fr->nWatMol = 0;
588     for (mb = 0; mb < mtop->nmolblock; mb++)
589     {
590         cgs  = &mtop->moltype[mtop->molblock[mb].type].cgs;
591         nmol = (mtop->molblock[mb].nmol*cgs->nr)/cginfo_mb[mb].cg_mod;
592         for (i = 0; i < cginfo_mb[mb].cg_mod; i++)
593         {
594             if (cg_sp[mb][i] == bestsp)
595             {
596                 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i], bestsol);
597                 fr->nWatMol += nmol;
598             }
599             else
600             {
601                 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[i], esolNO);
602             }
603         }
604         sfree(cg_sp[mb]);
605     }
606     sfree(cg_sp);
607
608     if (bestsol != esolNO && fp != NULL)
609     {
610         fprintf(fp, "\nEnabling %s-like water optimization for %d molecules.\n\n",
611                 esol_names[bestsol],
612                 solvent_parameters[bestsp].count);
613     }
614
615     sfree(solvent_parameters);
616     fr->solvent_opt = bestsol;
617 }
618
619 enum {
620     acNONE = 0, acCONSTRAINT, acSETTLE
621 };
622
623 static cginfo_mb_t *init_cginfo_mb(FILE *fplog, const gmx_mtop_t *mtop,
624                                    t_forcerec *fr, gmx_bool bNoSolvOpt,
625                                    gmx_bool *bFEP_NonBonded,
626                                    gmx_bool *bExcl_IntraCGAll_InterCGNone)
627 {
628     const t_block        *cgs;
629     const t_blocka       *excl;
630     const gmx_moltype_t  *molt;
631     const gmx_molblock_t *molb;
632     cginfo_mb_t          *cginfo_mb;
633     gmx_bool             *type_VDW;
634     int                  *cginfo;
635     int                   cg_offset, a_offset;
636     int                   mb, m, cg, a0, a1, gid, ai, j, aj, excl_nalloc;
637     int                  *a_con;
638     int                   ftype;
639     int                   ia;
640     gmx_bool              bId, *bExcl, bExclIntraAll, bExclInter, bHaveVDW, bHaveQ, bHavePerturbedAtoms;
641
642     snew(cginfo_mb, mtop->nmolblock);
643
644     snew(type_VDW, fr->ntype);
645     for (ai = 0; ai < fr->ntype; ai++)
646     {
647         type_VDW[ai] = FALSE;
648         for (j = 0; j < fr->ntype; j++)
649         {
650             type_VDW[ai] = type_VDW[ai] ||
651                 fr->bBHAM ||
652                 C6(fr->nbfp, fr->ntype, ai, j) != 0 ||
653                 C12(fr->nbfp, fr->ntype, ai, j) != 0;
654         }
655     }
656
657     *bFEP_NonBonded               = FALSE;
658     *bExcl_IntraCGAll_InterCGNone = TRUE;
659
660     excl_nalloc = 10;
661     snew(bExcl, excl_nalloc);
662     cg_offset = 0;
663     a_offset  = 0;
664     for (mb = 0; mb < mtop->nmolblock; mb++)
665     {
666         molb = &mtop->molblock[mb];
667         molt = &mtop->moltype[molb->type];
668         cgs  = &molt->cgs;
669         excl = &molt->excls;
670
671         /* Check if the cginfo is identical for all molecules in this block.
672          * If so, we only need an array of the size of one molecule.
673          * Otherwise we make an array of #mol times #cgs per molecule.
674          */
675         bId = TRUE;
676         for (m = 0; m < molb->nmol; m++)
677         {
678             int am = m*cgs->index[cgs->nr];
679             for (cg = 0; cg < cgs->nr; cg++)
680             {
681                 a0 = cgs->index[cg];
682                 a1 = cgs->index[cg+1];
683                 if (ggrpnr(&mtop->groups, egcENER, a_offset+am+a0) !=
684                     ggrpnr(&mtop->groups, egcENER, a_offset   +a0))
685                 {
686                     bId = FALSE;
687                 }
688                 if (mtop->groups.grpnr[egcQMMM] != NULL)
689                 {
690                     for (ai = a0; ai < a1; ai++)
691                     {
692                         if (mtop->groups.grpnr[egcQMMM][a_offset+am+ai] !=
693                             mtop->groups.grpnr[egcQMMM][a_offset   +ai])
694                         {
695                             bId = FALSE;
696                         }
697                     }
698                 }
699             }
700         }
701
702         cginfo_mb[mb].cg_start = cg_offset;
703         cginfo_mb[mb].cg_end   = cg_offset + molb->nmol*cgs->nr;
704         cginfo_mb[mb].cg_mod   = (bId ? 1 : molb->nmol)*cgs->nr;
705         snew(cginfo_mb[mb].cginfo, cginfo_mb[mb].cg_mod);
706         cginfo = cginfo_mb[mb].cginfo;
707
708         /* Set constraints flags for constrained atoms */
709         snew(a_con, molt->atoms.nr);
710         for (ftype = 0; ftype < F_NRE; ftype++)
711         {
712             if (interaction_function[ftype].flags & IF_CONSTRAINT)
713             {
714                 int nral;
715
716                 nral = NRAL(ftype);
717                 for (ia = 0; ia < molt->ilist[ftype].nr; ia += 1+nral)
718                 {
719                     int a;
720
721                     for (a = 0; a < nral; a++)
722                     {
723                         a_con[molt->ilist[ftype].iatoms[ia+1+a]] =
724                             (ftype == F_SETTLE ? acSETTLE : acCONSTRAINT);
725                     }
726                 }
727             }
728         }
729
730         for (m = 0; m < (bId ? 1 : molb->nmol); m++)
731         {
732             int cgm = m*cgs->nr;
733             int am  = m*cgs->index[cgs->nr];
734             for (cg = 0; cg < cgs->nr; cg++)
735             {
736                 a0 = cgs->index[cg];
737                 a1 = cgs->index[cg+1];
738
739                 /* Store the energy group in cginfo */
740                 gid = ggrpnr(&mtop->groups, egcENER, a_offset+am+a0);
741                 SET_CGINFO_GID(cginfo[cgm+cg], gid);
742
743                 /* Check the intra/inter charge group exclusions */
744                 if (a1-a0 > excl_nalloc)
745                 {
746                     excl_nalloc = a1 - a0;
747                     srenew(bExcl, excl_nalloc);
748                 }
749                 /* bExclIntraAll: all intra cg interactions excluded
750                  * bExclInter:    any inter cg interactions excluded
751                  */
752                 bExclIntraAll       = TRUE;
753                 bExclInter          = FALSE;
754                 bHaveVDW            = FALSE;
755                 bHaveQ              = FALSE;
756                 bHavePerturbedAtoms = FALSE;
757                 for (ai = a0; ai < a1; ai++)
758                 {
759                     /* Check VDW and electrostatic interactions */
760                     bHaveVDW = bHaveVDW || (type_VDW[molt->atoms.atom[ai].type] ||
761                                             type_VDW[molt->atoms.atom[ai].typeB]);
762                     bHaveQ  = bHaveQ    || (molt->atoms.atom[ai].q != 0 ||
763                                             molt->atoms.atom[ai].qB != 0);
764
765                     bHavePerturbedAtoms = bHavePerturbedAtoms || (PERTURBED(molt->atoms.atom[ai]) != 0);
766
767                     /* Clear the exclusion list for atom ai */
768                     for (aj = a0; aj < a1; aj++)
769                     {
770                         bExcl[aj-a0] = FALSE;
771                     }
772                     /* Loop over all the exclusions of atom ai */
773                     for (j = excl->index[ai]; j < excl->index[ai+1]; j++)
774                     {
775                         aj = excl->a[j];
776                         if (aj < a0 || aj >= a1)
777                         {
778                             bExclInter = TRUE;
779                         }
780                         else
781                         {
782                             bExcl[aj-a0] = TRUE;
783                         }
784                     }
785                     /* Check if ai excludes a0 to a1 */
786                     for (aj = a0; aj < a1; aj++)
787                     {
788                         if (!bExcl[aj-a0])
789                         {
790                             bExclIntraAll = FALSE;
791                         }
792                     }
793
794                     switch (a_con[ai])
795                     {
796                         case acCONSTRAINT:
797                             SET_CGINFO_CONSTR(cginfo[cgm+cg]);
798                             break;
799                         case acSETTLE:
800                             SET_CGINFO_SETTLE(cginfo[cgm+cg]);
801                             break;
802                         default:
803                             break;
804                     }
805                 }
806                 if (bExclIntraAll)
807                 {
808                     SET_CGINFO_EXCL_INTRA(cginfo[cgm+cg]);
809                 }
810                 if (bExclInter)
811                 {
812                     SET_CGINFO_EXCL_INTER(cginfo[cgm+cg]);
813                 }
814                 if (a1 - a0 > MAX_CHARGEGROUP_SIZE)
815                 {
816                     /* The size in cginfo is currently only read with DD */
817                     gmx_fatal(FARGS, "A charge group has size %d which is larger than the limit of %d atoms", a1-a0, MAX_CHARGEGROUP_SIZE);
818                 }
819                 if (bHaveVDW)
820                 {
821                     SET_CGINFO_HAS_VDW(cginfo[cgm+cg]);
822                 }
823                 if (bHaveQ)
824                 {
825                     SET_CGINFO_HAS_Q(cginfo[cgm+cg]);
826                 }
827                 if (bHavePerturbedAtoms && fr->efep != efepNO)
828                 {
829                     SET_CGINFO_FEP(cginfo[cgm+cg]);
830                     *bFEP_NonBonded = TRUE;
831                 }
832                 /* Store the charge group size */
833                 SET_CGINFO_NATOMS(cginfo[cgm+cg], a1-a0);
834
835                 if (!bExclIntraAll || bExclInter)
836                 {
837                     *bExcl_IntraCGAll_InterCGNone = FALSE;
838                 }
839             }
840         }
841
842         sfree(a_con);
843
844         cg_offset += molb->nmol*cgs->nr;
845         a_offset  += molb->nmol*cgs->index[cgs->nr];
846     }
847     sfree(bExcl);
848
849     /* the solvent optimizer is called after the QM is initialized,
850      * because we don't want to have the QM subsystemto become an
851      * optimized solvent
852      */
853
854     check_solvent(fplog, mtop, fr, cginfo_mb);
855
856     if (getenv("GMX_NO_SOLV_OPT"))
857     {
858         if (fplog)
859         {
860             fprintf(fplog, "Found environment variable GMX_NO_SOLV_OPT.\n"
861                     "Disabling all solvent optimization\n");
862         }
863         fr->solvent_opt = esolNO;
864     }
865     if (bNoSolvOpt)
866     {
867         fr->solvent_opt = esolNO;
868     }
869     if (!fr->solvent_opt)
870     {
871         for (mb = 0; mb < mtop->nmolblock; mb++)
872         {
873             for (cg = 0; cg < cginfo_mb[mb].cg_mod; cg++)
874             {
875                 SET_CGINFO_SOLOPT(cginfo_mb[mb].cginfo[cg], esolNO);
876             }
877         }
878     }
879
880     return cginfo_mb;
881 }
882
883 static int *cginfo_expand(int nmb, cginfo_mb_t *cgi_mb)
884 {
885     int  ncg, mb, cg;
886     int *cginfo;
887
888     ncg = cgi_mb[nmb-1].cg_end;
889     snew(cginfo, ncg);
890     mb = 0;
891     for (cg = 0; cg < ncg; cg++)
892     {
893         while (cg >= cgi_mb[mb].cg_end)
894         {
895             mb++;
896         }
897         cginfo[cg] =
898             cgi_mb[mb].cginfo[(cg - cgi_mb[mb].cg_start) % cgi_mb[mb].cg_mod];
899     }
900
901     return cginfo;
902 }
903
904 static void set_chargesum(FILE *log, t_forcerec *fr, const gmx_mtop_t *mtop)
905 {
906     /*This now calculates sum for q and c6*/
907     double         qsum, q2sum, q, c6sum, c6;
908     int            mb, nmol, i;
909     const t_atoms *atoms;
910
911     qsum   = 0;
912     q2sum  = 0;
913     c6sum  = 0;
914     for (mb = 0; mb < mtop->nmolblock; mb++)
915     {
916         nmol  = mtop->molblock[mb].nmol;
917         atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
918         for (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 (mb = 0; mb < mtop->nmolblock; mb++)
937         {
938             nmol  = mtop->molblock[mb].nmol;
939             atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
940             for (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
973 void update_forcerec(t_forcerec *fr, matrix box)
974 {
975     if (fr->eeltype == eelGRF)
976     {
977         calc_rffac(NULL, fr->eeltype, fr->epsilon_r, fr->epsilon_rf,
978                    fr->rcoulomb, fr->temp, fr->zsquare, box,
979                    &fr->kappa, &fr->k_rf, &fr->c_rf);
980     }
981 }
982
983 void set_avcsixtwelve(FILE *fplog, t_forcerec *fr, const gmx_mtop_t *mtop)
984 {
985     const t_atoms  *atoms, *atoms_tpi;
986     const t_blocka *excl;
987     int             mb, nmol, nmolc, i, j, tpi, tpj, j1, j2, k, nexcl, q;
988     gmx_int64_t     npair, npair_ij, tmpi, tmpj;
989     double          csix, ctwelve;
990     int             ntp, *typecount;
991     gmx_bool        bBHAM;
992     real           *nbfp;
993     real           *nbfp_comb = NULL;
994
995     ntp   = fr->ntype;
996     bBHAM = fr->bBHAM;
997     nbfp  = fr->nbfp;
998
999     /* For LJ-PME, we want to correct for the difference between the
1000      * actual C6 values and the C6 values used by the LJ-PME based on
1001      * combination rules. */
1002
1003     if (EVDW_PME(fr->vdwtype))
1004     {
1005         nbfp_comb = mk_nbfp_combination_rule(&mtop->ffparams,
1006                                              (fr->ljpme_combination_rule == eljpmeLB) ? eCOMB_ARITHMETIC : eCOMB_GEOMETRIC);
1007         for (tpi = 0; tpi < ntp; ++tpi)
1008         {
1009             for (tpj = 0; tpj < ntp; ++tpj)
1010             {
1011                 C6(nbfp_comb, ntp, tpi, tpj) =
1012                     C6(nbfp, ntp, tpi, tpj) - C6(nbfp_comb, ntp, tpi, tpj);
1013                 C12(nbfp_comb, ntp, tpi, tpj) = C12(nbfp, ntp, tpi, tpj);
1014             }
1015         }
1016         nbfp = nbfp_comb;
1017     }
1018     for (q = 0; q < (fr->efep == efepNO ? 1 : 2); q++)
1019     {
1020         csix    = 0;
1021         ctwelve = 0;
1022         npair   = 0;
1023         nexcl   = 0;
1024         if (!fr->n_tpi)
1025         {
1026             /* Count the types so we avoid natoms^2 operations */
1027             snew(typecount, ntp);
1028             gmx_mtop_count_atomtypes(mtop, q, typecount);
1029
1030             for (tpi = 0; tpi < ntp; tpi++)
1031             {
1032                 for (tpj = tpi; tpj < ntp; tpj++)
1033                 {
1034                     tmpi = typecount[tpi];
1035                     tmpj = typecount[tpj];
1036                     if (tpi != tpj)
1037                     {
1038                         npair_ij = tmpi*tmpj;
1039                     }
1040                     else
1041                     {
1042                         npair_ij = tmpi*(tmpi - 1)/2;
1043                     }
1044                     if (bBHAM)
1045                     {
1046                         /* nbfp now includes the 6.0 derivative prefactor */
1047                         csix    += npair_ij*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1048                     }
1049                     else
1050                     {
1051                         /* nbfp now includes the 6.0/12.0 derivative prefactors */
1052                         csix    += npair_ij*   C6(nbfp, ntp, tpi, tpj)/6.0;
1053                         ctwelve += npair_ij*  C12(nbfp, ntp, tpi, tpj)/12.0;
1054                     }
1055                     npair += npair_ij;
1056                 }
1057             }
1058             sfree(typecount);
1059             /* Subtract the excluded pairs.
1060              * The main reason for substracting exclusions is that in some cases
1061              * some combinations might never occur and the parameters could have
1062              * any value. These unused values should not influence the dispersion
1063              * correction.
1064              */
1065             for (mb = 0; mb < mtop->nmolblock; mb++)
1066             {
1067                 nmol  = mtop->molblock[mb].nmol;
1068                 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
1069                 excl  = &mtop->moltype[mtop->molblock[mb].type].excls;
1070                 for (i = 0; (i < atoms->nr); i++)
1071                 {
1072                     if (q == 0)
1073                     {
1074                         tpi = atoms->atom[i].type;
1075                     }
1076                     else
1077                     {
1078                         tpi = atoms->atom[i].typeB;
1079                     }
1080                     j1  = excl->index[i];
1081                     j2  = excl->index[i+1];
1082                     for (j = j1; j < j2; j++)
1083                     {
1084                         k = excl->a[j];
1085                         if (k > i)
1086                         {
1087                             if (q == 0)
1088                             {
1089                                 tpj = atoms->atom[k].type;
1090                             }
1091                             else
1092                             {
1093                                 tpj = atoms->atom[k].typeB;
1094                             }
1095                             if (bBHAM)
1096                             {
1097                                 /* nbfp now includes the 6.0 derivative prefactor */
1098                                 csix -= nmol*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1099                             }
1100                             else
1101                             {
1102                                 /* nbfp now includes the 6.0/12.0 derivative prefactors */
1103                                 csix    -= nmol*C6 (nbfp, ntp, tpi, tpj)/6.0;
1104                                 ctwelve -= nmol*C12(nbfp, ntp, tpi, tpj)/12.0;
1105                             }
1106                             nexcl += nmol;
1107                         }
1108                     }
1109                 }
1110             }
1111         }
1112         else
1113         {
1114             /* Only correct for the interaction of the test particle
1115              * with the rest of the system.
1116              */
1117             atoms_tpi =
1118                 &mtop->moltype[mtop->molblock[mtop->nmolblock-1].type].atoms;
1119
1120             npair = 0;
1121             for (mb = 0; mb < mtop->nmolblock; mb++)
1122             {
1123                 nmol  = mtop->molblock[mb].nmol;
1124                 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
1125                 for (j = 0; j < atoms->nr; j++)
1126                 {
1127                     nmolc = nmol;
1128                     /* Remove the interaction of the test charge group
1129                      * with itself.
1130                      */
1131                     if (mb == mtop->nmolblock-1)
1132                     {
1133                         nmolc--;
1134
1135                         if (mb == 0 && nmol == 1)
1136                         {
1137                             gmx_fatal(FARGS, "Old format tpr with TPI, please generate a new tpr file");
1138                         }
1139                     }
1140                     if (q == 0)
1141                     {
1142                         tpj = atoms->atom[j].type;
1143                     }
1144                     else
1145                     {
1146                         tpj = atoms->atom[j].typeB;
1147                     }
1148                     for (i = 0; i < fr->n_tpi; i++)
1149                     {
1150                         if (q == 0)
1151                         {
1152                             tpi = atoms_tpi->atom[i].type;
1153                         }
1154                         else
1155                         {
1156                             tpi = atoms_tpi->atom[i].typeB;
1157                         }
1158                         if (bBHAM)
1159                         {
1160                             /* nbfp now includes the 6.0 derivative prefactor */
1161                             csix    += nmolc*BHAMC(nbfp, ntp, tpi, tpj)/6.0;
1162                         }
1163                         else
1164                         {
1165                             /* nbfp now includes the 6.0/12.0 derivative prefactors */
1166                             csix    += nmolc*C6 (nbfp, ntp, tpi, tpj)/6.0;
1167                             ctwelve += nmolc*C12(nbfp, ntp, tpi, tpj)/12.0;
1168                         }
1169                         npair += nmolc;
1170                     }
1171                 }
1172             }
1173         }
1174         if (npair - nexcl <= 0 && fplog)
1175         {
1176             fprintf(fplog, "\nWARNING: There are no atom pairs for dispersion correction\n\n");
1177             csix     = 0;
1178             ctwelve  = 0;
1179         }
1180         else
1181         {
1182             csix    /= npair - nexcl;
1183             ctwelve /= npair - nexcl;
1184         }
1185         if (debug)
1186         {
1187             fprintf(debug, "Counted %d exclusions\n", nexcl);
1188             fprintf(debug, "Average C6 parameter is: %10g\n", (double)csix);
1189             fprintf(debug, "Average C12 parameter is: %10g\n", (double)ctwelve);
1190         }
1191         fr->avcsix[q]    = csix;
1192         fr->avctwelve[q] = ctwelve;
1193     }
1194
1195     if (EVDW_PME(fr->vdwtype))
1196     {
1197         sfree(nbfp_comb);
1198     }
1199
1200     if (fplog != NULL)
1201     {
1202         if (fr->eDispCorr == edispcAllEner ||
1203             fr->eDispCorr == edispcAllEnerPres)
1204         {
1205             fprintf(fplog, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
1206                     fr->avcsix[0], fr->avctwelve[0]);
1207         }
1208         else
1209         {
1210             fprintf(fplog, "Long Range LJ corr.: <C6> %10.4e\n", fr->avcsix[0]);
1211         }
1212     }
1213 }
1214
1215
1216 static void set_bham_b_max(FILE *fplog, t_forcerec *fr,
1217                            const gmx_mtop_t *mtop)
1218 {
1219     const t_atoms *at1, *at2;
1220     int            mt1, mt2, i, j, tpi, tpj, ntypes;
1221     real           b, bmin;
1222     real          *nbfp;
1223
1224     if (fplog)
1225     {
1226         fprintf(fplog, "Determining largest Buckingham b parameter for table\n");
1227     }
1228     nbfp   = fr->nbfp;
1229     ntypes = fr->ntype;
1230
1231     bmin           = -1;
1232     fr->bham_b_max = 0;
1233     for (mt1 = 0; mt1 < mtop->nmoltype; mt1++)
1234     {
1235         at1 = &mtop->moltype[mt1].atoms;
1236         for (i = 0; (i < at1->nr); i++)
1237         {
1238             tpi = at1->atom[i].type;
1239             if (tpi >= ntypes)
1240             {
1241                 gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", i, tpi, ntypes);
1242             }
1243
1244             for (mt2 = mt1; mt2 < mtop->nmoltype; mt2++)
1245             {
1246                 at2 = &mtop->moltype[mt2].atoms;
1247                 for (j = 0; (j < at2->nr); j++)
1248                 {
1249                     tpj = at2->atom[j].type;
1250                     if (tpj >= ntypes)
1251                     {
1252                         gmx_fatal(FARGS, "Atomtype[%d] = %d, maximum = %d", j, tpj, ntypes);
1253                     }
1254                     b = BHAMB(nbfp, ntypes, tpi, tpj);
1255                     if (b > fr->bham_b_max)
1256                     {
1257                         fr->bham_b_max = b;
1258                     }
1259                     if ((b < bmin) || (bmin == -1))
1260                     {
1261                         bmin = b;
1262                     }
1263                 }
1264             }
1265         }
1266     }
1267     if (fplog)
1268     {
1269         fprintf(fplog, "Buckingham b parameters, min: %g, max: %g\n",
1270                 bmin, fr->bham_b_max);
1271     }
1272 }
1273
1274 static void make_nbf_tables(FILE *fp, const output_env_t oenv,
1275                             t_forcerec *fr, real rtab,
1276                             const t_commrec *cr,
1277                             const char *tabfn, char *eg1, char *eg2,
1278                             t_nblists *nbl)
1279 {
1280     char buf[STRLEN];
1281     int  i, j;
1282
1283     if (tabfn == NULL)
1284     {
1285         if (debug)
1286         {
1287             fprintf(debug, "No table file name passed, can not read table, can not do non-bonded interactions\n");
1288         }
1289         return;
1290     }
1291
1292     sprintf(buf, "%s", tabfn);
1293     if (eg1 && eg2)
1294     {
1295         /* Append the two energy group names */
1296         sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "_%s_%s.%s",
1297                 eg1, eg2, ftp2ext(efXVG));
1298     }
1299     nbl->table_elec_vdw = make_tables(fp, oenv, fr, MASTER(cr), buf, rtab, 0);
1300     /* Copy the contents of the table to separate coulomb and LJ tables too,
1301      * to improve cache performance.
1302      */
1303     /* For performance reasons we want
1304      * the table data to be aligned to 16-byte. The pointers could be freed
1305      * but currently aren't.
1306      */
1307     nbl->table_elec.interaction   = GMX_TABLE_INTERACTION_ELEC;
1308     nbl->table_elec.format        = nbl->table_elec_vdw.format;
1309     nbl->table_elec.r             = nbl->table_elec_vdw.r;
1310     nbl->table_elec.n             = nbl->table_elec_vdw.n;
1311     nbl->table_elec.scale         = nbl->table_elec_vdw.scale;
1312     nbl->table_elec.scale_exp     = nbl->table_elec_vdw.scale_exp;
1313     nbl->table_elec.formatsize    = nbl->table_elec_vdw.formatsize;
1314     nbl->table_elec.ninteractions = 1;
1315     nbl->table_elec.stride        = nbl->table_elec.formatsize * nbl->table_elec.ninteractions;
1316     snew_aligned(nbl->table_elec.data, nbl->table_elec.stride*(nbl->table_elec.n+1), 32);
1317
1318     nbl->table_vdw.interaction   = GMX_TABLE_INTERACTION_VDWREP_VDWDISP;
1319     nbl->table_vdw.format        = nbl->table_elec_vdw.format;
1320     nbl->table_vdw.r             = nbl->table_elec_vdw.r;
1321     nbl->table_vdw.n             = nbl->table_elec_vdw.n;
1322     nbl->table_vdw.scale         = nbl->table_elec_vdw.scale;
1323     nbl->table_vdw.scale_exp     = nbl->table_elec_vdw.scale_exp;
1324     nbl->table_vdw.formatsize    = nbl->table_elec_vdw.formatsize;
1325     nbl->table_vdw.ninteractions = 2;
1326     nbl->table_vdw.stride        = nbl->table_vdw.formatsize * nbl->table_vdw.ninteractions;
1327     snew_aligned(nbl->table_vdw.data, nbl->table_vdw.stride*(nbl->table_vdw.n+1), 32);
1328
1329     for (i = 0; i <= nbl->table_elec_vdw.n; i++)
1330     {
1331         for (j = 0; j < 4; j++)
1332         {
1333             nbl->table_elec.data[4*i+j] = nbl->table_elec_vdw.data[12*i+j];
1334         }
1335         for (j = 0; j < 8; j++)
1336         {
1337             nbl->table_vdw.data[8*i+j] = nbl->table_elec_vdw.data[12*i+4+j];
1338         }
1339     }
1340 }
1341
1342 static void count_tables(int ftype1, int ftype2, const gmx_mtop_t *mtop,
1343                          int *ncount, int **count)
1344 {
1345     const gmx_moltype_t *molt;
1346     const t_ilist       *il;
1347     int                  mt, ftype, stride, i, j, tabnr;
1348
1349     for (mt = 0; mt < mtop->nmoltype; mt++)
1350     {
1351         molt = &mtop->moltype[mt];
1352         for (ftype = 0; ftype < F_NRE; ftype++)
1353         {
1354             if (ftype == ftype1 || ftype == ftype2)
1355             {
1356                 il     = &molt->ilist[ftype];
1357                 stride = 1 + NRAL(ftype);
1358                 for (i = 0; i < il->nr; i += stride)
1359                 {
1360                     tabnr = mtop->ffparams.iparams[il->iatoms[i]].tab.table;
1361                     if (tabnr < 0)
1362                     {
1363                         gmx_fatal(FARGS, "A bonded table number is smaller than 0: %d\n", tabnr);
1364                     }
1365                     if (tabnr >= *ncount)
1366                     {
1367                         srenew(*count, tabnr+1);
1368                         for (j = *ncount; j < tabnr+1; j++)
1369                         {
1370                             (*count)[j] = 0;
1371                         }
1372                         *ncount = tabnr+1;
1373                     }
1374                     (*count)[tabnr]++;
1375                 }
1376             }
1377         }
1378     }
1379 }
1380
1381 static bondedtable_t *make_bonded_tables(FILE *fplog,
1382                                          int ftype1, int ftype2,
1383                                          const gmx_mtop_t *mtop,
1384                                          const char *basefn, const char *tabext)
1385 {
1386     int            i, ncount, *count;
1387     char           tabfn[STRLEN];
1388     bondedtable_t *tab;
1389
1390     tab = NULL;
1391
1392     ncount = 0;
1393     count  = NULL;
1394     count_tables(ftype1, ftype2, mtop, &ncount, &count);
1395
1396     if (ncount > 0)
1397     {
1398         snew(tab, ncount);
1399         for (i = 0; i < ncount; i++)
1400         {
1401             if (count[i] > 0)
1402             {
1403                 sprintf(tabfn, "%s", basefn);
1404                 sprintf(tabfn + strlen(basefn) - strlen(ftp2ext(efXVG)) - 1, "_%s%d.%s",
1405                         tabext, i, ftp2ext(efXVG));
1406                 tab[i] = make_bonded_table(fplog, tabfn, NRAL(ftype1)-2);
1407             }
1408         }
1409         sfree(count);
1410     }
1411
1412     return tab;
1413 }
1414
1415 void forcerec_set_ranges(t_forcerec *fr,
1416                          int ncg_home, int ncg_force,
1417                          int natoms_force,
1418                          int natoms_force_constr, int natoms_f_novirsum)
1419 {
1420     fr->cg0 = 0;
1421     fr->hcg = ncg_home;
1422
1423     /* fr->ncg_force is unused in the standard code,
1424      * but it can be useful for modified code dealing with charge groups.
1425      */
1426     fr->ncg_force           = ncg_force;
1427     fr->natoms_force        = natoms_force;
1428     fr->natoms_force_constr = natoms_force_constr;
1429
1430     if (fr->natoms_force_constr > fr->nalloc_force)
1431     {
1432         fr->nalloc_force = over_alloc_dd(fr->natoms_force_constr);
1433
1434         if (fr->bTwinRange)
1435         {
1436             srenew(fr->f_twin, fr->nalloc_force);
1437         }
1438     }
1439
1440     if (fr->bF_NoVirSum)
1441     {
1442         fr->f_novirsum_n = natoms_f_novirsum;
1443         if (fr->f_novirsum_n > fr->f_novirsum_nalloc)
1444         {
1445             fr->f_novirsum_nalloc = over_alloc_dd(fr->f_novirsum_n);
1446             srenew(fr->f_novirsum_alloc, fr->f_novirsum_nalloc);
1447         }
1448     }
1449     else
1450     {
1451         fr->f_novirsum_n = 0;
1452     }
1453 }
1454
1455 static real cutoff_inf(real cutoff)
1456 {
1457     if (cutoff == 0)
1458     {
1459         cutoff = GMX_CUTOFF_INF;
1460     }
1461
1462     return cutoff;
1463 }
1464
1465 static void make_adress_tf_tables(FILE *fp, const output_env_t oenv,
1466                                   t_forcerec *fr, const t_inputrec *ir,
1467                                   const char *tabfn, const gmx_mtop_t *mtop,
1468                                   matrix     box)
1469 {
1470     char buf[STRLEN];
1471     int  i, j;
1472
1473     if (tabfn == NULL)
1474     {
1475         gmx_fatal(FARGS, "No thermoforce table file given. Use -tabletf to specify a file\n");
1476         return;
1477     }
1478
1479     snew(fr->atf_tabs, ir->adress->n_tf_grps);
1480
1481     sprintf(buf, "%s", tabfn);
1482     for (i = 0; i < ir->adress->n_tf_grps; i++)
1483     {
1484         j = ir->adress->tf_table_index[i]; /* get energy group index */
1485         sprintf(buf + strlen(tabfn) - strlen(ftp2ext(efXVG)) - 1, "tf_%s.%s",
1486                 *(mtop->groups.grpname[mtop->groups.grps[egcENER].nm_ind[j]]), ftp2ext(efXVG));
1487         if (fp)
1488         {
1489             fprintf(fp, "loading tf table for energygrp index %d from %s\n", ir->adress->tf_table_index[i], buf);
1490         }
1491         fr->atf_tabs[i] = make_atf_table(fp, oenv, fr, buf, box);
1492     }
1493
1494 }
1495
1496 gmx_bool can_use_allvsall(const t_inputrec *ir, gmx_bool bPrintNote, t_commrec *cr, FILE *fp)
1497 {
1498     gmx_bool bAllvsAll;
1499
1500     bAllvsAll =
1501         (
1502             ir->rlist == 0            &&
1503             ir->rcoulomb == 0         &&
1504             ir->rvdw == 0             &&
1505             ir->ePBC == epbcNONE      &&
1506             ir->vdwtype == evdwCUT    &&
1507             ir->coulombtype == eelCUT &&
1508             ir->efep == efepNO        &&
1509             (ir->implicit_solvent == eisNO ||
1510              (ir->implicit_solvent == eisGBSA && (ir->gb_algorithm == egbSTILL ||
1511                                                   ir->gb_algorithm == egbHCT   ||
1512                                                   ir->gb_algorithm == egbOBC))) &&
1513             getenv("GMX_NO_ALLVSALL") == NULL
1514         );
1515
1516     if (bAllvsAll && ir->opts.ngener > 1)
1517     {
1518         const char *note = "NOTE: Can not use all-vs-all force loops, because there are multiple energy monitor groups; you might get significantly higher performance when using only a single energy monitor group.\n";
1519
1520         if (bPrintNote)
1521         {
1522             if (MASTER(cr))
1523             {
1524                 fprintf(stderr, "\n%s\n", note);
1525             }
1526             if (fp != NULL)
1527             {
1528                 fprintf(fp, "\n%s\n", note);
1529             }
1530         }
1531         bAllvsAll = FALSE;
1532     }
1533
1534     if (bAllvsAll && fp && MASTER(cr))
1535     {
1536         fprintf(fp, "\nUsing SIMD all-vs-all kernels.\n\n");
1537     }
1538
1539     return bAllvsAll;
1540 }
1541
1542
1543 gmx_bool nbnxn_acceleration_supported(FILE             *fplog,
1544                                       const t_commrec  *cr,
1545                                       const t_inputrec *ir,
1546                                       gmx_bool          bGPU)
1547 {
1548     if (!bGPU && (ir->vdwtype == evdwPME && ir->ljpme_combination_rule == eljpmeLB))
1549     {
1550         md_print_warn(cr, fplog, "LJ-PME with Lorentz-Berthelot is not supported with %s, falling back to %s\n",
1551                       bGPU ? "GPUs" : "SIMD kernels",
1552                       bGPU ? "CPU only" : "plain-C kernels");
1553         return FALSE;
1554     }
1555
1556     return TRUE;
1557 }
1558
1559
1560 static void pick_nbnxn_kernel_cpu(const t_inputrec gmx_unused *ir,
1561                                   int                         *kernel_type,
1562                                   int                         *ewald_excl)
1563 {
1564     *kernel_type = nbnxnk4x4_PlainC;
1565     *ewald_excl  = ewaldexclTable;
1566
1567 #ifdef GMX_NBNXN_SIMD
1568     {
1569 #ifdef GMX_NBNXN_SIMD_4XN
1570         *kernel_type = nbnxnk4xN_SIMD_4xN;
1571 #endif
1572 #ifdef GMX_NBNXN_SIMD_2XNN
1573         *kernel_type = nbnxnk4xN_SIMD_2xNN;
1574 #endif
1575
1576 #if defined GMX_NBNXN_SIMD_2XNN && defined GMX_NBNXN_SIMD_4XN
1577         /* We need to choose if we want 2x(N+N) or 4xN kernels.
1578          * Currently this is based on the SIMD acceleration choice,
1579          * but it might be better to decide this at runtime based on CPU.
1580          *
1581          * 4xN calculates more (zero) interactions, but has less pair-search
1582          * work and much better kernel instruction scheduling.
1583          *
1584          * Up till now we have only seen that on Intel Sandy/Ivy Bridge,
1585          * which doesn't have FMA, both the analytical and tabulated Ewald
1586          * kernels have similar pair rates for 4x8 and 2x(4+4), so we choose
1587          * 2x(4+4) because it results in significantly fewer pairs.
1588          * For RF, the raw pair rate of the 4x8 kernel is higher than 2x(4+4),
1589          * 10% with HT, 50% without HT. As we currently don't detect the actual
1590          * use of HT, use 4x8 to avoid a potential performance hit.
1591          * On Intel Haswell 4x8 is always faster.
1592          */
1593         *kernel_type = nbnxnk4xN_SIMD_4xN;
1594
1595 #ifndef GMX_SIMD_HAVE_FMA
1596         if (EEL_PME_EWALD(ir->coulombtype) ||
1597             EVDW_PME(ir->vdwtype))
1598         {
1599             /* We have Ewald kernels without FMA (Intel Sandy/Ivy Bridge).
1600              * There are enough instructions to make 2x(4+4) efficient.
1601              */
1602             *kernel_type = nbnxnk4xN_SIMD_2xNN;
1603         }
1604 #endif
1605 #endif  /* GMX_NBNXN_SIMD_2XNN && GMX_NBNXN_SIMD_4XN */
1606
1607
1608         if (getenv("GMX_NBNXN_SIMD_4XN") != NULL)
1609         {
1610 #ifdef GMX_NBNXN_SIMD_4XN
1611             *kernel_type = nbnxnk4xN_SIMD_4xN;
1612 #else
1613             gmx_fatal(FARGS, "SIMD 4xN kernels requested, but GROMACS has been compiled without support for these kernels");
1614 #endif
1615         }
1616         if (getenv("GMX_NBNXN_SIMD_2XNN") != NULL)
1617         {
1618 #ifdef GMX_NBNXN_SIMD_2XNN
1619             *kernel_type = nbnxnk4xN_SIMD_2xNN;
1620 #else
1621             gmx_fatal(FARGS, "SIMD 2x(N+N) kernels requested, but GROMACS has been compiled without support for these kernels");
1622 #endif
1623         }
1624
1625         /* Analytical Ewald exclusion correction is only an option in
1626          * the SIMD kernel.
1627          * Since table lookup's don't parallelize with SIMD, analytical
1628          * will probably always be faster for a SIMD width of 8 or more.
1629          * With FMA analytical is sometimes faster for a width if 4 as well.
1630          * On BlueGene/Q, this is faster regardless of precision.
1631          * In single precision, this is faster on Bulldozer.
1632          */
1633 #if GMX_SIMD_REAL_WIDTH >= 8 || \
1634         (GMX_SIMD_REAL_WIDTH >= 4 && defined GMX_SIMD_HAVE_FMA && !defined GMX_DOUBLE) || \
1635         defined GMX_SIMD_IBM_QPX
1636         *ewald_excl = ewaldexclAnalytical;
1637 #endif
1638         if (getenv("GMX_NBNXN_EWALD_TABLE") != NULL)
1639         {
1640             *ewald_excl = ewaldexclTable;
1641         }
1642         if (getenv("GMX_NBNXN_EWALD_ANALYTICAL") != NULL)
1643         {
1644             *ewald_excl = ewaldexclAnalytical;
1645         }
1646
1647     }
1648 #endif /* GMX_NBNXN_SIMD */
1649 }
1650
1651
1652 const char *lookup_nbnxn_kernel_name(int kernel_type)
1653 {
1654     const char *returnvalue = NULL;
1655     switch (kernel_type)
1656     {
1657         case nbnxnkNotSet:
1658             returnvalue = "not set";
1659             break;
1660         case nbnxnk4x4_PlainC:
1661             returnvalue = "plain C";
1662             break;
1663         case nbnxnk4xN_SIMD_4xN:
1664         case nbnxnk4xN_SIMD_2xNN:
1665 #ifdef GMX_NBNXN_SIMD
1666 #if defined GMX_SIMD_X86_SSE2
1667             returnvalue = "SSE2";
1668 #elif defined GMX_SIMD_X86_SSE4_1
1669             returnvalue = "SSE4.1";
1670 #elif defined GMX_SIMD_X86_AVX_128_FMA
1671             returnvalue = "AVX_128_FMA";
1672 #elif defined GMX_SIMD_X86_AVX_256
1673             returnvalue = "AVX_256";
1674 #elif defined GMX_SIMD_X86_AVX2_256
1675             returnvalue = "AVX2_256";
1676 #else
1677             returnvalue = "SIMD";
1678 #endif
1679 #else  /* GMX_NBNXN_SIMD */
1680             returnvalue = "not available";
1681 #endif /* GMX_NBNXN_SIMD */
1682             break;
1683         case nbnxnk8x8x8_GPU: returnvalue    = "GPU"; break;
1684         case nbnxnk8x8x8_PlainC: returnvalue = "plain C"; break;
1685
1686         case nbnxnkNR:
1687         default:
1688             gmx_fatal(FARGS, "Illegal kernel type selected");
1689             returnvalue = NULL;
1690             break;
1691     }
1692     return returnvalue;
1693 };
1694
1695 static void pick_nbnxn_kernel(FILE                *fp,
1696                               const t_commrec     *cr,
1697                               gmx_bool             use_simd_kernels,
1698                               gmx_bool             bUseGPU,
1699                               gmx_bool             bEmulateGPU,
1700                               const t_inputrec    *ir,
1701                               int                 *kernel_type,
1702                               int                 *ewald_excl,
1703                               gmx_bool             bDoNonbonded)
1704 {
1705     assert(kernel_type);
1706
1707     *kernel_type = nbnxnkNotSet;
1708     *ewald_excl  = ewaldexclTable;
1709
1710     if (bEmulateGPU)
1711     {
1712         *kernel_type = nbnxnk8x8x8_PlainC;
1713
1714         if (bDoNonbonded)
1715         {
1716             md_print_warn(cr, fp, "Emulating a GPU run on the CPU (slow)");
1717         }
1718     }
1719     else if (bUseGPU)
1720     {
1721         *kernel_type = nbnxnk8x8x8_GPU;
1722     }
1723
1724     if (*kernel_type == nbnxnkNotSet)
1725     {
1726         /* LJ PME with LB combination rule does 7 mesh operations.
1727          * This so slow that we don't compile SIMD non-bonded kernels for that.
1728          */
1729         if (use_simd_kernels &&
1730             nbnxn_acceleration_supported(fp, cr, ir, FALSE))
1731         {
1732             pick_nbnxn_kernel_cpu(ir, kernel_type, ewald_excl);
1733         }
1734         else
1735         {
1736             *kernel_type = nbnxnk4x4_PlainC;
1737         }
1738     }
1739
1740     if (bDoNonbonded && fp != NULL)
1741     {
1742         fprintf(fp, "\nUsing %s %dx%d non-bonded kernels\n\n",
1743                 lookup_nbnxn_kernel_name(*kernel_type),
1744                 nbnxn_kernel_to_ci_size(*kernel_type),
1745                 nbnxn_kernel_to_cj_size(*kernel_type));
1746
1747         if (nbnxnk4x4_PlainC == *kernel_type ||
1748             nbnxnk8x8x8_PlainC == *kernel_type)
1749         {
1750             md_print_warn(cr, fp,
1751                           "WARNING: Using the slow %s kernels. This should\n"
1752                           "not happen during routine usage on supported platforms.\n\n",
1753                           lookup_nbnxn_kernel_name(*kernel_type));
1754         }
1755     }
1756 }
1757
1758 static void pick_nbnxn_resources(FILE                *fp,
1759                                  const t_commrec     *cr,
1760                                  const gmx_hw_info_t *hwinfo,
1761                                  gmx_bool             bDoNonbonded,
1762                                  gmx_bool            *bUseGPU,
1763                                  gmx_bool            *bEmulateGPU,
1764                                  const gmx_gpu_opt_t *gpu_opt)
1765 {
1766     gmx_bool bEmulateGPUEnvVarSet;
1767     char     gpu_err_str[STRLEN];
1768
1769     *bUseGPU = FALSE;
1770
1771     bEmulateGPUEnvVarSet = (getenv("GMX_EMULATE_GPU") != NULL);
1772
1773     /* Run GPU emulation mode if GMX_EMULATE_GPU is defined. Because
1774      * GPUs (currently) only handle non-bonded calculations, we will
1775      * automatically switch to emulation if non-bonded calculations are
1776      * turned off via GMX_NO_NONBONDED - this is the simple and elegant
1777      * way to turn off GPU initialization, data movement, and cleanup.
1778      *
1779      * GPU emulation can be useful to assess the performance one can expect by
1780      * adding GPU(s) to the machine. The conditional below allows this even
1781      * if mdrun is compiled without GPU acceleration support.
1782      * Note that you should freezing the system as otherwise it will explode.
1783      */
1784     *bEmulateGPU = (bEmulateGPUEnvVarSet ||
1785                     (!bDoNonbonded && gpu_opt->n_dev_use > 0));
1786
1787     /* Enable GPU mode when GPUs are available or no GPU emulation is requested.
1788      */
1789     if (gpu_opt->n_dev_use > 0 && !(*bEmulateGPU))
1790     {
1791         /* Each PP node will use the intra-node id-th device from the
1792          * list of detected/selected GPUs. */
1793         if (!init_gpu(fp, cr->rank_pp_intranode, gpu_err_str,
1794                       &hwinfo->gpu_info, gpu_opt))
1795         {
1796             /* At this point the init should never fail as we made sure that
1797              * we have all the GPUs we need. If it still does, we'll bail. */
1798             /* TODO the decorating of gpu_err_str is nicer if it
1799                happens inside init_gpu. Out here, the decorating with
1800                the MPI rank makes sense. */
1801             gmx_fatal(FARGS, "On rank %d failed to initialize GPU #%d: %s",
1802                       cr->nodeid,
1803                       get_cuda_gpu_device_id(&hwinfo->gpu_info, gpu_opt,
1804                                              cr->rank_pp_intranode),
1805                       gpu_err_str);
1806         }
1807
1808         /* Here we actually turn on hardware GPU acceleration */
1809         *bUseGPU = TRUE;
1810     }
1811 }
1812
1813 gmx_bool uses_simple_tables(int                 cutoff_scheme,
1814                             nonbonded_verlet_t *nbv,
1815                             int                 group)
1816 {
1817     gmx_bool bUsesSimpleTables = TRUE;
1818     int      grp_index;
1819
1820     switch (cutoff_scheme)
1821     {
1822         case ecutsGROUP:
1823             bUsesSimpleTables = TRUE;
1824             break;
1825         case ecutsVERLET:
1826             assert(NULL != nbv && NULL != nbv->grp);
1827             grp_index         = (group < 0) ? 0 : (nbv->ngrp - 1);
1828             bUsesSimpleTables = nbnxn_kernel_pairlist_simple(nbv->grp[grp_index].kernel_type);
1829             break;
1830         default:
1831             gmx_incons("unimplemented");
1832     }
1833     return bUsesSimpleTables;
1834 }
1835
1836 static void init_ewald_f_table(interaction_const_t *ic,
1837                                real                 rtab)
1838 {
1839     real maxr;
1840
1841     /* Get the Ewald table spacing based on Coulomb and/or LJ
1842      * Ewald coefficients and rtol.
1843      */
1844     ic->tabq_scale = ewald_spline3_table_scale(ic);
1845
1846     if (ic->cutoff_scheme == ecutsVERLET)
1847     {
1848         maxr = ic->rcoulomb;
1849     }
1850     else
1851     {
1852         maxr = std::max(ic->rcoulomb, rtab);
1853     }
1854     ic->tabq_size  = static_cast<int>(maxr*ic->tabq_scale) + 2;
1855
1856     sfree_aligned(ic->tabq_coul_FDV0);
1857     sfree_aligned(ic->tabq_coul_F);
1858     sfree_aligned(ic->tabq_coul_V);
1859
1860     sfree_aligned(ic->tabq_vdw_FDV0);
1861     sfree_aligned(ic->tabq_vdw_F);
1862     sfree_aligned(ic->tabq_vdw_V);
1863
1864     if (ic->eeltype == eelEWALD || EEL_PME(ic->eeltype))
1865     {
1866         /* Create the original table data in FDV0 */
1867         snew_aligned(ic->tabq_coul_FDV0, ic->tabq_size*4, 32);
1868         snew_aligned(ic->tabq_coul_F, ic->tabq_size, 32);
1869         snew_aligned(ic->tabq_coul_V, ic->tabq_size, 32);
1870         table_spline3_fill_ewald_lr(ic->tabq_coul_F, ic->tabq_coul_V, ic->tabq_coul_FDV0,
1871                                     ic->tabq_size, 1/ic->tabq_scale, ic->ewaldcoeff_q, v_q_ewald_lr);
1872     }
1873
1874     if (EVDW_PME(ic->vdwtype))
1875     {
1876         snew_aligned(ic->tabq_vdw_FDV0, ic->tabq_size*4, 32);
1877         snew_aligned(ic->tabq_vdw_F, ic->tabq_size, 32);
1878         snew_aligned(ic->tabq_vdw_V, ic->tabq_size, 32);
1879         table_spline3_fill_ewald_lr(ic->tabq_vdw_F, ic->tabq_vdw_V, ic->tabq_vdw_FDV0,
1880                                     ic->tabq_size, 1/ic->tabq_scale, ic->ewaldcoeff_lj, v_lj_ewald_lr);
1881     }
1882 }
1883
1884 void init_interaction_const_tables(FILE                *fp,
1885                                    interaction_const_t *ic,
1886                                    real                 rtab)
1887 {
1888     if (ic->eeltype == eelEWALD || EEL_PME(ic->eeltype) || EVDW_PME(ic->vdwtype))
1889     {
1890         init_ewald_f_table(ic, rtab);
1891
1892         if (fp != NULL)
1893         {
1894             fprintf(fp, "Initialized non-bonded Ewald correction tables, spacing: %.2e size: %d\n\n",
1895                     1/ic->tabq_scale, ic->tabq_size);
1896         }
1897     }
1898 }
1899
1900 static void clear_force_switch_constants(shift_consts_t *sc)
1901 {
1902     sc->c2   = 0;
1903     sc->c3   = 0;
1904     sc->cpot = 0;
1905 }
1906
1907 static void force_switch_constants(real p,
1908                                    real rsw, real rc,
1909                                    shift_consts_t *sc)
1910 {
1911     /* Here we determine the coefficient for shifting the force to zero
1912      * between distance rsw and the cut-off rc.
1913      * For a potential of r^-p, we have force p*r^-(p+1).
1914      * But to save flops we absorb p in the coefficient.
1915      * Thus we get:
1916      * force/p   = r^-(p+1) + c2*r^2 + c3*r^3
1917      * potential = r^-p + c2/3*r^3 + c3/4*r^4 + cpot
1918      */
1919     sc->c2   =  ((p + 1)*rsw - (p + 4)*rc)/(pow(rc, p + 2)*pow(rc - rsw, 2));
1920     sc->c3   = -((p + 1)*rsw - (p + 3)*rc)/(pow(rc, p + 2)*pow(rc - rsw, 3));
1921     sc->cpot = -pow(rc, -p) + p*sc->c2/3*pow(rc - rsw, 3) + p*sc->c3/4*pow(rc - rsw, 4);
1922 }
1923
1924 static void potential_switch_constants(real rsw, real rc,
1925                                        switch_consts_t *sc)
1926 {
1927     /* The switch function is 1 at rsw and 0 at rc.
1928      * The derivative and second derivate are zero at both ends.
1929      * rsw        = max(r - r_switch, 0)
1930      * sw         = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
1931      * dsw        = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
1932      * force      = force*dsw - potential*sw
1933      * potential *= sw
1934      */
1935     sc->c3 = -10*pow(rc - rsw, -3);
1936     sc->c4 =  15*pow(rc - rsw, -4);
1937     sc->c5 =  -6*pow(rc - rsw, -5);
1938 }
1939
1940 /*! \brief Construct interaction constants
1941  *
1942  * This data is used (particularly) by search and force code for
1943  * short-range interactions. Many of these are constant for the whole
1944  * simulation; some are constant only after PME tuning completes.
1945  */
1946 static void
1947 init_interaction_const(FILE                       *fp,
1948                        interaction_const_t       **interaction_const,
1949                        const t_forcerec           *fr)
1950 {
1951     interaction_const_t *ic;
1952     const real           minusSix          = -6.0;
1953     const real           minusTwelve       = -12.0;
1954
1955     snew(ic, 1);
1956
1957     ic->cutoff_scheme   = fr->cutoff_scheme;
1958
1959     /* Just allocate something so we can free it */
1960     snew_aligned(ic->tabq_coul_FDV0, 16, 32);
1961     snew_aligned(ic->tabq_coul_F, 16, 32);
1962     snew_aligned(ic->tabq_coul_V, 16, 32);
1963
1964     ic->rlist           = fr->rlist;
1965     ic->rlistlong       = fr->rlistlong;
1966
1967     /* Lennard-Jones */
1968     ic->vdwtype         = fr->vdwtype;
1969     ic->vdw_modifier    = fr->vdw_modifier;
1970     ic->rvdw            = fr->rvdw;
1971     ic->rvdw_switch     = fr->rvdw_switch;
1972     ic->ewaldcoeff_lj   = fr->ewaldcoeff_lj;
1973     ic->ljpme_comb_rule = fr->ljpme_combination_rule;
1974     ic->sh_lj_ewald     = 0;
1975     clear_force_switch_constants(&ic->dispersion_shift);
1976     clear_force_switch_constants(&ic->repulsion_shift);
1977
1978     switch (ic->vdw_modifier)
1979     {
1980         case eintmodPOTSHIFT:
1981             /* Only shift the potential, don't touch the force */
1982             ic->dispersion_shift.cpot = -pow(ic->rvdw, minusSix);
1983             ic->repulsion_shift.cpot  = -pow(ic->rvdw, minusTwelve);
1984             if (EVDW_PME(ic->vdwtype))
1985             {
1986                 real crc2;
1987
1988                 crc2            = sqr(ic->ewaldcoeff_lj*ic->rvdw);
1989                 ic->sh_lj_ewald = (exp(-crc2)*(1 + crc2 + 0.5*crc2*crc2) - 1)*pow(ic->rvdw, minusSix);
1990             }
1991             break;
1992         case eintmodFORCESWITCH:
1993             /* Switch the force, switch and shift the potential */
1994             force_switch_constants(6.0, ic->rvdw_switch, ic->rvdw,
1995                                    &ic->dispersion_shift);
1996             force_switch_constants(12.0, ic->rvdw_switch, ic->rvdw,
1997                                    &ic->repulsion_shift);
1998             break;
1999         case eintmodPOTSWITCH:
2000             /* Switch the potential and force */
2001             potential_switch_constants(ic->rvdw_switch, ic->rvdw,
2002                                        &ic->vdw_switch);
2003             break;
2004         case eintmodNONE:
2005         case eintmodEXACTCUTOFF:
2006             /* Nothing to do here */
2007             break;
2008         default:
2009             gmx_incons("unimplemented potential modifier");
2010     }
2011
2012     ic->sh_invrc6 = -ic->dispersion_shift.cpot;
2013
2014     /* Electrostatics */
2015     ic->eeltype          = fr->eeltype;
2016     ic->coulomb_modifier = fr->coulomb_modifier;
2017     ic->rcoulomb         = fr->rcoulomb;
2018     ic->epsilon_r        = fr->epsilon_r;
2019     ic->epsfac           = fr->epsfac;
2020     ic->ewaldcoeff_q     = fr->ewaldcoeff_q;
2021
2022     if (fr->coulomb_modifier == eintmodPOTSHIFT)
2023     {
2024         ic->sh_ewald = gmx_erfc(ic->ewaldcoeff_q*ic->rcoulomb);
2025     }
2026     else
2027     {
2028         ic->sh_ewald = 0;
2029     }
2030
2031     /* Reaction-field */
2032     if (EEL_RF(ic->eeltype))
2033     {
2034         ic->epsilon_rf = fr->epsilon_rf;
2035         ic->k_rf       = fr->k_rf;
2036         ic->c_rf       = fr->c_rf;
2037     }
2038     else
2039     {
2040         /* For plain cut-off we might use the reaction-field kernels */
2041         ic->epsilon_rf = ic->epsilon_r;
2042         ic->k_rf       = 0;
2043         if (fr->coulomb_modifier == eintmodPOTSHIFT)
2044         {
2045             ic->c_rf   = 1/ic->rcoulomb;
2046         }
2047         else
2048         {
2049             ic->c_rf   = 0;
2050         }
2051     }
2052
2053     if (fp != NULL)
2054     {
2055         real dispersion_shift;
2056
2057         dispersion_shift = ic->dispersion_shift.cpot;
2058         if (EVDW_PME(ic->vdwtype))
2059         {
2060             dispersion_shift -= ic->sh_lj_ewald;
2061         }
2062         fprintf(fp, "Potential shift: LJ r^-12: %.3e r^-6: %.3e",
2063                 ic->repulsion_shift.cpot, dispersion_shift);
2064
2065         if (ic->eeltype == eelCUT)
2066         {
2067             fprintf(fp, ", Coulomb %.e", -ic->c_rf);
2068         }
2069         else if (EEL_PME(ic->eeltype))
2070         {
2071             fprintf(fp, ", Ewald %.3e", -ic->sh_ewald);
2072         }
2073         fprintf(fp, "\n");
2074     }
2075
2076     *interaction_const = ic;
2077 }
2078
2079 /*! \brief Manage initialization within the NBNXN module of
2080  * run-time constants.
2081  */
2082 static void
2083 initialize_gpu_constants(const t_commrec gmx_unused      *cr,
2084                          interaction_const_t             *interaction_const,
2085                          const struct nonbonded_verlet_t *nbv)
2086 {
2087     if (nbv != NULL && nbv->bUseGPU)
2088     {
2089         nbnxn_gpu_init_const(nbv->gpu_nbv, interaction_const, nbv->grp);
2090
2091         /* With tMPI + GPUs some ranks may be sharing GPU(s) and therefore
2092          * also sharing texture references. To keep the code simple, we don't
2093          * treat texture references as shared resources, but this means that
2094          * the coulomb_tab and nbfp texture refs will get updated by multiple threads.
2095          * Hence, to ensure that the non-bonded kernels don't start before all
2096          * texture binding operations are finished, we need to wait for all ranks
2097          * to arrive here before continuing.
2098          *
2099          * Note that we could omit this barrier if GPUs are not shared (or
2100          * texture objects are used), but as this is initialization code, there
2101          * is no point in complicating things.
2102          */
2103 #ifdef GMX_THREAD_MPI
2104         if (PAR(cr))
2105         {
2106             gmx_barrier(cr);
2107         }
2108 #endif  /* GMX_THREAD_MPI */
2109     }
2110
2111 }
2112
2113 static void init_nb_verlet(FILE                *fp,
2114                            nonbonded_verlet_t **nb_verlet,
2115                            gmx_bool             bFEP_NonBonded,
2116                            const t_inputrec    *ir,
2117                            const t_forcerec    *fr,
2118                            const t_commrec     *cr,
2119                            const char          *nbpu_opt)
2120 {
2121     nonbonded_verlet_t *nbv;
2122     int                 i;
2123     char               *env;
2124     gmx_bool            bEmulateGPU, bHybridGPURun = FALSE;
2125
2126     nbnxn_alloc_t      *nb_alloc;
2127     nbnxn_free_t       *nb_free;
2128
2129     snew(nbv, 1);
2130
2131     pick_nbnxn_resources(fp, cr, fr->hwinfo,
2132                          fr->bNonbonded,
2133                          &nbv->bUseGPU,
2134                          &bEmulateGPU,
2135                          fr->gpu_opt);
2136
2137     nbv->nbs = NULL;
2138
2139     nbv->ngrp = (DOMAINDECOMP(cr) ? 2 : 1);
2140     for (i = 0; i < nbv->ngrp; i++)
2141     {
2142         nbv->grp[i].nbl_lists.nnbl = 0;
2143         nbv->grp[i].nbat           = NULL;
2144         nbv->grp[i].kernel_type    = nbnxnkNotSet;
2145
2146         if (i == 0) /* local */
2147         {
2148             pick_nbnxn_kernel(fp, cr, fr->use_simd_kernels,
2149                               nbv->bUseGPU, bEmulateGPU, ir,
2150                               &nbv->grp[i].kernel_type,
2151                               &nbv->grp[i].ewald_excl,
2152                               fr->bNonbonded);
2153         }
2154         else /* non-local */
2155         {
2156             if (nbpu_opt != NULL && strcmp(nbpu_opt, "gpu_cpu") == 0)
2157             {
2158                 /* Use GPU for local, select a CPU kernel for non-local */
2159                 pick_nbnxn_kernel(fp, cr, fr->use_simd_kernels,
2160                                   FALSE, FALSE, ir,
2161                                   &nbv->grp[i].kernel_type,
2162                                   &nbv->grp[i].ewald_excl,
2163                                   fr->bNonbonded);
2164
2165                 bHybridGPURun = TRUE;
2166             }
2167             else
2168             {
2169                 /* Use the same kernel for local and non-local interactions */
2170                 nbv->grp[i].kernel_type = nbv->grp[0].kernel_type;
2171                 nbv->grp[i].ewald_excl  = nbv->grp[0].ewald_excl;
2172             }
2173         }
2174     }
2175
2176     if (nbv->bUseGPU)
2177     {
2178         nbnxn_gpu_compile_kernels(cr->rank_pp_intranode, cr->nodeid, &fr->hwinfo->gpu_info, fr->gpu_opt, fr->ic);
2179
2180         /* init the NxN GPU data; the last argument tells whether we'll have
2181          * both local and non-local NB calculation on GPU */
2182         nbnxn_gpu_init(fp, &nbv->gpu_nbv,
2183                        &fr->hwinfo->gpu_info, fr->gpu_opt,
2184                        cr->rank_pp_intranode,
2185                        (nbv->ngrp > 1) && !bHybridGPURun);
2186
2187         if ((env = getenv("GMX_NB_MIN_CI")) != NULL)
2188         {
2189             char *end;
2190
2191             nbv->min_ci_balanced = strtol(env, &end, 10);
2192             if (!end || (*end != 0) || nbv->min_ci_balanced <= 0)
2193             {
2194                 gmx_fatal(FARGS, "Invalid value passed in GMX_NB_MIN_CI=%s, positive integer required", env);
2195             }
2196
2197             if (debug)
2198             {
2199                 fprintf(debug, "Neighbor-list balancing parameter: %d (passed as env. var.)\n",
2200                         nbv->min_ci_balanced);
2201             }
2202         }
2203         else
2204         {
2205             nbv->min_ci_balanced = nbnxn_gpu_min_ci_balanced(nbv->gpu_nbv);
2206             if (debug)
2207             {
2208                 fprintf(debug, "Neighbor-list balancing parameter: %d (auto-adjusted to the number of GPU multi-processors)\n",
2209                         nbv->min_ci_balanced);
2210             }
2211         }
2212     }
2213     else
2214     {
2215         nbv->min_ci_balanced = 0;
2216     }
2217
2218     *nb_verlet = nbv;
2219
2220     nbnxn_init_search(&nbv->nbs,
2221                       DOMAINDECOMP(cr) ? &cr->dd->nc : NULL,
2222                       DOMAINDECOMP(cr) ? domdec_zones(cr->dd) : NULL,
2223                       bFEP_NonBonded,
2224                       gmx_omp_nthreads_get(emntPairsearch));
2225
2226     for (i = 0; i < nbv->ngrp; i++)
2227     {
2228         gpu_set_host_malloc_and_free(nbv->grp[0].kernel_type == nbnxnk8x8x8_GPU,
2229                                      &nb_alloc, &nb_free);
2230
2231         nbnxn_init_pairlist_set(&nbv->grp[i].nbl_lists,
2232                                 nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type),
2233                                 /* 8x8x8 "non-simple" lists are ATM always combined */
2234                                 !nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type),
2235                                 nb_alloc, nb_free);
2236
2237         if (i == 0 ||
2238             nbv->grp[0].kernel_type != nbv->grp[i].kernel_type)
2239         {
2240             gmx_bool bSimpleList;
2241             int      enbnxninitcombrule;
2242
2243             bSimpleList = nbnxn_kernel_pairlist_simple(nbv->grp[i].kernel_type);
2244
2245             if (bSimpleList && (fr->vdwtype == evdwCUT && (fr->vdw_modifier == eintmodNONE || fr->vdw_modifier == eintmodPOTSHIFT)))
2246             {
2247                 /* Plain LJ cut-off: we can optimize with combination rules */
2248                 enbnxninitcombrule = enbnxninitcombruleDETECT;
2249             }
2250             else if (fr->vdwtype == evdwPME)
2251             {
2252                 /* LJ-PME: we need to use a combination rule for the grid */
2253                 if (fr->ljpme_combination_rule == eljpmeGEOM)
2254                 {
2255                     enbnxninitcombrule = enbnxninitcombruleGEOM;
2256                 }
2257                 else
2258                 {
2259                     enbnxninitcombrule = enbnxninitcombruleLB;
2260                 }
2261             }
2262             else
2263             {
2264                 /* We use a full combination matrix: no rule required */
2265                 enbnxninitcombrule = enbnxninitcombruleNONE;
2266             }
2267
2268
2269             snew(nbv->grp[i].nbat, 1);
2270             nbnxn_atomdata_init(fp,
2271                                 nbv->grp[i].nbat,
2272                                 nbv->grp[i].kernel_type,
2273                                 enbnxninitcombrule,
2274                                 fr->ntype, fr->nbfp,
2275                                 ir->opts.ngener,
2276                                 bSimpleList ? gmx_omp_nthreads_get(emntNonbonded) : 1,
2277                                 nb_alloc, nb_free);
2278         }
2279         else
2280         {
2281             nbv->grp[i].nbat = nbv->grp[0].nbat;
2282         }
2283     }
2284 }
2285
2286 gmx_bool usingGpu(nonbonded_verlet_t *nbv)
2287 {
2288     return nbv != NULL && nbv->bUseGPU;
2289 }
2290
2291 void init_forcerec(FILE              *fp,
2292                    const output_env_t oenv,
2293                    t_forcerec        *fr,
2294                    t_fcdata          *fcd,
2295                    const t_inputrec  *ir,
2296                    const gmx_mtop_t  *mtop,
2297                    const t_commrec   *cr,
2298                    matrix             box,
2299                    const char        *tabfn,
2300                    const char        *tabafn,
2301                    const char        *tabpfn,
2302                    const char        *tabbfn,
2303                    const char        *nbpu_opt,
2304                    gmx_bool           bNoSolvOpt,
2305                    real               print_force)
2306 {
2307     int            i, m, negp_pp, negptable, egi, egj;
2308     real           rtab;
2309     char          *env;
2310     double         dbl;
2311     const t_block *cgs;
2312     gmx_bool       bGenericKernelOnly;
2313     gmx_bool       bMakeTables, bMakeSeparate14Table, bSomeNormalNbListsAreInUse;
2314     gmx_bool       bFEP_NonBonded;
2315     int           *nm_ind, egp_flags;
2316
2317     if (fr->hwinfo == NULL)
2318     {
2319         /* Detect hardware, gather information.
2320          * In mdrun, hwinfo has already been set before calling init_forcerec.
2321          * Here we ignore GPUs, as tools will not use them anyhow.
2322          */
2323         fr->hwinfo = gmx_detect_hardware(fp, cr, FALSE);
2324     }
2325
2326     /* By default we turn SIMD kernels on, but it might be turned off further down... */
2327     fr->use_simd_kernels = TRUE;
2328
2329     fr->bDomDec = DOMAINDECOMP(cr);
2330
2331     if (check_box(ir->ePBC, box))
2332     {
2333         gmx_fatal(FARGS, check_box(ir->ePBC, box));
2334     }
2335
2336     /* Test particle insertion ? */
2337     if (EI_TPI(ir->eI))
2338     {
2339         /* Set to the size of the molecule to be inserted (the last one) */
2340         /* Because of old style topologies, we have to use the last cg
2341          * instead of the last molecule type.
2342          */
2343         cgs       = &mtop->moltype[mtop->molblock[mtop->nmolblock-1].type].cgs;
2344         fr->n_tpi = cgs->index[cgs->nr] - cgs->index[cgs->nr-1];
2345         if (fr->n_tpi != mtop->mols.index[mtop->mols.nr] - mtop->mols.index[mtop->mols.nr-1])
2346         {
2347             gmx_fatal(FARGS, "The molecule to insert can not consist of multiple charge groups.\nMake it a single charge group.");
2348         }
2349     }
2350     else
2351     {
2352         fr->n_tpi = 0;
2353     }
2354
2355     /* Copy AdResS parameters */
2356     if (ir->bAdress)
2357     {
2358         fr->adress_type           = ir->adress->type;
2359         fr->adress_const_wf       = ir->adress->const_wf;
2360         fr->adress_ex_width       = ir->adress->ex_width;
2361         fr->adress_hy_width       = ir->adress->hy_width;
2362         fr->adress_icor           = ir->adress->icor;
2363         fr->adress_site           = ir->adress->site;
2364         fr->adress_ex_forcecap    = ir->adress->ex_forcecap;
2365         fr->adress_do_hybridpairs = ir->adress->do_hybridpairs;
2366
2367
2368         snew(fr->adress_group_explicit, ir->adress->n_energy_grps);
2369         for (i = 0; i < ir->adress->n_energy_grps; i++)
2370         {
2371             fr->adress_group_explicit[i] = ir->adress->group_explicit[i];
2372         }
2373
2374         fr->n_adress_tf_grps = ir->adress->n_tf_grps;
2375         snew(fr->adress_tf_table_index, fr->n_adress_tf_grps);
2376         for (i = 0; i < fr->n_adress_tf_grps; i++)
2377         {
2378             fr->adress_tf_table_index[i] = ir->adress->tf_table_index[i];
2379         }
2380         copy_rvec(ir->adress->refs, fr->adress_refs);
2381     }
2382     else
2383     {
2384         fr->adress_type           = eAdressOff;
2385         fr->adress_do_hybridpairs = FALSE;
2386     }
2387
2388     /* Copy the user determined parameters */
2389     fr->userint1  = ir->userint1;
2390     fr->userint2  = ir->userint2;
2391     fr->userint3  = ir->userint3;
2392     fr->userint4  = ir->userint4;
2393     fr->userreal1 = ir->userreal1;
2394     fr->userreal2 = ir->userreal2;
2395     fr->userreal3 = ir->userreal3;
2396     fr->userreal4 = ir->userreal4;
2397
2398     /* Shell stuff */
2399     fr->fc_stepsize = ir->fc_stepsize;
2400
2401     /* Free energy */
2402     fr->efep        = ir->efep;
2403     fr->sc_alphavdw = ir->fepvals->sc_alpha;
2404     if (ir->fepvals->bScCoul)
2405     {
2406         fr->sc_alphacoul  = ir->fepvals->sc_alpha;
2407         fr->sc_sigma6_min = pow(ir->fepvals->sc_sigma_min, 6);
2408     }
2409     else
2410     {
2411         fr->sc_alphacoul  = 0;
2412         fr->sc_sigma6_min = 0; /* only needed when bScCoul is on */
2413     }
2414     fr->sc_power      = ir->fepvals->sc_power;
2415     fr->sc_r_power    = ir->fepvals->sc_r_power;
2416     fr->sc_sigma6_def = pow(ir->fepvals->sc_sigma, 6);
2417
2418     env = getenv("GMX_SCSIGMA_MIN");
2419     if (env != NULL)
2420     {
2421         dbl = 0;
2422         sscanf(env, "%20lf", &dbl);
2423         fr->sc_sigma6_min = pow(dbl, 6);
2424         if (fp)
2425         {
2426             fprintf(fp, "Setting the minimum soft core sigma to %g nm\n", dbl);
2427         }
2428     }
2429
2430     fr->bNonbonded = TRUE;
2431     if (getenv("GMX_NO_NONBONDED") != NULL)
2432     {
2433         /* turn off non-bonded calculations */
2434         fr->bNonbonded = FALSE;
2435         md_print_warn(cr, fp,
2436                       "Found environment variable GMX_NO_NONBONDED.\n"
2437                       "Disabling nonbonded calculations.\n");
2438     }
2439
2440     bGenericKernelOnly = FALSE;
2441
2442     /* We now check in the NS code whether a particular combination of interactions
2443      * can be used with water optimization, and disable it if that is not the case.
2444      */
2445
2446     if (getenv("GMX_NB_GENERIC") != NULL)
2447     {
2448         if (fp != NULL)
2449         {
2450             fprintf(fp,
2451                     "Found environment variable GMX_NB_GENERIC.\n"
2452                     "Disabling all interaction-specific nonbonded kernels, will only\n"
2453                     "use the slow generic ones in src/gmxlib/nonbonded/nb_generic.c\n\n");
2454         }
2455         bGenericKernelOnly = TRUE;
2456     }
2457
2458     if (bGenericKernelOnly == TRUE)
2459     {
2460         bNoSolvOpt         = TRUE;
2461     }
2462
2463     if ( (getenv("GMX_DISABLE_SIMD_KERNELS") != NULL) || (getenv("GMX_NOOPTIMIZEDKERNELS") != NULL) )
2464     {
2465         fr->use_simd_kernels = FALSE;
2466         if (fp != NULL)
2467         {
2468             fprintf(fp,
2469                     "\nFound environment variable GMX_DISABLE_SIMD_KERNELS.\n"
2470                     "Disabling the usage of any SIMD-specific non-bonded & bonded kernel routines\n"
2471                     "(e.g. SSE2/SSE4.1/AVX).\n\n");
2472         }
2473     }
2474
2475     fr->bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
2476
2477     /* Check if we can/should do all-vs-all kernels */
2478     fr->bAllvsAll       = can_use_allvsall(ir, FALSE, NULL, NULL);
2479     fr->AllvsAll_work   = NULL;
2480     fr->AllvsAll_workgb = NULL;
2481
2482     /* All-vs-all kernels have not been implemented in 4.6 and later.
2483      * See Redmine #1249. */
2484     if (fr->bAllvsAll)
2485     {
2486         fr->bAllvsAll            = FALSE;
2487         if (fp != NULL)
2488         {
2489             fprintf(fp,
2490                     "\nYour simulation settings would have triggered the efficient all-vs-all\n"
2491                     "kernels in GROMACS 4.5, but these have not been implemented in GROMACS\n"
2492                     "4.6 and 5.x. If performance is important, please use GROMACS 4.5.7\n"
2493                     "or try cutoff-scheme = Verlet.\n\n");
2494         }
2495     }
2496
2497     /* Neighbour searching stuff */
2498     fr->cutoff_scheme = ir->cutoff_scheme;
2499     fr->bGrid         = (ir->ns_type == ensGRID);
2500     fr->ePBC          = ir->ePBC;
2501
2502     if (fr->cutoff_scheme == ecutsGROUP)
2503     {
2504         const char *note = "NOTE: This file uses the deprecated 'group' cutoff_scheme. This will be\n"
2505             "removed in a future release when 'verlet' supports all interaction forms.\n";
2506
2507         if (MASTER(cr))
2508         {
2509             fprintf(stderr, "\n%s\n", note);
2510         }
2511         if (fp != NULL)
2512         {
2513             fprintf(fp, "\n%s\n", note);
2514         }
2515     }
2516
2517     /* Determine if we will do PBC for distances in bonded interactions */
2518     if (fr->ePBC == epbcNONE)
2519     {
2520         fr->bMolPBC = FALSE;
2521     }
2522     else
2523     {
2524         if (!DOMAINDECOMP(cr))
2525         {
2526             gmx_bool bSHAKE;
2527
2528             bSHAKE = (ir->eConstrAlg == econtSHAKE &&
2529                       (gmx_mtop_ftype_count(mtop, F_CONSTR) > 0 ||
2530                        gmx_mtop_ftype_count(mtop, F_CONSTRNC) > 0));
2531
2532             /* The group cut-off scheme and SHAKE assume charge groups
2533              * are whole, but not using molpbc is faster in most cases.
2534              * With intermolecular interactions we need PBC for calculating
2535              * distances between atoms in different molecules.
2536              */
2537             if ((fr->cutoff_scheme == ecutsGROUP || bSHAKE) &&
2538                 !mtop->bIntermolecularInteractions)
2539             {
2540                 fr->bMolPBC = ir->bPeriodicMols;
2541
2542                 if (bSHAKE && fr->bMolPBC)
2543                 {
2544                     gmx_fatal(FARGS, "SHAKE is not supported with periodic molecules");
2545                 }
2546             }
2547             else
2548             {
2549                 fr->bMolPBC = TRUE;
2550
2551                 if (getenv("GMX_USE_GRAPH") != NULL)
2552                 {
2553                     fr->bMolPBC = FALSE;
2554                     if (fp)
2555                     {
2556                         md_print_warn(cr, fp, "GMX_USE_GRAPH is set, using the graph for bonded interactions\n");
2557                     }
2558
2559                     if (mtop->bIntermolecularInteractions)
2560                     {
2561                         md_print_warn(cr, fp, "WARNING: Molecules linked by intermolecular interactions have to reside in the same periodic image, otherwise artifacts will occur!\n");
2562                     }
2563                 }
2564
2565                 if (bSHAKE && fr->bMolPBC)
2566                 {
2567                     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");
2568                 }
2569             }
2570         }
2571         else
2572         {
2573             fr->bMolPBC = dd_bonded_molpbc(cr->dd, fr->ePBC);
2574         }
2575     }
2576     fr->bGB = (ir->implicit_solvent == eisGBSA);
2577
2578     fr->rc_scaling = ir->refcoord_scaling;
2579     copy_rvec(ir->posres_com, fr->posres_com);
2580     copy_rvec(ir->posres_comB, fr->posres_comB);
2581     fr->rlist                    = cutoff_inf(ir->rlist);
2582     fr->rlistlong                = cutoff_inf(ir->rlistlong);
2583     fr->eeltype                  = ir->coulombtype;
2584     fr->vdwtype                  = ir->vdwtype;
2585     fr->ljpme_combination_rule   = ir->ljpme_combination_rule;
2586
2587     fr->coulomb_modifier = ir->coulomb_modifier;
2588     fr->vdw_modifier     = ir->vdw_modifier;
2589
2590     /* Electrostatics: Translate from interaction-setting-in-mdp-file to kernel interaction format */
2591     switch (fr->eeltype)
2592     {
2593         case eelCUT:
2594             fr->nbkernel_elec_interaction = (fr->bGB) ? GMX_NBKERNEL_ELEC_GENERALIZEDBORN : GMX_NBKERNEL_ELEC_COULOMB;
2595             break;
2596
2597         case eelRF:
2598         case eelGRF:
2599         case eelRF_NEC:
2600             fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
2601             break;
2602
2603         case eelRF_ZERO:
2604             fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_REACTIONFIELD;
2605             fr->coulomb_modifier          = eintmodEXACTCUTOFF;
2606             break;
2607
2608         case eelSWITCH:
2609         case eelSHIFT:
2610         case eelUSER:
2611         case eelENCADSHIFT:
2612         case eelPMESWITCH:
2613         case eelPMEUSER:
2614         case eelPMEUSERSWITCH:
2615             fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
2616             break;
2617
2618         case eelPME:
2619         case eelEWALD:
2620             fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_EWALD;
2621             break;
2622
2623         default:
2624             gmx_fatal(FARGS, "Unsupported electrostatic interaction: %s", eel_names[fr->eeltype]);
2625             break;
2626     }
2627
2628     /* Vdw: Translate from mdp settings to kernel format */
2629     switch (fr->vdwtype)
2630     {
2631         case evdwCUT:
2632             if (fr->bBHAM)
2633             {
2634                 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_BUCKINGHAM;
2635             }
2636             else
2637             {
2638                 fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LENNARDJONES;
2639             }
2640             break;
2641         case evdwPME:
2642             fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_LJEWALD;
2643             break;
2644
2645         case evdwSWITCH:
2646         case evdwSHIFT:
2647         case evdwUSER:
2648         case evdwENCADSHIFT:
2649             fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
2650             break;
2651
2652         default:
2653             gmx_fatal(FARGS, "Unsupported vdw interaction: %s", evdw_names[fr->vdwtype]);
2654             break;
2655     }
2656
2657     /* These start out identical to ir, but might be altered if we e.g. tabulate the interaction in the kernel */
2658     fr->nbkernel_elec_modifier    = fr->coulomb_modifier;
2659     fr->nbkernel_vdw_modifier     = fr->vdw_modifier;
2660
2661     fr->rvdw             = cutoff_inf(ir->rvdw);
2662     fr->rvdw_switch      = ir->rvdw_switch;
2663     fr->rcoulomb         = cutoff_inf(ir->rcoulomb);
2664     fr->rcoulomb_switch  = ir->rcoulomb_switch;
2665
2666     fr->bTwinRange = fr->rlistlong > fr->rlist;
2667     fr->bEwald     = (EEL_PME(fr->eeltype) || fr->eeltype == eelEWALD);
2668
2669     fr->reppow     = mtop->ffparams.reppow;
2670
2671     if (ir->cutoff_scheme == ecutsGROUP)
2672     {
2673         fr->bvdwtab    = ((fr->vdwtype != evdwCUT || !gmx_within_tol(fr->reppow, 12.0, 10*GMX_DOUBLE_EPS))
2674                           && !EVDW_PME(fr->vdwtype));
2675         /* We have special kernels for standard Ewald and PME, but the pme-switch ones are tabulated above */
2676         fr->bcoultab   = !(fr->eeltype == eelCUT ||
2677                            fr->eeltype == eelEWALD ||
2678                            fr->eeltype == eelPME ||
2679                            fr->eeltype == eelRF ||
2680                            fr->eeltype == eelRF_ZERO);
2681
2682         /* If the user absolutely wants different switch/shift settings for coul/vdw, it is likely
2683          * going to be faster to tabulate the interaction than calling the generic kernel.
2684          * However, if generic kernels have been requested we keep things analytically.
2685          */
2686         if (fr->nbkernel_elec_modifier == eintmodPOTSWITCH &&
2687             fr->nbkernel_vdw_modifier == eintmodPOTSWITCH &&
2688             bGenericKernelOnly == FALSE)
2689         {
2690             if ((fr->rcoulomb_switch != fr->rvdw_switch) || (fr->rcoulomb != fr->rvdw))
2691             {
2692                 fr->bcoultab = TRUE;
2693                 /* Once we tabulate electrostatics, we can use the switch function for LJ,
2694                  * which would otherwise need two tables.
2695                  */
2696             }
2697         }
2698         else if ((fr->nbkernel_elec_modifier == eintmodPOTSHIFT && fr->nbkernel_vdw_modifier == eintmodPOTSHIFT) ||
2699                  ((fr->nbkernel_elec_interaction == GMX_NBKERNEL_ELEC_REACTIONFIELD &&
2700                    fr->nbkernel_elec_modifier == eintmodEXACTCUTOFF &&
2701                    (fr->nbkernel_vdw_modifier == eintmodPOTSWITCH || fr->nbkernel_vdw_modifier == eintmodPOTSHIFT))))
2702         {
2703             if ((fr->rcoulomb != fr->rvdw) && (bGenericKernelOnly == FALSE))
2704             {
2705                 fr->bcoultab = TRUE;
2706             }
2707         }
2708
2709         if (fr->nbkernel_elec_modifier == eintmodFORCESWITCH)
2710         {
2711             fr->bcoultab = TRUE;
2712         }
2713         if (fr->nbkernel_vdw_modifier == eintmodFORCESWITCH)
2714         {
2715             fr->bvdwtab = TRUE;
2716         }
2717
2718         if (getenv("GMX_REQUIRE_TABLES"))
2719         {
2720             fr->bvdwtab  = TRUE;
2721             fr->bcoultab = TRUE;
2722         }
2723
2724         if (fp)
2725         {
2726             fprintf(fp, "Table routines are used for coulomb: %s\n", bool_names[fr->bcoultab]);
2727             fprintf(fp, "Table routines are used for vdw:     %s\n", bool_names[fr->bvdwtab ]);
2728         }
2729
2730         if (fr->bvdwtab == TRUE)
2731         {
2732             fr->nbkernel_vdw_interaction = GMX_NBKERNEL_VDW_CUBICSPLINETABLE;
2733             fr->nbkernel_vdw_modifier    = eintmodNONE;
2734         }
2735         if (fr->bcoultab == TRUE)
2736         {
2737             fr->nbkernel_elec_interaction = GMX_NBKERNEL_ELEC_CUBICSPLINETABLE;
2738             fr->nbkernel_elec_modifier    = eintmodNONE;
2739         }
2740     }
2741
2742     if (ir->cutoff_scheme == ecutsVERLET)
2743     {
2744         if (!gmx_within_tol(fr->reppow, 12.0, 10*GMX_DOUBLE_EPS))
2745         {
2746             gmx_fatal(FARGS, "Cut-off scheme %S only supports LJ repulsion power 12", ecutscheme_names[ir->cutoff_scheme]);
2747         }
2748         fr->bvdwtab  = FALSE;
2749         fr->bcoultab = FALSE;
2750     }
2751
2752     /* Tables are used for direct ewald sum */
2753     if (fr->bEwald)
2754     {
2755         if (EEL_PME(ir->coulombtype))
2756         {
2757             if (fp)
2758             {
2759                 fprintf(fp, "Will do PME sum in reciprocal space for electrostatic interactions.\n");
2760             }
2761             if (ir->coulombtype == eelP3M_AD)
2762             {
2763                 please_cite(fp, "Hockney1988");
2764                 please_cite(fp, "Ballenegger2012");
2765             }
2766             else
2767             {
2768                 please_cite(fp, "Essmann95a");
2769             }
2770
2771             if (ir->ewald_geometry == eewg3DC)
2772             {
2773                 if (fp)
2774                 {
2775                     fprintf(fp, "Using the Ewald3DC correction for systems with a slab geometry.\n");
2776                 }
2777                 please_cite(fp, "In-Chul99a");
2778             }
2779         }
2780         fr->ewaldcoeff_q = calc_ewaldcoeff_q(ir->rcoulomb, ir->ewald_rtol);
2781         init_ewald_tab(&(fr->ewald_table), ir, fp);
2782         if (fp)
2783         {
2784             fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for Ewald\n",
2785                     1/fr->ewaldcoeff_q);
2786         }
2787     }
2788
2789     if (EVDW_PME(ir->vdwtype))
2790     {
2791         if (fp)
2792         {
2793             fprintf(fp, "Will do PME sum in reciprocal space for LJ dispersion interactions.\n");
2794         }
2795         please_cite(fp, "Essmann95a");
2796         fr->ewaldcoeff_lj = calc_ewaldcoeff_lj(ir->rvdw, ir->ewald_rtol_lj);
2797         if (fp)
2798         {
2799             fprintf(fp, "Using a Gaussian width (1/beta) of %g nm for LJ Ewald\n",
2800                     1/fr->ewaldcoeff_lj);
2801         }
2802     }
2803
2804     /* Electrostatics */
2805     fr->epsilon_r       = ir->epsilon_r;
2806     fr->epsilon_rf      = ir->epsilon_rf;
2807     fr->fudgeQQ         = mtop->ffparams.fudgeQQ;
2808
2809     /* Parameters for generalized RF */
2810     fr->zsquare = 0.0;
2811     fr->temp    = 0.0;
2812
2813     if (fr->eeltype == eelGRF)
2814     {
2815         init_generalized_rf(fp, mtop, ir, fr);
2816     }
2817
2818     fr->bF_NoVirSum = (EEL_FULL(fr->eeltype) || EVDW_PME(fr->vdwtype) ||
2819                        gmx_mtop_ftype_count(mtop, F_POSRES) > 0 ||
2820                        gmx_mtop_ftype_count(mtop, F_FBPOSRES) > 0 ||
2821                        IR_ELEC_FIELD(*ir) ||
2822                        (fr->adress_icor != eAdressICOff)
2823                        );
2824
2825     if (fr->cutoff_scheme == ecutsGROUP &&
2826         ncg_mtop(mtop) > fr->cg_nalloc && !DOMAINDECOMP(cr))
2827     {
2828         /* Count the total number of charge groups */
2829         fr->cg_nalloc = ncg_mtop(mtop);
2830         srenew(fr->cg_cm, fr->cg_nalloc);
2831     }
2832     if (fr->shift_vec == NULL)
2833     {
2834         snew(fr->shift_vec, SHIFTS);
2835     }
2836
2837     if (fr->fshift == NULL)
2838     {
2839         snew(fr->fshift, SHIFTS);
2840     }
2841
2842     if (fr->nbfp == NULL)
2843     {
2844         fr->ntype = mtop->ffparams.atnr;
2845         fr->nbfp  = mk_nbfp(&mtop->ffparams, fr->bBHAM);
2846         if (EVDW_PME(fr->vdwtype))
2847         {
2848             fr->ljpme_c6grid  = make_ljpme_c6grid(&mtop->ffparams, fr);
2849         }
2850     }
2851
2852     /* Copy the energy group exclusions */
2853     fr->egp_flags = ir->opts.egp_flags;
2854
2855     /* Van der Waals stuff */
2856     if ((fr->vdwtype != evdwCUT) && (fr->vdwtype != evdwUSER) && !fr->bBHAM)
2857     {
2858         if (fr->rvdw_switch >= fr->rvdw)
2859         {
2860             gmx_fatal(FARGS, "rvdw_switch (%f) must be < rvdw (%f)",
2861                       fr->rvdw_switch, fr->rvdw);
2862         }
2863         if (fp)
2864         {
2865             fprintf(fp, "Using %s Lennard-Jones, switch between %g and %g nm\n",
2866                     (fr->eeltype == eelSWITCH) ? "switched" : "shifted",
2867                     fr->rvdw_switch, fr->rvdw);
2868         }
2869     }
2870
2871     if (fr->bBHAM && EVDW_PME(fr->vdwtype))
2872     {
2873         gmx_fatal(FARGS, "LJ PME not supported with Buckingham");
2874     }
2875
2876     if (fr->bBHAM && (fr->vdwtype == evdwSHIFT || fr->vdwtype == evdwSWITCH))
2877     {
2878         gmx_fatal(FARGS, "Switch/shift interaction not supported with Buckingham");
2879     }
2880
2881     if (fr->bBHAM && fr->cutoff_scheme == ecutsVERLET)
2882     {
2883         gmx_fatal(FARGS, "Verlet cutoff-scheme is not supported with Buckingham");
2884     }
2885
2886     if (fp)
2887     {
2888         fprintf(fp, "Cut-off's:   NS: %g   Coulomb: %g   %s: %g\n",
2889                 fr->rlist, fr->rcoulomb, fr->bBHAM ? "BHAM" : "LJ", fr->rvdw);
2890     }
2891
2892     fr->eDispCorr = ir->eDispCorr;
2893     if (ir->eDispCorr != edispcNO)
2894     {
2895         set_avcsixtwelve(fp, fr, mtop);
2896     }
2897
2898     if (fr->bBHAM)
2899     {
2900         set_bham_b_max(fp, fr, mtop);
2901     }
2902
2903     fr->gb_epsilon_solvent = ir->gb_epsilon_solvent;
2904
2905     /* Copy the GBSA data (radius, volume and surftens for each
2906      * atomtype) from the topology atomtype section to forcerec.
2907      */
2908     snew(fr->atype_radius, fr->ntype);
2909     snew(fr->atype_vol, fr->ntype);
2910     snew(fr->atype_surftens, fr->ntype);
2911     snew(fr->atype_gb_radius, fr->ntype);
2912     snew(fr->atype_S_hct, fr->ntype);
2913
2914     if (mtop->atomtypes.nr > 0)
2915     {
2916         for (i = 0; i < fr->ntype; i++)
2917         {
2918             fr->atype_radius[i] = mtop->atomtypes.radius[i];
2919         }
2920         for (i = 0; i < fr->ntype; i++)
2921         {
2922             fr->atype_vol[i] = mtop->atomtypes.vol[i];
2923         }
2924         for (i = 0; i < fr->ntype; i++)
2925         {
2926             fr->atype_surftens[i] = mtop->atomtypes.surftens[i];
2927         }
2928         for (i = 0; i < fr->ntype; i++)
2929         {
2930             fr->atype_gb_radius[i] = mtop->atomtypes.gb_radius[i];
2931         }
2932         for (i = 0; i < fr->ntype; i++)
2933         {
2934             fr->atype_S_hct[i] = mtop->atomtypes.S_hct[i];
2935         }
2936     }
2937
2938     /* Generate the GB table if needed */
2939     if (fr->bGB)
2940     {
2941 #ifdef GMX_DOUBLE
2942         fr->gbtabscale = 2000;
2943 #else
2944         fr->gbtabscale = 500;
2945 #endif
2946
2947         fr->gbtabr = 100;
2948         fr->gbtab  = make_gb_table(oenv, fr);
2949
2950         init_gb(&fr->born, fr, ir, mtop, ir->gb_algorithm);
2951
2952         /* Copy local gb data (for dd, this is done in dd_partition_system) */
2953         if (!DOMAINDECOMP(cr))
2954         {
2955             make_local_gb(cr, fr->born, ir->gb_algorithm);
2956         }
2957     }
2958
2959     /* Set the charge scaling */
2960     if (fr->epsilon_r != 0)
2961     {
2962         fr->epsfac = ONE_4PI_EPS0/fr->epsilon_r;
2963     }
2964     else
2965     {
2966         /* eps = 0 is infinite dieletric: no coulomb interactions */
2967         fr->epsfac = 0;
2968     }
2969
2970     /* Reaction field constants */
2971     if (EEL_RF(fr->eeltype))
2972     {
2973         calc_rffac(fp, fr->eeltype, fr->epsilon_r, fr->epsilon_rf,
2974                    fr->rcoulomb, fr->temp, fr->zsquare, box,
2975                    &fr->kappa, &fr->k_rf, &fr->c_rf);
2976     }
2977
2978     /*This now calculates sum for q and c6*/
2979     set_chargesum(fp, fr, mtop);
2980
2981     /* if we are using LR electrostatics, and they are tabulated,
2982      * the tables will contain modified coulomb interactions.
2983      * Since we want to use the non-shifted ones for 1-4
2984      * coulombic interactions, we must have an extra set of tables.
2985      */
2986
2987     /* Construct tables.
2988      * A little unnecessary to make both vdw and coul tables sometimes,
2989      * but what the heck... */
2990
2991     bMakeTables = fr->bcoultab || fr->bvdwtab || fr->bEwald ||
2992         (ir->eDispCorr != edispcNO && ir_vdw_switched(ir));
2993
2994     bMakeSeparate14Table = ((!bMakeTables || fr->eeltype != eelCUT || fr->vdwtype != evdwCUT ||
2995                              fr->coulomb_modifier != eintmodNONE ||
2996                              fr->vdw_modifier != eintmodNONE ||
2997                              fr->bBHAM || fr->bEwald) &&
2998                             (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 ||
2999                              gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0 ||
3000                              gmx_mtop_ftype_count(mtop, F_LJC_PAIRS_NB) > 0));
3001
3002     negp_pp   = ir->opts.ngener - ir->nwall;
3003     negptable = 0;
3004     if (!bMakeTables)
3005     {
3006         bSomeNormalNbListsAreInUse = TRUE;
3007         fr->nnblists               = 1;
3008     }
3009     else
3010     {
3011         bSomeNormalNbListsAreInUse = (ir->eDispCorr != edispcNO);
3012         for (egi = 0; egi < negp_pp; egi++)
3013         {
3014             for (egj = egi; egj < negp_pp; egj++)
3015             {
3016                 egp_flags = ir->opts.egp_flags[GID(egi, egj, ir->opts.ngener)];
3017                 if (!(egp_flags & EGP_EXCL))
3018                 {
3019                     if (egp_flags & EGP_TABLE)
3020                     {
3021                         negptable++;
3022                     }
3023                     else
3024                     {
3025                         bSomeNormalNbListsAreInUse = TRUE;
3026                     }
3027                 }
3028             }
3029         }
3030         if (bSomeNormalNbListsAreInUse)
3031         {
3032             fr->nnblists = negptable + 1;
3033         }
3034         else
3035         {
3036             fr->nnblists = negptable;
3037         }
3038         if (fr->nnblists > 1)
3039         {
3040             snew(fr->gid2nblists, ir->opts.ngener*ir->opts.ngener);
3041         }
3042     }
3043
3044     if (ir->adress)
3045     {
3046         fr->nnblists *= 2;
3047     }
3048
3049     snew(fr->nblists, fr->nnblists);
3050
3051     /* This code automatically gives table length tabext without cut-off's,
3052      * in that case grompp should already have checked that we do not need
3053      * normal tables and we only generate tables for 1-4 interactions.
3054      */
3055     rtab = ir->rlistlong + ir->tabext;
3056
3057     if (bMakeTables)
3058     {
3059         /* make tables for ordinary interactions */
3060         if (bSomeNormalNbListsAreInUse)
3061         {
3062             make_nbf_tables(fp, oenv, fr, rtab, cr, tabfn, NULL, NULL, &fr->nblists[0]);
3063             if (ir->adress)
3064             {
3065                 make_nbf_tables(fp, oenv, fr, rtab, cr, tabfn, NULL, NULL, &fr->nblists[fr->nnblists/2]);
3066             }
3067             if (!bMakeSeparate14Table)
3068             {
3069                 fr->tab14 = fr->nblists[0].table_elec_vdw;
3070             }
3071             m = 1;
3072         }
3073         else
3074         {
3075             m = 0;
3076         }
3077         if (negptable > 0)
3078         {
3079             /* Read the special tables for certain energy group pairs */
3080             nm_ind = mtop->groups.grps[egcENER].nm_ind;
3081             for (egi = 0; egi < negp_pp; egi++)
3082             {
3083                 for (egj = egi; egj < negp_pp; egj++)
3084                 {
3085                     egp_flags = ir->opts.egp_flags[GID(egi, egj, ir->opts.ngener)];
3086                     if ((egp_flags & EGP_TABLE) && !(egp_flags & EGP_EXCL))
3087                     {
3088                         if (fr->nnblists > 1)
3089                         {
3090                             fr->gid2nblists[GID(egi, egj, ir->opts.ngener)] = m;
3091                         }
3092                         /* Read the table file with the two energy groups names appended */
3093                         make_nbf_tables(fp, oenv, fr, rtab, cr, tabfn,
3094                                         *mtop->groups.grpname[nm_ind[egi]],
3095                                         *mtop->groups.grpname[nm_ind[egj]],
3096                                         &fr->nblists[m]);
3097                         if (ir->adress)
3098                         {
3099                             make_nbf_tables(fp, oenv, fr, rtab, cr, tabfn,
3100                                             *mtop->groups.grpname[nm_ind[egi]],
3101                                             *mtop->groups.grpname[nm_ind[egj]],
3102                                             &fr->nblists[fr->nnblists/2+m]);
3103                         }
3104                         m++;
3105                     }
3106                     else if (fr->nnblists > 1)
3107                     {
3108                         fr->gid2nblists[GID(egi, egj, ir->opts.ngener)] = 0;
3109                     }
3110                 }
3111             }
3112         }
3113     }
3114     else if ((fr->eDispCorr != edispcNO) &&
3115              ((fr->vdw_modifier == eintmodPOTSWITCH) ||
3116               (fr->vdw_modifier == eintmodFORCESWITCH) ||
3117               (fr->vdw_modifier == eintmodPOTSHIFT)))
3118     {
3119         /* Tables might not be used for the potential modifier interactions per se, but
3120          * we still need them to evaluate switch/shift dispersion corrections in this case.
3121          */
3122         make_nbf_tables(fp, oenv, fr, rtab, cr, tabfn, NULL, NULL, &fr->nblists[0]);
3123     }
3124
3125     if (bMakeSeparate14Table)
3126     {
3127         /* generate extra tables with plain Coulomb for 1-4 interactions only */
3128         fr->tab14 = make_tables(fp, oenv, fr, MASTER(cr), tabpfn, rtab,
3129                                 GMX_MAKETABLES_14ONLY);
3130     }
3131
3132     /* Read AdResS Thermo Force table if needed */
3133     if (fr->adress_icor == eAdressICThermoForce)
3134     {
3135         /* old todo replace */
3136
3137         if (ir->adress->n_tf_grps > 0)
3138         {
3139             make_adress_tf_tables(fp, oenv, fr, ir, tabfn, mtop, box);
3140
3141         }
3142         else
3143         {
3144             /* load the default table */
3145             snew(fr->atf_tabs, 1);
3146             fr->atf_tabs[DEFAULT_TF_TABLE] = make_atf_table(fp, oenv, fr, tabafn, box);
3147         }
3148     }
3149
3150     /* Wall stuff */
3151     fr->nwall = ir->nwall;
3152     if (ir->nwall && ir->wall_type == ewtTABLE)
3153     {
3154         make_wall_tables(fp, oenv, ir, tabfn, &mtop->groups, fr);
3155     }
3156
3157     if (fcd && tabbfn)
3158     {
3159         fcd->bondtab  = make_bonded_tables(fp,
3160                                            F_TABBONDS, F_TABBONDSNC,
3161                                            mtop, tabbfn, "b");
3162         fcd->angletab = make_bonded_tables(fp,
3163                                            F_TABANGLES, -1,
3164                                            mtop, tabbfn, "a");
3165         fcd->dihtab   = make_bonded_tables(fp,
3166                                            F_TABDIHS, -1,
3167                                            mtop, tabbfn, "d");
3168     }
3169     else
3170     {
3171         if (debug)
3172         {
3173             fprintf(debug, "No fcdata or table file name passed, can not read table, can not do bonded interactions\n");
3174         }
3175     }
3176
3177     /* QM/MM initialization if requested
3178      */
3179     if (ir->bQMMM)
3180     {
3181         fprintf(stderr, "QM/MM calculation requested.\n");
3182     }
3183
3184     fr->bQMMM      = ir->bQMMM;
3185     fr->qr         = mk_QMMMrec();
3186
3187     /* Set all the static charge group info */
3188     fr->cginfo_mb = init_cginfo_mb(fp, mtop, fr, bNoSolvOpt,
3189                                    &bFEP_NonBonded,
3190                                    &fr->bExcl_IntraCGAll_InterCGNone);
3191     if (DOMAINDECOMP(cr))
3192     {
3193         fr->cginfo = NULL;
3194     }
3195     else
3196     {
3197         fr->cginfo = cginfo_expand(mtop->nmolblock, fr->cginfo_mb);
3198     }
3199
3200     if (!DOMAINDECOMP(cr))
3201     {
3202         forcerec_set_ranges(fr, ncg_mtop(mtop), ncg_mtop(mtop),
3203                             mtop->natoms, mtop->natoms, mtop->natoms);
3204     }
3205
3206     fr->print_force = print_force;
3207
3208
3209     /* coarse load balancing vars */
3210     fr->t_fnbf    = 0.;
3211     fr->t_wait    = 0.;
3212     fr->timesteps = 0;
3213
3214     /* Initialize neighbor search */
3215     init_ns(fp, cr, &fr->ns, fr, mtop);
3216
3217     if (cr->duty & DUTY_PP)
3218     {
3219         gmx_nonbonded_setup(fr, bGenericKernelOnly);
3220         /*
3221            if (ir->bAdress)
3222             {
3223                 gmx_setup_adress_kernels(fp,bGenericKernelOnly);
3224             }
3225          */
3226     }
3227
3228     /* Initialize the thread working data for bonded interactions */
3229     init_bonded_threading(fp, fr, mtop->groups.grps[egcENER].nr);
3230
3231     snew(fr->excl_load, fr->nthreads+1);
3232
3233     /* fr->ic is used both by verlet and group kernels (to some extent) now */
3234     init_interaction_const(fp, &fr->ic, fr);
3235
3236     if (fr->cutoff_scheme == ecutsVERLET)
3237     {
3238         if (ir->rcoulomb != ir->rvdw)
3239         {
3240             gmx_fatal(FARGS, "With Verlet lists rcoulomb and rvdw should be identical");
3241         }
3242
3243         init_nb_verlet(fp, &fr->nbv, bFEP_NonBonded, ir, fr, cr, nbpu_opt);
3244     }
3245
3246     init_interaction_const_tables(fp, fr->ic, rtab);
3247
3248     initialize_gpu_constants(cr, fr->ic, fr->nbv);
3249
3250     if (ir->eDispCorr != edispcNO)
3251     {
3252         calc_enervirdiff(fp, ir->eDispCorr, fr);
3253     }
3254 }
3255
3256 #define pr_real(fp, r) fprintf(fp, "%s: %e\n",#r, r)
3257 #define pr_int(fp, i)  fprintf((fp), "%s: %d\n",#i, i)
3258 #define pr_bool(fp, b) fprintf((fp), "%s: %s\n",#b, bool_names[b])
3259
3260 void pr_forcerec(FILE *fp, t_forcerec *fr)
3261 {
3262     int i;
3263
3264     pr_real(fp, fr->rlist);
3265     pr_real(fp, fr->rcoulomb);
3266     pr_real(fp, fr->fudgeQQ);
3267     pr_bool(fp, fr->bGrid);
3268     pr_bool(fp, fr->bTwinRange);
3269     /*pr_int(fp,fr->cg0);
3270        pr_int(fp,fr->hcg);*/
3271     for (i = 0; i < fr->nnblists; i++)
3272     {
3273         pr_int(fp, fr->nblists[i].table_elec_vdw.n);
3274     }
3275     pr_real(fp, fr->rcoulomb_switch);
3276     pr_real(fp, fr->rcoulomb);
3277
3278     fflush(fp);
3279 }
3280
3281 void forcerec_set_excl_load(t_forcerec           *fr,
3282                             const gmx_localtop_t *top)
3283 {
3284     const int *ind, *a;
3285     int        t, i, j, ntot, n, ntarget;
3286
3287     ind = top->excls.index;
3288     a   = top->excls.a;
3289
3290     ntot = 0;
3291     for (i = 0; i < top->excls.nr; i++)
3292     {
3293         for (j = ind[i]; j < ind[i+1]; j++)
3294         {
3295             if (a[j] > i)
3296             {
3297                 ntot++;
3298             }
3299         }
3300     }
3301
3302     fr->excl_load[0] = 0;
3303     n                = 0;
3304     i                = 0;
3305     for (t = 1; t <= fr->nthreads; t++)
3306     {
3307         ntarget = (ntot*t)/fr->nthreads;
3308         while (i < top->excls.nr && n < ntarget)
3309         {
3310             for (j = ind[i]; j < ind[i+1]; j++)
3311             {
3312                 if (a[j] > i)
3313                 {
3314                     n++;
3315                 }
3316             }
3317             i++;
3318         }
3319         fr->excl_load[t] = i;
3320     }
3321 }
3322
3323 /* Frees GPU memory and destroys the GPU context.
3324  *
3325  * Note that this function needs to be called even if GPUs are not used
3326  * in this run because the PME ranks have no knowledge of whether GPUs
3327  * are used or not, but all ranks need to enter the barrier below.
3328  */
3329 void free_gpu_resources(const t_forcerec     *fr,
3330                         const t_commrec      *cr,
3331                         const gmx_gpu_info_t *gpu_info,
3332                         const gmx_gpu_opt_t  *gpu_opt)
3333 {
3334     gmx_bool bIsPPrankUsingGPU;
3335     char     gpu_err_str[STRLEN];
3336
3337     bIsPPrankUsingGPU = (cr->duty & DUTY_PP) && fr && fr->nbv && fr->nbv->bUseGPU;
3338
3339     if (bIsPPrankUsingGPU)
3340     {
3341         /* free nbnxn data in GPU memory */
3342         nbnxn_gpu_free(fr->nbv->gpu_nbv);
3343
3344         /* With tMPI we need to wait for all ranks to finish deallocation before
3345          * destroying the context in free_gpu() as some ranks may be sharing
3346          * GPU and context.
3347          * Note: as only PP ranks need to free GPU resources, so it is safe to
3348          * not call the barrier on PME ranks.
3349          */
3350 #ifdef GMX_THREAD_MPI
3351         if (PAR(cr))
3352         {
3353             gmx_barrier(cr);
3354         }
3355 #endif  /* GMX_THREAD_MPI */
3356
3357         /* uninitialize GPU (by destroying the context) */
3358         if (!free_cuda_gpu(cr->rank_pp_intranode, gpu_err_str, gpu_info, gpu_opt))
3359         {
3360             gmx_warning("On rank %d failed to free GPU #%d: %s",
3361                         cr->nodeid, get_current_cuda_gpu_device_id(), gpu_err_str);
3362         }
3363     }
3364 }