2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
37 /* This file is completely threadsafe - keep it that way! */
40 #include "gromacs/legacyheaders/nsgrid.h"
45 #include "gromacs/fileio/pdbio.h"
46 #include "gromacs/legacyheaders/domdec.h"
47 #include "gromacs/legacyheaders/macros.h"
48 #include "gromacs/legacyheaders/network.h"
49 #include "gromacs/legacyheaders/typedefs.h"
50 #include "gromacs/legacyheaders/types/commrec.h"
51 #include "gromacs/math/vec.h"
52 #include "gromacs/pbcutil/pbc.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/futil.h"
55 #include "gromacs/utility/smalloc.h"
57 /***********************************
59 ***********************************/
61 const char *range_warn =
62 "Explanation: During neighborsearching, we assign each particle to a grid\n"
63 "based on its coordinates. If your system contains collisions or parameter\n"
64 "errors that give particles very high velocities you might end up with some\n"
65 "coordinates being +-Infinity or NaN (not-a-number). Obviously, we cannot\n"
66 "put these on a grid, so this is usually where we detect those errors.\n"
67 "Make sure your system is properly energy-minimized and that the potential\n"
68 "energy seems reasonable before trying again.";
70 static void calc_x_av_stddev(int n, rvec *x, rvec av, rvec stddev)
78 for (i = 0; i < n; i++)
80 for (d = 0; d < DIM; d++)
83 s2[d] += x[i][d]*x[i][d];
87 dsvmul(1.0/n, s1, s1);
88 dsvmul(1.0/n, s2, s2);
90 for (d = 0; d < DIM; d++)
93 stddev[d] = sqrt(s2[d] - s1[d]*s1[d]);
97 void get_nsgrid_boundaries_vac(real av, real stddev,
98 real *bound0, real *bound1,
99 real *bdens0, real *bdens1)
101 /* Set the grid to 2 times the standard deviation of
102 * the charge group centers in both directions.
103 * For a uniform bounded distribution the width is sqrt(3)*stddev,
104 * so all charge groups fall within the width.
105 * For a sphere stddev is r/sqrt(5): 99.2% falls within the width.
106 * For a Gaussian distribution 98% fall within the width.
108 *bound0 = av - NSGRID_STDDEV_FAC*stddev;
109 *bound1 = av + NSGRID_STDDEV_FAC*stddev;
111 *bdens0 = av - GRID_STDDEV_FAC*stddev;
112 *bdens1 = av + GRID_STDDEV_FAC*stddev;
115 static void dd_box_bounds_to_ns_bounds(real box0, real box_size,
116 real *gr0, real *gr1)
120 /* Redetermine av and stddev from the DD box boundaries */
121 av = box0 + 0.5*box_size;
122 stddev = 0.5*box_size/GRID_STDDEV_FAC;
124 *gr0 = av - NSGRID_STDDEV_FAC*stddev;
125 *gr1 = av + NSGRID_STDDEV_FAC*stddev;
128 void get_nsgrid_boundaries(int nboundeddim, matrix box,
130 gmx_ddbox_t *ddbox, rvec *gr0, rvec *gr1,
132 rvec grid_x0, rvec grid_x1,
136 real vol, bdens0, bdens1;
139 if (nboundeddim < DIM)
141 calc_x_av_stddev(ncg, cgcm, av, stddev);
145 for (d = 0; d < DIM; d++)
149 grid_x0[d] = (gr0 != NULL ? (*gr0)[d] : 0);
150 grid_x1[d] = (gr1 != NULL ? (*gr1)[d] : box[d][d]);
151 vol *= (grid_x1[d] - grid_x0[d]);
157 get_nsgrid_boundaries_vac(av[d], stddev[d],
158 &grid_x0[d], &grid_x1[d],
163 /* Temporary fix which uses global ddbox boundaries
164 * for unbounded dimensions.
165 * Should be replaced by local boundaries, which makes
166 * the ns grid smaller and does not require global comm.
168 dd_box_bounds_to_ns_bounds(ddbox->box0[d], ddbox->box_size[d],
169 &grid_x0[d], &grid_x1[d]);
173 /* Check for a DD cell not at a lower edge */
174 if (dd != NULL && gr0 != NULL && dd->ci[d] > 0)
176 grid_x0[d] = (*gr0)[d];
179 /* Check for a DD cell not at a higher edge */
180 if (dd != NULL && gr1 != NULL && dd->ci[d] < dd->nc[d]-1)
182 grid_x1[d] = (*gr1)[d];
185 vol *= (bdens1 - bdens0);
190 fprintf(debug, "Set grid boundaries dim %d: %f %f\n",
191 d, grid_x0[d], grid_x1[d]);
195 *grid_density = ncg/vol;
198 static void set_grid_sizes(matrix box, rvec izones_x0, rvec izones_x1, real rlist,
199 const gmx_domdec_t *dd, const gmx_ddbox_t *ddbox,
204 gmx_bool bDD, bDDRect;
207 real inv_r_ideal, size, add_tric, radd;
209 for (i = 0; (i < DIM); i++)
214 "set_grid_sizes, i-zone bounds for dim %d: %6.3f %6.3f\n",
215 i, izones_x0[i], izones_x1[i]);
217 izones_size[i] = izones_x1[i] - izones_x0[i];
220 /* Use the ideal number of cg's per cell to set the ideal cell size */
221 inv_r_ideal = pow(grid_density/grid->ncg_ideal, 1.0/3.0);
222 if (rlist > 0 && inv_r_ideal*rlist < 1)
224 inv_r_ideal = 1/rlist;
228 fprintf(debug, "CG density %f ideal ns cell size %f\n",
229 grid_density, 1/inv_r_ideal);
232 clear_rvec(grid->cell_offset);
233 for (i = 0; (i < DIM); i++)
235 /* Initial settings, for DD might change below */
236 grid->cell_offset[i] = izones_x0[i];
237 size = izones_size[i];
239 bDD = dd && (dd->nc[i] > 1);
246 /* With DD grid cell jumps only the first decomposition
247 * direction has uniform DD cell boundaries.
249 bDDRect = !(ddbox->tric_dir[i] ||
250 (dd->bGridJump && i != dd->dim[0]));
253 if (i >= ddbox->npbcdim &&
255 izones_x1[i] + radd > ddbox->box0[i] + ddbox->box_size[i]))
257 radd = ddbox->box0[i] + ddbox->box_size[i] - izones_x1[i];
264 /* With DD we only need a grid of one DD cell size + rlist */
271 size += radd/ddbox->skew_fac[i];
274 /* Check if the cell boundary in this direction is
275 * perpendicular to the Cartesian axis.
277 for (j = i+1; j < grid->npbcdim; j++)
281 /* Correct the offset for the home cell location */
282 grid->cell_offset[i] += izones_x0[j]*box[j][i]/box[j][j];
284 /* Correct the offset and size for the off-diagonal
285 * displacement of opposing DD cell corners.
287 /* Without rouding we would need to add:
288 * box[j][i]*rlist/(dd->skew_fac[i]*box[j][j])
290 /* Determine the shift for the corners of the triclinic box */
291 add_tric = izones_size[j]*box[j][i]/box[j][j];
292 if (dd->ndim == 1 && j == ZZ)
294 /* With 1D domain decomposition the cg's are not in
295 * the triclinic box, but trilinic x-y and rectangular y-z.
296 * Therefore we need to add the shift from the trilinic
297 * corner to the corner at y=0.
299 add_tric += -box[YY][XX]*box[ZZ][YY]/box[YY][YY];
303 grid->cell_offset[i] += add_tric;
315 /* No DD or the box is triclinic is this direction:
316 * we will use the normal grid ns that checks all cells
317 * that are within cut-off distance of the i-particle.
319 grid->n[i] = (int)(size*inv_r_ideal + 0.5);
324 grid->cell_size[i] = size/grid->n[i];
329 /* We use grid->ncpddc[i] such that all particles
330 * in one ns cell belong to a single DD cell only.
331 * We can then beforehand exclude certain ns grid cells
332 * for non-home i-particles.
334 grid->ncpddc[i] = (int)(izones_size[i]*inv_r_ideal + 0.5);
335 if (grid->ncpddc[i] < 2)
339 grid->cell_size[i] = izones_size[i]/grid->ncpddc[i];
340 grid->n[i] = grid->ncpddc[i] + (int)(radd/grid->cell_size[i]) + 1;
344 fprintf(debug, "grid dim %d size %d x %f: %f - %f\n",
345 i, grid->n[i], grid->cell_size[i],
346 grid->cell_offset[i],
347 grid->cell_offset[i]+grid->n[i]*grid->cell_size[i]);
353 fprintf(debug, "CG ncg ideal %d, actual density %.1f\n",
354 grid->ncg_ideal, grid_density*grid->cell_size[XX]*grid->cell_size[YY]*grid->cell_size[ZZ]);
358 t_grid *init_grid(FILE *fplog, t_forcerec *fr)
366 grid->npbcdim = ePBC2npbcdim(fr->ePBC);
368 if (fr->ePBC == epbcXY && fr->nwall == 2)
370 grid->nboundeddim = 3;
374 grid->nboundeddim = grid->npbcdim;
379 fprintf(debug, "The coordinates are bounded in %d dimensions\n",
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");
388 sscanf(ptr, "%d", &grid->ncg_ideal);
391 fprintf(fplog, "Set ncg_ideal to %d\n", grid->ncg_ideal);
393 if (grid->ncg_ideal <= 0)
395 gmx_fatal(FARGS, "The number of cg's per cell should be > 0");
400 fprintf(debug, "Set ncg_ideal to %d\n", grid->ncg_ideal);
406 void done_grid(t_grid *grid)
411 sfree(grid->cell_index);
415 grid->cells_nalloc = 0;
423 fprintf(debug, "Successfully freed memory for grid pointers.");
427 int xyz2ci_(int nry, int nrz, int x, int y, int z)
428 /* Return the cell index */
430 return (nry*nrz*x+nrz*y+z);
433 void ci2xyz(t_grid *grid, int i, int *x, int *y, int *z)
434 /* Return x,y and z from the cell index */
438 range_check_mesg(i, 0, grid->nr, range_warn);
440 ci = grid->cell_index[i];
441 *x = ci / (grid->n[YY]*grid->n[ZZ]);
442 ci -= (*x)*grid->n[YY]*grid->n[ZZ];
443 *y = ci / grid->n[ZZ];
444 ci -= (*y)*grid->n[ZZ];
448 static int ci_not_used(ivec n)
450 /* Return an improbable value */
451 return xyz2ci(n[YY], n[ZZ], 3*n[XX], 3*n[YY], 3*n[ZZ]);
454 static void set_grid_ncg(t_grid *grid, int ncg)
459 if (grid->nr+1 > grid->nr_alloc)
461 nr_old = grid->nr_alloc;
462 grid->nr_alloc = over_alloc_dd(grid->nr) + 1;
463 srenew(grid->cell_index, grid->nr_alloc);
464 for (i = nr_old; i < grid->nr_alloc; i++)
466 grid->cell_index[i] = 0;
468 srenew(grid->a, grid->nr_alloc);
472 void grid_first(FILE *fplog, t_grid *grid,
473 gmx_domdec_t *dd, const gmx_ddbox_t *ddbox,
474 matrix box, rvec izones_x0, rvec izones_x1,
475 real rlistlong, real grid_density)
480 set_grid_sizes(box, izones_x0, izones_x1, rlistlong, dd, ddbox, grid, grid_density);
482 grid->ncells = grid->n[XX]*grid->n[YY]*grid->n[ZZ];
484 if (grid->ncells+1 > grid->cells_nalloc)
486 /* Allocate double the size so we have some headroom */
487 grid->cells_nalloc = 2*grid->ncells;
488 srenew(grid->nra, grid->cells_nalloc+1);
489 srenew(grid->index, grid->cells_nalloc+1);
493 fprintf(fplog, "Grid: %d x %d x %d cells\n",
494 grid->n[XX], grid->n[YY], grid->n[ZZ]);
498 m = max(grid->n[XX], max(grid->n[YY], grid->n[ZZ]));
499 if (m > grid->dc_nalloc)
501 /* Allocate with double the initial size for box scaling */
502 grid->dc_nalloc = 2*m;
503 srenew(grid->dcx2, grid->dc_nalloc);
504 srenew(grid->dcy2, grid->dc_nalloc);
505 srenew(grid->dcz2, grid->dc_nalloc);
509 for (i = 0; (i < grid->ncells); i++)
515 static void calc_bor(int cg0, int cg1, int ncg, int CG0[2], int CG1[2])
535 fprintf(debug, "calc_bor: cg0=%d, cg1=%d, ncg=%d\n", cg0, cg1, ncg);
536 for (m = 0; (m < 2); m++)
538 fprintf(debug, "CG0[%d]=%d, CG1[%d]=%d\n", m, CG0[m], m, CG1[m]);
544 void calc_elemnr(t_grid *grid, int cg0, int cg1, int ncg)
547 int *cell_index = grid->cell_index;
548 int *nra = grid->nra;
552 ncells = grid->ncells;
555 gmx_fatal(FARGS, "Number of grid cells is zero. Probably the system and box collapsed.\n");
558 not_used = ci_not_used(grid->n);
560 calc_bor(cg0, cg1, ncg, CG0, CG1);
561 for (m = 0; (m < 2); m++)
563 for (i = CG0[m]; (i < CG1[m]); i++)
568 range_check_mesg(ci, 0, ncells, range_warn);
575 void calc_ptrs(t_grid *grid)
577 int *index = grid->index;
578 int *nra = grid->nra;
579 int ix, iy, iz, ci, nr;
582 ncells = grid->ncells;
585 gmx_fatal(FARGS, "Number of grid cells is zero. Probably the system and box collapsed.\n");
589 for (ix = 0; (ix < grid->n[XX]); ix++)
591 for (iy = 0; (iy < grid->n[YY]); iy++)
593 for (iz = 0; (iz < grid->n[ZZ]); iz++, ci++)
595 range_check_mesg(ci, 0, ncells, range_warn);
605 void grid_last(t_grid *grid, int cg0, int cg1, int ncg)
609 int ci, not_used, ind, ncells;
610 int *cell_index = grid->cell_index;
611 int *nra = grid->nra;
612 int *index = grid->index;
615 ncells = grid->ncells;
618 gmx_fatal(FARGS, "Number of grid cells is zero. Probably the system and box collapsed.\n");
621 not_used = ci_not_used(grid->n);
623 calc_bor(cg0, cg1, ncg, CG0, CG1);
624 for (m = 0; (m < 2); m++)
626 for (i = CG0[m]; (i < CG1[m]); i++)
631 range_check_mesg(ci, 0, ncells, range_warn);
632 ind = index[ci]+nra[ci]++;
633 range_check_mesg(ind, 0, grid->nr, range_warn);
640 void fill_grid(gmx_domdec_zones_t *dd_zones,
641 t_grid *grid, int ncg_tot,
642 int cg0, int cg1, rvec cg_cm[])
647 int zone, ccg0, ccg1, cg, d, not_used;
648 ivec shift0, useall, b0, b1, ind;
653 /* We have already filled the grid up to grid->ncg,
654 * continue from there.
659 set_grid_ncg(grid, ncg_tot);
661 cell_index = grid->cell_index;
663 /* Initiate cell borders */
667 for (d = 0; d < DIM; d++)
669 if (grid->cell_size[d] > 0)
671 n_box[d] = 1/grid->cell_size[d];
678 copy_rvec(grid->cell_offset, offset);
682 fprintf(debug, "Filling grid from %d to %d\n", cg0, cg1);
686 if (dd_zones == NULL)
688 for (cg = cg0; cg < cg1; cg++)
690 for (d = 0; d < DIM; d++)
692 ind[d] = (cg_cm[cg][d] - offset[d])*n_box[d];
693 /* With pbc we should be done here.
694 * Without pbc cg's outside the grid
695 * should be assigned to the closest grid cell.
701 else if (ind[d] >= grid->n[d])
703 ind[d] = grid->n[d] - 1;
706 cell_index[cg] = xyz2ci(nry, nrz, ind[XX], ind[YY], ind[ZZ]);
711 for (zone = 0; zone < dd_zones->n; zone++)
713 ccg0 = dd_zones->cg_range[zone];
714 ccg1 = dd_zones->cg_range[zone+1];
715 if (ccg1 <= cg0 || ccg0 >= cg1)
720 /* Determine the ns grid cell limits for this DD zone */
721 for (d = 0; d < DIM; d++)
723 shift0[d] = dd_zones->shift[zone][d];
724 useall[d] = (shift0[d] == 0 || d >= grid->npbcdim);
725 /* Check if we need to do normal or optimized grid assignments.
726 * Normal is required for dims without DD or triclinic dims.
727 * DD edge cell on dims without pbc will be automatically
728 * be correct, since the shift=0 zones with have b0 and b1
729 * set to the grid boundaries and there are no shift=1 zones.
731 if (grid->ncpddc[d] == 0)
741 b1[d] = grid->ncpddc[d];
746 b0[d] = grid->ncpddc[d];
752 not_used = ci_not_used(grid->n);
754 /* Put all the charge groups of this DD zone on the grid */
755 for (cg = ccg0; cg < ccg1; cg++)
757 if (cell_index[cg] == -1)
759 /* This cg has moved to another node */
760 cell_index[cg] = NSGRID_SIGNAL_MOVED_FAC*grid->ncells;
765 for (d = 0; d < DIM; d++)
767 ind[d] = (cg_cm[cg][d] - offset[d])*n_box[d];
768 /* Here we have to correct for rounding problems,
769 * as this cg_cm to cell index operation is not necessarily
770 * binary identical to the operation for the DD zone assignment
771 * and therefore a cg could end up in an unused grid cell.
772 * For dimensions without pbc we need to check
773 * for cells on the edge if charge groups are beyond
774 * the grid and if so, store them in the closest cell.
780 else if (ind[d] >= b1[d])
788 /* Charge groups in this DD zone further away than the cut-off
789 * in direction do not participate in non-bonded interactions.
795 if (cg > grid->nr_alloc)
797 fprintf(stderr, "WARNING: nra_alloc %d cg0 %d cg1 %d cg %d\n",
798 grid->nr_alloc, cg0, cg1, cg);
802 cell_index[cg] = xyz2ci(nry, nrz, ind[XX], ind[YY], ind[ZZ]);
806 cell_index[cg] = not_used;
815 void check_grid(t_grid *grid)
817 int ix, iy, iz, ci, cci, nra;
819 if (grid->ncells <= 0)
821 gmx_fatal(FARGS, "Number of grid cells is zero. Probably the system and box collapsed.\n");
826 for (ix = 0; (ix < grid->n[XX]); ix++)
828 for (iy = 0; (iy < grid->n[YY]); iy++)
830 for (iz = 0; (iz < grid->n[ZZ]); iz++, ci++)
834 nra = grid->index[ci]-grid->index[cci];
835 if (nra != grid->nra[cci])
837 gmx_fatal(FARGS, "nra=%d, grid->nra=%d, cci=%d",
838 nra, grid->nra[cci], cci);
841 cci = xyz2ci(grid->n[YY], grid->n[ZZ], ix, iy, iz);
842 range_check_mesg(cci, 0, grid->ncells, range_warn);
846 gmx_fatal(FARGS, "ci = %d, cci = %d", ci, cci);
853 void print_grid(FILE *log, t_grid *grid)
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++)
866 fprintf(log, "%5d %5d\n", i, grid->cell_index[i]);
868 fprintf(log, "cells\n");
869 fprintf(log, " ix iy iz nr index cgs...\n");
871 for (ix = 0; (ix < grid->n[XX]); ix++)
873 for (iy = 0; (iy < grid->n[YY]); iy++)
875 for (iz = 0; (iz < grid->n[ZZ]); iz++, ci++)
877 index = grid->index[ci];
879 fprintf(log, "%3d%3d%3d%5d%5d", ix, iy, iz, nra, index);
880 for (i = 0; (i < nra); i++)
882 fprintf(log, "%5d", grid->a[index+i]);