PLUTO
al_exchange.c
Go to the documentation of this file.
1 /* ///////////////////////////////////////////////////////////////////// */
2 /*!
3  \file
4  \brief Fill the ghost boundaries
5 
6  Fill the ghost boundaries
7 
8  \author A. Malagoli (University of Chicago)
9  \author A. Mignone (mignone@to.astro.it)
10 
11  \date Jun 13, 2007
12 */
13 /* ///////////////////////////////////////////////////////////////////// */
14 #include "al_hidden.h" /*I "al_hidden.h" I*/
15 
16 /*
17  The SZ structure stack is defined and maintained
18  in al_szptr_.c
19  Here we include an external reference to it in
20  order to be able to make internal references to it.
21 */
22 extern SZ *sz_stack[AL_MAX_ARRAYS];
23 extern int stack_ptr[AL_MAX_ARRAYS];
24 
25 /* ********************************************************************* */
26 int AL_Exchange(void *vbuf, int sz_ptr)
27 /*!
28  * Fill the ghost boundaries
29  *
30  * \param [in] vbuf pointer to buffer
31  * \param [in] sz_ptr integer pointer to the distributed array descriptor
32  *********************************************************************** */
33 {
34  char *buf;
35  register int nd;
36  int myrank, nproc;
37  int ndim, gp, nleft, nright, tag1, tag2;
38  int sendb, recvb;
39  MPI_Datatype itype;
40  MPI_Comm comm;
41  MPI_Status status;
42  SZ *s;
43 
44  buf = (char *) vbuf;
45 
46  /* DIAGNOSTICS
47  Check that sz_ptr points to an allocated SZ
48  */
49  if( stack_ptr[sz_ptr] == AL_STACK_FREE){
50  printf("AL_Decompose: wrong SZ pointer\n");
51  }
52 
53  s = sz_stack[sz_ptr];
54 
55  myrank = s->rank;
56  nproc = s->size;
57  comm = s->comm;
58  ndim = s->ndim;
59 
60  for(nd=0;nd<ndim;nd++){
61  gp = s->bg[nd];
62  /* If gp=0, do nothing */
63  if( gp > 0 ){
64  nleft = s->left[nd];
65  nright = s->right[nd];
66  itype = s->type_rl[nd];
67  tag1 = s->tag1[nd];
68  sendb = s->sendb1[nd];
69  recvb = s->recvb1[nd];
70 
71  MPI_Sendrecv(&buf[sendb], 1, itype, nleft, tag1,
72  &buf[recvb], 1, itype, nright,tag1,
73  comm, &status);
74 
75  nleft = s->left[nd];
76  nright = s->right[nd];
77  itype = s->type_lr[nd];
78  tag2 = s->tag2[nd];
79 
80  sendb = s->sendb2[nd];
81  recvb = s->recvb2[nd];
82 
83  MPI_Sendrecv(&buf[sendb], 1, itype, nright, tag2,
84  &buf[recvb], 1, itype, nleft,tag2,
85  comm, &status);
86  }
87  }
88 
89  /* DIAGNOSTICS */
90 #ifdef DEBUG
91  if(myrank==0) printf("AL_Exchange: filled ghost regions\n");
92 #endif
93 
94  return (int) AL_SUCCESS;
95 }
96 
97 /* ********************************************************************* */
98 int AL_Exchange_periods (void *vbuf, int *periods, int sz_ptr)
99 /*!
100  * Same as AL_Exchange, but exchanges periodic
101  * boundaries at physical domain in the dim direction
102  * only if periods[dim] = 1.
103  * If a dimension is not periodic and periods[dim] = 1
104  * nothing changes.
105  *
106  * \param [in] vbuf pointer to buffer
107  * \param [in] periods
108  * \param [in] sz_ptr integer pointer to the distributed array descriptor
109  *********************************************************************** */
110 {
111  char *buf;
112  register int nd;
113  int myrank, nproc;
114  int ndim, gp, nleft, nright, tag1, tag2;
115  int sendb, recvb;
116  MPI_Datatype itype;
117  MPI_Comm comm;
118  MPI_Status status;
119  SZ *s;
120  int is_beg[3], is_end[3];
121 
122  buf = (char *) vbuf;
123 
124  /* -- DIAGNOSTICS
125  Check that sz_ptr points to an allocated SZ
126  -- */
127  if( stack_ptr[sz_ptr] == AL_STACK_FREE){
128  printf("AL_Decompose: wrong SZ pointer\n");
129  }
130 
131  s = sz_stack[sz_ptr];
132 
133  myrank = s->rank;
134  nproc = s->size;
135  comm = s->comm;
136  ndim = s->ndim;
137 
138  AL_Is_boundary (sz_ptr, is_beg, is_end);
139 
140  for(nd=0;nd<ndim;nd++){
141  gp = s->bg[nd];
142 
143  /* If gp=0, do nothing */
144 
145  if( gp > 0 ){
146  nleft = s->left[nd];
147  nright = s->right[nd];
148 
149  if (is_beg[nd] && periods[nd] == 0) nleft = MPI_PROC_NULL;
150  if (is_end[nd] && periods[nd] == 0) nright = MPI_PROC_NULL;
151 
152  itype = s->type_rl[nd];
153  tag1 = s->tag1[nd];
154 
155  sendb = s->sendb1[nd];
156  recvb = s->recvb1[nd];
157 
158  MPI_Sendrecv(&buf[sendb], 1, itype, nleft, tag1,
159  &buf[recvb], 1, itype, nright,tag1,
160  comm, &status);
161 
162  nleft = s->left[nd];
163  nright = s->right[nd];
164 
165  if (is_beg[nd] && periods[nd] == 0) nleft = MPI_PROC_NULL;
166  if (is_end[nd] && periods[nd] == 0) nright = MPI_PROC_NULL;
167 
168  itype = s->type_lr[nd];
169  tag2 = s->tag2[nd];
170 
171  sendb = s->sendb2[nd];
172  recvb = s->recvb2[nd];
173 
174  MPI_Sendrecv(&buf[sendb], 1, itype, nright, tag2,
175  &buf[recvb], 1, itype, nleft,tag2,
176  comm, &status);
177  }
178  }
179 
180  /* DIAGNOSTICS */
181 #ifdef DEBUG
182  printf("AL_Exchange: filled ghost regions\n");
183 #endif
184 
185  return (int) AL_SUCCESS;
186 }
187 
188 
189 
int recvb2[AL_MAX_DIM]
Definition: al_hidden.h:95
int AL_Is_boundary(int sz_ptr, int *is_gbeg, int *is_gend)
Definition: al_boundary.c:28
int AL_Exchange_periods(void *vbuf, int *periods, int sz_ptr)
Definition: al_exchange.c:98
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
int tag1[AL_MAX_DIM]
Definition: al_hidden.h:96
#define AL_SUCCESS
Definition: al_codes.h:32
int right[AL_MAX_DIM]
Definition: al_hidden.h:86
#define AL_STACK_FREE
Definition: al_codes.h:24
int recvb1[AL_MAX_DIM]
Definition: al_hidden.h:94
int sendb1[AL_MAX_DIM]
Definition: al_hidden.h:92
int rank
Definition: al_hidden.h:44
MPI_Datatype type_rl[AL_MAX_DIM]
Definition: al_hidden.h:106
int AL_Exchange(void *vbuf, int sz_ptr)
Definition: al_exchange.c:26
int bg[AL_MAX_DIM]
Definition: al_hidden.h:72
int left[AL_MAX_DIM]
Definition: al_hidden.h:84
int sendb2[AL_MAX_DIM]
Definition: al_hidden.h:93
#define s
int size
Definition: al_hidden.h:45
int tag2[AL_MAX_DIM]
Definition: al_hidden.h:96
MPI_Comm comm
Definition: al_hidden.h:47
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
Definition: al_hidden.h:38
Internal include file for the ArrayLib.
MPI_Datatype type_lr[AL_MAX_DIM]
Definition: al_hidden.h:104
#define AL_MAX_ARRAYS
Definition: al_codes.h:21
int ndim
Definition: al_hidden.h:43