Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / mdlib / ns.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 "sysstuff.h"
43 #include "smalloc.h"
44 #include "macros.h"
45 #include "maths.h"
46 #include "vec.h"
47 #include "network.h"
48 #include "nsgrid.h"
49 #include "force.h"
50 #include "nonbonded.h"
51 #include "ns.h"
52 #include "pbc.h"
53 #include "names.h"
54 #include "gmx_fatal.h"
55 #include "nrnb.h"
56 #include "txtdump.h"
57 #include "mtop_util.h"
58
59 #include "domdec.h"
60 #include "adress.h"
61
62
63 /* 
64  *    E X C L U S I O N   H A N D L I N G
65  */
66
67 #ifdef DEBUG
68 static void SETEXCL_(t_excl e[],atom_id i,atom_id j)
69 {   e[j] = e[j] | (1<<i); }
70 static void RMEXCL_(t_excl e[],atom_id i,atom_id j) 
71 { e[j]=e[j] & ~(1<<i); }
72 static gmx_bool ISEXCL_(t_excl e[],atom_id i,atom_id j) 
73 { return (gmx_bool)(e[j] & (1<<i)); }
74 static gmx_bool NOTEXCL_(t_excl e[],atom_id i,atom_id j)
75 {  return !(ISEXCL(e,i,j)); }
76 #else
77 #define SETEXCL(e,i,j) (e)[((atom_id) (j))] |= (1<<((atom_id) (i)))
78 #define RMEXCL(e,i,j)  (e)[((atom_id) (j))] &= (~(1<<((atom_id) (i))))
79 #define ISEXCL(e,i,j)  (gmx_bool) ((e)[((atom_id) (j))] & (1<<((atom_id) (i))))
80 #define NOTEXCL(e,i,j) !(ISEXCL(e,i,j))
81 #endif
82
83 static int
84 round_up_to_simd_width(int length, int simd_width)
85 {
86     int offset,newlength;
87     
88     offset = (simd_width>0) ? length % simd_width : 0;
89
90     return (offset==0) ? length : length-offset+simd_width;
91 }
92 /************************************************
93  *
94  *  U T I L I T I E S    F O R    N S
95  *
96  ************************************************/
97
98 static void reallocate_nblist(t_nblist *nl)
99 {
100     if (gmx_debug_at)
101     {
102         fprintf(debug,"reallocating neigborlist (ielec=%d, ivdw=%d, igeometry=%d, type=%d), maxnri=%d\n",
103                 nl->ielec,nl->ivdw,nl->igeometry,nl->type,nl->maxnri);
104     }
105     srenew(nl->iinr,   nl->maxnri);
106     if (nl->igeometry == GMX_NBLIST_GEOMETRY_CG_CG)
107     {
108         srenew(nl->iinr_end,nl->maxnri);
109     }
110     srenew(nl->gid,    nl->maxnri);
111     srenew(nl->shift,  nl->maxnri);
112     srenew(nl->jindex, nl->maxnri+1);
113 }
114
115
116 static void init_nblist(FILE *log, t_nblist *nl_sr,t_nblist *nl_lr,
117                         int maxsr,int maxlr,
118                         int ivdw, int ivdwmod,
119                         int ielec, int ielecmod,
120                         int igeometry, int type)
121 {
122     t_nblist *nl;
123     int      homenr;
124     int      i,nn;
125     
126     for(i=0; (i<2); i++)
127     {
128         nl     = (i == 0) ? nl_sr : nl_lr;
129         homenr = (i == 0) ? maxsr : maxlr;
130
131         if (nl == NULL)
132         {
133             continue;
134         }
135
136
137         /* Set coul/vdw in neighborlist, and for the normal loops we determine
138          * an index of which one to call.
139          */
140         nl->ivdw        = ivdw;
141         nl->ivdwmod     = ivdwmod;
142         nl->ielec       = ielec;
143         nl->ielecmod    = ielecmod;
144         nl->type        = type;
145         nl->igeometry   = igeometry;
146
147         if (nl->type==GMX_NBLIST_INTERACTION_FREE_ENERGY)
148         {
149             nl->igeometry  = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
150         }
151
152         /* This will also set the simd_padding_width field */
153         gmx_nonbonded_set_kernel_pointers( (i==0) ? log : NULL,nl);
154         
155         /* maxnri is influenced by the number of shifts (maximum is 8)
156          * and the number of energy groups.
157          * If it is not enough, nl memory will be reallocated during the run.
158          * 4 seems to be a reasonable factor, which only causes reallocation
159          * during runs with tiny and many energygroups.
160          */
161         nl->maxnri      = homenr*4;
162         nl->maxnrj      = 0;
163         nl->maxlen      = 0;
164         nl->nri         = -1;
165         nl->nrj         = 0;
166         nl->iinr        = NULL;
167         nl->gid         = NULL;
168         nl->shift       = NULL;
169         nl->jindex      = NULL;
170         reallocate_nblist(nl);
171         nl->jindex[0] = 0;
172
173         if(debug)
174         {
175             fprintf(debug,"Initiating neighbourlist (ielec=%d, ivdw=%d, type=%d) for %s interactions,\nwith %d SR, %d LR atoms.\n",
176                     nl->ielec,nl->ivdw,nl->type,gmx_nblist_geometry_names[nl->igeometry],maxsr,maxlr);
177         }
178     }
179 }
180
181 void init_neighbor_list(FILE *log,t_forcerec *fr,int homenr)
182 {
183    /* Make maxlr tunable! (does not seem to be a big difference though) 
184     * This parameter determines the number of i particles in a long range 
185     * neighbourlist. Too few means many function calls, too many means
186     * cache trashing.
187     */
188    int maxsr,maxsr_wat,maxlr,maxlr_wat;
189    int ielec,ielecf,ivdw,ielecmod,ielecmodf,ivdwmod,type;
190    int solvent;
191    int igeometry_def,igeometry_w,igeometry_ww;
192    int i;
193    t_nblists *nbl;
194
195    /* maxsr     = homenr-fr->nWatMol*3; */
196    maxsr     = homenr;
197
198    if (maxsr < 0)
199    {
200      gmx_fatal(FARGS,"%s, %d: Negative number of short range atoms.\n"
201                  "Call your Gromacs dealer for assistance.",__FILE__,__LINE__);
202    }
203    /* This is just for initial allocation, so we do not reallocate
204     * all the nlist arrays many times in a row.
205     * The numbers seem very accurate, but they are uncritical.
206     */
207    maxsr_wat = min(fr->nWatMol,(homenr+2)/3); 
208    if (fr->bTwinRange) 
209    {
210        maxlr     = 50;
211        maxlr_wat = min(maxsr_wat,maxlr);
212    }
213    else
214    {
215      maxlr = maxlr_wat = 0;
216    }  
217
218    /* Determine the values for ielec/ivdw. */
219    ielec = fr->nbkernel_elec_interaction;
220    ivdw  = fr->nbkernel_vdw_interaction;
221    ielecmod = fr->nbkernel_elec_modifier;
222    ivdwmod  = fr->nbkernel_vdw_modifier;
223    type     = GMX_NBLIST_INTERACTION_STANDARD;
224
225    fr->ns.bCGlist = (getenv("GMX_NBLISTCG") != 0);
226    if (!fr->ns.bCGlist)
227    {
228        igeometry_def = GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE;
229    }
230    else
231    {
232        igeometry_def = GMX_NBLIST_GEOMETRY_CG_CG;
233        if (log != NULL)
234        {
235            fprintf(log,"\nUsing charge-group - charge-group neighbor lists and kernels\n\n");
236        }
237    }
238    
239    if (fr->solvent_opt == esolTIP4P) {
240        igeometry_w  = GMX_NBLIST_GEOMETRY_WATER4_PARTICLE;
241        igeometry_ww = GMX_NBLIST_GEOMETRY_WATER4_WATER4;
242    } else {
243        igeometry_w  = GMX_NBLIST_GEOMETRY_WATER3_PARTICLE;
244        igeometry_ww = GMX_NBLIST_GEOMETRY_WATER3_WATER3;
245    }
246
247    for(i=0; i<fr->nnblists; i++) 
248    {
249        nbl = &(fr->nblists[i]);
250
251        if ((fr->adress_type!=eAdressOff) && (i>=fr->nnblists/2)){
252            type=GMX_NBLIST_INTERACTION_ADRESS;
253        }
254        init_nblist(log,&nbl->nlist_sr[eNL_VDWQQ],&nbl->nlist_lr[eNL_VDWQQ],
255                    maxsr,maxlr,ivdw,ivdwmod,ielec,ielecmod,igeometry_def, type);
256        init_nblist(log,&nbl->nlist_sr[eNL_VDW],&nbl->nlist_lr[eNL_VDW],
257                    maxsr,maxlr,ivdw,ivdwmod,GMX_NBKERNEL_ELEC_NONE,eintmodNONE,igeometry_def, type);
258        init_nblist(log,&nbl->nlist_sr[eNL_QQ],&nbl->nlist_lr[eNL_QQ],
259                    maxsr,maxlr,GMX_NBKERNEL_VDW_NONE,eintmodNONE,ielec,ielecmod,igeometry_def, type);
260        init_nblist(log,&nbl->nlist_sr[eNL_VDWQQ_WATER],&nbl->nlist_lr[eNL_VDWQQ_WATER],
261                    maxsr_wat,maxlr_wat,ivdw,ivdwmod,ielec,ielecmod, igeometry_w, type);
262        init_nblist(log,&nbl->nlist_sr[eNL_QQ_WATER],&nbl->nlist_lr[eNL_QQ_WATER],
263                    maxsr_wat,maxlr_wat,GMX_NBKERNEL_VDW_NONE,eintmodNONE,ielec,ielecmod, igeometry_w, type);
264        init_nblist(log,&nbl->nlist_sr[eNL_VDWQQ_WATERWATER],&nbl->nlist_lr[eNL_VDWQQ_WATERWATER],
265                    maxsr_wat,maxlr_wat,ivdw,ivdwmod,ielec,ielecmod, igeometry_ww, type);
266        init_nblist(log,&nbl->nlist_sr[eNL_QQ_WATERWATER],&nbl->nlist_lr[eNL_QQ_WATERWATER],
267                    maxsr_wat,maxlr_wat,GMX_NBKERNEL_VDW_NONE,eintmodNONE,ielec,ielecmod, igeometry_ww, type);
268
269        /* Did we get the solvent loops so we can use optimized water kernels? */
270        if(nbl->nlist_sr[eNL_VDWQQ_WATER].kernelptr_vf==NULL
271           || nbl->nlist_sr[eNL_QQ_WATER].kernelptr_vf==NULL
272 #ifndef DISABLE_WATERWATER_NLIST
273           || nbl->nlist_sr[eNL_VDWQQ_WATERWATER].kernelptr_vf==NULL
274           || nbl->nlist_sr[eNL_QQ_WATERWATER].kernelptr_vf==NULL
275 #endif
276           )
277        {
278            fr->solvent_opt = esolNO;
279            fprintf(log,"Note: The available nonbonded kernels do not support water optimization - disabling.\n");
280        }
281        
282        if (fr->efep != efepNO) 
283        {
284            if ((fr->bEwald) && (fr->sc_alphacoul > 0)) /* need to handle long range differently if using softcore */
285            {
286                ielecf = GMX_NBKERNEL_ELEC_EWALD;
287                ielecmodf = eintmodNONE;
288            }
289            else
290            {
291                ielecf = ielec;
292                ielecmodf = ielecmod;
293            }
294
295            init_nblist(log,&nbl->nlist_sr[eNL_VDWQQ_FREE],&nbl->nlist_lr[eNL_VDWQQ_FREE],
296                        maxsr,maxlr,ivdw,ivdwmod,ielecf,ielecmod,GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
297            init_nblist(log,&nbl->nlist_sr[eNL_VDW_FREE],&nbl->nlist_lr[eNL_VDW_FREE],
298                        maxsr,maxlr,ivdw,ivdwmod,GMX_NBKERNEL_ELEC_NONE,eintmodNONE,GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
299            init_nblist(log,&nbl->nlist_sr[eNL_QQ_FREE],&nbl->nlist_lr[eNL_QQ_FREE],
300                        maxsr,maxlr,GMX_NBKERNEL_VDW_NONE,eintmodNONE,ielecf,ielecmod,GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_FREE_ENERGY);
301        }  
302    }
303    /* QMMM MM list */
304    if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
305    {
306        init_nblist(log,&fr->QMMMlist,NULL,
307                    maxsr,maxlr,0,0,ielec,ielecmod,GMX_NBLIST_GEOMETRY_PARTICLE_PARTICLE, GMX_NBLIST_INTERACTION_STANDARD);
308    }
309
310    if(log!=NULL)
311    {
312        fprintf(log,"\n");
313    }
314
315    fr->ns.nblist_initialized=TRUE;
316 }
317
318 static void reset_nblist(t_nblist *nl)
319 {
320      nl->nri       = -1;
321      nl->nrj       = 0;
322      nl->maxlen    = 0;
323      if (nl->jindex)
324      {
325          nl->jindex[0] = 0;
326      }
327 }
328
329 static void reset_neighbor_lists(t_forcerec *fr,gmx_bool bResetSR, gmx_bool bResetLR)
330 {
331     int n,i;
332   
333     if (fr->bQMMM)
334     {
335         /* only reset the short-range nblist */
336         reset_nblist(&(fr->QMMMlist));
337     }
338
339     for(n=0; n<fr->nnblists; n++)
340     {
341         for(i=0; i<eNL_NR; i++)
342         {
343             if(bResetSR)
344             {
345                 reset_nblist( &(fr->nblists[n].nlist_sr[i]) );
346             }
347             if(bResetLR)
348             {
349                 reset_nblist( &(fr->nblists[n].nlist_lr[i]) );
350             }
351         }
352     }
353 }
354
355
356
357
358 static inline void new_i_nblist(t_nblist *nlist,
359                                 gmx_bool bLR,atom_id i_atom,int shift,int gid)
360 {
361     int    i,k,nri,nshift;
362     
363     nri = nlist->nri;
364     
365     /* Check whether we have to increase the i counter */
366     if ((nri == -1) ||
367         (nlist->iinr[nri]  != i_atom) || 
368         (nlist->shift[nri] != shift) || 
369         (nlist->gid[nri]   != gid))
370     {
371         /* This is something else. Now see if any entries have 
372          * been added in the list of the previous atom.
373          */
374         if ((nri == -1) ||
375             ((nlist->jindex[nri+1] > nlist->jindex[nri]) && 
376              (nlist->gid[nri] != -1)))
377         {
378             /* If so increase the counter */
379             nlist->nri++;
380             nri++;
381             if (nlist->nri >= nlist->maxnri)
382             {
383                 nlist->maxnri += over_alloc_large(nlist->nri);
384                 reallocate_nblist(nlist);
385             }
386         }
387         /* Set the number of neighbours and the atom number */
388         nlist->jindex[nri+1] = nlist->jindex[nri];
389         nlist->iinr[nri]     = i_atom;
390         nlist->gid[nri]      = gid;
391         nlist->shift[nri]    = shift;
392     }
393 }
394
395 static inline void close_i_nblist(t_nblist *nlist) 
396 {
397     int nri = nlist->nri;
398     int len;
399     
400     if (nri >= 0)
401     {
402         /* Add elements up to padding. Since we allocate memory in units
403          * of the simd_padding width, we do not have to check for possible
404          * list reallocation here.
405          */
406         while((nlist->nrj % nlist->simd_padding_width)!=0)
407         {
408             /* Use -4 here, so we can write forces for 4 atoms before real data */
409             nlist->jjnr[nlist->nrj++]=-4;
410         }
411         nlist->jindex[nri+1] = nlist->nrj;
412         
413         len=nlist->nrj -  nlist->jindex[nri];
414         
415         /* nlist length for water i molecules is treated statically 
416          * in the innerloops 
417          */
418         if (len > nlist->maxlen)
419         {
420             nlist->maxlen = len;
421         }
422     }
423 }
424
425 static inline void close_nblist(t_nblist *nlist)
426 {
427     /* Only close this nblist when it has been initialized.
428      * Avoid the creation of i-lists with no j-particles.
429      */
430     if (nlist->nrj == 0)
431     {
432         /* Some assembly kernels do not support empty lists,
433          * make sure here that we don't generate any empty lists.
434          * With the current ns code this branch is taken in two cases:
435          * No i-particles at all: nri=-1 here
436          * There are i-particles, but no j-particles; nri=0 here
437          */
438         nlist->nri = 0;
439     }
440     else
441     {
442         /* Close list number nri by incrementing the count */
443         nlist->nri++;
444     }
445 }
446
447 static inline void close_neighbor_lists(t_forcerec *fr,gmx_bool bMakeQMMMnblist)
448 {
449     int n,i;
450     
451     if (bMakeQMMMnblist)
452     {
453             close_nblist(&(fr->QMMMlist));
454     }
455
456     for(n=0; n<fr->nnblists; n++)
457     {
458         for(i=0; (i<eNL_NR); i++)
459         {
460             close_nblist(&(fr->nblists[n].nlist_sr[i]));
461             close_nblist(&(fr->nblists[n].nlist_lr[i]));
462         }
463     }
464 }
465
466
467 static inline void add_j_to_nblist(t_nblist *nlist,atom_id j_atom,gmx_bool bLR)
468 {
469     int nrj=nlist->nrj;
470     
471     if (nlist->nrj >= nlist->maxnrj)
472     {
473         nlist->maxnrj = round_up_to_simd_width(over_alloc_small(nlist->nrj + 1),nlist->simd_padding_width);
474         
475         if (gmx_debug_at)
476             fprintf(debug,"Increasing %s nblist (ielec=%d,ivdw=%d,type=%d,igeometry=%d) j size to %d\n",
477                     bLR ? "LR" : "SR",nlist->ielec,nlist->ivdw,nlist->type,nlist->igeometry,nlist->maxnrj);
478         
479         srenew(nlist->jjnr,nlist->maxnrj);
480     }
481
482     nlist->jjnr[nrj] = j_atom;
483     nlist->nrj ++;
484 }
485
486 static inline void add_j_to_nblist_cg(t_nblist *nlist,
487                                       atom_id j_start,int j_end,
488                                       t_excl *bexcl,gmx_bool i_is_j,
489                                       gmx_bool bLR)
490 {
491     int nrj=nlist->nrj;
492     int j;
493
494     if (nlist->nrj >= nlist->maxnrj)
495     {
496         nlist->maxnrj = over_alloc_small(nlist->nrj + 1);
497         if (gmx_debug_at)
498             fprintf(debug,"Increasing %s nblist (ielec=%d,ivdw=%d,type=%d,igeometry=%d) j size to %d\n",
499                     bLR ? "LR" : "SR",nlist->ielec,nlist->ivdw,nlist->type,nlist->igeometry,nlist->maxnrj);
500         
501         srenew(nlist->jjnr    ,nlist->maxnrj);
502         srenew(nlist->jjnr_end,nlist->maxnrj);
503         srenew(nlist->excl    ,nlist->maxnrj*MAX_CGCGSIZE);
504     }
505
506     nlist->jjnr[nrj]     = j_start;
507     nlist->jjnr_end[nrj] = j_end;
508
509     if (j_end - j_start > MAX_CGCGSIZE)
510     {
511         gmx_fatal(FARGS,"The charge-group - charge-group neighborlist do not support charge groups larger than %d, found a charge group of size %d",MAX_CGCGSIZE,j_end-j_start);
512     }
513
514     /* Set the exclusions */
515     for(j=j_start; j<j_end; j++)
516     {
517         nlist->excl[nrj*MAX_CGCGSIZE + j - j_start] = bexcl[j];
518     }
519     if (i_is_j)
520     {
521         /* Avoid double counting of intra-cg interactions */
522         for(j=1; j<j_end-j_start; j++)
523         {
524             nlist->excl[nrj*MAX_CGCGSIZE + j] |= (1<<j) - 1;
525         }
526     }
527
528     nlist->nrj ++;
529 }
530
531 typedef void
532 put_in_list_t(gmx_bool              bHaveVdW[],
533               int               ngid,
534               t_mdatoms *       md,
535               int               icg,
536               int               jgid,
537               int               nj,
538               atom_id           jjcg[],
539               atom_id           index[],
540               t_excl            bExcl[],
541               int               shift,
542               t_forcerec *      fr,
543               gmx_bool          bLR,
544               gmx_bool          bDoVdW,
545               gmx_bool          bDoCoul,
546               int               solvent_opt);
547
548 static void 
549 put_in_list_at(gmx_bool              bHaveVdW[],
550                int               ngid,
551                t_mdatoms *       md,
552                int               icg,
553                int               jgid,
554                int               nj,
555                atom_id           jjcg[],
556                atom_id           index[],
557                t_excl            bExcl[],
558                int               shift,
559                t_forcerec *      fr,
560                gmx_bool          bLR,
561                gmx_bool          bDoVdW,
562                gmx_bool          bDoCoul,
563                int               solvent_opt)
564 {
565     /* The a[] index has been removed,
566      * to put it back in i_atom should be a[i0] and jj should be a[jj].
567      */
568     t_nblist *   vdwc;
569     t_nblist *   vdw;
570     t_nblist *   coul;
571     t_nblist *   vdwc_free  = NULL;
572     t_nblist *   vdw_free   = NULL;
573     t_nblist *   coul_free  = NULL;
574     t_nblist *   vdwc_ww    = NULL;
575     t_nblist *   coul_ww    = NULL;
576     
577     int             i,j,jcg,igid,gid,nbl_ind,ind_ij;
578     atom_id   jj,jj0,jj1,i_atom;
579     int       i0,nicg,len;
580     
581     int       *cginfo;
582     int       *type,*typeB;
583     real      *charge,*chargeB;
584     real      qi,qiB,qq,rlj;
585     gmx_bool      bFreeEnergy,bFree,bFreeJ,bNotEx,*bPert;
586     gmx_bool      bDoVdW_i,bDoCoul_i,bDoCoul_i_sol;
587     int       iwater,jwater;
588     t_nblist  *nlist;
589     
590     /* Copy some pointers */
591     cginfo  = fr->cginfo;
592     charge  = md->chargeA;
593     chargeB = md->chargeB;
594     type    = md->typeA;
595     typeB   = md->typeB;
596     bPert   = md->bPerturbed;
597     
598     /* Get atom range */
599     i0     = index[icg];
600     nicg   = index[icg+1]-i0;
601     
602     /* Get the i charge group info */
603     igid   = GET_CGINFO_GID(cginfo[icg]);
604
605     iwater = (solvent_opt!=esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
606     
607     bFreeEnergy = FALSE;
608     if (md->nPerturbed) 
609     {
610         /* Check if any of the particles involved are perturbed. 
611          * If not we can do the cheaper normal put_in_list
612          * and use more solvent optimization.
613          */
614         for(i=0; i<nicg; i++)
615         {
616             bFreeEnergy |= bPert[i0+i];
617         }
618         /* Loop over the j charge groups */
619         for(j=0; (j<nj && !bFreeEnergy); j++) 
620         {
621             jcg = jjcg[j];
622             jj0 = index[jcg];
623             jj1 = index[jcg+1];
624             /* Finally loop over the atoms in the j-charge group */     
625             for(jj=jj0; jj<jj1; jj++)
626             {
627                 bFreeEnergy |= bPert[jj];
628             }
629         }
630     }
631     
632     /* Unpack pointers to neighbourlist structs */
633     if (fr->nnblists == 1)
634     {
635         nbl_ind = 0;
636     }
637     else
638     {
639         nbl_ind = fr->gid2nblists[GID(igid,jgid,ngid)];
640     }
641     if (bLR)
642     {
643         nlist = fr->nblists[nbl_ind].nlist_lr;
644     }
645     else
646     {
647         nlist = fr->nblists[nbl_ind].nlist_sr;
648     }
649     
650     if (iwater != esolNO)
651     {
652         vdwc = &nlist[eNL_VDWQQ_WATER];
653         vdw  = &nlist[eNL_VDW];
654         coul = &nlist[eNL_QQ_WATER];
655 #ifndef DISABLE_WATERWATER_NLIST
656         vdwc_ww = &nlist[eNL_VDWQQ_WATERWATER];
657         coul_ww = &nlist[eNL_QQ_WATERWATER];
658 #endif
659     } 
660     else 
661     {
662         vdwc = &nlist[eNL_VDWQQ];
663         vdw  = &nlist[eNL_VDW];
664         coul = &nlist[eNL_QQ];
665     }
666     
667     if (!bFreeEnergy) 
668     {
669         if (iwater != esolNO) 
670         {
671             /* Loop over the atoms in the i charge group */    
672             i_atom  = i0;
673             gid     = GID(igid,jgid,ngid);
674             /* Create new i_atom for each energy group */
675             if (bDoCoul && bDoVdW)
676             {
677                 new_i_nblist(vdwc,bLR,i_atom,shift,gid);
678 #ifndef DISABLE_WATERWATER_NLIST
679                 new_i_nblist(vdwc_ww,bLR,i_atom,shift,gid);
680 #endif
681             }
682             if (bDoVdW)
683             {
684                 new_i_nblist(vdw,bLR,i_atom,shift,gid);
685             }
686             if (bDoCoul) 
687             {
688                 new_i_nblist(coul,bLR,i_atom,shift,gid);
689 #ifndef DISABLE_WATERWATER_NLIST
690                 new_i_nblist(coul_ww,bLR,i_atom,shift,gid);
691 #endif
692             }      
693           /* Loop over the j charge groups */
694             for(j=0; (j<nj); j++) 
695             {
696                 jcg=jjcg[j];
697                 
698                 if (jcg == icg)
699                 {
700                     continue;
701                 }
702
703                 jj0 = index[jcg];
704                 jwater = GET_CGINFO_SOLOPT(cginfo[jcg]);
705                 
706                 if (iwater == esolSPC && jwater == esolSPC)
707                 {
708                     /* Interaction between two SPC molecules */
709                     if (!bDoCoul)
710                     {
711                         /* VdW only - only first atoms in each water interact */
712                         add_j_to_nblist(vdw,jj0,bLR);
713                     }
714                     else 
715                     {
716 #ifdef DISABLE_WATERWATER_NLIST 
717                         /* Add entries for the three atoms - only do VdW if we need to */
718                         if (!bDoVdW)
719                         {
720                             add_j_to_nblist(coul,jj0,bLR);
721                         }
722                         else
723                         {
724                             add_j_to_nblist(vdwc,jj0,bLR);
725                         }
726                         add_j_to_nblist(coul,jj0+1,bLR);
727                         add_j_to_nblist(coul,jj0+2,bLR);            
728 #else
729                         /* One entry for the entire water-water interaction */
730                         if (!bDoVdW)
731                         {
732                             add_j_to_nblist(coul_ww,jj0,bLR);
733                         }
734                         else
735                         {
736                             add_j_to_nblist(vdwc_ww,jj0,bLR);
737                         }
738 #endif
739                     }  
740                 } 
741                 else if (iwater == esolTIP4P && jwater == esolTIP4P)
742                 {
743                     /* Interaction between two TIP4p molecules */
744                     if (!bDoCoul)
745                     {
746                         /* VdW only - only first atoms in each water interact */
747                         add_j_to_nblist(vdw,jj0,bLR);
748                     }
749                     else 
750                     {
751 #ifdef DISABLE_WATERWATER_NLIST 
752                         /* Add entries for the four atoms - only do VdW if we need to */
753                         if (bDoVdW)
754                         {
755                             add_j_to_nblist(vdw,jj0,bLR);
756                         }
757                         add_j_to_nblist(coul,jj0+1,bLR);
758                         add_j_to_nblist(coul,jj0+2,bLR);            
759                         add_j_to_nblist(coul,jj0+3,bLR);            
760 #else
761                         /* One entry for the entire water-water interaction */
762                         if (!bDoVdW)
763                         {
764                             add_j_to_nblist(coul_ww,jj0,bLR);
765                         }
766                         else
767                         {
768                             add_j_to_nblist(vdwc_ww,jj0,bLR);
769                         }
770 #endif
771                     }                                   
772                 }
773                 else 
774                 {
775                     /* j charge group is not water, but i is.
776                      * Add entries to the water-other_atom lists; the geometry of the water
777                      * molecule doesn't matter - that is taken care of in the nonbonded kernel,
778                      * so we don't care if it is SPC or TIP4P...
779                      */
780                     
781                     jj1 = index[jcg+1];
782                     
783                     if (!bDoVdW) 
784                     {
785                         for(jj=jj0; (jj<jj1); jj++) 
786                         {
787                             if (charge[jj] != 0)
788                             {
789                                 add_j_to_nblist(coul,jj,bLR);
790                             }
791                         }
792                     }
793                     else if (!bDoCoul)
794                     {
795                         for(jj=jj0; (jj<jj1); jj++)
796                         {
797                             if (bHaveVdW[type[jj]])
798                             {
799                                 add_j_to_nblist(vdw,jj,bLR);
800                             }
801                         }
802                     }
803                     else 
804                     {
805                         /* _charge_ _groups_ interact with both coulomb and LJ */
806                         /* Check which atoms we should add to the lists!       */
807                         for(jj=jj0; (jj<jj1); jj++) 
808                         {
809                             if (bHaveVdW[type[jj]]) 
810                             {
811                                 if (charge[jj] != 0)
812                                 {
813                                     add_j_to_nblist(vdwc,jj,bLR);
814                                 }
815                                 else
816                                 {
817                                     add_j_to_nblist(vdw,jj,bLR);
818                                 }
819                             }
820                             else if (charge[jj] != 0)
821                             {
822                                 add_j_to_nblist(coul,jj,bLR);
823                             }
824                         }
825                     }
826                 }
827             }
828             close_i_nblist(vdw); 
829             close_i_nblist(coul); 
830             close_i_nblist(vdwc);  
831 #ifndef DISABLE_WATERWATER_NLIST
832             close_i_nblist(coul_ww);
833             close_i_nblist(vdwc_ww); 
834 #endif
835         } 
836         else
837         { 
838             /* no solvent as i charge group */
839             /* Loop over the atoms in the i charge group */    
840             for(i=0; i<nicg; i++) 
841             {
842                 i_atom  = i0+i;
843                 gid     = GID(igid,jgid,ngid);
844                 qi      = charge[i_atom];
845                 
846                 /* Create new i_atom for each energy group */
847                 if (bDoVdW && bDoCoul)
848                 {
849                     new_i_nblist(vdwc,bLR,i_atom,shift,gid);
850                 }
851                 if (bDoVdW)
852                 {
853                     new_i_nblist(vdw,bLR,i_atom,shift,gid);
854                 }
855                 if (bDoCoul)
856                 {
857                     new_i_nblist(coul,bLR,i_atom,shift,gid);
858                 }
859                 bDoVdW_i  = (bDoVdW  && bHaveVdW[type[i_atom]]);
860                 bDoCoul_i = (bDoCoul && qi!=0);
861                 
862                 if (bDoVdW_i || bDoCoul_i) 
863                 {
864                     /* Loop over the j charge groups */
865                     for(j=0; (j<nj); j++) 
866                     {
867                         jcg=jjcg[j];
868                         
869                         /* Check for large charge groups */
870                         if (jcg == icg)
871                         {
872                             jj0 = i0 + i + 1;
873                         }
874                         else
875                         {
876                             jj0 = index[jcg];
877                         }
878                         
879                         jj1=index[jcg+1];
880                         /* Finally loop over the atoms in the j-charge group */ 
881                         for(jj=jj0; jj<jj1; jj++) 
882                         {
883                             bNotEx = NOTEXCL(bExcl,i,jj);
884                             
885                             if (bNotEx) 
886                             {
887                                 if (!bDoVdW_i) 
888                                 { 
889                                     if (charge[jj] != 0)
890                                     {
891                                         add_j_to_nblist(coul,jj,bLR);
892                                     }
893                                 }
894                                 else if (!bDoCoul_i) 
895                                 {
896                                     if (bHaveVdW[type[jj]])
897                                     {
898                                         add_j_to_nblist(vdw,jj,bLR);
899                                     }
900                                 }
901                                 else 
902                                 {
903                                     if (bHaveVdW[type[jj]]) 
904                                     {
905                                         if (charge[jj] != 0)
906                                         {
907                                             add_j_to_nblist(vdwc,jj,bLR);
908                                         }
909                                         else
910                                         {
911                                             add_j_to_nblist(vdw,jj,bLR);
912                                         }
913                                     } 
914                                     else if (charge[jj] != 0)
915                                     {
916                                         add_j_to_nblist(coul,jj,bLR);
917                                     }
918                                 }
919                             }
920                         }
921                     }
922                 }
923                 close_i_nblist(vdw);
924                 close_i_nblist(coul);
925                 close_i_nblist(vdwc);
926             }
927         }
928     }
929     else
930     {
931         /* we are doing free energy */
932         vdwc_free = &nlist[eNL_VDWQQ_FREE];
933         vdw_free  = &nlist[eNL_VDW_FREE];
934         coul_free = &nlist[eNL_QQ_FREE];
935         /* Loop over the atoms in the i charge group */    
936         for(i=0; i<nicg; i++) 
937         {
938             i_atom  = i0+i;
939             gid     = GID(igid,jgid,ngid);
940             qi      = charge[i_atom];
941             qiB     = chargeB[i_atom];
942             
943             /* Create new i_atom for each energy group */
944             if (bDoVdW && bDoCoul) 
945                 new_i_nblist(vdwc,bLR,i_atom,shift,gid);
946             if (bDoVdW)   
947                 new_i_nblist(vdw,bLR,i_atom,shift,gid);
948             if (bDoCoul) 
949                 new_i_nblist(coul,bLR,i_atom,shift,gid);
950             
951             new_i_nblist(vdw_free,bLR,i_atom,shift,gid);
952             new_i_nblist(coul_free,bLR,i_atom,shift,gid);
953             new_i_nblist(vdwc_free,bLR,i_atom,shift,gid);
954             
955             bDoVdW_i  = (bDoVdW  &&
956                          (bHaveVdW[type[i_atom]] || bHaveVdW[typeB[i_atom]]));
957             bDoCoul_i = (bDoCoul && (qi!=0 || qiB!=0));
958             /* For TIP4P the first atom does not have a charge,
959              * but the last three do. So we should still put an atom
960              * without LJ but with charge in the water-atom neighborlist
961              * for a TIP4p i charge group.
962              * For SPC type water the first atom has LJ and charge,
963              * so there is no such problem.
964              */
965             if (iwater == esolNO)
966             {
967                 bDoCoul_i_sol = bDoCoul_i;
968             }
969             else
970             {
971                 bDoCoul_i_sol = bDoCoul;
972             }
973             
974             if (bDoVdW_i || bDoCoul_i_sol) 
975             {
976                 /* Loop over the j charge groups */
977                 for(j=0; (j<nj); j++)
978                 {
979                     jcg=jjcg[j];
980                     
981                     /* Check for large charge groups */
982                     if (jcg == icg)
983                     {
984                         jj0 = i0 + i + 1;
985                     }
986                     else
987                     {
988                         jj0 = index[jcg];
989                     }
990                     
991                     jj1=index[jcg+1];
992                     /* Finally loop over the atoms in the j-charge group */     
993                     bFree = bPert[i_atom];
994                     for(jj=jj0; (jj<jj1); jj++) 
995                     {
996                         bFreeJ = bFree || bPert[jj];
997                         /* Complicated if, because the water H's should also
998                          * see perturbed j-particles
999                          */
1000                         if (iwater==esolNO || i==0 || bFreeJ) 
1001                         {
1002                             bNotEx = NOTEXCL(bExcl,i,jj);
1003                             
1004                             if (bNotEx) 
1005                             {
1006                                 if (bFreeJ)
1007                                 {
1008                                     if (!bDoVdW_i) 
1009                                     {
1010                                         if (charge[jj]!=0 || chargeB[jj]!=0)
1011                                         {
1012                                             add_j_to_nblist(coul_free,jj,bLR);
1013                                         }
1014                                     }
1015                                     else if (!bDoCoul_i) 
1016                                     {
1017                                         if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]])
1018                                         {
1019                                             add_j_to_nblist(vdw_free,jj,bLR);
1020                                         }
1021                                     }
1022                                     else 
1023                                     {
1024                                         if (bHaveVdW[type[jj]] || bHaveVdW[typeB[jj]]) 
1025                                         {
1026                                             if (charge[jj]!=0 || chargeB[jj]!=0)
1027                                             {
1028                                                 add_j_to_nblist(vdwc_free,jj,bLR);
1029                                             }
1030                                             else
1031                                             {
1032                                                 add_j_to_nblist(vdw_free,jj,bLR);
1033                                             }
1034                                         }
1035                                         else if (charge[jj]!=0 || chargeB[jj]!=0)
1036                                             add_j_to_nblist(coul_free,jj,bLR);
1037                                     }
1038                                 }
1039                                 else if (!bDoVdW_i) 
1040                                 { 
1041                                     /* This is done whether or not bWater is set */
1042                                     if (charge[jj] != 0)
1043                                     {
1044                                         add_j_to_nblist(coul,jj,bLR);
1045                                     }
1046                                 }
1047                                 else if (!bDoCoul_i_sol) 
1048                                 { 
1049                                     if (bHaveVdW[type[jj]])
1050                                     {
1051                                         add_j_to_nblist(vdw,jj,bLR);
1052                                     }
1053                                 }
1054                                 else 
1055                                 {
1056                                     if (bHaveVdW[type[jj]]) 
1057                                     {
1058                                         if (charge[jj] != 0)
1059                                         {
1060                                             add_j_to_nblist(vdwc,jj,bLR);
1061                                         }
1062                                         else
1063                                         {
1064                                             add_j_to_nblist(vdw,jj,bLR);
1065                                         }
1066                                     } 
1067                                     else if (charge[jj] != 0)
1068                                     {
1069                                         add_j_to_nblist(coul,jj,bLR);
1070                                     }
1071                                 }
1072                             }
1073                         }
1074                     }
1075                 }
1076             }
1077             close_i_nblist(vdw);
1078             close_i_nblist(coul);
1079             close_i_nblist(vdwc);
1080             close_i_nblist(vdw_free);
1081             close_i_nblist(coul_free);
1082             close_i_nblist(vdwc_free);
1083         }
1084     }
1085 }
1086
1087 static void
1088 put_in_list_adress(gmx_bool              bHaveVdW[],
1089                int               ngid,
1090                t_mdatoms *       md,
1091                int               icg,
1092                int               jgid,
1093                int               nj,
1094                atom_id           jjcg[],
1095                atom_id           index[],
1096                t_excl            bExcl[],
1097                int               shift,
1098                t_forcerec *      fr,
1099                gmx_bool          bLR,
1100                gmx_bool          bDoVdW,
1101                gmx_bool          bDoCoul,
1102                int               solvent_opt)
1103 {
1104     /* The a[] index has been removed,
1105      * to put it back in i_atom should be a[i0] and jj should be a[jj].
1106      */
1107     t_nblist *   vdwc;
1108     t_nblist *   vdw;
1109     t_nblist *   coul;
1110     t_nblist *   vdwc_adress  = NULL;
1111     t_nblist *   vdw_adress   = NULL;
1112     t_nblist *   coul_adress  = NULL;
1113     t_nblist *   vdwc_ww    = NULL;
1114     t_nblist *   coul_ww    = NULL;
1115
1116     int             i,j,jcg,igid,gid,nbl_ind,nbl_ind_adress;
1117     atom_id   jj,jj0,jj1,i_atom;
1118     int       i0,nicg,len;
1119
1120     int       *cginfo;
1121     int       *type,*typeB;
1122     real      *charge,*chargeB;
1123     real      *wf;
1124     real      qi,qiB,qq,rlj;
1125     gmx_bool      bFreeEnergy,bFree,bFreeJ,bNotEx,*bPert;
1126     gmx_bool      bDoVdW_i,bDoCoul_i,bDoCoul_i_sol;
1127     gmx_bool      b_hybrid;
1128     gmx_bool      j_all_atom;
1129     int       iwater,jwater;
1130     t_nblist  *nlist, *nlist_adress;
1131
1132     /* Copy some pointers */
1133     cginfo  = fr->cginfo;
1134     charge  = md->chargeA;
1135     chargeB = md->chargeB;
1136     type    = md->typeA;
1137     typeB   = md->typeB;
1138     bPert   = md->bPerturbed;
1139     wf      = md->wf;
1140
1141     /* Get atom range */
1142     i0     = index[icg];
1143     nicg   = index[icg+1]-i0;
1144
1145     /* Get the i charge group info */
1146     igid   = GET_CGINFO_GID(cginfo[icg]);
1147
1148     iwater = (solvent_opt!=esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
1149
1150     if (md->nPerturbed)
1151     {
1152         gmx_fatal(FARGS,"AdResS does not support free energy pertubation\n");
1153     }
1154
1155     /* Unpack pointers to neighbourlist structs */
1156     if (fr->nnblists == 2)
1157     {
1158         nbl_ind = 0;
1159         nbl_ind_adress = 1;
1160     }
1161     else
1162     {
1163         nbl_ind = fr->gid2nblists[GID(igid,jgid,ngid)];
1164         nbl_ind_adress = nbl_ind+fr->nnblists/2;
1165     }
1166     if (bLR)
1167     {
1168         nlist = fr->nblists[nbl_ind].nlist_lr;
1169         nlist_adress= fr->nblists[nbl_ind_adress].nlist_lr;
1170     }
1171     else
1172     {
1173         nlist = fr->nblists[nbl_ind].nlist_sr;
1174         nlist_adress = fr->nblists[nbl_ind_adress].nlist_sr;
1175     }
1176
1177
1178         vdwc = &nlist[eNL_VDWQQ];
1179         vdw  = &nlist[eNL_VDW];
1180         coul = &nlist[eNL_QQ];
1181     
1182         vdwc_adress = &nlist_adress[eNL_VDWQQ];
1183         vdw_adress  = &nlist_adress[eNL_VDW];
1184         coul_adress = &nlist_adress[eNL_QQ];
1185
1186     /* We do not support solvent optimization with AdResS for now.
1187       For this we would need hybrid solvent-other kernels */
1188     
1189             /* no solvent as i charge group */
1190             /* Loop over the atoms in the i charge group */
1191             for(i=0; i<nicg; i++)
1192             {
1193                 i_atom  = i0+i;
1194                 gid     = GID(igid,jgid,ngid);
1195                 qi      = charge[i_atom];
1196
1197                 /* Create new i_atom for each energy group */
1198                 if (bDoVdW && bDoCoul)
1199                 {
1200                     new_i_nblist(vdwc,bLR,i_atom,shift,gid);
1201                     new_i_nblist(vdwc_adress,bLR,i_atom,shift,gid);
1202
1203                 }
1204                 if (bDoVdW)
1205                 {
1206                     new_i_nblist(vdw,bLR,i_atom,shift,gid);
1207                     new_i_nblist(vdw_adress,bLR,i_atom,shift,gid);
1208
1209                 }
1210                 if (bDoCoul)
1211                 {
1212                     new_i_nblist(coul,bLR,i_atom,shift,gid);
1213                     new_i_nblist(coul_adress,bLR,i_atom,shift,gid);
1214                 }
1215                 bDoVdW_i  = (bDoVdW  && bHaveVdW[type[i_atom]]);
1216                 bDoCoul_i = (bDoCoul && qi!=0);
1217
1218                 if (bDoVdW_i || bDoCoul_i)
1219                 {
1220                     /* Loop over the j charge groups */
1221                     for(j=0; (j<nj); j++)
1222                     {
1223                         jcg=jjcg[j];
1224
1225                         /* Check for large charge groups */
1226                         if (jcg == icg)
1227                         {
1228                             jj0 = i0 + i + 1;
1229                         }
1230                         else
1231                         {
1232                             jj0 = index[jcg];
1233                         }
1234
1235                         jj1=index[jcg+1];
1236                         /* Finally loop over the atoms in the j-charge group */
1237                         for(jj=jj0; jj<jj1; jj++)
1238                         {
1239                             bNotEx = NOTEXCL(bExcl,i,jj);
1240
1241                             b_hybrid=!((wf[i_atom]==1&&wf[jj]==1)||(wf[i_atom] ==0 && wf[jj]==0));
1242
1243                             if (bNotEx)
1244                             {
1245                                 if (!bDoVdW_i)
1246                                 {
1247                                     if (charge[jj] != 0)
1248                                     {
1249                                         if(!b_hybrid){
1250                                             add_j_to_nblist(coul,jj,bLR);
1251                                         }else{
1252                                             add_j_to_nblist(coul_adress,jj,bLR);
1253                                         }
1254                                     }
1255                                 }
1256                                 else if (!bDoCoul_i)
1257                                 {
1258                                     if (bHaveVdW[type[jj]])
1259                                     {
1260                                         if(!b_hybrid){
1261                                             add_j_to_nblist(vdw,jj,bLR);
1262                                         }else{
1263                                             add_j_to_nblist(vdw_adress,jj,bLR);
1264                                         }
1265                                     }
1266                                 }
1267                                 else
1268                                 {
1269                                     if (bHaveVdW[type[jj]])
1270                                     {
1271                                         if (charge[jj] != 0)
1272                                         {
1273                                             if(!b_hybrid){
1274                                                 add_j_to_nblist(vdwc,jj,bLR);
1275                                             }else{
1276                                                 add_j_to_nblist(vdwc_adress,jj,bLR);
1277                                             }
1278                                         }
1279                                         else
1280                                         {
1281                                             if(!b_hybrid){
1282                                                 add_j_to_nblist(vdw,jj,bLR);
1283                                             }else{
1284                                                 add_j_to_nblist(vdw_adress,jj,bLR);
1285                                             }
1286
1287                                         }
1288                                     }
1289                                     else if (charge[jj] != 0)
1290                                     {
1291                                         if(!b_hybrid){
1292                                             add_j_to_nblist(coul,jj,bLR);
1293                                         }else{
1294                                             add_j_to_nblist(coul_adress,jj,bLR);
1295                                         }
1296
1297                                     }
1298                                 }
1299                             }
1300                         }
1301                     }
1302                 
1303                 close_i_nblist(vdw);
1304                 close_i_nblist(coul);
1305                 close_i_nblist(vdwc);
1306                 close_i_nblist(vdw_adress);
1307                 close_i_nblist(coul_adress);
1308                 close_i_nblist(vdwc_adress);
1309             }
1310     }
1311 }
1312
1313 static void 
1314 put_in_list_qmmm(gmx_bool              bHaveVdW[],
1315                  int               ngid,
1316                  t_mdatoms *       md,
1317                  int               icg,
1318                  int               jgid,
1319                  int               nj,
1320                  atom_id           jjcg[],
1321                  atom_id           index[],
1322                  t_excl            bExcl[],
1323                  int               shift,
1324                  t_forcerec *      fr,
1325                  gmx_bool          bLR,
1326                  gmx_bool          bDoVdW,
1327                  gmx_bool          bDoCoul,
1328                  int               solvent_opt)
1329 {
1330     t_nblist *   coul;
1331     int           i,j,jcg,igid,gid;
1332     atom_id   jj,jj0,jj1,i_atom;
1333     int       i0,nicg;
1334     gmx_bool      bNotEx;
1335     
1336     /* Get atom range */
1337     i0     = index[icg];
1338     nicg   = index[icg+1]-i0;
1339     
1340     /* Get the i charge group info */
1341     igid   = GET_CGINFO_GID(fr->cginfo[icg]);
1342     
1343     coul = &fr->QMMMlist;
1344     
1345     /* Loop over atoms in the ith charge group */
1346     for (i=0;i<nicg;i++)
1347     {
1348         i_atom = i0+i;
1349         gid    = GID(igid,jgid,ngid);
1350         /* Create new i_atom for each energy group */
1351         new_i_nblist(coul,bLR,i_atom,shift,gid);
1352         
1353         /* Loop over the j charge groups */
1354         for (j=0;j<nj;j++)
1355         {
1356             jcg=jjcg[j];
1357             
1358             /* Charge groups cannot have QM and MM atoms simultaneously */
1359             if (jcg!=icg)
1360             {
1361                 jj0 = index[jcg];
1362                 jj1 = index[jcg+1];
1363                 /* Finally loop over the atoms in the j-charge group */
1364                 for(jj=jj0; jj<jj1; jj++)
1365                 {
1366                     bNotEx = NOTEXCL(bExcl,i,jj);
1367                     if(bNotEx)
1368                         add_j_to_nblist(coul,jj,bLR);
1369                 }
1370             }
1371         }
1372         close_i_nblist(coul);
1373     }
1374 }
1375
1376 static void 
1377 put_in_list_cg(gmx_bool              bHaveVdW[],
1378                int               ngid,
1379                t_mdatoms *       md,
1380                int               icg,
1381                int               jgid,
1382                int               nj,
1383                atom_id           jjcg[],
1384                atom_id           index[],
1385                t_excl            bExcl[],
1386                int               shift,
1387                t_forcerec *      fr,
1388                gmx_bool          bLR,
1389                gmx_bool          bDoVdW,
1390                gmx_bool          bDoCoul,
1391                int               solvent_opt)
1392 {
1393     int          cginfo;
1394     int          igid,gid,nbl_ind;
1395     t_nblist *   vdwc;
1396     int          j,jcg;
1397
1398     cginfo = fr->cginfo[icg];
1399
1400     igid = GET_CGINFO_GID(cginfo);
1401     gid  = GID(igid,jgid,ngid);
1402
1403     /* Unpack pointers to neighbourlist structs */
1404     if (fr->nnblists == 1)
1405     {
1406         nbl_ind = 0;
1407     }
1408     else
1409     {
1410         nbl_ind = fr->gid2nblists[gid];
1411     }
1412     if (bLR)
1413     {
1414         vdwc = &fr->nblists[nbl_ind].nlist_lr[eNL_VDWQQ];
1415     }
1416     else
1417     {
1418         vdwc = &fr->nblists[nbl_ind].nlist_sr[eNL_VDWQQ];
1419     }
1420
1421     /* Make a new neighbor list for charge group icg.
1422      * Currently simply one neighbor list is made with LJ and Coulomb.
1423      * If required, zero interactions could be removed here
1424      * or in the force loop.
1425      */
1426     new_i_nblist(vdwc,bLR,index[icg],shift,gid);
1427     vdwc->iinr_end[vdwc->nri] = index[icg+1];
1428
1429     for(j=0; (j<nj); j++) 
1430     {
1431         jcg = jjcg[j];
1432         /* Skip the icg-icg pairs if all self interactions are excluded */
1433         if (!(jcg == icg && GET_CGINFO_EXCL_INTRA(cginfo)))
1434         {
1435             /* Here we add the j charge group jcg to the list,
1436              * exclusions are also added to the list.
1437              */
1438             add_j_to_nblist_cg(vdwc,index[jcg],index[jcg+1],bExcl,icg==jcg,bLR);
1439         }
1440     }
1441
1442     close_i_nblist(vdwc);  
1443 }
1444
1445 static void setexcl(atom_id start,atom_id end,t_blocka *excl,gmx_bool b,
1446                     t_excl bexcl[])
1447 {
1448     atom_id i,k;
1449     
1450     if (b)
1451     {
1452         for(i=start; i<end; i++)
1453         {
1454             for(k=excl->index[i]; k<excl->index[i+1]; k++)
1455             {
1456                 SETEXCL(bexcl,i-start,excl->a[k]);
1457             }
1458         }
1459     }
1460     else
1461     {
1462         for(i=start; i<end; i++)
1463         {
1464             for(k=excl->index[i]; k<excl->index[i+1]; k++)
1465             {
1466                 RMEXCL(bexcl,i-start,excl->a[k]);
1467             }
1468         }
1469     }
1470 }
1471
1472 int calc_naaj(int icg,int cgtot)
1473 {
1474     int naaj;
1475     
1476     if ((cgtot % 2) == 1)
1477     {
1478         /* Odd number of charge groups, easy */
1479         naaj = 1 + (cgtot/2);
1480     }
1481     else if ((cgtot % 4) == 0)
1482     {
1483     /* Multiple of four is hard */
1484         if (icg < cgtot/2)
1485         {
1486             if ((icg % 2) == 0)
1487             {
1488                 naaj=1+(cgtot/2);
1489             }
1490             else
1491             {
1492                 naaj=cgtot/2;
1493             }
1494         }
1495         else
1496         {
1497             if ((icg % 2) == 1)
1498             {
1499                 naaj=1+(cgtot/2);
1500             }
1501             else
1502             {
1503                 naaj=cgtot/2;
1504             }
1505         }
1506     }
1507     else
1508     {
1509         /* cgtot/2 = odd */
1510         if ((icg % 2) == 0)
1511         {
1512             naaj=1+(cgtot/2);
1513         }
1514         else
1515         {
1516             naaj=cgtot/2;
1517         }
1518     }
1519 #ifdef DEBUG
1520     fprintf(log,"naaj=%d\n",naaj);
1521 #endif
1522
1523     return naaj;
1524 }
1525
1526 /************************************************
1527  *
1528  *  S I M P L E      C O R E     S T U F F
1529  *
1530  ************************************************/
1531
1532 static real calc_image_tric(rvec xi,rvec xj,matrix box,
1533                             rvec b_inv,int *shift)
1534 {
1535     /* This code assumes that the cut-off is smaller than
1536      * a half times the smallest diagonal element of the box.
1537      */
1538     const real h25=2.5;
1539     real dx,dy,dz;
1540     real r2;
1541     int  tx,ty,tz;
1542     
1543     /* Compute diff vector */
1544     dz = xj[ZZ] - xi[ZZ];
1545     dy = xj[YY] - xi[YY];
1546     dx = xj[XX] - xi[XX];
1547     
1548   /* Perform NINT operation, using trunc operation, therefore
1549    * we first add 2.5 then subtract 2 again
1550    */
1551     tz = dz*b_inv[ZZ] + h25;
1552     tz -= 2;
1553     dz -= tz*box[ZZ][ZZ];
1554     dy -= tz*box[ZZ][YY];
1555     dx -= tz*box[ZZ][XX];
1556
1557     ty = dy*b_inv[YY] + h25;
1558     ty -= 2;
1559     dy -= ty*box[YY][YY];
1560     dx -= ty*box[YY][XX];
1561     
1562     tx = dx*b_inv[XX]+h25;
1563     tx -= 2;
1564     dx -= tx*box[XX][XX];
1565   
1566     /* Distance squared */
1567     r2 = (dx*dx) + (dy*dy) + (dz*dz);
1568
1569     *shift = XYZ2IS(tx,ty,tz);
1570
1571     return r2;
1572 }
1573
1574 static real calc_image_rect(rvec xi,rvec xj,rvec box_size,
1575                             rvec b_inv,int *shift)
1576 {
1577     const real h15=1.5;
1578     real ddx,ddy,ddz;
1579     real dx,dy,dz;
1580     real r2;
1581     int  tx,ty,tz;
1582     
1583     /* Compute diff vector */
1584     dx = xj[XX] - xi[XX];
1585     dy = xj[YY] - xi[YY];
1586     dz = xj[ZZ] - xi[ZZ];
1587   
1588     /* Perform NINT operation, using trunc operation, therefore
1589      * we first add 1.5 then subtract 1 again
1590      */
1591     tx = dx*b_inv[XX] + h15;
1592     ty = dy*b_inv[YY] + h15;
1593     tz = dz*b_inv[ZZ] + h15;
1594     tx--;
1595     ty--;
1596     tz--;
1597     
1598     /* Correct diff vector for translation */
1599     ddx = tx*box_size[XX] - dx;
1600     ddy = ty*box_size[YY] - dy;
1601     ddz = tz*box_size[ZZ] - dz;
1602     
1603     /* Distance squared */
1604     r2 = (ddx*ddx) + (ddy*ddy) + (ddz*ddz);
1605     
1606     *shift = XYZ2IS(tx,ty,tz);
1607     
1608     return r2;
1609 }
1610
1611 static void add_simple(t_ns_buf *nsbuf,int nrj,atom_id cg_j,
1612                        gmx_bool bHaveVdW[],int ngid,t_mdatoms *md,
1613                        int icg,int jgid,t_block *cgs,t_excl bexcl[],
1614                        int shift,t_forcerec *fr,put_in_list_t *put_in_list)
1615 {
1616     if (nsbuf->nj + nrj > MAX_CG)
1617     {
1618         put_in_list(bHaveVdW,ngid,md,icg,jgid,nsbuf->ncg,nsbuf->jcg,
1619                     cgs->index,bexcl,shift,fr,FALSE,TRUE,TRUE,fr->solvent_opt);
1620         /* Reset buffer contents */
1621         nsbuf->ncg = nsbuf->nj = 0;
1622     }
1623     nsbuf->jcg[nsbuf->ncg++] = cg_j;
1624     nsbuf->nj += nrj;
1625 }
1626
1627 static void ns_inner_tric(rvec x[],int icg,int *i_egp_flags,
1628                           int njcg,atom_id jcg[],
1629                           matrix box,rvec b_inv,real rcut2,
1630                           t_block *cgs,t_ns_buf **ns_buf,
1631                           gmx_bool bHaveVdW[],int ngid,t_mdatoms *md,
1632                           t_excl bexcl[],t_forcerec *fr,
1633                           put_in_list_t *put_in_list)
1634 {
1635     int      shift;
1636     int      j,nrj,jgid;
1637     int      *cginfo=fr->cginfo;
1638     atom_id  cg_j,*cgindex;
1639     t_ns_buf *nsbuf;
1640     
1641     cgindex = cgs->index;
1642     shift   = CENTRAL;
1643     for(j=0; (j<njcg); j++)
1644     {
1645         cg_j   = jcg[j];
1646         nrj    = cgindex[cg_j+1]-cgindex[cg_j];
1647         if (calc_image_tric(x[icg],x[cg_j],box,b_inv,&shift) < rcut2)
1648         {
1649             jgid  = GET_CGINFO_GID(cginfo[cg_j]);
1650             if (!(i_egp_flags[jgid] & EGP_EXCL))
1651             {
1652                 add_simple(&ns_buf[jgid][shift],nrj,cg_j,
1653                            bHaveVdW,ngid,md,icg,jgid,cgs,bexcl,shift,fr,
1654                            put_in_list);
1655             }
1656         }
1657     }
1658 }
1659
1660 static void ns_inner_rect(rvec x[],int icg,int *i_egp_flags,
1661                           int njcg,atom_id jcg[],
1662                           gmx_bool bBox,rvec box_size,rvec b_inv,real rcut2,
1663                           t_block *cgs,t_ns_buf **ns_buf,
1664                           gmx_bool bHaveVdW[],int ngid,t_mdatoms *md,
1665                           t_excl bexcl[],t_forcerec *fr,
1666                           put_in_list_t *put_in_list)
1667 {
1668     int      shift;
1669     int      j,nrj,jgid;
1670     int      *cginfo=fr->cginfo;
1671     atom_id  cg_j,*cgindex;
1672     t_ns_buf *nsbuf;
1673
1674     cgindex = cgs->index;
1675     if (bBox)
1676     {
1677         shift = CENTRAL;
1678         for(j=0; (j<njcg); j++)
1679         {
1680             cg_j   = jcg[j];
1681             nrj    = cgindex[cg_j+1]-cgindex[cg_j];
1682             if (calc_image_rect(x[icg],x[cg_j],box_size,b_inv,&shift) < rcut2)
1683             {
1684                 jgid  = GET_CGINFO_GID(cginfo[cg_j]);
1685                 if (!(i_egp_flags[jgid] & EGP_EXCL))
1686                 {
1687                     add_simple(&ns_buf[jgid][shift],nrj,cg_j,
1688                                bHaveVdW,ngid,md,icg,jgid,cgs,bexcl,shift,fr,
1689                                put_in_list);
1690                 }
1691             }
1692         }
1693     }
1694     else
1695     {
1696         for(j=0; (j<njcg); j++)
1697         {
1698             cg_j   = jcg[j];
1699             nrj    = cgindex[cg_j+1]-cgindex[cg_j];
1700             if ((rcut2 == 0) || (distance2(x[icg],x[cg_j]) < rcut2)) {
1701                 jgid  = GET_CGINFO_GID(cginfo[cg_j]);
1702                 if (!(i_egp_flags[jgid] & EGP_EXCL))
1703                 {
1704                     add_simple(&ns_buf[jgid][CENTRAL],nrj,cg_j,
1705                                bHaveVdW,ngid,md,icg,jgid,cgs,bexcl,CENTRAL,fr,
1706                                put_in_list);
1707                 }
1708             }
1709         }
1710     }
1711 }
1712
1713 /* ns_simple_core needs to be adapted for QMMM still 2005 */
1714
1715 static int ns_simple_core(t_forcerec *fr,
1716                           gmx_localtop_t *top,
1717                           t_mdatoms *md,
1718                           matrix box,rvec box_size,
1719                           t_excl bexcl[],atom_id *aaj,
1720                           int ngid,t_ns_buf **ns_buf,
1721                           put_in_list_t *put_in_list,gmx_bool bHaveVdW[])
1722 {
1723     int      naaj,k;
1724     real     rlist2;
1725     int      nsearch,icg,jcg,igid,i0,nri,nn;
1726     int      *cginfo;
1727     t_ns_buf *nsbuf;
1728     /* atom_id  *i_atoms; */
1729     t_block  *cgs=&(top->cgs);
1730     t_blocka *excl=&(top->excls);
1731     rvec     b_inv;
1732     int      m;
1733     gmx_bool     bBox,bTriclinic;
1734     int      *i_egp_flags;
1735     
1736     rlist2 = sqr(fr->rlist);
1737     
1738     bBox = (fr->ePBC != epbcNONE);
1739     if (bBox)
1740     {
1741         for(m=0; (m<DIM); m++)
1742         {
1743             b_inv[m] = divide_err(1.0,box_size[m]);
1744         }
1745         bTriclinic = TRICLINIC(box);
1746     }
1747     else
1748     {
1749         bTriclinic = FALSE;
1750     }
1751     
1752     cginfo = fr->cginfo;
1753     
1754     nsearch=0;
1755     for (icg=fr->cg0; (icg<fr->hcg); icg++)
1756     {
1757         /*
1758           i0        = cgs->index[icg];
1759           nri       = cgs->index[icg+1]-i0;
1760           i_atoms   = &(cgs->a[i0]);
1761           i_eg_excl = fr->eg_excl + ngid*md->cENER[*i_atoms];
1762           setexcl(nri,i_atoms,excl,TRUE,bexcl);
1763         */
1764         igid = GET_CGINFO_GID(cginfo[icg]);
1765         i_egp_flags = fr->egp_flags + ngid*igid;
1766         setexcl(cgs->index[icg],cgs->index[icg+1],excl,TRUE,bexcl);
1767         
1768         naaj=calc_naaj(icg,cgs->nr);
1769         if (bTriclinic)
1770         {
1771             ns_inner_tric(fr->cg_cm,icg,i_egp_flags,naaj,&(aaj[icg]),
1772                           box,b_inv,rlist2,cgs,ns_buf,
1773                           bHaveVdW,ngid,md,bexcl,fr,put_in_list);
1774         }
1775         else
1776         {
1777             ns_inner_rect(fr->cg_cm,icg,i_egp_flags,naaj,&(aaj[icg]),
1778                           bBox,box_size,b_inv,rlist2,cgs,ns_buf,
1779                           bHaveVdW,ngid,md,bexcl,fr,put_in_list);
1780         }
1781         nsearch += naaj;
1782         
1783         for(nn=0; (nn<ngid); nn++)
1784         {
1785             for(k=0; (k<SHIFTS); k++)
1786             {
1787                 nsbuf = &(ns_buf[nn][k]);
1788                 if (nsbuf->ncg > 0)
1789                 {
1790                     put_in_list(bHaveVdW,ngid,md,icg,nn,nsbuf->ncg,nsbuf->jcg,
1791                                 cgs->index,bexcl,k,fr,FALSE,TRUE,TRUE,fr->solvent_opt);
1792                     nsbuf->ncg=nsbuf->nj=0;
1793                 }
1794             }
1795         }
1796         /* setexcl(nri,i_atoms,excl,FALSE,bexcl); */
1797         setexcl(cgs->index[icg],cgs->index[icg+1],excl,FALSE,bexcl);
1798     }
1799     close_neighbor_lists(fr,FALSE);
1800     
1801     return nsearch;
1802 }
1803
1804 /************************************************
1805  *
1806  *    N S 5     G R I D     S T U F F
1807  *
1808  ************************************************/
1809
1810 static inline void get_dx(int Nx,real gridx,real rc2,int xgi,real x,
1811                           int *dx0,int *dx1,real *dcx2)
1812 {
1813     real dcx,tmp;
1814     int  xgi0,xgi1,i;
1815     
1816     if (xgi < 0)
1817     {
1818         *dx0 = 0;
1819         xgi0 = -1;
1820         *dx1 = -1;
1821         xgi1 = 0;
1822     }
1823     else if (xgi >= Nx)
1824     {
1825         *dx0 = Nx;
1826         xgi0 = Nx-1;
1827         *dx1 = Nx-1;
1828         xgi1 = Nx;
1829     }
1830     else
1831     {
1832         dcx2[xgi] = 0;
1833         *dx0 = xgi;
1834         xgi0 = xgi-1;
1835         *dx1 = xgi;
1836         xgi1 = xgi+1;
1837     }
1838     
1839     for(i=xgi0; i>=0; i--)
1840     {
1841         dcx = (i+1)*gridx-x;
1842         tmp = dcx*dcx;
1843         if (tmp >= rc2)
1844             break;
1845         *dx0 = i;
1846         dcx2[i] = tmp;
1847     }
1848     for(i=xgi1; i<Nx; i++)
1849     {
1850         dcx = i*gridx-x;
1851         tmp = dcx*dcx;
1852         if (tmp >= rc2)
1853         {
1854             break;
1855         }
1856         *dx1 = i;
1857         dcx2[i] = tmp;
1858     }
1859 }
1860
1861 static inline void get_dx_dd(int Nx,real gridx,real rc2,int xgi,real x,
1862                              int ncpddc,int shift_min,int shift_max,
1863                              int *g0,int *g1,real *dcx2)
1864 {
1865     real dcx,tmp;
1866     int  g_min,g_max,shift_home;
1867     
1868     if (xgi < 0)
1869     {
1870         g_min = 0;
1871         g_max = Nx - 1;
1872         *g0   = 0;
1873         *g1   = -1;
1874     }
1875     else if (xgi >= Nx)
1876     {
1877         g_min = 0;
1878         g_max = Nx - 1;
1879         *g0   = Nx;
1880         *g1   = Nx - 1;
1881     }
1882     else
1883     {
1884         if (ncpddc == 0)
1885         {
1886             g_min = 0;
1887             g_max = Nx - 1;
1888         }
1889         else
1890         {
1891             if (xgi < ncpddc)
1892             {
1893                 shift_home = 0;
1894             }
1895             else
1896             {
1897                 shift_home = -1;
1898             }
1899             g_min = (shift_min == shift_home ? 0          : ncpddc);
1900             g_max = (shift_max == shift_home ? ncpddc - 1 : Nx - 1);
1901         }
1902         if (shift_min > 0)
1903         {
1904             *g0 = g_min;
1905             *g1 = g_min - 1;
1906         }
1907         else if (shift_max < 0)
1908         {
1909             *g0 = g_max + 1;
1910             *g1 = g_max;
1911         }
1912         else
1913         {
1914             *g0 = xgi;
1915             *g1 = xgi;
1916             dcx2[xgi] = 0;
1917         }
1918     }
1919     
1920     while (*g0 > g_min)
1921     {
1922         /* Check one grid cell down */
1923         dcx = ((*g0 - 1) + 1)*gridx - x;
1924         tmp = dcx*dcx;
1925         if (tmp >= rc2)
1926         {
1927             break;
1928         }
1929         (*g0)--;
1930         dcx2[*g0] = tmp;
1931     }
1932     
1933     while (*g1 < g_max)
1934     {
1935         /* Check one grid cell up */
1936         dcx = (*g1 + 1)*gridx - x;
1937         tmp = dcx*dcx;
1938         if (tmp >= rc2)
1939         {
1940             break;
1941         }
1942         (*g1)++;
1943         dcx2[*g1] = tmp;
1944     }
1945 }
1946
1947
1948 #define sqr(x) ((x)*(x))
1949 #define calc_dx2(XI,YI,ZI,y) (sqr(XI-y[XX]) + sqr(YI-y[YY]) + sqr(ZI-y[ZZ]))
1950 #define calc_cyl_dx2(XI,YI,y) (sqr(XI-y[XX]) + sqr(YI-y[YY]))
1951 /****************************************************
1952  *
1953  *    F A S T   N E I G H B O R  S E A R C H I N G
1954  *
1955  *    Optimized neighboursearching routine using grid 
1956  *    at least 1x1x1, see GROMACS manual
1957  *
1958  ****************************************************/
1959
1960
1961 static void get_cutoff2(t_forcerec *fr,gmx_bool bDoLongRange,
1962                         real *rvdw2,real *rcoul2,
1963                         real *rs2,real *rm2,real *rl2)
1964 {
1965     *rs2 = sqr(fr->rlist);
1966
1967     if (bDoLongRange && fr->bTwinRange)
1968     {
1969         /* The VdW and elec. LR cut-off's could be different,
1970          * so we can not simply set them to rlistlong.
1971          */
1972         if (EVDW_MIGHT_BE_ZERO_AT_CUTOFF(fr->vdwtype) &&
1973             fr->rvdw > fr->rlist)
1974         {
1975             *rvdw2  = sqr(fr->rlistlong);
1976         }
1977         else
1978         {
1979             *rvdw2  = sqr(fr->rvdw);
1980         }
1981         if (EEL_MIGHT_BE_ZERO_AT_CUTOFF(fr->eeltype) &&
1982             fr->rcoulomb > fr->rlist)
1983         {
1984             *rcoul2 = sqr(fr->rlistlong);
1985         }
1986         else
1987         {
1988             *rcoul2 = sqr(fr->rcoulomb);
1989         }
1990     }
1991     else
1992     {
1993         /* Workaround for a gcc -O3 or -ffast-math problem */
1994         *rvdw2  = *rs2;
1995         *rcoul2 = *rs2;
1996     }
1997     *rm2 = min(*rvdw2,*rcoul2);
1998     *rl2 = max(*rvdw2,*rcoul2);
1999 }
2000
2001 static void init_nsgrid_lists(t_forcerec *fr,int ngid,gmx_ns_t *ns)
2002 {
2003     real rvdw2,rcoul2,rs2,rm2,rl2;
2004     int j;
2005
2006     get_cutoff2(fr,TRUE,&rvdw2,&rcoul2,&rs2,&rm2,&rl2);
2007
2008     /* Short range buffers */
2009     snew(ns->nl_sr,ngid);
2010     /* Counters */
2011     snew(ns->nsr,ngid);
2012     snew(ns->nlr_ljc,ngid);
2013     snew(ns->nlr_one,ngid);
2014     
2015     /* Always allocate both list types, since rcoulomb might now change with PME load balancing */
2016     /* Long range VdW and Coul buffers */
2017     snew(ns->nl_lr_ljc,ngid);
2018     /* Long range VdW or Coul only buffers */
2019     snew(ns->nl_lr_one,ngid);
2020
2021     for(j=0; (j<ngid); j++) {
2022         snew(ns->nl_sr[j],MAX_CG);
2023         snew(ns->nl_lr_ljc[j],MAX_CG);
2024         snew(ns->nl_lr_one[j],MAX_CG);
2025     }
2026     if (debug)
2027     {
2028         fprintf(debug,
2029                 "ns5_core: rs2 = %g, rm2 = %g, rl2 = %g (nm^2)\n",
2030                 rs2,rm2,rl2);
2031     }
2032 }
2033
2034 static int nsgrid_core(FILE *log,t_commrec *cr,t_forcerec *fr,
2035                        matrix box,rvec box_size,int ngid,
2036                        gmx_localtop_t *top,
2037                        t_grid *grid,rvec x[],
2038                        t_excl bexcl[],gmx_bool *bExcludeAlleg,
2039                        t_nrnb *nrnb,t_mdatoms *md,
2040                        real *lambda,real *dvdlambda,
2041                        gmx_grppairener_t *grppener,
2042                        put_in_list_t *put_in_list,
2043                        gmx_bool bHaveVdW[],
2044                        gmx_bool bDoLongRange,gmx_bool bMakeQMMMnblist)
2045 {
2046     gmx_ns_t *ns;
2047     atom_id **nl_lr_ljc,**nl_lr_one,**nl_sr;
2048     int     *nlr_ljc,*nlr_one,*nsr;
2049     gmx_domdec_t *dd=NULL;
2050     t_block *cgs=&(top->cgs);
2051     int     *cginfo=fr->cginfo;
2052     /* atom_id *i_atoms,*cgsindex=cgs->index; */
2053     ivec    sh0,sh1,shp;
2054     int     cell_x,cell_y,cell_z;
2055     int     d,tx,ty,tz,dx,dy,dz,cj;
2056 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2057     int     zsh_ty,zsh_tx,ysh_tx;
2058 #endif
2059     int     dx0,dx1,dy0,dy1,dz0,dz1;
2060     int     Nx,Ny,Nz,shift=-1,j,nrj,nns,nn=-1;
2061     real    gridx,gridy,gridz,grid_x,grid_y,grid_z;
2062     real    *dcx2,*dcy2,*dcz2;
2063     int     zgi,ygi,xgi;
2064     int     cg0,cg1,icg=-1,cgsnr,i0,igid,nri,naaj,max_jcg;
2065     int     jcg0,jcg1,jjcg,cgj0,jgid;
2066     int     *grida,*gridnra,*gridind;
2067     gmx_bool    rvdw_lt_rcoul,rcoul_lt_rvdw;
2068     rvec    xi,*cgcm,grid_offset;
2069     real    r2,rs2,rvdw2,rcoul2,rm2,rl2,XI,YI,ZI,dcx,dcy,dcz,tmp1,tmp2;
2070     int     *i_egp_flags;
2071     gmx_bool    bDomDec,bTriclinicX,bTriclinicY;
2072     ivec    ncpddc;
2073     
2074     ns = &fr->ns;
2075     
2076     bDomDec = DOMAINDECOMP(cr);
2077     if (bDomDec)
2078     {
2079         dd = cr->dd;
2080     }
2081     
2082     bTriclinicX = ((YY < grid->npbcdim &&
2083                     (!bDomDec || dd->nc[YY]==1) && box[YY][XX] != 0) ||
2084                    (ZZ < grid->npbcdim &&
2085                     (!bDomDec || dd->nc[ZZ]==1) && box[ZZ][XX] != 0));
2086     bTriclinicY =  (ZZ < grid->npbcdim &&
2087                     (!bDomDec || dd->nc[ZZ]==1) && box[ZZ][YY] != 0);
2088     
2089     cgsnr    = cgs->nr;
2090
2091     get_cutoff2(fr,bDoLongRange,&rvdw2,&rcoul2,&rs2,&rm2,&rl2);
2092
2093     rvdw_lt_rcoul = (rvdw2 >= rcoul2);
2094     rcoul_lt_rvdw = (rcoul2 >= rvdw2);
2095     
2096     if (bMakeQMMMnblist)
2097     {
2098         rm2 = rl2;
2099         rs2 = rl2;
2100     }
2101
2102     nl_sr     = ns->nl_sr;
2103     nsr       = ns->nsr;
2104     nl_lr_ljc = ns->nl_lr_ljc;
2105     nl_lr_one = ns->nl_lr_one;
2106     nlr_ljc   = ns->nlr_ljc;
2107     nlr_one   = ns->nlr_one;
2108     
2109     /* Unpack arrays */
2110     cgcm    = fr->cg_cm;
2111     Nx      = grid->n[XX];
2112     Ny      = grid->n[YY];
2113     Nz      = grid->n[ZZ];
2114     grida   = grid->a;
2115     gridind = grid->index;
2116     gridnra = grid->nra;
2117     nns     = 0;
2118     
2119     gridx      = grid->cell_size[XX];
2120     gridy      = grid->cell_size[YY];
2121     gridz      = grid->cell_size[ZZ];
2122     grid_x     = 1/gridx;
2123     grid_y     = 1/gridy;
2124     grid_z     = 1/gridz;
2125     copy_rvec(grid->cell_offset,grid_offset);
2126     copy_ivec(grid->ncpddc,ncpddc);
2127     dcx2       = grid->dcx2;
2128     dcy2       = grid->dcy2;
2129     dcz2       = grid->dcz2;
2130     
2131 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2132     zsh_ty = floor(-box[ZZ][YY]/box[YY][YY]+0.5);
2133     zsh_tx = floor(-box[ZZ][XX]/box[XX][XX]+0.5);
2134     ysh_tx = floor(-box[YY][XX]/box[XX][XX]+0.5);
2135     if (zsh_tx!=0 && ysh_tx!=0)
2136     {
2137         /* This could happen due to rounding, when both ratios are 0.5 */
2138         ysh_tx = 0;
2139     }
2140 #endif
2141     
2142     debug_gmx();
2143
2144     if (fr->n_tpi)
2145     {
2146         /* We only want a list for the test particle */
2147         cg0 = cgsnr - 1;
2148     }
2149     else
2150     {
2151         cg0 = grid->icg0;
2152     }
2153     cg1 = grid->icg1;
2154
2155     /* Set the shift range */
2156     for(d=0; d<DIM; d++)
2157     {
2158         sh0[d] = -1;
2159         sh1[d] = 1;
2160         /* Check if we need periodicity shifts.
2161          * Without PBC or with domain decomposition we don't need them.
2162          */
2163         if (d >= ePBC2npbcdim(fr->ePBC) || (bDomDec && dd->nc[d] > 1))
2164         {
2165             shp[d] = 0;
2166         }
2167         else
2168         {
2169             if (d == XX &&
2170                 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < sqrt(rl2))
2171             {
2172                 shp[d] = 2;
2173             }
2174             else
2175             {
2176                 shp[d] = 1;
2177             }
2178         }
2179     }
2180     
2181     /* Loop over charge groups */
2182     for(icg=cg0; (icg < cg1); icg++)
2183     {
2184         igid = GET_CGINFO_GID(cginfo[icg]);
2185         /* Skip this charge group if all energy groups are excluded! */
2186         if (bExcludeAlleg[igid])
2187         {
2188             continue;
2189         }
2190         
2191         i0   = cgs->index[icg];
2192         
2193         if (bMakeQMMMnblist)
2194         { 
2195             /* Skip this charge group if it is not a QM atom while making a
2196              * QM/MM neighbourlist
2197              */
2198             if (md->bQM[i0]==FALSE)
2199             {
2200                 continue; /* MM particle, go to next particle */ 
2201             }
2202             
2203             /* Compute the number of charge groups that fall within the control
2204              * of this one (icg)
2205              */
2206             naaj    = calc_naaj(icg,cgsnr);
2207             jcg0    = icg;
2208             jcg1    = icg + naaj;
2209             max_jcg = cgsnr;       
2210         } 
2211         else
2212         { 
2213             /* make a normal neighbourlist */
2214             
2215             if (bDomDec)
2216             {
2217                 /* Get the j charge-group and dd cell shift ranges */
2218                 dd_get_ns_ranges(cr->dd,icg,&jcg0,&jcg1,sh0,sh1);
2219                 max_jcg = 0;
2220             }
2221             else
2222             {
2223                 /* Compute the number of charge groups that fall within the control
2224                  * of this one (icg)
2225                  */
2226                 naaj = calc_naaj(icg,cgsnr);
2227                 jcg0 = icg;
2228                 jcg1 = icg + naaj;
2229                 
2230                 if (fr->n_tpi)
2231                 {
2232                     /* The i-particle is awlways the test particle,
2233                      * so we want all j-particles
2234                      */
2235                     max_jcg = cgsnr - 1;
2236                 }
2237                 else
2238                 {
2239                     max_jcg  = jcg1 - cgsnr;
2240                 }
2241             }
2242         }
2243         
2244         i_egp_flags = fr->egp_flags + igid*ngid;
2245         
2246         /* Set the exclusions for the atoms in charge group icg using a bitmask */
2247         setexcl(i0,cgs->index[icg+1],&top->excls,TRUE,bexcl);
2248         
2249         ci2xyz(grid,icg,&cell_x,&cell_y,&cell_z);
2250         
2251         /* Changed iicg to icg, DvdS 990115 
2252          * (but see consistency check above, DvdS 990330) 
2253          */
2254 #ifdef NS5DB
2255         fprintf(log,"icg=%5d, naaj=%5d, cell %d %d %d\n",
2256                 icg,naaj,cell_x,cell_y,cell_z);
2257 #endif
2258         /* Loop over shift vectors in three dimensions */
2259         for (tz=-shp[ZZ]; tz<=shp[ZZ]; tz++)
2260         {
2261             ZI = cgcm[icg][ZZ]+tz*box[ZZ][ZZ];
2262             /* Calculate range of cells in Z direction that have the shift tz */
2263             zgi = cell_z + tz*Nz;
2264 #define FAST_DD_NS
2265 #ifndef FAST_DD_NS
2266             get_dx(Nz,gridz,rl2,zgi,ZI,&dz0,&dz1,dcz2);
2267 #else
2268             get_dx_dd(Nz,gridz,rl2,zgi,ZI-grid_offset[ZZ],
2269                       ncpddc[ZZ],sh0[ZZ],sh1[ZZ],&dz0,&dz1,dcz2);
2270 #endif
2271             if (dz0 > dz1)
2272             {
2273                 continue;
2274             }
2275             for (ty=-shp[YY]; ty<=shp[YY]; ty++)
2276             {
2277                 YI = cgcm[icg][YY]+ty*box[YY][YY]+tz*box[ZZ][YY];
2278                 /* Calculate range of cells in Y direction that have the shift ty */
2279                 if (bTriclinicY)
2280                 {
2281                     ygi = (int)(Ny + (YI - grid_offset[YY])*grid_y) - Ny;
2282                 }
2283                 else
2284                 {
2285                     ygi = cell_y + ty*Ny;
2286                 }
2287 #ifndef FAST_DD_NS
2288                 get_dx(Ny,gridy,rl2,ygi,YI,&dy0,&dy1,dcy2);
2289 #else
2290                 get_dx_dd(Ny,gridy,rl2,ygi,YI-grid_offset[YY],
2291                           ncpddc[YY],sh0[YY],sh1[YY],&dy0,&dy1,dcy2);
2292 #endif
2293                 if (dy0 > dy1)
2294                 {
2295                     continue;
2296                 }
2297                 for (tx=-shp[XX]; tx<=shp[XX]; tx++)
2298                 {
2299                     XI = cgcm[icg][XX]+tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
2300                     /* Calculate range of cells in X direction that have the shift tx */
2301                     if (bTriclinicX)
2302                     {
2303                         xgi = (int)(Nx + (XI - grid_offset[XX])*grid_x) - Nx;
2304                     }
2305                     else
2306                     {
2307                         xgi = cell_x + tx*Nx;
2308                     }
2309 #ifndef FAST_DD_NS
2310                     get_dx(Nx,gridx,rl2,xgi*Nx,XI,&dx0,&dx1,dcx2);
2311 #else
2312                     get_dx_dd(Nx,gridx,rl2,xgi,XI-grid_offset[XX],
2313                               ncpddc[XX],sh0[XX],sh1[XX],&dx0,&dx1,dcx2);
2314 #endif
2315                     if (dx0 > dx1)
2316                     {
2317                         continue;
2318                     }
2319                     /* Adress: an explicit cg that has a weigthing function of 0 is excluded
2320                      *  from the neigbour list as it will not interact  */
2321                     if (fr->adress_type != eAdressOff){
2322                         if (md->wf[cgs->index[icg]]==0 && egp_explicit(fr, igid)){
2323                             continue;
2324                         }
2325                     }
2326                     /* Get shift vector */        
2327                     shift=XYZ2IS(tx,ty,tz);
2328 #ifdef NS5DB
2329                     range_check(shift,0,SHIFTS);
2330 #endif
2331                     for(nn=0; (nn<ngid); nn++)
2332                     {
2333                         nsr[nn]      = 0;
2334                         nlr_ljc[nn]  = 0;
2335                         nlr_one[nn] = 0;
2336                     }
2337 #ifdef NS5DB
2338                     fprintf(log,"shift: %2d, dx0,1: %2d,%2d, dy0,1: %2d,%2d, dz0,1: %2d,%2d\n",
2339                             shift,dx0,dx1,dy0,dy1,dz0,dz1);
2340                     fprintf(log,"cgcm: %8.3f  %8.3f  %8.3f\n",cgcm[icg][XX],
2341                             cgcm[icg][YY],cgcm[icg][ZZ]);
2342                     fprintf(log,"xi:   %8.3f  %8.3f  %8.3f\n",XI,YI,ZI);
2343 #endif
2344                     for (dx=dx0; (dx<=dx1); dx++)
2345                     {
2346                         tmp1 = rl2 - dcx2[dx];
2347                         for (dy=dy0; (dy<=dy1); dy++)
2348                         {
2349                             tmp2 = tmp1 - dcy2[dy];
2350                             if (tmp2 > 0)
2351                             {
2352                                 for (dz=dz0; (dz<=dz1); dz++) {
2353                                     if (tmp2 > dcz2[dz]) {
2354                                         /* Find grid-cell cj in which possible neighbours are */
2355                                         cj   = xyz2ci(Ny,Nz,dx,dy,dz);
2356                                         
2357                                         /* Check out how many cgs (nrj) there in this cell */
2358                                         nrj  = gridnra[cj];
2359                                         
2360                                         /* Find the offset in the cg list */
2361                                         cgj0 = gridind[cj];
2362                                         
2363                                         /* Check if all j's are out of range so we
2364                                          * can skip the whole cell.
2365                                          * Should save some time, especially with DD.
2366                                          */
2367                                         if (nrj == 0 ||
2368                                             (grida[cgj0] >= max_jcg &&
2369                                              (grida[cgj0] >= jcg1 || grida[cgj0+nrj-1] < jcg0)))
2370                                         {
2371                                             continue;
2372                                         }
2373                                         
2374                                         /* Loop over cgs */
2375                                         for (j=0; (j<nrj); j++)
2376                                         {
2377                                             jjcg = grida[cgj0+j];
2378                                             
2379                                             /* check whether this guy is in range! */
2380                                             if ((jjcg >= jcg0 && jjcg < jcg1) ||
2381                                                 (jjcg < max_jcg))
2382                                             {
2383                                                 r2=calc_dx2(XI,YI,ZI,cgcm[jjcg]);
2384                                                 if (r2 < rl2) {
2385                                                     /* jgid = gid[cgsatoms[cgsindex[jjcg]]]; */
2386                                                     jgid = GET_CGINFO_GID(cginfo[jjcg]);
2387                                                     /* check energy group exclusions */
2388                                                     if (!(i_egp_flags[jgid] & EGP_EXCL))
2389                                                     {
2390                                                         if (r2 < rs2)
2391                                                         {
2392                                                             if (nsr[jgid] >= MAX_CG)
2393                                                             {
2394                                                                 /* Add to short-range list */
2395                                                                 put_in_list(bHaveVdW,ngid,md,icg,jgid,
2396                                                                             nsr[jgid],nl_sr[jgid],
2397                                                                             cgs->index,/* cgsatoms, */ bexcl,
2398                                                                             shift,fr,FALSE,TRUE,TRUE,fr->solvent_opt);
2399                                                                 nsr[jgid]=0;
2400                                                             }
2401                                                             nl_sr[jgid][nsr[jgid]++]=jjcg;
2402                                                         } 
2403                                                         else if (r2 < rm2)
2404                                                         {
2405                                                             if (nlr_ljc[jgid] >= MAX_CG)
2406                                                             {
2407                                                                 /* Add to LJ+coulomb long-range list */
2408                                                                 put_in_list(bHaveVdW,ngid,md,icg,jgid,
2409                                                                             nlr_ljc[jgid],nl_lr_ljc[jgid],top->cgs.index,
2410                                                                             bexcl,shift,fr,TRUE,TRUE,TRUE,fr->solvent_opt);
2411                                                                 nlr_ljc[jgid]=0;
2412                                                             }
2413                                                             nl_lr_ljc[jgid][nlr_ljc[jgid]++]=jjcg;
2414                                                         }
2415                                                         else
2416                                                         {
2417                                                             if (nlr_one[jgid] >= MAX_CG)
2418                                                             {
2419                                                                 /* Add to long-range list with only coul, or only LJ */
2420                                                                 put_in_list(bHaveVdW,ngid,md,icg,jgid,
2421                                                                             nlr_one[jgid],nl_lr_one[jgid],top->cgs.index,
2422                                                                             bexcl,shift,fr,TRUE,rvdw_lt_rcoul,rcoul_lt_rvdw,fr->solvent_opt);
2423                                                                 nlr_one[jgid]=0;
2424                                                             }
2425                                                             nl_lr_one[jgid][nlr_one[jgid]++]=jjcg;
2426                                                         }
2427                                                     }
2428                                                 }
2429                                                 nns++;
2430                                             }
2431                                         }
2432                                     }
2433                                 }
2434                             }
2435                         }
2436                     }
2437                     /* CHECK whether there is anything left in the buffers */
2438                     for(nn=0; (nn<ngid); nn++)
2439                     {
2440                         if (nsr[nn] > 0)
2441                         {
2442                             put_in_list(bHaveVdW,ngid,md,icg,nn,nsr[nn],nl_sr[nn],
2443                                         cgs->index, /* cgsatoms, */ bexcl,
2444                                         shift,fr,FALSE,TRUE,TRUE,fr->solvent_opt);
2445                         }
2446                         
2447                         if (nlr_ljc[nn] > 0)
2448                         {
2449                             put_in_list(bHaveVdW,ngid,md,icg,nn,nlr_ljc[nn],
2450                                         nl_lr_ljc[nn],top->cgs.index,
2451                                         bexcl,shift,fr,TRUE,TRUE,TRUE,fr->solvent_opt);
2452                         }
2453                         
2454                         if (nlr_one[nn] > 0)
2455                         {
2456                             put_in_list(bHaveVdW,ngid,md,icg,nn,nlr_one[nn],
2457                                         nl_lr_one[nn],top->cgs.index,
2458                                         bexcl,shift,fr,TRUE,rvdw_lt_rcoul,rcoul_lt_rvdw,fr->solvent_opt);
2459                         }
2460                     }
2461                 }
2462             }
2463         }
2464         /* setexcl(nri,i_atoms,&top->atoms.excl,FALSE,bexcl); */
2465         setexcl(cgs->index[icg],cgs->index[icg+1],&top->excls,FALSE,bexcl);
2466     }
2467     /* No need to perform any left-over force calculations anymore (as we used to do here)
2468      * since we now save the proper long-range lists for later evaluation.
2469      */
2470
2471     debug_gmx();
2472      
2473     /* Close neighbourlists */
2474     close_neighbor_lists(fr,bMakeQMMMnblist);
2475     
2476     return nns;
2477 }
2478
2479 void ns_realloc_natoms(gmx_ns_t *ns,int natoms)
2480 {
2481     int i;
2482     
2483     if (natoms > ns->nra_alloc)
2484     {
2485         ns->nra_alloc = over_alloc_dd(natoms);
2486         srenew(ns->bexcl,ns->nra_alloc);
2487         for(i=0; i<ns->nra_alloc; i++)
2488         {
2489             ns->bexcl[i] = 0;
2490         }
2491     }
2492 }
2493
2494 void init_ns(FILE *fplog,const t_commrec *cr,
2495              gmx_ns_t *ns,t_forcerec *fr,
2496              const gmx_mtop_t *mtop,
2497              matrix box)
2498 {
2499     int  mt,icg,nr_in_cg,maxcg,i,j,jcg,ngid,ncg;
2500     t_block *cgs;
2501     char *ptr;
2502     
2503     /* Compute largest charge groups size (# atoms) */
2504     nr_in_cg=1;
2505     for(mt=0; mt<mtop->nmoltype; mt++) {
2506         cgs = &mtop->moltype[mt].cgs;
2507         for (icg=0; (icg < cgs->nr); icg++)
2508         {
2509             nr_in_cg=max(nr_in_cg,(int)(cgs->index[icg+1]-cgs->index[icg]));
2510         }
2511     }
2512
2513     /* Verify whether largest charge group is <= max cg.
2514      * This is determined by the type of the local exclusion type 
2515      * Exclusions are stored in bits. (If the type is not large
2516      * enough, enlarge it, unsigned char -> unsigned short -> unsigned long)
2517      */
2518     maxcg = sizeof(t_excl)*8;
2519     if (nr_in_cg > maxcg)
2520     {
2521         gmx_fatal(FARGS,"Max #atoms in a charge group: %d > %d\n",
2522                   nr_in_cg,maxcg);
2523     }
2524     
2525     ngid = mtop->groups.grps[egcENER].nr;
2526     snew(ns->bExcludeAlleg,ngid);
2527     for(i=0; i<ngid; i++) {
2528         ns->bExcludeAlleg[i] = TRUE;
2529         for(j=0; j<ngid; j++)
2530         {
2531             if (!(fr->egp_flags[i*ngid+j] & EGP_EXCL))
2532             {
2533                 ns->bExcludeAlleg[i] = FALSE;
2534             }
2535         }
2536     }
2537     
2538     if (fr->bGrid) {
2539         /* Grid search */
2540         ns->grid = init_grid(fplog,fr);
2541         init_nsgrid_lists(fr,ngid,ns);
2542     }
2543     else
2544     {
2545         /* Simple search */
2546         snew(ns->ns_buf,ngid);
2547         for(i=0; (i<ngid); i++)
2548         {
2549             snew(ns->ns_buf[i],SHIFTS);
2550         }
2551         ncg = ncg_mtop(mtop);
2552         snew(ns->simple_aaj,2*ncg);
2553         for(jcg=0; (jcg<ncg); jcg++)
2554         {
2555             ns->simple_aaj[jcg]     = jcg;
2556             ns->simple_aaj[jcg+ncg] = jcg;
2557         }
2558     }
2559     
2560     /* Create array that determines whether or not atoms have VdW */
2561     snew(ns->bHaveVdW,fr->ntype);
2562     for(i=0; (i<fr->ntype); i++)
2563     {
2564         for(j=0; (j<fr->ntype); j++)
2565         {
2566             ns->bHaveVdW[i] = (ns->bHaveVdW[i] || 
2567                                (fr->bBHAM ? 
2568                                 ((BHAMA(fr->nbfp,fr->ntype,i,j) != 0) ||
2569                                  (BHAMB(fr->nbfp,fr->ntype,i,j) != 0) ||
2570                                  (BHAMC(fr->nbfp,fr->ntype,i,j) != 0)) :
2571                                 ((C6(fr->nbfp,fr->ntype,i,j) != 0) ||
2572                                  (C12(fr->nbfp,fr->ntype,i,j) != 0))));
2573         }
2574     }
2575     if (debug) 
2576         pr_bvec(debug,0,"bHaveVdW",ns->bHaveVdW,fr->ntype,TRUE);
2577     
2578     ns->nra_alloc = 0;
2579     ns->bexcl = NULL;
2580     if (!DOMAINDECOMP(cr))
2581     {
2582         /* This could be reduced with particle decomposition */
2583         ns_realloc_natoms(ns,mtop->natoms);
2584     }
2585
2586     ns->nblist_initialized=FALSE;
2587
2588     /* nbr list debug dump */
2589     {
2590         char *ptr=getenv("GMX_DUMP_NL");
2591         if (ptr)
2592         {
2593             ns->dump_nl=strtol(ptr,NULL,10);
2594             if (fplog)
2595             {
2596                 fprintf(fplog, "GMX_DUMP_NL = %d", ns->dump_nl);
2597             }
2598         }
2599         else
2600         {
2601             ns->dump_nl=0;
2602         }
2603     }
2604 }
2605
2606                          
2607 int search_neighbours(FILE *log,t_forcerec *fr,
2608                       rvec x[],matrix box,
2609                       gmx_localtop_t *top,
2610                       gmx_groups_t *groups,
2611                       t_commrec *cr,
2612                       t_nrnb *nrnb,t_mdatoms *md,
2613                       real *lambda,real *dvdlambda,
2614                       gmx_grppairener_t *grppener,
2615                       gmx_bool bFillGrid,
2616                       gmx_bool bDoLongRangeNS,
2617                       gmx_bool bPadListsForKernels)
2618 {
2619     t_block  *cgs=&(top->cgs);
2620     rvec     box_size,grid_x0,grid_x1;
2621     int      i,j,m,ngid;
2622     real     min_size,grid_dens;
2623     int      nsearch;
2624     gmx_bool     bGrid;
2625     char     *ptr;
2626     gmx_bool     *i_egp_flags;
2627     int      cg_start,cg_end,start,end;
2628     gmx_ns_t *ns;
2629     t_grid   *grid;
2630     gmx_domdec_zones_t *dd_zones;
2631     put_in_list_t *put_in_list;
2632
2633     ns = &fr->ns;
2634
2635     /* Set some local variables */
2636     bGrid = fr->bGrid;
2637     ngid = groups->grps[egcENER].nr;
2638     
2639     for(m=0; (m<DIM); m++)
2640     {
2641         box_size[m] = box[m][m];
2642     }
2643   
2644     if (fr->ePBC != epbcNONE)
2645     {
2646         if (sqr(fr->rlistlong) >= max_cutoff2(fr->ePBC,box))
2647         {
2648             gmx_fatal(FARGS,"One of the box vectors has become shorter than twice the cut-off length or box_yy-|box_zy| or box_zz has become smaller than the cut-off.");
2649         }
2650         if (!bGrid)
2651         {
2652             min_size = min(box_size[XX],min(box_size[YY],box_size[ZZ]));
2653             if (2*fr->rlistlong >= min_size)
2654                 gmx_fatal(FARGS,"One of the box diagonal elements has become smaller than twice the cut-off length.");
2655         }
2656     }
2657     
2658     if (DOMAINDECOMP(cr))
2659     {
2660         ns_realloc_natoms(ns,cgs->index[cgs->nr]);
2661     }
2662     debug_gmx();
2663     
2664     /* Reset the neighbourlists */
2665     reset_neighbor_lists(fr,TRUE,TRUE);
2666     
2667     if (bGrid && bFillGrid)
2668     {
2669                 
2670         grid = ns->grid;
2671         if (DOMAINDECOMP(cr))
2672         {
2673             dd_zones = domdec_zones(cr->dd);
2674         }
2675         else
2676         {
2677             dd_zones = NULL;
2678
2679             get_nsgrid_boundaries(grid->nboundeddim,box,NULL,NULL,NULL,NULL,
2680                                   cgs->nr,fr->cg_cm,grid_x0,grid_x1,&grid_dens);
2681
2682             grid_first(log,grid,NULL,NULL,fr->ePBC,box,grid_x0,grid_x1,
2683                        fr->rlistlong,grid_dens);
2684         }
2685         debug_gmx();
2686         
2687         /* Don't know why this all is... (DvdS 3/99) */
2688 #ifndef SEGV
2689         start = 0;
2690         end   = cgs->nr;
2691 #else
2692         start = fr->cg0;
2693         end   = (cgs->nr+1)/2;
2694 #endif
2695         
2696         if (DOMAINDECOMP(cr))
2697         {
2698             end = cgs->nr;
2699             fill_grid(log,dd_zones,grid,end,-1,end,fr->cg_cm);
2700             grid->icg0 = 0;
2701             grid->icg1 = dd_zones->izone[dd_zones->nizone-1].cg1;
2702         }
2703         else
2704         {
2705             fill_grid(log,NULL,grid,cgs->nr,fr->cg0,fr->hcg,fr->cg_cm);
2706             grid->icg0 = fr->cg0;
2707             grid->icg1 = fr->hcg;
2708             debug_gmx();
2709             
2710             if (PARTDECOMP(cr))
2711                 mv_grid(cr,grid);
2712             debug_gmx();
2713         }
2714         
2715         calc_elemnr(log,grid,start,end,cgs->nr);
2716         calc_ptrs(grid);
2717         grid_last(log,grid,start,end,cgs->nr);
2718         
2719         if (gmx_debug_at)
2720         {
2721             check_grid(debug,grid);
2722             print_grid(debug,grid);
2723         }
2724     }
2725     else if (fr->n_tpi)
2726     {
2727         /* Set the grid cell index for the test particle only.
2728          * The cell to cg index is not corrected, but that does not matter.
2729          */
2730         fill_grid(log,NULL,ns->grid,fr->hcg,fr->hcg-1,fr->hcg,fr->cg_cm);
2731     }
2732     debug_gmx();
2733
2734     if (fr->adress_type == eAdressOff){
2735         if (!fr->ns.bCGlist)
2736         {
2737             put_in_list = put_in_list_at;
2738         }
2739         else
2740         {
2741             put_in_list = put_in_list_cg;
2742         }
2743     }else{
2744          put_in_list = put_in_list_adress;
2745     }
2746
2747     /* Do the core! */
2748     if (bGrid)
2749     {
2750         grid = ns->grid;
2751         nsearch = nsgrid_core(log,cr,fr,box,box_size,ngid,top,
2752                               grid,x,ns->bexcl,ns->bExcludeAlleg,
2753                               nrnb,md,lambda,dvdlambda,grppener,
2754                               put_in_list,ns->bHaveVdW,
2755                               bDoLongRangeNS,FALSE);
2756         
2757         /* neighbour searching withouth QMMM! QM atoms have zero charge in
2758          * the classical calculation. The charge-charge interaction
2759          * between QM and MM atoms is handled in the QMMM core calculation
2760          * (see QMMM.c). The VDW however, we'd like to compute classically
2761          * and the QM MM atom pairs have just been put in the
2762          * corresponding neighbourlists. in case of QMMM we still need to
2763          * fill a special QMMM neighbourlist that contains all neighbours
2764          * of the QM atoms. If bQMMM is true, this list will now be made: 
2765          */
2766         if (fr->bQMMM && fr->qr->QMMMscheme!=eQMMMschemeoniom)
2767         {
2768             nsearch += nsgrid_core(log,cr,fr,box,box_size,ngid,top,
2769                                    grid,x,ns->bexcl,ns->bExcludeAlleg,
2770                                    nrnb,md,lambda,dvdlambda,grppener,
2771                                    put_in_list_qmmm,ns->bHaveVdW,
2772                                    bDoLongRangeNS,TRUE);
2773         }
2774     }
2775     else 
2776     {
2777         nsearch = ns_simple_core(fr,top,md,box,box_size,
2778                                  ns->bexcl,ns->simple_aaj,
2779                                  ngid,ns->ns_buf,put_in_list,ns->bHaveVdW);
2780     }
2781     debug_gmx();
2782
2783 #ifdef DEBUG
2784     pr_nsblock(log);
2785 #endif
2786     
2787     inc_nrnb(nrnb,eNR_NS,nsearch);
2788     /* inc_nrnb(nrnb,eNR_LR,fr->nlr); */
2789     
2790     return nsearch;
2791 }
2792
2793 int natoms_beyond_ns_buffer(t_inputrec *ir,t_forcerec *fr,t_block *cgs,
2794                             matrix scale_tot,rvec *x)
2795 {
2796     int  cg0,cg1,cg,a0,a1,a,i,j;
2797     real rint,hbuf2,scale;
2798     rvec *cg_cm,cgsc;
2799     gmx_bool bIsotropic;
2800     int  nBeyond;
2801     
2802     nBeyond = 0;
2803     
2804     rint = max(ir->rcoulomb,ir->rvdw);
2805     if (ir->rlist < rint)
2806     {
2807         gmx_fatal(FARGS,"The neighbor search buffer has negative size: %f nm",
2808                   ir->rlist - rint);
2809     }
2810     cg_cm = fr->cg_cm;
2811     
2812     cg0 = fr->cg0;
2813     cg1 = fr->hcg;
2814     
2815     if (!EI_DYNAMICS(ir->eI) || !DYNAMIC_BOX(*ir))
2816     {
2817         hbuf2 = sqr(0.5*(ir->rlist - rint));
2818         for(cg=cg0; cg<cg1; cg++)
2819         {
2820             a0 = cgs->index[cg];
2821             a1 = cgs->index[cg+1];
2822             for(a=a0; a<a1; a++)
2823             {
2824                 if (distance2(cg_cm[cg],x[a]) > hbuf2)
2825                 {
2826                     nBeyond++;
2827                 }
2828             }
2829         }
2830     }
2831     else
2832     {
2833         bIsotropic = TRUE;
2834         scale = scale_tot[0][0];
2835         for(i=1; i<DIM; i++)
2836         {
2837             /* With anisotropic scaling, the original spherical ns volumes become
2838              * ellipsoids. To avoid costly transformations we use the minimum
2839              * eigenvalue of the scaling matrix for determining the buffer size.
2840              * Since the lower half is 0, the eigenvalues are the diagonal elements.
2841              */
2842             scale = min(scale,scale_tot[i][i]);
2843             if (scale_tot[i][i] != scale_tot[i-1][i-1])
2844             {
2845                 bIsotropic = FALSE;
2846             }
2847             for(j=0; j<i; j++)
2848             {
2849                 if (scale_tot[i][j] != 0)
2850                 {
2851                     bIsotropic = FALSE;
2852                 }
2853             }
2854         }
2855         hbuf2 = sqr(0.5*(scale*ir->rlist - rint));
2856         if (bIsotropic)
2857         {
2858             for(cg=cg0; cg<cg1; cg++)
2859             {
2860                 svmul(scale,cg_cm[cg],cgsc);
2861                 a0 = cgs->index[cg];
2862                 a1 = cgs->index[cg+1];
2863                 for(a=a0; a<a1; a++)
2864                 {
2865                     if (distance2(cgsc,x[a]) > hbuf2)
2866                     {                    
2867                         nBeyond++;
2868                     }
2869                 }
2870             }
2871         }
2872         else
2873         {
2874             /* Anistropic scaling */
2875             for(cg=cg0; cg<cg1; cg++)
2876             {
2877                 /* Since scale_tot contains the transpose of the scaling matrix,
2878                  * we need to multiply with the transpose.
2879                  */
2880                 tmvmul_ur0(scale_tot,cg_cm[cg],cgsc);
2881                 a0 = cgs->index[cg];
2882                 a1 = cgs->index[cg+1];
2883                 for(a=a0; a<a1; a++)
2884                 {
2885                     if (distance2(cgsc,x[a]) > hbuf2)
2886                     {
2887                         nBeyond++;
2888                     }
2889                 }
2890             }
2891         }
2892     }
2893     
2894     return nBeyond;
2895 }