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