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