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