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