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