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