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