Merge release-4-6 into master
[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     gmx_bool      bEnergyGroupCG;
1156
1157     /* Copy some pointers */
1158     cginfo  = fr->cginfo;
1159     charge  = md->chargeA;
1160     chargeB = md->chargeB;
1161     type    = md->typeA;
1162     typeB   = md->typeB;
1163     bPert   = md->bPerturbed;
1164     wf      = md->wf;
1165
1166     /* Get atom range */
1167     i0     = index[icg];
1168     nicg   = index[icg+1]-i0;
1169
1170     /* Get the i charge group info */
1171     igid   = GET_CGINFO_GID(cginfo[icg]);
1172
1173     iwater = (solvent_opt != esolNO) ? GET_CGINFO_SOLOPT(cginfo[icg]) : esolNO;
1174
1175     if (md->nPerturbed)
1176     {
1177         gmx_fatal(FARGS, "AdResS does not support free energy pertubation\n");
1178     }
1179
1180     /* Unpack pointers to neighbourlist structs */
1181     if (fr->nnblists == 2)
1182     {
1183         nbl_ind        = 0;
1184         nbl_ind_adress = 1;
1185     }
1186     else
1187     {
1188         nbl_ind        = fr->gid2nblists[GID(igid, jgid, ngid)];
1189         nbl_ind_adress = nbl_ind+fr->nnblists/2;
1190     }
1191     if (bLR)
1192     {
1193         nlist        = fr->nblists[nbl_ind].nlist_lr;
1194         nlist_adress = fr->nblists[nbl_ind_adress].nlist_lr;
1195     }
1196     else
1197     {
1198         nlist        = fr->nblists[nbl_ind].nlist_sr;
1199         nlist_adress = fr->nblists[nbl_ind_adress].nlist_sr;
1200     }
1201
1202
1203     vdwc = &nlist[eNL_VDWQQ];
1204     vdw  = &nlist[eNL_VDW];
1205     coul = &nlist[eNL_QQ];
1206
1207     vdwc_adress = &nlist_adress[eNL_VDWQQ];
1208     vdw_adress  = &nlist_adress[eNL_VDW];
1209     coul_adress = &nlist_adress[eNL_QQ];
1210
1211     /* We do not support solvent optimization with AdResS for now.
1212        For this we would need hybrid solvent-other kernels */
1213
1214     /* no solvent as i charge group */
1215     /* Loop over the atoms in the i charge group */
1216     for (i = 0; i < nicg; i++)
1217     {
1218         i_atom  = i0+i;
1219         gid     = GID(igid, jgid, ngid);
1220         qi      = charge[i_atom];
1221
1222         /* Create new i_atom for each energy group */
1223         if (bDoVdW && bDoCoul)
1224         {
1225             new_i_nblist(vdwc, bLR, i_atom, shift, gid);
1226             new_i_nblist(vdwc_adress, bLR, i_atom, shift, gid);
1227
1228         }
1229         if (bDoVdW)
1230         {
1231             new_i_nblist(vdw, bLR, i_atom, shift, gid);
1232             new_i_nblist(vdw_adress, bLR, i_atom, shift, gid);
1233
1234         }
1235         if (bDoCoul)
1236         {
1237             new_i_nblist(coul, bLR, i_atom, shift, gid);
1238             new_i_nblist(coul_adress, bLR, i_atom, shift, gid);
1239         }
1240         bDoVdW_i  = (bDoVdW  && bHaveVdW[type[i_atom]]);
1241         bDoCoul_i = (bDoCoul && qi != 0);
1242
1243         /* Here we find out whether the energy groups interaction belong to a
1244          * coarse-grained (vsite) or atomistic interaction. Note that, beacuse
1245          * interactions between coarse-grained and other (atomistic) energygroups
1246          * are excluded automatically by grompp, it is sufficient to check for
1247          * the group id of atom i (igid) */
1248         bEnergyGroupCG = !egp_explicit(fr, igid);
1249
1250         if (bDoVdW_i || bDoCoul_i)
1251         {
1252             /* Loop over the j charge groups */
1253             for (j = 0; (j < nj); j++)
1254             {
1255                 jcg = jjcg[j];
1256
1257                 /* Check for large charge groups */
1258                 if (jcg == icg)
1259                 {
1260                     jj0 = i0 + i + 1;
1261                 }
1262                 else
1263                 {
1264                     jj0 = index[jcg];
1265                 }
1266
1267                 jj1 = index[jcg+1];
1268                 /* Finally loop over the atoms in the j-charge group */
1269                 for (jj = jj0; jj < jj1; jj++)
1270                 {
1271                     bNotEx = NOTEXCL(bExcl, i, jj);
1272
1273                     /* Now we have to exclude interactions which will be zero
1274                      * anyway due to the AdResS weights (in previous implementations
1275                      * this was done in the force kernel). This is necessary as
1276                      * pure interactions (those with b_hybrid=false, i.e. w_i*w_j==1 or 0)
1277                      * are put into neighbour lists which will be passed to the
1278                      * standard (optimized) kernels for speed. The interactions with
1279                      * b_hybrid=true are placed into the _adress neighbour lists and
1280                      * processed by the generic AdResS kernel.
1281                      */
1282                     if ( (bEnergyGroupCG &&
1283                          wf[i_atom] >= 1-GMX_REAL_EPS && wf[jj] >= 1-GMX_REAL_EPS ) ||
1284                            ( !bEnergyGroupCG && wf[jj] <= GMX_REAL_EPS ) )
1285                     {
1286                         continue;
1287                     }
1288
1289                     b_hybrid = !((wf[i_atom] >= 1-GMX_REAL_EPS && wf[jj] >= 1-GMX_REAL_EPS) ||
1290                         (wf[i_atom] <= GMX_REAL_EPS && wf[jj] <= GMX_REAL_EPS));
1291
1292                     if (bNotEx)
1293                     {
1294                         if (!bDoVdW_i)
1295                         {
1296                             if (charge[jj] != 0)
1297                             {
1298                                 if (!b_hybrid)
1299                                 {
1300                                     add_j_to_nblist(coul, jj, bLR);
1301                                 }
1302                                 else
1303                                 {
1304                                     add_j_to_nblist(coul_adress, jj, bLR);
1305                                 }
1306                             }
1307                         }
1308                         else if (!bDoCoul_i)
1309                         {
1310                             if (bHaveVdW[type[jj]])
1311                             {
1312                                 if (!b_hybrid)
1313                                 {
1314                                     add_j_to_nblist(vdw, jj, bLR);
1315                                 }
1316                                 else
1317                                 {
1318                                     add_j_to_nblist(vdw_adress, jj, bLR);
1319                                 }
1320                             }
1321                         }
1322                         else
1323                         {
1324                             if (bHaveVdW[type[jj]])
1325                             {
1326                                 if (charge[jj] != 0)
1327                                 {
1328                                     if (!b_hybrid)
1329                                     {
1330                                         add_j_to_nblist(vdwc, jj, bLR);
1331                                     }
1332                                     else
1333                                     {
1334                                         add_j_to_nblist(vdwc_adress, jj, bLR);
1335                                     }
1336                                 }
1337                                 else
1338                                 {
1339                                     if (!b_hybrid)
1340                                     {
1341                                         add_j_to_nblist(vdw, jj, bLR);
1342                                     }
1343                                     else
1344                                     {
1345                                         add_j_to_nblist(vdw_adress, jj, bLR);
1346                                     }
1347
1348                                 }
1349                             }
1350                             else if (charge[jj] != 0)
1351                             {
1352                                 if (!b_hybrid)
1353                                 {
1354                                     add_j_to_nblist(coul, jj, bLR);
1355                                 }
1356                                 else
1357                                 {
1358                                     add_j_to_nblist(coul_adress, jj, bLR);
1359                                 }
1360
1361                             }
1362                         }
1363                     }
1364                 }
1365             }
1366
1367             close_i_nblist(vdw);
1368             close_i_nblist(coul);
1369             close_i_nblist(vdwc);
1370             close_i_nblist(vdw_adress);
1371             close_i_nblist(coul_adress);
1372             close_i_nblist(vdwc_adress);
1373         }
1374     }
1375 }
1376
1377 static void
1378 put_in_list_qmmm(gmx_bool              bHaveVdW[],
1379                  int                   ngid,
1380                  t_mdatoms     *       md,
1381                  int                   icg,
1382                  int                   jgid,
1383                  int                   nj,
1384                  atom_id               jjcg[],
1385                  atom_id               index[],
1386                  t_excl                bExcl[],
1387                  int                   shift,
1388                  t_forcerec     *      fr,
1389                  gmx_bool              bLR,
1390                  gmx_bool              bDoVdW,
1391                  gmx_bool              bDoCoul,
1392                  int                   solvent_opt)
1393 {
1394     t_nblist  *   coul;
1395     int           i, j, jcg, igid, gid;
1396     atom_id       jj, jj0, jj1, i_atom;
1397     int           i0, nicg;
1398     gmx_bool      bNotEx;
1399
1400     /* Get atom range */
1401     i0     = index[icg];
1402     nicg   = index[icg+1]-i0;
1403
1404     /* Get the i charge group info */
1405     igid   = GET_CGINFO_GID(fr->cginfo[icg]);
1406
1407     coul = &fr->QMMMlist;
1408
1409     /* Loop over atoms in the ith charge group */
1410     for (i = 0; i < nicg; i++)
1411     {
1412         i_atom = i0+i;
1413         gid    = GID(igid, jgid, ngid);
1414         /* Create new i_atom for each energy group */
1415         new_i_nblist(coul, bLR, i_atom, shift, gid);
1416
1417         /* Loop over the j charge groups */
1418         for (j = 0; j < nj; j++)
1419         {
1420             jcg = jjcg[j];
1421
1422             /* Charge groups cannot have QM and MM atoms simultaneously */
1423             if (jcg != icg)
1424             {
1425                 jj0 = index[jcg];
1426                 jj1 = index[jcg+1];
1427                 /* Finally loop over the atoms in the j-charge group */
1428                 for (jj = jj0; jj < jj1; jj++)
1429                 {
1430                     bNotEx = NOTEXCL(bExcl, i, jj);
1431                     if (bNotEx)
1432                     {
1433                         add_j_to_nblist(coul, jj, bLR);
1434                     }
1435                 }
1436             }
1437         }
1438         close_i_nblist(coul);
1439     }
1440 }
1441
1442 static void
1443 put_in_list_cg(gmx_bool              bHaveVdW[],
1444                int                   ngid,
1445                t_mdatoms     *       md,
1446                int                   icg,
1447                int                   jgid,
1448                int                   nj,
1449                atom_id               jjcg[],
1450                atom_id               index[],
1451                t_excl                bExcl[],
1452                int                   shift,
1453                t_forcerec     *      fr,
1454                gmx_bool              bLR,
1455                gmx_bool              bDoVdW,
1456                gmx_bool              bDoCoul,
1457                int                   solvent_opt)
1458 {
1459     int          cginfo;
1460     int          igid, gid, nbl_ind;
1461     t_nblist *   vdwc;
1462     int          j, jcg;
1463
1464     cginfo = fr->cginfo[icg];
1465
1466     igid = GET_CGINFO_GID(cginfo);
1467     gid  = GID(igid, jgid, ngid);
1468
1469     /* Unpack pointers to neighbourlist structs */
1470     if (fr->nnblists == 1)
1471     {
1472         nbl_ind = 0;
1473     }
1474     else
1475     {
1476         nbl_ind = fr->gid2nblists[gid];
1477     }
1478     if (bLR)
1479     {
1480         vdwc = &fr->nblists[nbl_ind].nlist_lr[eNL_VDWQQ];
1481     }
1482     else
1483     {
1484         vdwc = &fr->nblists[nbl_ind].nlist_sr[eNL_VDWQQ];
1485     }
1486
1487     /* Make a new neighbor list for charge group icg.
1488      * Currently simply one neighbor list is made with LJ and Coulomb.
1489      * If required, zero interactions could be removed here
1490      * or in the force loop.
1491      */
1492     new_i_nblist(vdwc, bLR, index[icg], shift, gid);
1493     vdwc->iinr_end[vdwc->nri] = index[icg+1];
1494
1495     for (j = 0; (j < nj); j++)
1496     {
1497         jcg = jjcg[j];
1498         /* Skip the icg-icg pairs if all self interactions are excluded */
1499         if (!(jcg == icg && GET_CGINFO_EXCL_INTRA(cginfo)))
1500         {
1501             /* Here we add the j charge group jcg to the list,
1502              * exclusions are also added to the list.
1503              */
1504             add_j_to_nblist_cg(vdwc, index[jcg], index[jcg+1], bExcl, icg == jcg, bLR);
1505         }
1506     }
1507
1508     close_i_nblist(vdwc);
1509 }
1510
1511 static void setexcl(atom_id start, atom_id end, t_blocka *excl, gmx_bool b,
1512                     t_excl bexcl[])
1513 {
1514     atom_id i, k;
1515
1516     if (b)
1517     {
1518         for (i = start; i < end; i++)
1519         {
1520             for (k = excl->index[i]; k < excl->index[i+1]; k++)
1521             {
1522                 SETEXCL(bexcl, i-start, excl->a[k]);
1523             }
1524         }
1525     }
1526     else
1527     {
1528         for (i = start; i < end; i++)
1529         {
1530             for (k = excl->index[i]; k < excl->index[i+1]; k++)
1531             {
1532                 RMEXCL(bexcl, i-start, excl->a[k]);
1533             }
1534         }
1535     }
1536 }
1537
1538 int calc_naaj(int icg, int cgtot)
1539 {
1540     int naaj;
1541
1542     if ((cgtot % 2) == 1)
1543     {
1544         /* Odd number of charge groups, easy */
1545         naaj = 1 + (cgtot/2);
1546     }
1547     else if ((cgtot % 4) == 0)
1548     {
1549         /* Multiple of four is hard */
1550         if (icg < cgtot/2)
1551         {
1552             if ((icg % 2) == 0)
1553             {
1554                 naaj = 1+(cgtot/2);
1555             }
1556             else
1557             {
1558                 naaj = cgtot/2;
1559             }
1560         }
1561         else
1562         {
1563             if ((icg % 2) == 1)
1564             {
1565                 naaj = 1+(cgtot/2);
1566             }
1567             else
1568             {
1569                 naaj = cgtot/2;
1570             }
1571         }
1572     }
1573     else
1574     {
1575         /* cgtot/2 = odd */
1576         if ((icg % 2) == 0)
1577         {
1578             naaj = 1+(cgtot/2);
1579         }
1580         else
1581         {
1582             naaj = cgtot/2;
1583         }
1584     }
1585 #ifdef DEBUG
1586     fprintf(log, "naaj=%d\n", naaj);
1587 #endif
1588
1589     return naaj;
1590 }
1591
1592 /************************************************
1593  *
1594  *  S I M P L E      C O R E     S T U F F
1595  *
1596  ************************************************/
1597
1598 static real calc_image_tric(rvec xi, rvec xj, matrix box,
1599                             rvec b_inv, int *shift)
1600 {
1601     /* This code assumes that the cut-off is smaller than
1602      * a half times the smallest diagonal element of the box.
1603      */
1604     const real h25 = 2.5;
1605     real       dx, dy, dz;
1606     real       r2;
1607     int        tx, ty, tz;
1608
1609     /* Compute diff vector */
1610     dz = xj[ZZ] - xi[ZZ];
1611     dy = xj[YY] - xi[YY];
1612     dx = xj[XX] - xi[XX];
1613
1614     /* Perform NINT operation, using trunc operation, therefore
1615      * we first add 2.5 then subtract 2 again
1616      */
1617     tz  = dz*b_inv[ZZ] + h25;
1618     tz -= 2;
1619     dz -= tz*box[ZZ][ZZ];
1620     dy -= tz*box[ZZ][YY];
1621     dx -= tz*box[ZZ][XX];
1622
1623     ty  = dy*b_inv[YY] + h25;
1624     ty -= 2;
1625     dy -= ty*box[YY][YY];
1626     dx -= ty*box[YY][XX];
1627
1628     tx  = dx*b_inv[XX]+h25;
1629     tx -= 2;
1630     dx -= tx*box[XX][XX];
1631
1632     /* Distance squared */
1633     r2 = (dx*dx) + (dy*dy) + (dz*dz);
1634
1635     *shift = XYZ2IS(tx, ty, tz);
1636
1637     return r2;
1638 }
1639
1640 static real calc_image_rect(rvec xi, rvec xj, rvec box_size,
1641                             rvec b_inv, int *shift)
1642 {
1643     const real h15 = 1.5;
1644     real       ddx, ddy, ddz;
1645     real       dx, dy, dz;
1646     real       r2;
1647     int        tx, ty, tz;
1648
1649     /* Compute diff vector */
1650     dx = xj[XX] - xi[XX];
1651     dy = xj[YY] - xi[YY];
1652     dz = xj[ZZ] - xi[ZZ];
1653
1654     /* Perform NINT operation, using trunc operation, therefore
1655      * we first add 1.5 then subtract 1 again
1656      */
1657     tx = dx*b_inv[XX] + h15;
1658     ty = dy*b_inv[YY] + h15;
1659     tz = dz*b_inv[ZZ] + h15;
1660     tx--;
1661     ty--;
1662     tz--;
1663
1664     /* Correct diff vector for translation */
1665     ddx = tx*box_size[XX] - dx;
1666     ddy = ty*box_size[YY] - dy;
1667     ddz = tz*box_size[ZZ] - dz;
1668
1669     /* Distance squared */
1670     r2 = (ddx*ddx) + (ddy*ddy) + (ddz*ddz);
1671
1672     *shift = XYZ2IS(tx, ty, tz);
1673
1674     return r2;
1675 }
1676
1677 static void add_simple(t_ns_buf *nsbuf, int nrj, atom_id cg_j,
1678                        gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1679                        int icg, int jgid, t_block *cgs, t_excl bexcl[],
1680                        int shift, t_forcerec *fr, put_in_list_t *put_in_list)
1681 {
1682     if (nsbuf->nj + nrj > MAX_CG)
1683     {
1684         put_in_list(bHaveVdW, ngid, md, icg, jgid, nsbuf->ncg, nsbuf->jcg,
1685                     cgs->index, bexcl, shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
1686         /* Reset buffer contents */
1687         nsbuf->ncg = nsbuf->nj = 0;
1688     }
1689     nsbuf->jcg[nsbuf->ncg++] = cg_j;
1690     nsbuf->nj               += nrj;
1691 }
1692
1693 static void ns_inner_tric(rvec x[], int icg, int *i_egp_flags,
1694                           int njcg, atom_id jcg[],
1695                           matrix box, rvec b_inv, real rcut2,
1696                           t_block *cgs, t_ns_buf **ns_buf,
1697                           gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1698                           t_excl bexcl[], t_forcerec *fr,
1699                           put_in_list_t *put_in_list)
1700 {
1701     int       shift;
1702     int       j, nrj, jgid;
1703     int      *cginfo = fr->cginfo;
1704     atom_id   cg_j, *cgindex;
1705     t_ns_buf *nsbuf;
1706
1707     cgindex = cgs->index;
1708     shift   = CENTRAL;
1709     for (j = 0; (j < njcg); j++)
1710     {
1711         cg_j   = jcg[j];
1712         nrj    = cgindex[cg_j+1]-cgindex[cg_j];
1713         if (calc_image_tric(x[icg], x[cg_j], box, b_inv, &shift) < rcut2)
1714         {
1715             jgid  = GET_CGINFO_GID(cginfo[cg_j]);
1716             if (!(i_egp_flags[jgid] & EGP_EXCL))
1717             {
1718                 add_simple(&ns_buf[jgid][shift], nrj, cg_j,
1719                            bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, shift, fr,
1720                            put_in_list);
1721             }
1722         }
1723     }
1724 }
1725
1726 static void ns_inner_rect(rvec x[], int icg, int *i_egp_flags,
1727                           int njcg, atom_id jcg[],
1728                           gmx_bool bBox, rvec box_size, rvec b_inv, real rcut2,
1729                           t_block *cgs, t_ns_buf **ns_buf,
1730                           gmx_bool bHaveVdW[], int ngid, t_mdatoms *md,
1731                           t_excl bexcl[], t_forcerec *fr,
1732                           put_in_list_t *put_in_list)
1733 {
1734     int       shift;
1735     int       j, nrj, jgid;
1736     int      *cginfo = fr->cginfo;
1737     atom_id   cg_j, *cgindex;
1738     t_ns_buf *nsbuf;
1739
1740     cgindex = cgs->index;
1741     if (bBox)
1742     {
1743         shift = CENTRAL;
1744         for (j = 0; (j < njcg); j++)
1745         {
1746             cg_j   = jcg[j];
1747             nrj    = cgindex[cg_j+1]-cgindex[cg_j];
1748             if (calc_image_rect(x[icg], x[cg_j], box_size, b_inv, &shift) < rcut2)
1749             {
1750                 jgid  = GET_CGINFO_GID(cginfo[cg_j]);
1751                 if (!(i_egp_flags[jgid] & EGP_EXCL))
1752                 {
1753                     add_simple(&ns_buf[jgid][shift], nrj, cg_j,
1754                                bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, shift, fr,
1755                                put_in_list);
1756                 }
1757             }
1758         }
1759     }
1760     else
1761     {
1762         for (j = 0; (j < njcg); j++)
1763         {
1764             cg_j   = jcg[j];
1765             nrj    = cgindex[cg_j+1]-cgindex[cg_j];
1766             if ((rcut2 == 0) || (distance2(x[icg], x[cg_j]) < rcut2))
1767             {
1768                 jgid  = GET_CGINFO_GID(cginfo[cg_j]);
1769                 if (!(i_egp_flags[jgid] & EGP_EXCL))
1770                 {
1771                     add_simple(&ns_buf[jgid][CENTRAL], nrj, cg_j,
1772                                bHaveVdW, ngid, md, icg, jgid, cgs, bexcl, CENTRAL, fr,
1773                                put_in_list);
1774                 }
1775             }
1776         }
1777     }
1778 }
1779
1780 /* ns_simple_core needs to be adapted for QMMM still 2005 */
1781
1782 static int ns_simple_core(t_forcerec *fr,
1783                           gmx_localtop_t *top,
1784                           t_mdatoms *md,
1785                           matrix box, rvec box_size,
1786                           t_excl bexcl[], atom_id *aaj,
1787                           int ngid, t_ns_buf **ns_buf,
1788                           put_in_list_t *put_in_list, gmx_bool bHaveVdW[])
1789 {
1790     int          naaj, k;
1791     real         rlist2;
1792     int          nsearch, icg, jcg, igid, i0, nri, nn;
1793     int         *cginfo;
1794     t_ns_buf    *nsbuf;
1795     /* atom_id  *i_atoms; */
1796     t_block     *cgs  = &(top->cgs);
1797     t_blocka    *excl = &(top->excls);
1798     rvec         b_inv;
1799     int          m;
1800     gmx_bool     bBox, bTriclinic;
1801     int         *i_egp_flags;
1802
1803     rlist2 = sqr(fr->rlist);
1804
1805     bBox = (fr->ePBC != epbcNONE);
1806     if (bBox)
1807     {
1808         for (m = 0; (m < DIM); m++)
1809         {
1810             b_inv[m] = divide_err(1.0, box_size[m]);
1811         }
1812         bTriclinic = TRICLINIC(box);
1813     }
1814     else
1815     {
1816         bTriclinic = FALSE;
1817     }
1818
1819     cginfo = fr->cginfo;
1820
1821     nsearch = 0;
1822     for (icg = fr->cg0; (icg < fr->hcg); icg++)
1823     {
1824         /*
1825            i0        = cgs->index[icg];
1826            nri       = cgs->index[icg+1]-i0;
1827            i_atoms   = &(cgs->a[i0]);
1828            i_eg_excl = fr->eg_excl + ngid*md->cENER[*i_atoms];
1829            setexcl(nri,i_atoms,excl,TRUE,bexcl);
1830          */
1831         igid        = GET_CGINFO_GID(cginfo[icg]);
1832         i_egp_flags = fr->egp_flags + ngid*igid;
1833         setexcl(cgs->index[icg], cgs->index[icg+1], excl, TRUE, bexcl);
1834
1835         naaj = calc_naaj(icg, cgs->nr);
1836         if (bTriclinic)
1837         {
1838             ns_inner_tric(fr->cg_cm, icg, i_egp_flags, naaj, &(aaj[icg]),
1839                           box, b_inv, rlist2, cgs, ns_buf,
1840                           bHaveVdW, ngid, md, bexcl, fr, put_in_list);
1841         }
1842         else
1843         {
1844             ns_inner_rect(fr->cg_cm, icg, i_egp_flags, naaj, &(aaj[icg]),
1845                           bBox, box_size, b_inv, rlist2, cgs, ns_buf,
1846                           bHaveVdW, ngid, md, bexcl, fr, put_in_list);
1847         }
1848         nsearch += naaj;
1849
1850         for (nn = 0; (nn < ngid); nn++)
1851         {
1852             for (k = 0; (k < SHIFTS); k++)
1853             {
1854                 nsbuf = &(ns_buf[nn][k]);
1855                 if (nsbuf->ncg > 0)
1856                 {
1857                     put_in_list(bHaveVdW, ngid, md, icg, nn, nsbuf->ncg, nsbuf->jcg,
1858                                 cgs->index, bexcl, k, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
1859                     nsbuf->ncg = nsbuf->nj = 0;
1860                 }
1861             }
1862         }
1863         /* setexcl(nri,i_atoms,excl,FALSE,bexcl); */
1864         setexcl(cgs->index[icg], cgs->index[icg+1], excl, FALSE, bexcl);
1865     }
1866     close_neighbor_lists(fr, FALSE);
1867
1868     return nsearch;
1869 }
1870
1871 /************************************************
1872  *
1873  *    N S 5     G R I D     S T U F F
1874  *
1875  ************************************************/
1876
1877 static inline void get_dx(int Nx, real gridx, real rc2, int xgi, real x,
1878                           int *dx0, int *dx1, real *dcx2)
1879 {
1880     real dcx, tmp;
1881     int  xgi0, xgi1, i;
1882
1883     if (xgi < 0)
1884     {
1885         *dx0 = 0;
1886         xgi0 = -1;
1887         *dx1 = -1;
1888         xgi1 = 0;
1889     }
1890     else if (xgi >= Nx)
1891     {
1892         *dx0 = Nx;
1893         xgi0 = Nx-1;
1894         *dx1 = Nx-1;
1895         xgi1 = Nx;
1896     }
1897     else
1898     {
1899         dcx2[xgi] = 0;
1900         *dx0      = xgi;
1901         xgi0      = xgi-1;
1902         *dx1      = xgi;
1903         xgi1      = xgi+1;
1904     }
1905
1906     for (i = xgi0; i >= 0; i--)
1907     {
1908         dcx = (i+1)*gridx-x;
1909         tmp = dcx*dcx;
1910         if (tmp >= rc2)
1911         {
1912             break;
1913         }
1914         *dx0    = i;
1915         dcx2[i] = tmp;
1916     }
1917     for (i = xgi1; i < Nx; i++)
1918     {
1919         dcx = i*gridx-x;
1920         tmp = dcx*dcx;
1921         if (tmp >= rc2)
1922         {
1923             break;
1924         }
1925         *dx1    = i;
1926         dcx2[i] = tmp;
1927     }
1928 }
1929
1930 static inline void get_dx_dd(int Nx, real gridx, real rc2, int xgi, real x,
1931                              int ncpddc, int shift_min, int shift_max,
1932                              int *g0, int *g1, real *dcx2)
1933 {
1934     real dcx, tmp;
1935     int  g_min, g_max, shift_home;
1936
1937     if (xgi < 0)
1938     {
1939         g_min = 0;
1940         g_max = Nx - 1;
1941         *g0   = 0;
1942         *g1   = -1;
1943     }
1944     else if (xgi >= Nx)
1945     {
1946         g_min = 0;
1947         g_max = Nx - 1;
1948         *g0   = Nx;
1949         *g1   = Nx - 1;
1950     }
1951     else
1952     {
1953         if (ncpddc == 0)
1954         {
1955             g_min = 0;
1956             g_max = Nx - 1;
1957         }
1958         else
1959         {
1960             if (xgi < ncpddc)
1961             {
1962                 shift_home = 0;
1963             }
1964             else
1965             {
1966                 shift_home = -1;
1967             }
1968             g_min = (shift_min == shift_home ? 0          : ncpddc);
1969             g_max = (shift_max == shift_home ? ncpddc - 1 : Nx - 1);
1970         }
1971         if (shift_min > 0)
1972         {
1973             *g0 = g_min;
1974             *g1 = g_min - 1;
1975         }
1976         else if (shift_max < 0)
1977         {
1978             *g0 = g_max + 1;
1979             *g1 = g_max;
1980         }
1981         else
1982         {
1983             *g0       = xgi;
1984             *g1       = xgi;
1985             dcx2[xgi] = 0;
1986         }
1987     }
1988
1989     while (*g0 > g_min)
1990     {
1991         /* Check one grid cell down */
1992         dcx = ((*g0 - 1) + 1)*gridx - x;
1993         tmp = dcx*dcx;
1994         if (tmp >= rc2)
1995         {
1996             break;
1997         }
1998         (*g0)--;
1999         dcx2[*g0] = tmp;
2000     }
2001
2002     while (*g1 < g_max)
2003     {
2004         /* Check one grid cell up */
2005         dcx = (*g1 + 1)*gridx - x;
2006         tmp = dcx*dcx;
2007         if (tmp >= rc2)
2008         {
2009             break;
2010         }
2011         (*g1)++;
2012         dcx2[*g1] = tmp;
2013     }
2014 }
2015
2016
2017 #define sqr(x) ((x)*(x))
2018 #define calc_dx2(XI, YI, ZI, y) (sqr(XI-y[XX]) + sqr(YI-y[YY]) + sqr(ZI-y[ZZ]))
2019 #define calc_cyl_dx2(XI, YI, y) (sqr(XI-y[XX]) + sqr(YI-y[YY]))
2020 /****************************************************
2021  *
2022  *    F A S T   N E I G H B O R  S E A R C H I N G
2023  *
2024  *    Optimized neighboursearching routine using grid
2025  *    at least 1x1x1, see GROMACS manual
2026  *
2027  ****************************************************/
2028
2029
2030 static void get_cutoff2(t_forcerec *fr, gmx_bool bDoLongRange,
2031                         real *rvdw2, real *rcoul2,
2032                         real *rs2, real *rm2, real *rl2)
2033 {
2034     *rs2 = sqr(fr->rlist);
2035
2036     if (bDoLongRange && fr->bTwinRange)
2037     {
2038         /* The VdW and elec. LR cut-off's could be different,
2039          * so we can not simply set them to rlistlong.
2040          */
2041         if (EVDW_MIGHT_BE_ZERO_AT_CUTOFF(fr->vdwtype) &&
2042             fr->rvdw > fr->rlist)
2043         {
2044             *rvdw2  = sqr(fr->rlistlong);
2045         }
2046         else
2047         {
2048             *rvdw2  = sqr(fr->rvdw);
2049         }
2050         if (EEL_MIGHT_BE_ZERO_AT_CUTOFF(fr->eeltype) &&
2051             fr->rcoulomb > fr->rlist)
2052         {
2053             *rcoul2 = sqr(fr->rlistlong);
2054         }
2055         else
2056         {
2057             *rcoul2 = sqr(fr->rcoulomb);
2058         }
2059     }
2060     else
2061     {
2062         /* Workaround for a gcc -O3 or -ffast-math problem */
2063         *rvdw2  = *rs2;
2064         *rcoul2 = *rs2;
2065     }
2066     *rm2 = min(*rvdw2, *rcoul2);
2067     *rl2 = max(*rvdw2, *rcoul2);
2068 }
2069
2070 static void init_nsgrid_lists(t_forcerec *fr, int ngid, gmx_ns_t *ns)
2071 {
2072     real rvdw2, rcoul2, rs2, rm2, rl2;
2073     int  j;
2074
2075     get_cutoff2(fr, TRUE, &rvdw2, &rcoul2, &rs2, &rm2, &rl2);
2076
2077     /* Short range buffers */
2078     snew(ns->nl_sr, ngid);
2079     /* Counters */
2080     snew(ns->nsr, ngid);
2081     snew(ns->nlr_ljc, ngid);
2082     snew(ns->nlr_one, ngid);
2083
2084     /* Always allocate both list types, since rcoulomb might now change with PME load balancing */
2085     /* Long range VdW and Coul buffers */
2086     snew(ns->nl_lr_ljc, ngid);
2087     /* Long range VdW or Coul only buffers */
2088     snew(ns->nl_lr_one, ngid);
2089
2090     for (j = 0; (j < ngid); j++)
2091     {
2092         snew(ns->nl_sr[j], MAX_CG);
2093         snew(ns->nl_lr_ljc[j], MAX_CG);
2094         snew(ns->nl_lr_one[j], MAX_CG);
2095     }
2096     if (debug)
2097     {
2098         fprintf(debug,
2099                 "ns5_core: rs2 = %g, rm2 = %g, rl2 = %g (nm^2)\n",
2100                 rs2, rm2, rl2);
2101     }
2102 }
2103
2104 static int nsgrid_core(FILE *log, t_commrec *cr, t_forcerec *fr,
2105                        matrix box, rvec box_size, int ngid,
2106                        gmx_localtop_t *top,
2107                        t_grid *grid, rvec x[],
2108                        t_excl bexcl[], gmx_bool *bExcludeAlleg,
2109                        t_nrnb *nrnb, t_mdatoms *md,
2110                        real *lambda, real *dvdlambda,
2111                        gmx_grppairener_t *grppener,
2112                        put_in_list_t *put_in_list,
2113                        gmx_bool bHaveVdW[],
2114                        gmx_bool bDoLongRange, gmx_bool bMakeQMMMnblist)
2115 {
2116     gmx_ns_t     *ns;
2117     atom_id     **nl_lr_ljc, **nl_lr_one, **nl_sr;
2118     int          *nlr_ljc, *nlr_one, *nsr;
2119     gmx_domdec_t *dd     = NULL;
2120     t_block      *cgs    = &(top->cgs);
2121     int          *cginfo = fr->cginfo;
2122     /* atom_id *i_atoms,*cgsindex=cgs->index; */
2123     ivec          sh0, sh1, shp;
2124     int           cell_x, cell_y, cell_z;
2125     int           d, tx, ty, tz, dx, dy, dz, cj;
2126 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2127     int           zsh_ty, zsh_tx, ysh_tx;
2128 #endif
2129     int           dx0, dx1, dy0, dy1, dz0, dz1;
2130     int           Nx, Ny, Nz, shift = -1, j, nrj, nns, nn = -1;
2131     real          gridx, gridy, gridz, grid_x, grid_y, grid_z;
2132     real         *dcx2, *dcy2, *dcz2;
2133     int           zgi, ygi, xgi;
2134     int           cg0, cg1, icg = -1, cgsnr, i0, igid, nri, naaj, max_jcg;
2135     int           jcg0, jcg1, jjcg, cgj0, jgid;
2136     int          *grida, *gridnra, *gridind;
2137     gmx_bool      rvdw_lt_rcoul, rcoul_lt_rvdw;
2138     rvec          xi, *cgcm, grid_offset;
2139     real          r2, rs2, rvdw2, rcoul2, rm2, rl2, XI, YI, ZI, dcx, dcy, dcz, tmp1, tmp2;
2140     int          *i_egp_flags;
2141     gmx_bool      bDomDec, bTriclinicX, bTriclinicY;
2142     ivec          ncpddc;
2143
2144     ns = &fr->ns;
2145
2146     bDomDec = DOMAINDECOMP(cr);
2147     if (bDomDec)
2148     {
2149         dd = cr->dd;
2150     }
2151
2152     bTriclinicX = ((YY < grid->npbcdim &&
2153                     (!bDomDec || dd->nc[YY] == 1) && box[YY][XX] != 0) ||
2154                    (ZZ < grid->npbcdim &&
2155                     (!bDomDec || dd->nc[ZZ] == 1) && box[ZZ][XX] != 0));
2156     bTriclinicY =  (ZZ < grid->npbcdim &&
2157                     (!bDomDec || dd->nc[ZZ] == 1) && box[ZZ][YY] != 0);
2158
2159     cgsnr    = cgs->nr;
2160
2161     get_cutoff2(fr, bDoLongRange, &rvdw2, &rcoul2, &rs2, &rm2, &rl2);
2162
2163     rvdw_lt_rcoul = (rvdw2 >= rcoul2);
2164     rcoul_lt_rvdw = (rcoul2 >= rvdw2);
2165
2166     if (bMakeQMMMnblist)
2167     {
2168         rm2 = rl2;
2169         rs2 = rl2;
2170     }
2171
2172     nl_sr     = ns->nl_sr;
2173     nsr       = ns->nsr;
2174     nl_lr_ljc = ns->nl_lr_ljc;
2175     nl_lr_one = ns->nl_lr_one;
2176     nlr_ljc   = ns->nlr_ljc;
2177     nlr_one   = ns->nlr_one;
2178
2179     /* Unpack arrays */
2180     cgcm    = fr->cg_cm;
2181     Nx      = grid->n[XX];
2182     Ny      = grid->n[YY];
2183     Nz      = grid->n[ZZ];
2184     grida   = grid->a;
2185     gridind = grid->index;
2186     gridnra = grid->nra;
2187     nns     = 0;
2188
2189     gridx      = grid->cell_size[XX];
2190     gridy      = grid->cell_size[YY];
2191     gridz      = grid->cell_size[ZZ];
2192     grid_x     = 1/gridx;
2193     grid_y     = 1/gridy;
2194     grid_z     = 1/gridz;
2195     copy_rvec(grid->cell_offset, grid_offset);
2196     copy_ivec(grid->ncpddc, ncpddc);
2197     dcx2       = grid->dcx2;
2198     dcy2       = grid->dcy2;
2199     dcz2       = grid->dcz2;
2200
2201 #ifdef ALLOW_OFFDIAG_LT_HALFDIAG
2202     zsh_ty = floor(-box[ZZ][YY]/box[YY][YY]+0.5);
2203     zsh_tx = floor(-box[ZZ][XX]/box[XX][XX]+0.5);
2204     ysh_tx = floor(-box[YY][XX]/box[XX][XX]+0.5);
2205     if (zsh_tx != 0 && ysh_tx != 0)
2206     {
2207         /* This could happen due to rounding, when both ratios are 0.5 */
2208         ysh_tx = 0;
2209     }
2210 #endif
2211
2212     debug_gmx();
2213
2214     if (fr->n_tpi)
2215     {
2216         /* We only want a list for the test particle */
2217         cg0 = cgsnr - 1;
2218     }
2219     else
2220     {
2221         cg0 = grid->icg0;
2222     }
2223     cg1 = grid->icg1;
2224
2225     /* Set the shift range */
2226     for (d = 0; d < DIM; d++)
2227     {
2228         sh0[d] = -1;
2229         sh1[d] = 1;
2230         /* Check if we need periodicity shifts.
2231          * Without PBC or with domain decomposition we don't need them.
2232          */
2233         if (d >= ePBC2npbcdim(fr->ePBC) || (bDomDec && dd->nc[d] > 1))
2234         {
2235             shp[d] = 0;
2236         }
2237         else
2238         {
2239             if (d == XX &&
2240                 box[XX][XX] - fabs(box[YY][XX]) - fabs(box[ZZ][XX]) < sqrt(rl2))
2241             {
2242                 shp[d] = 2;
2243             }
2244             else
2245             {
2246                 shp[d] = 1;
2247             }
2248         }
2249     }
2250
2251     /* Loop over charge groups */
2252     for (icg = cg0; (icg < cg1); icg++)
2253     {
2254         igid = GET_CGINFO_GID(cginfo[icg]);
2255         /* Skip this charge group if all energy groups are excluded! */
2256         if (bExcludeAlleg[igid])
2257         {
2258             continue;
2259         }
2260
2261         i0   = cgs->index[icg];
2262
2263         if (bMakeQMMMnblist)
2264         {
2265             /* Skip this charge group if it is not a QM atom while making a
2266              * QM/MM neighbourlist
2267              */
2268             if (md->bQM[i0] == FALSE)
2269             {
2270                 continue; /* MM particle, go to next particle */
2271             }
2272
2273             /* Compute the number of charge groups that fall within the control
2274              * of this one (icg)
2275              */
2276             naaj    = calc_naaj(icg, cgsnr);
2277             jcg0    = icg;
2278             jcg1    = icg + naaj;
2279             max_jcg = cgsnr;
2280         }
2281         else
2282         {
2283             /* make a normal neighbourlist */
2284
2285             if (bDomDec)
2286             {
2287                 /* Get the j charge-group and dd cell shift ranges */
2288                 dd_get_ns_ranges(cr->dd, icg, &jcg0, &jcg1, sh0, sh1);
2289                 max_jcg = 0;
2290             }
2291             else
2292             {
2293                 /* Compute the number of charge groups that fall within the control
2294                  * of this one (icg)
2295                  */
2296                 naaj = calc_naaj(icg, cgsnr);
2297                 jcg0 = icg;
2298                 jcg1 = icg + naaj;
2299
2300                 if (fr->n_tpi)
2301                 {
2302                     /* The i-particle is awlways the test particle,
2303                      * so we want all j-particles
2304                      */
2305                     max_jcg = cgsnr - 1;
2306                 }
2307                 else
2308                 {
2309                     max_jcg  = jcg1 - cgsnr;
2310                 }
2311             }
2312         }
2313
2314         i_egp_flags = fr->egp_flags + igid*ngid;
2315
2316         /* Set the exclusions for the atoms in charge group icg using a bitmask */
2317         setexcl(i0, cgs->index[icg+1], &top->excls, TRUE, bexcl);
2318
2319         ci2xyz(grid, icg, &cell_x, &cell_y, &cell_z);
2320
2321         /* Changed iicg to icg, DvdS 990115
2322          * (but see consistency check above, DvdS 990330)
2323          */
2324 #ifdef NS5DB
2325         fprintf(log, "icg=%5d, naaj=%5d, cell %d %d %d\n",
2326                 icg, naaj, cell_x, cell_y, cell_z);
2327 #endif
2328         /* Loop over shift vectors in three dimensions */
2329         for (tz = -shp[ZZ]; tz <= shp[ZZ]; tz++)
2330         {
2331             ZI = cgcm[icg][ZZ]+tz*box[ZZ][ZZ];
2332             /* Calculate range of cells in Z direction that have the shift tz */
2333             zgi = cell_z + tz*Nz;
2334 #define FAST_DD_NS
2335 #ifndef FAST_DD_NS
2336             get_dx(Nz, gridz, rl2, zgi, ZI, &dz0, &dz1, dcz2);
2337 #else
2338             get_dx_dd(Nz, gridz, rl2, zgi, ZI-grid_offset[ZZ],
2339                       ncpddc[ZZ], sh0[ZZ], sh1[ZZ], &dz0, &dz1, dcz2);
2340 #endif
2341             if (dz0 > dz1)
2342             {
2343                 continue;
2344             }
2345             for (ty = -shp[YY]; ty <= shp[YY]; ty++)
2346             {
2347                 YI = cgcm[icg][YY]+ty*box[YY][YY]+tz*box[ZZ][YY];
2348                 /* Calculate range of cells in Y direction that have the shift ty */
2349                 if (bTriclinicY)
2350                 {
2351                     ygi = (int)(Ny + (YI - grid_offset[YY])*grid_y) - Ny;
2352                 }
2353                 else
2354                 {
2355                     ygi = cell_y + ty*Ny;
2356                 }
2357 #ifndef FAST_DD_NS
2358                 get_dx(Ny, gridy, rl2, ygi, YI, &dy0, &dy1, dcy2);
2359 #else
2360                 get_dx_dd(Ny, gridy, rl2, ygi, YI-grid_offset[YY],
2361                           ncpddc[YY], sh0[YY], sh1[YY], &dy0, &dy1, dcy2);
2362 #endif
2363                 if (dy0 > dy1)
2364                 {
2365                     continue;
2366                 }
2367                 for (tx = -shp[XX]; tx <= shp[XX]; tx++)
2368                 {
2369                     XI = cgcm[icg][XX]+tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
2370                     /* Calculate range of cells in X direction that have the shift tx */
2371                     if (bTriclinicX)
2372                     {
2373                         xgi = (int)(Nx + (XI - grid_offset[XX])*grid_x) - Nx;
2374                     }
2375                     else
2376                     {
2377                         xgi = cell_x + tx*Nx;
2378                     }
2379 #ifndef FAST_DD_NS
2380                     get_dx(Nx, gridx, rl2, xgi*Nx, XI, &dx0, &dx1, dcx2);
2381 #else
2382                     get_dx_dd(Nx, gridx, rl2, xgi, XI-grid_offset[XX],
2383                               ncpddc[XX], sh0[XX], sh1[XX], &dx0, &dx1, dcx2);
2384 #endif
2385                     if (dx0 > dx1)
2386                     {
2387                         continue;
2388                     }
2389                     /* Adress: an explicit cg that has a weigthing function of 0 is excluded
2390                      *  from the neigbour list as it will not interact  */
2391                     if (fr->adress_type != eAdressOff)
2392                     {
2393                         if (md->wf[cgs->index[icg]] <= GMX_REAL_EPS && egp_explicit(fr, igid))
2394                         {
2395                             continue;
2396                         }
2397                     }
2398                     /* Get shift vector */
2399                     shift = XYZ2IS(tx, ty, tz);
2400 #ifdef NS5DB
2401                     range_check(shift, 0, SHIFTS);
2402 #endif
2403                     for (nn = 0; (nn < ngid); nn++)
2404                     {
2405                         nsr[nn]      = 0;
2406                         nlr_ljc[nn]  = 0;
2407                         nlr_one[nn]  = 0;
2408                     }
2409 #ifdef NS5DB
2410                     fprintf(log, "shift: %2d, dx0,1: %2d,%2d, dy0,1: %2d,%2d, dz0,1: %2d,%2d\n",
2411                             shift, dx0, dx1, dy0, dy1, dz0, dz1);
2412                     fprintf(log, "cgcm: %8.3f  %8.3f  %8.3f\n", cgcm[icg][XX],
2413                             cgcm[icg][YY], cgcm[icg][ZZ]);
2414                     fprintf(log, "xi:   %8.3f  %8.3f  %8.3f\n", XI, YI, ZI);
2415 #endif
2416                     for (dx = dx0; (dx <= dx1); dx++)
2417                     {
2418                         tmp1 = rl2 - dcx2[dx];
2419                         for (dy = dy0; (dy <= dy1); dy++)
2420                         {
2421                             tmp2 = tmp1 - dcy2[dy];
2422                             if (tmp2 > 0)
2423                             {
2424                                 for (dz = dz0; (dz <= dz1); dz++)
2425                                 {
2426                                     if (tmp2 > dcz2[dz])
2427                                     {
2428                                         /* Find grid-cell cj in which possible neighbours are */
2429                                         cj   = xyz2ci(Ny, Nz, dx, dy, dz);
2430
2431                                         /* Check out how many cgs (nrj) there in this cell */
2432                                         nrj  = gridnra[cj];
2433
2434                                         /* Find the offset in the cg list */
2435                                         cgj0 = gridind[cj];
2436
2437                                         /* Check if all j's are out of range so we
2438                                          * can skip the whole cell.
2439                                          * Should save some time, especially with DD.
2440                                          */
2441                                         if (nrj == 0 ||
2442                                             (grida[cgj0] >= max_jcg &&
2443                                              (grida[cgj0] >= jcg1 || grida[cgj0+nrj-1] < jcg0)))
2444                                         {
2445                                             continue;
2446                                         }
2447
2448                                         /* Loop over cgs */
2449                                         for (j = 0; (j < nrj); j++)
2450                                         {
2451                                             jjcg = grida[cgj0+j];
2452
2453                                             /* check whether this guy is in range! */
2454                                             if ((jjcg >= jcg0 && jjcg < jcg1) ||
2455                                                 (jjcg < max_jcg))
2456                                             {
2457                                                 r2 = calc_dx2(XI, YI, ZI, cgcm[jjcg]);
2458                                                 if (r2 < rl2)
2459                                                 {
2460                                                     /* jgid = gid[cgsatoms[cgsindex[jjcg]]]; */
2461                                                     jgid = GET_CGINFO_GID(cginfo[jjcg]);
2462                                                     /* check energy group exclusions */
2463                                                     if (!(i_egp_flags[jgid] & EGP_EXCL))
2464                                                     {
2465                                                         if (r2 < rs2)
2466                                                         {
2467                                                             if (nsr[jgid] >= MAX_CG)
2468                                                             {
2469                                                                 /* Add to short-range list */
2470                                                                 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2471                                                                             nsr[jgid], nl_sr[jgid],
2472                                                                             cgs->index, /* cgsatoms, */ bexcl,
2473                                                                             shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
2474                                                                 nsr[jgid] = 0;
2475                                                             }
2476                                                             nl_sr[jgid][nsr[jgid]++] = jjcg;
2477                                                         }
2478                                                         else if (r2 < rm2)
2479                                                         {
2480                                                             if (nlr_ljc[jgid] >= MAX_CG)
2481                                                             {
2482                                                                 /* Add to LJ+coulomb long-range list */
2483                                                                 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2484                                                                             nlr_ljc[jgid], nl_lr_ljc[jgid], top->cgs.index,
2485                                                                             bexcl, shift, fr, TRUE, TRUE, TRUE, fr->solvent_opt);
2486                                                                 nlr_ljc[jgid] = 0;
2487                                                             }
2488                                                             nl_lr_ljc[jgid][nlr_ljc[jgid]++] = jjcg;
2489                                                         }
2490                                                         else
2491                                                         {
2492                                                             if (nlr_one[jgid] >= MAX_CG)
2493                                                             {
2494                                                                 /* Add to long-range list with only coul, or only LJ */
2495                                                                 put_in_list(bHaveVdW, ngid, md, icg, jgid,
2496                                                                             nlr_one[jgid], nl_lr_one[jgid], top->cgs.index,
2497                                                                             bexcl, shift, fr, TRUE, rvdw_lt_rcoul, rcoul_lt_rvdw, fr->solvent_opt);
2498                                                                 nlr_one[jgid] = 0;
2499                                                             }
2500                                                             nl_lr_one[jgid][nlr_one[jgid]++] = jjcg;
2501                                                         }
2502                                                     }
2503                                                 }
2504                                                 nns++;
2505                                             }
2506                                         }
2507                                     }
2508                                 }
2509                             }
2510                         }
2511                     }
2512                     /* CHECK whether there is anything left in the buffers */
2513                     for (nn = 0; (nn < ngid); nn++)
2514                     {
2515                         if (nsr[nn] > 0)
2516                         {
2517                             put_in_list(bHaveVdW, ngid, md, icg, nn, nsr[nn], nl_sr[nn],
2518                                         cgs->index, /* cgsatoms, */ bexcl,
2519                                         shift, fr, FALSE, TRUE, TRUE, fr->solvent_opt);
2520                         }
2521
2522                         if (nlr_ljc[nn] > 0)
2523                         {
2524                             put_in_list(bHaveVdW, ngid, md, icg, nn, nlr_ljc[nn],
2525                                         nl_lr_ljc[nn], top->cgs.index,
2526                                         bexcl, shift, fr, TRUE, TRUE, TRUE, fr->solvent_opt);
2527                         }
2528
2529                         if (nlr_one[nn] > 0)
2530                         {
2531                             put_in_list(bHaveVdW, ngid, md, icg, nn, nlr_one[nn],
2532                                         nl_lr_one[nn], top->cgs.index,
2533                                         bexcl, shift, fr, TRUE, rvdw_lt_rcoul, rcoul_lt_rvdw, fr->solvent_opt);
2534                         }
2535                     }
2536                 }
2537             }
2538         }
2539         /* setexcl(nri,i_atoms,&top->atoms.excl,FALSE,bexcl); */
2540         setexcl(cgs->index[icg], cgs->index[icg+1], &top->excls, FALSE, bexcl);
2541     }
2542     /* No need to perform any left-over force calculations anymore (as we used to do here)
2543      * since we now save the proper long-range lists for later evaluation.
2544      */
2545
2546     debug_gmx();
2547
2548     /* Close neighbourlists */
2549     close_neighbor_lists(fr, bMakeQMMMnblist);
2550
2551     return nns;
2552 }
2553
2554 void ns_realloc_natoms(gmx_ns_t *ns, int natoms)
2555 {
2556     int i;
2557
2558     if (natoms > ns->nra_alloc)
2559     {
2560         ns->nra_alloc = over_alloc_dd(natoms);
2561         srenew(ns->bexcl, ns->nra_alloc);
2562         for (i = 0; i < ns->nra_alloc; i++)
2563         {
2564             ns->bexcl[i] = 0;
2565         }
2566     }
2567 }
2568
2569 void init_ns(FILE *fplog, const t_commrec *cr,
2570              gmx_ns_t *ns, t_forcerec *fr,
2571              const gmx_mtop_t *mtop,
2572              matrix box)
2573 {
2574     int  mt, icg, nr_in_cg, maxcg, i, j, jcg, ngid, ncg;
2575     t_block *cgs;
2576     char *ptr;
2577
2578     /* Compute largest charge groups size (# atoms) */
2579     nr_in_cg = 1;
2580     for (mt = 0; mt < mtop->nmoltype; mt++)
2581     {
2582         cgs = &mtop->moltype[mt].cgs;
2583         for (icg = 0; (icg < cgs->nr); icg++)
2584         {
2585             nr_in_cg = max(nr_in_cg, (int)(cgs->index[icg+1]-cgs->index[icg]));
2586         }
2587     }
2588
2589     /* Verify whether largest charge group is <= max cg.
2590      * This is determined by the type of the local exclusion type
2591      * Exclusions are stored in bits. (If the type is not large
2592      * enough, enlarge it, unsigned char -> unsigned short -> unsigned long)
2593      */
2594     maxcg = sizeof(t_excl)*8;
2595     if (nr_in_cg > maxcg)
2596     {
2597         gmx_fatal(FARGS, "Max #atoms in a charge group: %d > %d\n",
2598                   nr_in_cg, maxcg);
2599     }
2600
2601     ngid = mtop->groups.grps[egcENER].nr;
2602     snew(ns->bExcludeAlleg, ngid);
2603     for (i = 0; i < ngid; i++)
2604     {
2605         ns->bExcludeAlleg[i] = TRUE;
2606         for (j = 0; j < ngid; j++)
2607         {
2608             if (!(fr->egp_flags[i*ngid+j] & EGP_EXCL))
2609             {
2610                 ns->bExcludeAlleg[i] = FALSE;
2611             }
2612         }
2613     }
2614
2615     if (fr->bGrid)
2616     {
2617         /* Grid search */
2618         ns->grid = init_grid(fplog, fr);
2619         init_nsgrid_lists(fr, ngid, ns);
2620     }
2621     else
2622     {
2623         /* Simple search */
2624         snew(ns->ns_buf, ngid);
2625         for (i = 0; (i < ngid); i++)
2626         {
2627             snew(ns->ns_buf[i], SHIFTS);
2628         }
2629         ncg = ncg_mtop(mtop);
2630         snew(ns->simple_aaj, 2*ncg);
2631         for (jcg = 0; (jcg < ncg); jcg++)
2632         {
2633             ns->simple_aaj[jcg]     = jcg;
2634             ns->simple_aaj[jcg+ncg] = jcg;
2635         }
2636     }
2637
2638     /* Create array that determines whether or not atoms have VdW */
2639     snew(ns->bHaveVdW, fr->ntype);
2640     for (i = 0; (i < fr->ntype); i++)
2641     {
2642         for (j = 0; (j < fr->ntype); j++)
2643         {
2644             ns->bHaveVdW[i] = (ns->bHaveVdW[i] ||
2645                                (fr->bBHAM ?
2646                                 ((BHAMA(fr->nbfp, fr->ntype, i, j) != 0) ||
2647                                  (BHAMB(fr->nbfp, fr->ntype, i, j) != 0) ||
2648                                  (BHAMC(fr->nbfp, fr->ntype, i, j) != 0)) :
2649                                 ((C6(fr->nbfp, fr->ntype, i, j) != 0) ||
2650                                  (C12(fr->nbfp, fr->ntype, i, j) != 0))));
2651         }
2652     }
2653     if (debug)
2654     {
2655         pr_bvec(debug, 0, "bHaveVdW", ns->bHaveVdW, fr->ntype, TRUE);
2656     }
2657
2658     ns->nra_alloc = 0;
2659     ns->bexcl     = NULL;
2660     if (!DOMAINDECOMP(cr))
2661     {
2662         /* This could be reduced with particle decomposition */
2663         ns_realloc_natoms(ns, mtop->natoms);
2664     }
2665
2666     ns->nblist_initialized = FALSE;
2667
2668     /* nbr list debug dump */
2669     {
2670         char *ptr = getenv("GMX_DUMP_NL");
2671         if (ptr)
2672         {
2673             ns->dump_nl = strtol(ptr, NULL, 10);
2674             if (fplog)
2675             {
2676                 fprintf(fplog, "GMX_DUMP_NL = %d", ns->dump_nl);
2677             }
2678         }
2679         else
2680         {
2681             ns->dump_nl = 0;
2682         }
2683     }
2684 }
2685
2686
2687 int search_neighbours(FILE *log, t_forcerec *fr,
2688                       rvec x[], matrix box,
2689                       gmx_localtop_t *top,
2690                       gmx_groups_t *groups,
2691                       t_commrec *cr,
2692                       t_nrnb *nrnb, t_mdatoms *md,
2693                       real *lambda, real *dvdlambda,
2694                       gmx_grppairener_t *grppener,
2695                       gmx_bool bFillGrid,
2696                       gmx_bool bDoLongRangeNS,
2697                       gmx_bool bPadListsForKernels)
2698 {
2699     t_block  *cgs = &(top->cgs);
2700     rvec     box_size, grid_x0, grid_x1;
2701     int      i, j, m, ngid;
2702     real     min_size, grid_dens;
2703     int      nsearch;
2704     gmx_bool     bGrid;
2705     char     *ptr;
2706     gmx_bool     *i_egp_flags;
2707     int      cg_start, cg_end, start, end;
2708     gmx_ns_t *ns;
2709     t_grid   *grid;
2710     gmx_domdec_zones_t *dd_zones;
2711     put_in_list_t *put_in_list;
2712
2713     ns = &fr->ns;
2714
2715     /* Set some local variables */
2716     bGrid = fr->bGrid;
2717     ngid  = groups->grps[egcENER].nr;
2718
2719     for (m = 0; (m < DIM); m++)
2720     {
2721         box_size[m] = box[m][m];
2722     }
2723
2724     if (fr->ePBC != epbcNONE)
2725     {
2726         if (sqr(fr->rlistlong) >= max_cutoff2(fr->ePBC, box))
2727         {
2728             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.");
2729         }
2730         if (!bGrid)
2731         {
2732             min_size = min(box_size[XX], min(box_size[YY], box_size[ZZ]));
2733             if (2*fr->rlistlong >= min_size)
2734             {
2735                 gmx_fatal(FARGS, "One of the box diagonal elements has become smaller than twice the cut-off length.");
2736             }
2737         }
2738     }
2739
2740     if (DOMAINDECOMP(cr))
2741     {
2742         ns_realloc_natoms(ns, cgs->index[cgs->nr]);
2743     }
2744     debug_gmx();
2745
2746     /* Reset the neighbourlists */
2747     reset_neighbor_lists(fr, TRUE, TRUE);
2748
2749     if (bGrid && bFillGrid)
2750     {
2751
2752         grid = ns->grid;
2753         if (DOMAINDECOMP(cr))
2754         {
2755             dd_zones = domdec_zones(cr->dd);
2756         }
2757         else
2758         {
2759             dd_zones = NULL;
2760
2761             get_nsgrid_boundaries(grid->nboundeddim, box, NULL, NULL, NULL, NULL,
2762                                   cgs->nr, fr->cg_cm, grid_x0, grid_x1, &grid_dens);
2763
2764             grid_first(log, grid, NULL, NULL, fr->ePBC, box, grid_x0, grid_x1,
2765                        fr->rlistlong, grid_dens);
2766         }
2767         debug_gmx();
2768
2769         /* Don't know why this all is... (DvdS 3/99) */
2770 #ifndef SEGV
2771         start = 0;
2772         end   = cgs->nr;
2773 #else
2774         start = fr->cg0;
2775         end   = (cgs->nr+1)/2;
2776 #endif
2777
2778         if (DOMAINDECOMP(cr))
2779         {
2780             end = cgs->nr;
2781             fill_grid(log, dd_zones, grid, end, -1, end, fr->cg_cm);
2782             grid->icg0 = 0;
2783             grid->icg1 = dd_zones->izone[dd_zones->nizone-1].cg1;
2784         }
2785         else
2786         {
2787             fill_grid(log, NULL, grid, cgs->nr, fr->cg0, fr->hcg, fr->cg_cm);
2788             grid->icg0 = fr->cg0;
2789             grid->icg1 = fr->hcg;
2790             debug_gmx();
2791
2792             if (PARTDECOMP(cr))
2793             {
2794                 mv_grid(cr, grid);
2795             }
2796             debug_gmx();
2797         }
2798
2799         calc_elemnr(log, grid, start, end, cgs->nr);
2800         calc_ptrs(grid);
2801         grid_last(log, grid, start, end, cgs->nr);
2802
2803         if (gmx_debug_at)
2804         {
2805             check_grid(debug, grid);
2806             print_grid(debug, grid);
2807         }
2808     }
2809     else if (fr->n_tpi)
2810     {
2811         /* Set the grid cell index for the test particle only.
2812          * The cell to cg index is not corrected, but that does not matter.
2813          */
2814         fill_grid(log, NULL, ns->grid, fr->hcg, fr->hcg-1, fr->hcg, fr->cg_cm);
2815     }
2816     debug_gmx();
2817
2818     if (fr->adress_type == eAdressOff)
2819     {
2820         if (!fr->ns.bCGlist)
2821         {
2822             put_in_list = put_in_list_at;
2823         }
2824         else
2825         {
2826             put_in_list = put_in_list_cg;
2827         }
2828     }
2829     else
2830     {
2831         put_in_list = put_in_list_adress;
2832     }
2833
2834     /* Do the core! */
2835     if (bGrid)
2836     {
2837         grid    = ns->grid;
2838         nsearch = nsgrid_core(log, cr, fr, box, box_size, ngid, top,
2839                               grid, x, ns->bexcl, ns->bExcludeAlleg,
2840                               nrnb, md, lambda, dvdlambda, grppener,
2841                               put_in_list, ns->bHaveVdW,
2842                               bDoLongRangeNS, FALSE);
2843
2844         /* neighbour searching withouth QMMM! QM atoms have zero charge in
2845          * the classical calculation. The charge-charge interaction
2846          * between QM and MM atoms is handled in the QMMM core calculation
2847          * (see QMMM.c). The VDW however, we'd like to compute classically
2848          * and the QM MM atom pairs have just been put in the
2849          * corresponding neighbourlists. in case of QMMM we still need to
2850          * fill a special QMMM neighbourlist that contains all neighbours
2851          * of the QM atoms. If bQMMM is true, this list will now be made:
2852          */
2853         if (fr->bQMMM && fr->qr->QMMMscheme != eQMMMschemeoniom)
2854         {
2855             nsearch += nsgrid_core(log, cr, fr, box, box_size, ngid, top,
2856                                    grid, x, ns->bexcl, ns->bExcludeAlleg,
2857                                    nrnb, md, lambda, dvdlambda, grppener,
2858                                    put_in_list_qmmm, ns->bHaveVdW,
2859                                    bDoLongRangeNS, TRUE);
2860         }
2861     }
2862     else
2863     {
2864         nsearch = ns_simple_core(fr, top, md, box, box_size,
2865                                  ns->bexcl, ns->simple_aaj,
2866                                  ngid, ns->ns_buf, put_in_list, ns->bHaveVdW);
2867     }
2868     debug_gmx();
2869
2870 #ifdef DEBUG
2871     pr_nsblock(log);
2872 #endif
2873
2874     inc_nrnb(nrnb, eNR_NS, nsearch);
2875     /* inc_nrnb(nrnb,eNR_LR,fr->nlr); */
2876
2877     return nsearch;
2878 }
2879
2880 int natoms_beyond_ns_buffer(t_inputrec *ir, t_forcerec *fr, t_block *cgs,
2881                             matrix scale_tot, rvec *x)
2882 {
2883     int  cg0, cg1, cg, a0, a1, a, i, j;
2884     real rint, hbuf2, scale;
2885     rvec *cg_cm, cgsc;
2886     gmx_bool bIsotropic;
2887     int  nBeyond;
2888
2889     nBeyond = 0;
2890
2891     rint = max(ir->rcoulomb, ir->rvdw);
2892     if (ir->rlist < rint)
2893     {
2894         gmx_fatal(FARGS, "The neighbor search buffer has negative size: %f nm",
2895                   ir->rlist - rint);
2896     }
2897     cg_cm = fr->cg_cm;
2898
2899     cg0 = fr->cg0;
2900     cg1 = fr->hcg;
2901
2902     if (!EI_DYNAMICS(ir->eI) || !DYNAMIC_BOX(*ir))
2903     {
2904         hbuf2 = sqr(0.5*(ir->rlist - rint));
2905         for (cg = cg0; cg < cg1; cg++)
2906         {
2907             a0 = cgs->index[cg];
2908             a1 = cgs->index[cg+1];
2909             for (a = a0; a < a1; a++)
2910             {
2911                 if (distance2(cg_cm[cg], x[a]) > hbuf2)
2912                 {
2913                     nBeyond++;
2914                 }
2915             }
2916         }
2917     }
2918     else
2919     {
2920         bIsotropic = TRUE;
2921         scale      = scale_tot[0][0];
2922         for (i = 1; i < DIM; i++)
2923         {
2924             /* With anisotropic scaling, the original spherical ns volumes become
2925              * ellipsoids. To avoid costly transformations we use the minimum
2926              * eigenvalue of the scaling matrix for determining the buffer size.
2927              * Since the lower half is 0, the eigenvalues are the diagonal elements.
2928              */
2929             scale = min(scale, scale_tot[i][i]);
2930             if (scale_tot[i][i] != scale_tot[i-1][i-1])
2931             {
2932                 bIsotropic = FALSE;
2933             }
2934             for (j = 0; j < i; j++)
2935             {
2936                 if (scale_tot[i][j] != 0)
2937                 {
2938                     bIsotropic = FALSE;
2939                 }
2940             }
2941         }
2942         hbuf2 = sqr(0.5*(scale*ir->rlist - rint));
2943         if (bIsotropic)
2944         {
2945             for (cg = cg0; cg < cg1; cg++)
2946             {
2947                 svmul(scale, cg_cm[cg], cgsc);
2948                 a0 = cgs->index[cg];
2949                 a1 = cgs->index[cg+1];
2950                 for (a = a0; a < a1; a++)
2951                 {
2952                     if (distance2(cgsc, x[a]) > hbuf2)
2953                     {
2954                         nBeyond++;
2955                     }
2956                 }
2957             }
2958         }
2959         else
2960         {
2961             /* Anistropic scaling */
2962             for (cg = cg0; cg < cg1; cg++)
2963             {
2964                 /* Since scale_tot contains the transpose of the scaling matrix,
2965                  * we need to multiply with the transpose.
2966                  */
2967                 tmvmul_ur0(scale_tot, cg_cm[cg], cgsc);
2968                 a0 = cgs->index[cg];
2969                 a1 = cgs->index[cg+1];
2970                 for (a = a0; a < a1; a++)
2971                 {
2972                     if (distance2(cgsc, x[a]) > hbuf2)
2973                     {
2974                         nBeyond++;
2975                     }
2976                 }
2977             }
2978         }
2979     }
2980
2981     return nBeyond;
2982 }