4fbb9bb8825dd2417409d7f8f5c34995a0ea256c
[alexxy/gromacs.git] / src / gromacs / legacyheaders / types / commrec.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, 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 _commrec_h
38 #define _commrec_h
39
40 #include "../../utility/gmxmpi.h"
41 #include "../typedefs.h"
42
43 #ifdef __cplusplus
44 extern "C" {
45 #endif
46
47
48 #define DD_MAXZONE  8
49 #define DD_MAXIZONE 4
50
51 typedef struct gmx_domdec_master *gmx_domdec_master_p_t;
52
53 typedef struct {
54     int  j0;     /* j-zone start               */
55     int  j1;     /* j-zone end                 */
56     int  cg1;    /* i-charge-group end         */
57     int  jcg0;   /* j-charge-group start       */
58     int  jcg1;   /* j-charge-group end         */
59     ivec shift0; /* Minimum shifts to consider */
60     ivec shift1; /* Maximum shifts to consider */
61 } gmx_domdec_ns_ranges_t;
62
63 typedef struct {
64     rvec x0;     /* Zone lower corner in triclinic coordinates         */
65     rvec x1;     /* Zone upper corner in triclinic coordinates         */
66     rvec bb_x0;  /* Zone bounding box lower corner in Cartesian coords */
67     rvec bb_x1;  /* Zone bounding box upper corner in Cartesian coords */
68 } gmx_domdec_zone_size_t;
69
70 struct gmx_domdec_zones_t {
71     /* The number of zones including the home zone */
72     int                    n;
73     /* The shift of the zones with respect to the home zone */
74     ivec                   shift[DD_MAXZONE];
75     /* The charge group boundaries for the zones */
76     int                    cg_range[DD_MAXZONE+1];
77     /* The number of neighbor search zones with i-particles */
78     int                    nizone;
79     /* The neighbor search charge group ranges for each i-zone */
80     gmx_domdec_ns_ranges_t izone[DD_MAXIZONE];
81     /* Boundaries of the zones */
82     gmx_domdec_zone_size_t size[DD_MAXZONE];
83     /* The cg density of the home zone */
84     real                   dens_zone0;
85 };
86
87 typedef struct gmx_ga2la *gmx_ga2la_t;
88
89 typedef struct gmx_hash *gmx_hash_t;
90
91 typedef struct gmx_reverse_top *gmx_reverse_top_p_t;
92
93 typedef struct gmx_domdec_constraints *gmx_domdec_constraints_p_t;
94
95 typedef struct gmx_domdec_specat_comm *gmx_domdec_specat_comm_p_t;
96
97 typedef struct gmx_domdec_comm *gmx_domdec_comm_p_t;
98
99 typedef struct gmx_pme_comm_n_box *gmx_pme_comm_n_box_p_t;
100
101 struct gmx_ddbox_t {
102     int  npbcdim;
103     int  nboundeddim;
104     rvec box0;
105     rvec box_size;
106     /* Tells if the box is skewed for each of the three cartesian directions */
107     ivec tric_dir;
108     rvec skew_fac;
109     /* Orthogonal vectors for triclinic cells, Cartesian index */
110     rvec v[DIM][DIM];
111     /* Normal vectors for the cells walls */
112     rvec normal[DIM];
113 };
114
115
116 typedef struct {
117     /* these buffers are used as destination buffers if MPI_IN_PLACE isn't
118        supported.*/
119     int             *ibuf; /* for ints */
120     int              ibuf_alloc;
121
122     gmx_int64_t     *libuf;
123     int              libuf_alloc;
124
125     float           *fbuf; /* for floats */
126     int              fbuf_alloc;
127
128     double          *dbuf; /* for doubles */
129     int              dbuf_alloc;
130 } mpi_in_place_buf_t;
131
132
133 struct gmx_domdec_t {
134     /* The DD particle-particle nodes only */
135     /* The communication setup within the communicator all
136      * defined in dd->comm in domdec.c
137      */
138     int                    nnodes;
139     MPI_Comm               mpi_comm_all;
140     /* Use MPI_Sendrecv communication instead of non-blocking calls */
141     gmx_bool               bSendRecv2;
142     /* The local DD cell index and rank */
143     ivec                   ci;
144     int                    rank;
145     ivec                   master_ci;
146     int                    masterrank;
147     /* Communication with the PME only nodes */
148     int                    pme_nodeid;
149     gmx_bool               pme_receive_vir_ener;
150     gmx_pme_comm_n_box_p_t cnb;
151     int                    nreq_pme;
152     MPI_Request            req_pme[8];
153
154
155     /* The communication setup, identical for each cell, cartesian index */
156     ivec     nc;
157     int      ndim;
158     ivec     dim; /* indexed by 0 to ndim */
159     gmx_bool bGridJump;
160
161     /* PBC from dim 0 to npbcdim */
162     int npbcdim;
163
164     /* Screw PBC? */
165     gmx_bool bScrewPBC;
166
167     /* Forward and backward neighboring cells, indexed by 0 to ndim */
168     int  neighbor[DIM][2];
169
170     /* Only available on the master node */
171     gmx_domdec_master_p_t ma;
172
173     /* Are there inter charge group constraints */
174     gmx_bool bInterCGcons;
175     gmx_bool bInterCGsettles;
176
177     /* Global atom number to interaction list */
178     gmx_reverse_top_p_t reverse_top;
179     int                 nbonded_global;
180     int                 nbonded_local;
181
182     /* The number of inter charge-group exclusions */
183     int  n_intercg_excl;
184
185     /* Vsite stuff */
186     gmx_hash_t                 ga2la_vsite;
187     gmx_domdec_specat_comm_p_t vsite_comm;
188
189     /* Constraint stuff */
190     gmx_domdec_constraints_p_t constraints;
191     gmx_domdec_specat_comm_p_t constraint_comm;
192
193     /* The local to gobal charge group index and local cg to local atom index */
194     int   ncg_home;
195     int   ncg_tot;
196     int  *index_gl;
197     int  *cgindex;
198     int   cg_nalloc;
199     /* Local atom to local cg index, only for special cases */
200     int  *la2lc;
201     int   la2lc_nalloc;
202
203     /* The number of home atoms */
204     int   nat_home;
205     /* The total number of atoms: home and received zones */
206     int   nat_tot;
207     /* Index from the local atoms to the global atoms */
208     int  *gatindex;
209     int   gatindex_nalloc;
210
211     /* Global atom number to local atom number list */
212     gmx_ga2la_t ga2la;
213
214     /* Communication stuff */
215     gmx_domdec_comm_p_t comm;
216
217     /* The partioning count, to keep track of the state */
218     gmx_int64_t ddp_count;
219
220
221     /* gmx_pme_recv_f buffer */
222     int   pme_recv_f_alloc;
223     rvec *pme_recv_f_buf;
224
225 };
226
227 struct gmx_multisim_t {
228     int       nsim;
229     int       sim;
230     MPI_Group mpi_group_masters;
231     MPI_Comm  mpi_comm_masters;
232     /* these buffers are used as destination buffers if MPI_IN_PLACE isn't
233        supported.*/
234     mpi_in_place_buf_t *mpb;
235 };
236
237 #define DUTY_PP  (1<<0)
238 #define DUTY_PME (1<<1)
239
240 typedef struct {
241     int      bUse;
242     MPI_Comm comm_intra;
243     int      rank_intra;
244     MPI_Comm comm_inter;
245
246 } gmx_nodecomm_t;
247
248 struct t_commrec {
249     /* The nodeids in one sim are numbered sequentially from 0.
250      * All communication within some simulation should happen
251      * in mpi_comm_mysim, or its subset mpi_comm_mygroup.
252      */
253     int sim_nodeid, nnodes, npmenodes;
254
255     /* thread numbers: */
256     /* Not used yet: int threadid, nthreads; */
257     /* The nodeid in the PP/PME, PP or PME group */
258     int      nodeid;
259     MPI_Comm mpi_comm_mysim;
260     MPI_Comm mpi_comm_mygroup;
261
262     /* MPI ranks within a physical node for hardware access */
263     int            nrank_intranode;    /* nr of ranks on this physical node */
264     int            rank_intranode;     /* our rank on this physical node */
265     int            nrank_pp_intranode; /* as nrank_intranode, for particle-particle only */
266     int            rank_pp_intranode;  /* as rank_intranode, for particle-particle only */
267
268     gmx_nodecomm_t nc;
269
270     /* For domain decomposition */
271     gmx_domdec_t *dd;
272
273     /* The duties of this node, see the defines above */
274     int             duty;
275
276     gmx_multisim_t *ms;
277
278     /* these buffers are used as destination buffers if MPI_IN_PLACE isn't
279        supported.*/
280     mpi_in_place_buf_t *mpb;
281 };
282
283 #define MASTERNODE(cr)     (((cr)->nodeid == 0) || !PAR(cr))
284 /* #define MASTERTHREAD(cr)   ((cr)->threadid == 0) */
285 /* #define MASTER(cr)         (MASTERNODE(cr) && MASTERTHREAD(cr)) */
286 #define MASTER(cr)         MASTERNODE(cr)
287 #define SIMMASTER(cr)      ((MASTER(cr) && ((cr)->duty & DUTY_PP)) || !PAR(cr))
288 #define NODEPAR(cr)        ((cr)->nnodes > 1)
289 /* #define THREADPAR(cr)      ((cr)->nthreads > 1) */
290 /* #define PAR(cr)            (NODEPAR(cr) || THREADPAR(cr)) */
291 #define PAR(cr)            NODEPAR(cr)
292 #define RANK(cr, nodeid)    (nodeid)
293 #define MASTERRANK(cr)     (0)
294
295 /* Note that even with particle decomposition removed, the use of
296  * non-DD parallelization in TPI, NM and multi-simulations means that
297  * PAR(cr) and DOMAINDECOMP(cr) are not universally synonymous. In
298  * particular, DOMAINDECOMP(cr) == true indicates that there is more
299  * than one domain, not just that the dd algorithm is active. */
300 #define DOMAINDECOMP(cr)   (((cr)->dd != NULL) && PAR(cr))
301 #define DDMASTER(dd)       ((dd)->rank == (dd)->masterrank)
302
303 #define MULTISIM(cr)       ((cr)->ms)
304 #define MSRANK(ms, nodeid)  (nodeid)
305 #define MASTERSIM(ms)      ((ms)->sim == 0)
306
307 /* The master of all (the node that prints the remaining run time etc.) */
308 #define MULTIMASTER(cr)    (SIMMASTER(cr) && (!MULTISIM(cr) || MASTERSIM((cr)->ms)))
309
310 #ifdef __cplusplus
311 }
312 #endif
313 #endif