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