Split lines with many copyright years
[alexxy/gromacs.git] / src / gromacs / mdlib / nsgrid.cpp
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,2015,2017,2018 by the GROMACS development team.
7  * Copyright (c) 2019,2020, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 /* This file is completely threadsafe - keep it that way! */
39 #include "gmxpre.h"
40
41 #include "nsgrid.h"
42
43 #include <cmath>
44 #include <cstdio>
45 #include <cstdlib>
46
47 #include <algorithm>
48
49 #include "gromacs/domdec/dlb.h"
50 #include "gromacs/domdec/domdec.h"
51 #include "gromacs/domdec/domdec_struct.h"
52 #include "gromacs/fileio/pdbio.h"
53 #include "gromacs/gmxlib/network.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/mdtypes/forcerec.h"
56 #include "gromacs/pbcutil/pbc.h"
57 #include "gromacs/utility/fatalerror.h"
58 #include "gromacs/utility/futil.h"
59 #include "gromacs/utility/smalloc.h"
60
61 /***********************************
62  *         Grid Routines
63  ***********************************/
64
65 static const char* range_warn =
66         "Explanation: During neighborsearching, we assign each particle to a grid\n"
67         "based on its coordinates. If your system contains collisions or parameter\n"
68         "errors that give particles very high velocities you might end up with some\n"
69         "coordinates being +-Infinity or NaN (not-a-number). Obviously, we cannot\n"
70         "put these on a grid, so this is usually where we detect those errors.\n"
71         "Make sure your system is properly energy-minimized and that the potential\n"
72         "energy seems reasonable before trying again.";
73
74 static void calc_x_av_stddev(int n, rvec* x, rvec av, rvec stddev)
75 {
76     dvec s1, s2;
77     int  i, d;
78
79     clear_dvec(s1);
80     clear_dvec(s2);
81
82     for (i = 0; i < n; i++)
83     {
84         for (d = 0; d < DIM; d++)
85         {
86             s1[d] += x[i][d];
87             s2[d] += x[i][d] * x[i][d];
88         }
89     }
90
91     dsvmul(1.0 / n, s1, s1);
92     dsvmul(1.0 / n, s2, s2);
93
94     for (d = 0; d < DIM; d++)
95     {
96         av[d]     = s1[d];
97         stddev[d] = std::sqrt(s2[d] - s1[d] * s1[d]);
98     }
99 }
100
101 static void get_nsgrid_boundaries_vac(real av, real stddev, real* bound0, real* bound1, real* bdens0, real* bdens1)
102 {
103     /* Set the grid to 2 times the standard deviation of
104      * the charge group centers in both directions.
105      * For a uniform bounded distribution the width is sqrt(3)*stddev,
106      * so all charge groups fall within the width.
107      * For a sphere stddev is r/sqrt(5): 99.2% falls within the width.
108      * For a Gaussian distribution 98% fall within the width.
109      */
110     *bound0 = av - NSGRID_STDDEV_FAC * stddev;
111     *bound1 = av + NSGRID_STDDEV_FAC * stddev;
112
113     *bdens0 = av - GRID_STDDEV_FAC * stddev;
114     *bdens1 = av + GRID_STDDEV_FAC * stddev;
115 }
116
117 static void dd_box_bounds_to_ns_bounds(real box0, real box_size, real* gr0, real* gr1)
118 {
119     real av, stddev;
120
121     /* Redetermine av and stddev from the DD box boundaries */
122     av     = box0 + 0.5 * box_size;
123     stddev = 0.5 * box_size / GRID_STDDEV_FAC;
124
125     *gr0 = av - NSGRID_STDDEV_FAC * stddev;
126     *gr1 = av + NSGRID_STDDEV_FAC * stddev;
127 }
128
129 void get_nsgrid_boundaries(int           nboundeddim,
130                            matrix        box,
131                            gmx_domdec_t* dd,
132                            gmx_ddbox_t*  ddbox,
133                            rvec*         gr0,
134                            rvec*         gr1,
135                            int           ncg,
136                            rvec*         cgcm,
137                            rvec          grid_x0,
138                            rvec          grid_x1,
139                            real*         grid_density)
140 {
141     rvec av, stddev;
142     real vol, bdens0, bdens1;
143     int  d;
144
145     if (nboundeddim < DIM)
146     {
147         calc_x_av_stddev(ncg, cgcm, av, stddev);
148     }
149
150     vol = 1;
151     for (d = 0; d < DIM; d++)
152     {
153         if (d < nboundeddim)
154         {
155             grid_x0[d] = (gr0 != nullptr ? (*gr0)[d] : 0);
156             grid_x1[d] = (gr1 != nullptr ? (*gr1)[d] : box[d][d]);
157             vol *= (grid_x1[d] - grid_x0[d]);
158         }
159         else
160         {
161             if (ddbox == nullptr)
162             {
163                 get_nsgrid_boundaries_vac(av[d], stddev[d], &grid_x0[d], &grid_x1[d], &bdens0, &bdens1);
164             }
165             else
166             {
167                 /* Temporary fix which uses global ddbox boundaries
168                  * for unbounded dimensions.
169                  * Should be replaced by local boundaries, which makes
170                  * the ns grid smaller and does not require global comm.
171                  */
172                 dd_box_bounds_to_ns_bounds(ddbox->box0[d], ddbox->box_size[d], &grid_x0[d], &grid_x1[d]);
173                 bdens0 = grid_x0[d];
174                 bdens1 = grid_x1[d];
175             }
176             /* Check for a DD cell not at a lower edge */
177             if (dd != nullptr && gr0 != nullptr && dd->ci[d] > 0)
178             {
179                 grid_x0[d] = (*gr0)[d];
180                 bdens0     = (*gr0)[d];
181             }
182             /* Check for a DD cell not at a higher edge */
183             if (dd != nullptr && gr1 != nullptr && dd->ci[d] < dd->numCells[d] - 1)
184             {
185                 grid_x1[d] = (*gr1)[d];
186                 bdens1     = (*gr1)[d];
187             }
188             vol *= (bdens1 - bdens0);
189         }
190
191         if (debug)
192         {
193             fprintf(debug, "Set grid boundaries dim %d: %f %f\n", d, grid_x0[d], grid_x1[d]);
194         }
195     }
196
197     *grid_density = ncg / vol;
198 }
199
200 static void set_grid_sizes(matrix              box,
201                            rvec                izones_x0,
202                            rvec                izones_x1,
203                            real                rlist,
204                            const gmx_domdec_t* dd,
205                            const gmx_ddbox_t*  ddbox,
206                            t_grid*             grid,
207                            real                grid_density)
208 {
209     int      i, j;
210     gmx_bool bDD, bDDRect;
211     rvec     izones_size;
212     real     inv_r_ideal, size, add_tric, radd;
213
214     for (i = 0; (i < DIM); i++)
215     {
216         if (debug)
217         {
218             fprintf(debug, "set_grid_sizes, i-zone bounds for dim %d: %6.3f %6.3f\n", i,
219                     izones_x0[i], izones_x1[i]);
220         }
221         izones_size[i] = izones_x1[i] - izones_x0[i];
222     }
223
224     /* Use the ideal number of cg's per cell to set the ideal cell size */
225     inv_r_ideal = std::cbrt(grid_density / grid->ncg_ideal);
226     if (rlist > 0 && inv_r_ideal * rlist < 1)
227     {
228         inv_r_ideal = 1 / rlist;
229     }
230     if (debug)
231     {
232         fprintf(debug, "CG density %f ideal ns cell size %f\n", grid_density, 1 / inv_r_ideal);
233     }
234
235     clear_rvec(grid->cell_offset);
236     for (i = 0; (i < DIM); i++)
237     {
238         /* Initial settings, for DD might change below */
239         grid->cell_offset[i] = izones_x0[i];
240         size                 = izones_size[i];
241
242         bDD = (dd != nullptr) && (dd->numCells[i] > 1);
243         if (!bDD)
244         {
245             bDDRect = FALSE;
246         }
247         else
248         {
249             /* With DD grid cell jumps only the first decomposition
250              * direction has uniform DD cell boundaries.
251              */
252             bDDRect = !((ddbox->tric_dir[i] != 0) || (dd_dlb_is_on(dd) && i != dd->dim[0]));
253
254             radd = rlist;
255             if (i >= ddbox->npbcdim
256                 && (rlist == 0 || izones_x1[i] + radd > ddbox->box0[i] + ddbox->box_size[i]))
257             {
258                 radd = ddbox->box0[i] + ddbox->box_size[i] - izones_x1[i];
259                 if (radd < 0)
260                 {
261                     radd = 0;
262                 }
263             }
264
265             /* With DD we only need a grid of one DD cell size + rlist */
266             if (bDDRect)
267             {
268                 size += radd;
269             }
270             else
271             {
272                 size += radd / ddbox->skew_fac[i];
273             }
274
275             /* Check if the cell boundary in this direction is
276              * perpendicular to the Cartesian axis.
277              * Since grid->npbcdim isan integer that in principle can take
278              * any value, we help the compiler avoid warnings and potentially
279              * optimize by ensuring that j < DIM here.
280              */
281             for (j = i + 1; j < grid->npbcdim && j < DIM; j++)
282             {
283                 if (box[j][i] != 0)
284                 {
285                     /* Correct the offset for the home cell location */
286                     grid->cell_offset[i] += izones_x0[j] * box[j][i] / box[j][j];
287
288                     /* Correct the offset and size for the off-diagonal
289                      * displacement of opposing DD cell corners.
290                      */
291                     /* Without rouding we would need to add:
292                      * box[j][i]*rlist/(dd->skew_fac[i]*box[j][j])
293                      */
294                     /* Determine the shift for the corners of the triclinic box */
295                     add_tric = izones_size[j] * box[j][i] / box[j][j];
296                     if (dd->ndim == 1 && j == ZZ)
297                     {
298                         /* With 1D domain decomposition the cg's are not in
299                          * the triclinic box, but trilinic x-y and rectangular y-z.
300                          * Therefore we need to add the shift from the trilinic
301                          * corner to the corner at y=0.
302                          */
303                         add_tric += -box[YY][XX] * box[ZZ][YY] / box[YY][YY];
304                     }
305                     if (box[j][i] < 0)
306                     {
307                         grid->cell_offset[i] += add_tric;
308                         size -= add_tric;
309                     }
310                     else
311                     {
312                         size += add_tric;
313                     }
314                 }
315             }
316         }
317         if (!bDDRect)
318         {
319             /* No DD or the box is triclinic is this direction:
320              * we will use the normal grid ns that checks all cells
321              * that are within cut-off distance of the i-particle.
322              */
323             grid->n[i] = gmx::roundToInt(size * inv_r_ideal);
324             if (grid->n[i] < 2)
325             {
326                 grid->n[i] = 2;
327             }
328             grid->cell_size[i] = size / grid->n[i];
329             grid->ncpddc[i]    = 0;
330         }
331         else
332         {
333             /* We use grid->ncpddc[i] such that all particles
334              * in one ns cell belong to a single DD cell only.
335              * We can then beforehand exclude certain ns grid cells
336              * for non-home i-particles.
337              */
338             grid->ncpddc[i] = gmx::roundToInt(izones_size[i] * inv_r_ideal);
339             if (grid->ncpddc[i] < 2)
340             {
341                 grid->ncpddc[i] = 2;
342             }
343             grid->cell_size[i] = izones_size[i] / grid->ncpddc[i];
344             grid->n[i]         = grid->ncpddc[i] + static_cast<int>(radd / grid->cell_size[i]) + 1;
345         }
346         if (debug)
347         {
348             fprintf(debug, "grid dim %d size %d x %f: %f - %f\n", i, grid->n[i], grid->cell_size[i],
349                     grid->cell_offset[i], grid->cell_offset[i] + grid->n[i] * grid->cell_size[i]);
350         }
351     }
352
353     if (debug)
354     {
355         fprintf(debug, "CG ncg ideal %d, actual density %.1f\n", grid->ncg_ideal,
356                 grid_density * grid->cell_size[XX] * grid->cell_size[YY] * grid->cell_size[ZZ]);
357     }
358 }
359
360 t_grid* init_grid(FILE* fplog, t_forcerec* fr)
361 {
362     char*   ptr;
363     t_grid* grid;
364
365     snew(grid, 1);
366
367     grid->npbcdim = ePBC2npbcdim(fr->ePBC);
368
369     if (fr->ePBC == epbcXY && fr->nwall == 2)
370     {
371         grid->nboundeddim = 3;
372     }
373     else
374     {
375         grid->nboundeddim = grid->npbcdim;
376     }
377
378     if (debug)
379     {
380         fprintf(debug, "The coordinates are bounded in %d dimensions\n", grid->nboundeddim);
381     }
382
383     /* The ideal number of cg's per ns grid cell seems to be 10 */
384     grid->ncg_ideal = 10;
385     ptr             = getenv("GMX_NSCELL_NCG");
386     if (ptr)
387     {
388         sscanf(ptr, "%d", &grid->ncg_ideal);
389         if (fplog)
390         {
391             fprintf(fplog, "Set ncg_ideal to %d\n", grid->ncg_ideal);
392         }
393         if (grid->ncg_ideal <= 0)
394         {
395             gmx_fatal(FARGS, "The number of cg's per cell should be > 0");
396         }
397     }
398     if (debug)
399     {
400         fprintf(debug, "Set ncg_ideal to %d\n", grid->ncg_ideal);
401     }
402
403     return grid;
404 }
405
406 void done_grid(t_grid* grid)
407 {
408     if (grid == nullptr)
409     {
410         return;
411     }
412     grid->nr = 0;
413     clear_ivec(grid->n);
414     grid->ncells = 0;
415     sfree(grid->cell_index);
416     sfree(grid->a);
417     sfree(grid->index);
418     sfree(grid->nra);
419     grid->cells_nalloc = 0;
420     sfree(grid->dcx2);
421     sfree(grid->dcy2);
422     sfree(grid->dcz2);
423     grid->dc_nalloc = 0;
424
425     if (debug)
426     {
427         fprintf(debug, "Successfully freed memory for grid pointers.");
428     }
429     sfree(grid);
430 }
431
432 int xyz2ci_(int nry, int nrz, int x, int y, int z)
433 /* Return the cell index */
434 {
435     return (nry * nrz * x + nrz * y + z);
436 }
437
438 void ci2xyz(t_grid* grid, int i, int* x, int* y, int* z)
439 /* Return x,y and z from the cell index */
440 {
441     int ci;
442
443     range_check_mesg(i, 0, grid->nr, range_warn);
444
445     ci = grid->cell_index[i];
446     *x = ci / (grid->n[YY] * grid->n[ZZ]);
447     ci -= (*x) * grid->n[YY] * grid->n[ZZ];
448     *y = ci / grid->n[ZZ];
449     ci -= (*y) * grid->n[ZZ];
450     *z = ci;
451 }
452
453 static int ci_not_used(const ivec n)
454 {
455     /* Return an improbable value */
456     return xyz2ci(n[YY], n[ZZ], 3 * n[XX], 3 * n[YY], 3 * n[ZZ]);
457 }
458
459 static void set_grid_ncg(t_grid* grid, int ncg)
460 {
461     int nr_old, i;
462
463     grid->nr = ncg;
464     if (grid->nr + 1 > grid->nr_alloc)
465     {
466         nr_old         = grid->nr_alloc;
467         grid->nr_alloc = over_alloc_dd(grid->nr) + 1;
468         srenew(grid->cell_index, grid->nr_alloc);
469         for (i = nr_old; i < grid->nr_alloc; i++)
470         {
471             grid->cell_index[i] = 0;
472         }
473         srenew(grid->a, grid->nr_alloc);
474     }
475 }
476
477 void grid_first(FILE*              fplog,
478                 t_grid*            grid,
479                 gmx_domdec_t*      dd,
480                 const gmx_ddbox_t* ddbox,
481                 matrix             box,
482                 rvec               izones_x0,
483                 rvec               izones_x1,
484                 real               rlist,
485                 real               grid_density)
486 {
487     int i, m;
488
489     set_grid_sizes(box, izones_x0, izones_x1, rlist, dd, ddbox, grid, grid_density);
490
491     grid->ncells = grid->n[XX] * grid->n[YY] * grid->n[ZZ];
492
493     if (grid->ncells + 1 > grid->cells_nalloc)
494     {
495         /* Allocate double the size so we have some headroom */
496         grid->cells_nalloc = 2 * grid->ncells;
497         srenew(grid->nra, grid->cells_nalloc + 1);
498         srenew(grid->index, grid->cells_nalloc + 1);
499
500         if (fplog)
501         {
502             fprintf(fplog, "Grid: %d x %d x %d cells\n", grid->n[XX], grid->n[YY], grid->n[ZZ]);
503         }
504     }
505
506     m = std::max(grid->n[XX], std::max(grid->n[YY], grid->n[ZZ]));
507     if (m > grid->dc_nalloc)
508     {
509         /* Allocate with double the initial size for box scaling */
510         grid->dc_nalloc = 2 * m;
511         srenew(grid->dcx2, grid->dc_nalloc);
512         srenew(grid->dcy2, grid->dc_nalloc);
513         srenew(grid->dcz2, grid->dc_nalloc);
514     }
515
516     grid->nr = 0;
517     for (i = 0; (i < grid->ncells); i++)
518     {
519         grid->nra[i] = 0;
520     }
521 }
522
523 static void calc_bor(int cg0, int cg1, int ncg, int CG0[2], int CG1[2])
524 {
525     if (cg1 > ncg)
526     {
527         CG0[0] = cg0;
528         CG1[0] = ncg;
529         CG0[1] = 0;
530         CG1[1] = cg1 - ncg;
531     }
532     else
533     {
534         CG0[0] = cg0;
535         CG1[0] = cg1;
536         CG0[1] = 0;
537         CG1[1] = 0;
538     }
539     if (debug)
540     {
541         int m;
542
543         fprintf(debug, "calc_bor: cg0=%d, cg1=%d, ncg=%d\n", cg0, cg1, ncg);
544         for (m = 0; (m < 2); m++)
545         {
546             fprintf(debug, "CG0[%d]=%d, CG1[%d]=%d\n", m, CG0[m], m, CG1[m]);
547         }
548     }
549 }
550
551 void calc_elemnr(t_grid* grid, int cg0, int cg1, int ncg)
552 {
553     int  CG0[2], CG1[2];
554     int* cell_index = grid->cell_index;
555     int* nra        = grid->nra;
556     int  i, m, ncells;
557     int  ci, not_used;
558
559     ncells = grid->ncells;
560     if (ncells <= 0)
561     {
562         gmx_fatal(FARGS, "Number of grid cells is zero. Probably the system and box collapsed.\n");
563     }
564
565     not_used = ci_not_used(grid->n);
566
567     calc_bor(cg0, cg1, ncg, CG0, CG1);
568     for (m = 0; (m < 2); m++)
569     {
570         for (i = CG0[m]; (i < CG1[m]); i++)
571         {
572             ci = cell_index[i];
573             if (ci != not_used)
574             {
575                 range_check_mesg(ci, 0, ncells, range_warn);
576                 nra[ci]++;
577             }
578         }
579     }
580 }
581
582 void calc_ptrs(t_grid* grid)
583 {
584     int* index = grid->index;
585     int* nra   = grid->nra;
586     int  ix, iy, iz, ci, nr;
587     int  nnra, ncells;
588
589     ncells = grid->ncells;
590     if (ncells <= 0)
591     {
592         gmx_fatal(FARGS, "Number of grid cells is zero. Probably the system and box collapsed.\n");
593     }
594
595     ci = nr = 0;
596     for (ix = 0; (ix < grid->n[XX]); ix++)
597     {
598         for (iy = 0; (iy < grid->n[YY]); iy++)
599         {
600             for (iz = 0; (iz < grid->n[ZZ]); iz++, ci++)
601             {
602                 range_check_mesg(ci, 0, ncells, range_warn);
603                 index[ci] = nr;
604                 nnra      = nra[ci];
605                 nr += nnra;
606                 nra[ci] = 0;
607             }
608         }
609     }
610 }
611
612 void grid_last(t_grid* grid, int cg0, int cg1, int ncg)
613 {
614     int  CG0[2], CG1[2];
615     int  i, m;
616     int  ci, not_used, ind, ncells;
617     int* cell_index = grid->cell_index;
618     int* nra        = grid->nra;
619     int* index      = grid->index;
620     int* a          = grid->a;
621
622     ncells = grid->ncells;
623     if (ncells <= 0)
624     {
625         gmx_fatal(FARGS, "Number of grid cells is zero. Probably the system and box collapsed.\n");
626     }
627
628     not_used = ci_not_used(grid->n);
629
630     calc_bor(cg0, cg1, ncg, CG0, CG1);
631     for (m = 0; (m < 2); m++)
632     {
633         for (i = CG0[m]; (i < CG1[m]); i++)
634         {
635             ci = cell_index[i];
636             if (ci != not_used)
637             {
638                 range_check_mesg(ci, 0, ncells, range_warn);
639                 ind = index[ci] + nra[ci]++;
640                 range_check_mesg(ind, 0, grid->nr, range_warn);
641                 a[ind] = i;
642             }
643         }
644     }
645 }
646
647 void fill_grid(gmx_domdec_zones_t* dd_zones, t_grid* grid, int ncg_tot, int cg0, int cg1, rvec cg_cm[])
648 {
649     int*     cell_index;
650     int      nry, nrz;
651     rvec     n_box, offset;
652     int      zone, ccg0, ccg1, cg, d, not_used;
653     ivec     shift0, useall, b0, b1, ind;
654     gmx_bool bUse;
655
656     if (cg0 == -1)
657     {
658         /* We have already filled the grid up to grid->ncg,
659          * continue from there.
660          */
661         cg0 = grid->nr;
662     }
663
664     set_grid_ncg(grid, ncg_tot);
665
666     cell_index = grid->cell_index;
667
668     /* Initiate cell borders */
669     nry = grid->n[YY];
670     nrz = grid->n[ZZ];
671     for (d = 0; d < DIM; d++)
672     {
673         if (grid->cell_size[d] > 0)
674         {
675             n_box[d] = 1 / grid->cell_size[d];
676         }
677         else
678         {
679             n_box[d] = 0;
680         }
681     }
682     copy_rvec(grid->cell_offset, offset);
683
684     if (debug)
685     {
686         fprintf(debug, "Filling grid from %d to %d\n", cg0, cg1);
687     }
688
689     if (dd_zones == nullptr)
690     {
691         for (cg = cg0; cg < cg1; cg++)
692         {
693             for (d = 0; d < DIM; d++)
694             {
695                 ind[d] = static_cast<int>((cg_cm[cg][d] - offset[d]) * n_box[d]);
696                 /* With pbc we should be done here.
697                  * Without pbc cg's outside the grid
698                  * should be assigned to the closest grid cell.
699                  */
700                 if (ind[d] < 0)
701                 {
702                     ind[d] = 0;
703                 }
704                 else if (ind[d] >= grid->n[d])
705                 {
706                     ind[d] = grid->n[d] - 1;
707                 }
708             }
709             cell_index[cg] = xyz2ci(nry, nrz, ind[XX], ind[YY], ind[ZZ]);
710         }
711     }
712     else
713     {
714         for (zone = 0; zone < dd_zones->n; zone++)
715         {
716             ccg0 = dd_zones->cg_range[zone];
717             ccg1 = dd_zones->cg_range[zone + 1];
718             if (ccg1 <= cg0 || ccg0 >= cg1)
719             {
720                 continue;
721             }
722
723             /* Determine the ns grid cell limits for this DD zone */
724             for (d = 0; d < DIM; d++)
725             {
726                 shift0[d] = dd_zones->shift[zone][d];
727                 useall[d] = static_cast<int>(shift0[d] == 0 || d >= grid->npbcdim);
728                 /* Check if we need to do normal or optimized grid assignments.
729                  * Normal is required for dims without DD or triclinic dims.
730                  * DD edge cell on dims without pbc will be automatically
731                  * be correct, since the shift=0 zones with have b0 and b1
732                  * set to the grid boundaries and there are no shift=1 zones.
733                  */
734                 if (grid->ncpddc[d] == 0)
735                 {
736                     b0[d] = 0;
737                     b1[d] = grid->n[d];
738                 }
739                 else
740                 {
741                     if (shift0[d] == 0)
742                     {
743                         b0[d] = 0;
744                         b1[d] = grid->ncpddc[d];
745                     }
746                     else
747                     {
748                         /* shift = 1 */
749                         b0[d] = grid->ncpddc[d];
750                         b1[d] = grid->n[d];
751                     }
752                 }
753             }
754
755             not_used = ci_not_used(grid->n);
756
757             /* Put all the charge groups of this DD zone on the grid */
758             for (cg = ccg0; cg < ccg1; cg++)
759             {
760                 if (cell_index[cg] == -1)
761                 {
762                     /* This cg has moved to another node */
763                     cell_index[cg] = NSGRID_SIGNAL_MOVED_FAC * grid->ncells;
764                     continue;
765                 }
766
767                 bUse = TRUE;
768                 for (d = 0; d < DIM; d++)
769                 {
770                     ind[d] = static_cast<int>((cg_cm[cg][d] - offset[d]) * n_box[d]);
771                     /* Here we have to correct for rounding problems,
772                      * as this cg_cm to cell index operation is not necessarily
773                      * binary identical to the operation for the DD zone assignment
774                      * and therefore a cg could end up in an unused grid cell.
775                      * For dimensions without pbc we need to check
776                      * for cells on the edge if charge groups are beyond
777                      * the grid and if so, store them in the closest cell.
778                      */
779                     if (ind[d] < b0[d])
780                     {
781                         ind[d] = b0[d];
782                     }
783                     else if (ind[d] >= b1[d])
784                     {
785                         if (useall[d])
786                         {
787                             ind[d] = b1[d] - 1;
788                         }
789                         else
790                         {
791                             /* Charge groups in this DD zone further away than the cut-off
792                              * in direction do not participate in non-bonded interactions.
793                              */
794                             bUse = FALSE;
795                         }
796                     }
797                 }
798                 if (cg > grid->nr_alloc)
799                 {
800                     fprintf(stderr, "WARNING: nra_alloc %d cg0 %d cg1 %d cg %d\n", grid->nr_alloc,
801                             cg0, cg1, cg);
802                 }
803                 if (bUse)
804                 {
805                     cell_index[cg] = xyz2ci(nry, nrz, ind[XX], ind[YY], ind[ZZ]);
806                 }
807                 else
808                 {
809                     cell_index[cg] = not_used;
810                 }
811             }
812         }
813     }
814 }
815
816 void check_grid(t_grid* grid)
817 {
818     int ix, iy, iz, ci, cci, nra;
819
820     if (grid->ncells <= 0)
821     {
822         gmx_fatal(FARGS, "Number of grid cells is zero. Probably the system and box collapsed.\n");
823     }
824
825     ci  = 0;
826     cci = 0;
827     for (ix = 0; (ix < grid->n[XX]); ix++)
828     {
829         for (iy = 0; (iy < grid->n[YY]); iy++)
830         {
831             for (iz = 0; (iz < grid->n[ZZ]); iz++, ci++)
832             {
833                 if (ci > 0)
834                 {
835                     nra = grid->index[ci] - grid->index[cci];
836                     if (nra != grid->nra[cci])
837                     {
838                         gmx_fatal(FARGS, "nra=%d, grid->nra=%d, cci=%d", nra, grid->nra[cci], cci);
839                     }
840                 }
841                 cci = xyz2ci(grid->n[YY], grid->n[ZZ], ix, iy, iz);
842                 range_check_mesg(cci, 0, grid->ncells, range_warn);
843
844                 if (cci != ci)
845                 {
846                     gmx_fatal(FARGS, "ci = %d, cci = %d", ci, cci);
847                 }
848             }
849         }
850     }
851 }
852
853 void print_grid(FILE* log, t_grid* grid)
854 {
855     int i, nra, index;
856     int ix, iy, iz, ci;
857
858     fprintf(log, "nr:        %d\n", grid->nr);
859     fprintf(log, "nrx:       %d\n", grid->n[XX]);
860     fprintf(log, "nry:       %d\n", grid->n[YY]);
861     fprintf(log, "nrz:       %d\n", grid->n[ZZ]);
862     fprintf(log, "ncg_ideal: %d\n", grid->ncg_ideal);
863     fprintf(log, "    i  cell_index\n");
864     for (i = 0; (i < grid->nr); i++)
865     {
866         fprintf(log, "%5d  %5d\n", i, grid->cell_index[i]);
867     }
868     fprintf(log, "cells\n");
869     fprintf(log, " ix iy iz   nr  index  cgs...\n");
870     ci = 0;
871     for (ix = 0; (ix < grid->n[XX]); ix++)
872     {
873         for (iy = 0; (iy < grid->n[YY]); iy++)
874         {
875             for (iz = 0; (iz < grid->n[ZZ]); iz++, ci++)
876             {
877                 index = grid->index[ci];
878                 nra   = grid->nra[ci];
879                 fprintf(log, "%3d%3d%3d%5d%5d", ix, iy, iz, nra, index);
880                 for (i = 0; (i < nra); i++)
881                 {
882                     fprintf(log, "%5d", grid->a[index + i]);
883                 }
884                 fprintf(log, "\n");
885             }
886         }
887     }
888     fflush(log);
889 }