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