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