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