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