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