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