3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Gromacs Runs On Most of All Computer Systems
42 #define GRID_STDDEV_FAC sqrt(3)
43 #define NSGRID_STDDEV_FAC 2.0
45 * GRID_STDDEV_FAC * stddev is used to estimate the interaction density.
46 * sqrt(3) gives a uniform load for a rectangular block of cg's.
47 * For a sphere it is not a bad approximation for 4x1x1 up to 4x2x2.
49 * The extent of the neighborsearch grid is a bit larger than sqrt(3)
50 * to account for less dense regions at the edges of the system.
53 #define NSGRID_SIGNAL_MOVED_FAC 4
54 /* A cell index of NSGRID_SIGNAL_MOVED_FAC*ncells signals
55 * that a charge group moved to another DD domain.
58 t_grid *init_grid(FILE *fplog, t_forcerec *fr);
60 void done_grid(t_grid *grid);
62 void get_nsgrid_boundaries(int nboundeddim, matrix box,
67 rvec grid_x0, rvec grid_x1,
69 /* Return the ns grid boundaries grid_x0 and grid_x1
70 * and the estimate for the grid density.
71 * For non-bounded dimensions the boundaries are determined
72 * from the average and std.dev. of cgcm.
73 * The are determined from box, unless gr0!=NULL or gr1!=NULL,
74 * then they are taken from gr0 or gr1.
75 * With dd and unbounded dimensions, the proper grid borders for cells
76 * on the edges are determined from cgcm.
79 void grid_first(FILE *log, t_grid *grid,
80 gmx_domdec_t *dd, const gmx_ddbox_t *ddbox, matrix box, rvec izones_x0, rvec izones_x1,
81 real rlong, real grid_density);
83 void fill_grid(gmx_domdec_zones_t *dd_zones,
84 t_grid *grid, int ncg_tot,
85 int cg0, int cg1, rvec cg_cm[]);
86 /* Allocates space on the grid for ncg_tot cg's.
87 * Fills the grid with cg's from cg0 to cg1.
88 * When cg0 is -1, contiues filling from grid->nr to cg1.
91 void calc_elemnr(t_grid *grid, int cg0, int cg1, int ncg);
93 void calc_ptrs(t_grid *grid);
95 void grid_last(t_grid *grid, int cg0, int cg1, int ncg);
97 int xyz2ci_(int nry, int nrz, int x, int y, int z);
98 #define xyz2ci(nry, nrz, x, y, z) ((nry)*(nrz)*(x)+(nrz)*(y)+(z))
99 /* Return the cell index */
101 void ci2xyz(t_grid *grid, int i, int *x, int *y, int *z);
103 void check_grid(t_grid *grid);
105 void print_grid(FILE *log, t_grid *grid);
107 void mv_grid(t_commrec *cr, t_grid *grid);
108 /* Move the grid over processors */