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,2015, 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 #ifndef GMX_MDLIB_NSGRID_H
38 #define GMX_MDLIB_NSGRID_H
42 #include "gromacs/math/vectypes.h"
43 #include "gromacs/utility/real.h"
46 struct gmx_domdec_zones_t
;
51 typedef struct t_grid
{
52 int nr
; // Total number of charge groups
53 int nboundeddim
; // The number of bounded dimensions
54 int npbcdim
; // The number of dimensions with pbc
55 int ncg_ideal
; // The ideal number of cg's per cell
56 ivec n
; // The dimension of the grid
57 int ncells
; // Total number of cells
58 int cells_nalloc
; // Allocation size of index and nra
59 ivec ncpddc
; // The number of cells per DD cell
60 rvec cell_size
; // The size of the cells
61 rvec cell_offset
; // The offset of the cell (0,0,0)
62 int *cell_index
; // The cell number of each cg
63 int *index
; // The index into a for each cell
64 // The location of the cell in the index
65 // array can be found by calling xyz2ci
66 int *nra
; // The number of entries in a cell
67 int icg0
; // The start of the i-cg range
68 int icg1
; // The end of the i-cg range
71 int *a
; // The grid of cgs
72 int nr_alloc
; // Allocation size of cell_index and a
73 real
*dcx2
; // Squared distance from atom to j-cell
74 real
*dcy2
; // Squared distance from atom to j-cell
75 real
*dcz2
; // Squared distance from atom to j-cell
76 int dc_nalloc
; // Allocation size of dcx2, dyc2, dcz2
79 /*! \brief Used when estimating the interaction density.
81 * GRID_STDDEV_FAC * stddev estimates the interaction density. The
82 * value sqrt(3) == 1.73205080757 gives a uniform load for a
83 * rectangular 3D block of charge groups. For a sphere, it is not a
84 * bad approximation for 4x1x1 up to 4x2x2.
86 * \todo It would be nicer to use sqrt(3) here, when all code that
87 * includes this file is in C++, which will let us cope with the
88 * std::sqrt<T> on Windows. */
89 static const real GRID_STDDEV_FAC
= 1.73205080757;
91 /*! \brief The extent of the neighborsearch grid is a bit larger than sqrt(3)
92 * to account for less dense regions at the edges of the system.
94 static const real NSGRID_STDDEV_FAC
= 2.0;
96 #define NSGRID_SIGNAL_MOVED_FAC 4
97 /* A cell index of NSGRID_SIGNAL_MOVED_FAC*ncells signals
98 * that a charge group moved to another DD domain.
101 t_grid
*init_grid(FILE *fplog
, t_forcerec
*fr
);
103 void done_grid(t_grid
*grid
);
105 void get_nsgrid_boundaries(int nboundeddim
, matrix box
,
106 struct gmx_domdec_t
*dd
,
108 rvec
*gr0
, rvec
*gr1
,
110 rvec grid_x0
, rvec grid_x1
,
112 /* Return the ns grid boundaries grid_x0 and grid_x1
113 * and the estimate for the grid density.
114 * For non-bounded dimensions the boundaries are determined
115 * from the average and std.dev. of cgcm.
116 * The are determined from box, unless gr0!=NULL or gr1!=NULL,
117 * then they are taken from gr0 or gr1.
118 * With dd and unbounded dimensions, the proper grid borders for cells
119 * on the edges are determined from cgcm.
122 void grid_first(FILE *log
, t_grid
*grid
,
123 struct gmx_domdec_t
*dd
, const gmx_ddbox_t
*ddbox
, matrix box
, rvec izones_x0
, rvec izones_x1
,
124 real rlong
, real grid_density
);
126 void fill_grid(struct gmx_domdec_zones_t
*dd_zones
,
127 t_grid
*grid
, int ncg_tot
,
128 int cg0
, int cg1
, rvec cg_cm
[]);
129 /* Allocates space on the grid for ncg_tot cg's.
130 * Fills the grid with cg's from cg0 to cg1.
131 * When cg0 is -1, contiues filling from grid->nr to cg1.
134 void calc_elemnr(t_grid
*grid
, int cg0
, int cg1
, int ncg
);
136 void calc_ptrs(t_grid
*grid
);
138 void grid_last(t_grid
*grid
, int cg0
, int cg1
, int ncg
);
140 int xyz2ci_(int nry
, int nrz
, int x
, int y
, int z
);
141 #define xyz2ci(nry, nrz, x, y, z) ((nry)*(nrz)*(x)+(nrz)*(y)+(z))
142 /* Return the cell index */
144 void ci2xyz(t_grid
*grid
, int i
, int *x
, int *y
, int *z
);
146 void check_grid(t_grid
*grid
);
148 void print_grid(FILE *log
, t_grid
*grid
);