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