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