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