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