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