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