PLUTO
al_proto.h File Reference

ArrayLib function prototypes header file. More...

This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Functions

int AL_Init (int *, char ***)
 
int AL_Finalize ()
 
int AL_Initialized ()
 
int AL_Sz_init (MPI_Comm, int *)
 
int AL_Free (int)
 
int AL_Sz_free (int)
 
int AL_Valid_ptr (int)
 
int AL_Set_comm (MPI_Comm, int)
 
int AL_Set_dimensions (int, int)
 
int AL_Set_type (AL_Datatype, int, int)
 
int AL_Set_global_dim (int *, int)
 
int AL_Set_local_dim (int *, int)
 
int AL_Set_parallel_dim (int *, int)
 
int AL_Set_periodic_dim (int *, int)
 
int AL_Set_staggered_dim (int *, int)
 
int AL_Set_ghosts (int *, int)
 
int AL_Get_size (int, int *)
 
int AL_Get_comm (int, MPI_Comm *)
 
int AL_Get_cart_comm (int, MPI_Comm *)
 
int AL_Get_dimensions (int, int *)
 
int AL_Get_type (int, AL_Datatype *)
 
int AL_Get_buffsize (int, int *)
 
int AL_Get_global_dim (int, int *)
 
int AL_Get_local_dim (int, int *)
 
int AL_Get_local_dim_gp (int, int *)
 
int AL_Get_parallel_dim (int, int *)
 
int AL_Get_periodic_dim (int, int *)
 
int AL_Get_staggered_dim (int, int *)
 
int AL_Get_ghosts (int, int *)
 
int AL_Get_offsets (int, int *)
 
int AL_Get_lbounds (int, int *, int *, int *, int)
 
int AL_Get_gbounds (int, int *, int *, int *, int)
 
int AL_Get_bounds (int, int *, int *, int *, int)
 
int AL_Is_boundary (int, int *, int *)
 
int AL_Get_stride (int, int *)
 
int AL_Decompose (int, int *, int)
 
int AL_Type_create_subarray (int, int *, int *, int *, int, MPI_Datatype, MPI_Datatype *)
 
void * AL_Allocate_array (int)
 
int AL_Exchange (void *, int)
 
int AL_Exchange_dim (char *, int *, int)
 
int AL_Exchange_periods (void *vbuf, int *periods, int sz_ptr)
 
int AL_File_open (char *, int)
 
long long AL_Get_offset (int)
 
int AL_Set_offset (int, long long)
 
int AL_Write_header (void *, int, AL_Datatype, int)
 
int AL_File_close (int)
 
int AL_Write_common (void *, int, AL_Datatype, int)
 
int AL_Read_common (void *, int, AL_Datatype, int)
 
int AL_Write_array (void *, int, int)
 
int AL_Read_array (void *, int, int)
 
int AL_Write_array_begin (void *, int, int *, int *, int)
 
int AL_Write_array_end (void *, int)
 
int AL_Init_stack_ ()
 
int AL_Allocate_sz_ ()
 
int AL_Deallocate_sz_ (int)
 
int AL_Auto_Decomp_ (int, int, int *, int *)
 
int AL_Sort_ (int, int *, int *)
 

Detailed Description

ArrayLib function prototypes header file.

Contains the function prototypes used in the ArrayLib routines.

Author
A. Malagoli (University of Chicago)
A. Mignone (migno.nosp@m.ne@p.nosp@m.h.uni.nosp@m.to.i.nosp@m.t)
G. Muscianisi (g.mus.nosp@m.cian.nosp@m.isi@c.nosp@m.inec.nosp@m.a.it)
Date
Aug 16, 2012

Definition in file al_proto.h.

Function Documentation

void* AL_Allocate_array ( int  sz_ptr)

Allocate the local buffer for a distributed array

Parameters
[in]sz_ptrpointer to SZ array (integer)

Definition at line 25 of file al_alloc.c.

31 {
32  char *a;
33  int myrank, nproc;
34  int size;
35  long int buffsize;
36  MPI_Comm comm;
37  SZ *s;
38 
39 
40 
41  /* DIAGNOSTICS
42  Check that sz_ptr points to an allocated SZ
43  */
44  if( stack_ptr[sz_ptr] == AL_STACK_FREE){
45  printf("AL_Allocate_array: wrong SZ pointer\n");
46  }
47 
48  s = sz_stack[sz_ptr];
49 
50  myrank = s->rank;
51  nproc = s->size;
52  comm = s->comm;
53 
54  buffsize = s->buffsize;
55 
56  MPI_Type_size(s->type, &size);
57 
58  if( !(a = (char *)AL_CALLOC_((int) buffsize,size)) ){
59  printf("[%d] AL_Allocate_array: allocation error\n",myrank);
60  return 0;
61  }
62 
63  /* DIAGNOSTICS */
64 #ifdef DEBUG
65  printf("AL_Allocate_array: allocated %ld bytes\n",buffsize*size);
66 #endif
67 
68  return (void *)a;
69 }
static double a
Definition: init.c:135
#define AL_CALLOC_(nelem, size)
Definition: al_defs.h:34
#define AL_STACK_FREE
Definition: al_codes.h:24
int rank
Definition: al_hidden.h:44
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
#define s
int size
Definition: al_hidden.h:45
MPI_Comm comm
Definition: al_hidden.h:47
long int buffsize
Definition: al_hidden.h:52
MPI_Datatype type
Definition: al_hidden.h:40
Definition: al_hidden.h:38
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
int AL_Allocate_sz_ ( )

Return an integer pointer to an SZ stack entry.

Definition at line 58 of file al_szptr_.c.

62 {
63  register int i;
64  int sz_ptr; /* The returned pointer to the SZ array */
65 
66  sz_ptr = -1;
67 
68  /*
69  First, search the stack for the first
70  available pointer
71  */
72  for( i=0; i<AL_MAX_ARRAYS; i++){
73  if( stack_ptr[i] == AL_STACK_FREE ){
74  sz_ptr = i;
76  break;
77  }
78  }
79 
80  /*
81  If the allocation did not fail, then proceed to
82  allocate a SZ structure
83  */
84  if( sz_ptr != -1 ){
85  if( !(sz_stack[i] = (SZ *)malloc(sizeof(SZ))) ){
86  printf("AL_Allocate_sz_: Failed to allocate SZ\n");
87  }
88  }
89 
90  /*
91  If we were successful, then set sz.compiled to AL_FALSE,
92  indicating the need for intialization.
93  */
94  sz_stack[sz_ptr]->compiled = AL_FALSE;
95 
96  return sz_ptr;
97 
98 }
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
#define AL_STACK_FREE
Definition: al_codes.h:24
#define AL_STACK_USED
Definition: al_codes.h:25
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
int i
Definition: analysis.c:2
#define AL_FALSE
Definition: al_codes.h:29
Definition: al_hidden.h:38
int compiled
Definition: al_hidden.h:39
#define AL_MAX_ARRAYS
Definition: al_codes.h:21

Here is the caller graph for this function:

int AL_Auto_Decomp_ ( int  nproc,
int  npdim,
int *  ldims,
int *  gdims 
)

Find a "maximally cubic" processors distribution.

Parameters
[in]nprocnumber of processors
[in]npdimnumber of parallel dimensions
[out]ldimsprocessor decomposition along directions
[in]gdimsglobal array sizes

Definition at line 19 of file al_decomp_.c.

28 {
29  /*
30  1D case. We do not need to do anything, really
31  */
32  if( npdim == 1 ) {
33  ldims[0] = nproc;
34  return (int) AL_SUCCESS;
35  }
36 
37  /*
38  2D case. Slightly more complicated
39  */
40  if( npdim == 2 ) {
41  if( AL_Decomp_2d_(nproc, npdim, ldims, gdims) == AL_FAILURE) {
42  return (int) AL_FAILURE;
43  } else {
44  return (int) AL_SUCCESS;
45  }
46  }
47 
48  /*
49  3D case. Even more complicated
50  */
51  if( npdim == 3 ) {
52  if( AL_Decomp_3d_(nproc, npdim, ldims, gdims) == AL_FAILURE) {
53  return (int) AL_FAILURE;
54  } else {
55  return (int) AL_SUCCESS;
56  }
57  }
58 
59 
60 return 0;
61 }
#define AL_SUCCESS
Definition: al_codes.h:32
int AL_Decomp_2d_(int nproc, int npdim, int *ldims, int *gdims)
Definition: al_decomp_.c:187
#define AL_FAILURE
Definition: al_codes.h:33
int AL_Decomp_3d_(int nproc, int npdim, int *ldims, int *gdims)
Definition: al_decomp_.c:64

Here is the call graph for this function:

Here is the caller graph for this function:

int AL_Deallocate_sz_ ( int  sz_ptr)

Deallocate an integer pointer to an SZ stack entry.

Parameters
[in]sz_ptrInteger pointer to an entry in the SZ stack

Definition at line 102 of file al_szptr_.c.

108 {
109  register int i;
110  int ndim;
111  SZ *s;
112 
113  if(sz_ptr<0) return AL_SUCCESS;
114  if(stack_ptr[sz_ptr] == AL_STACK_FREE) return AL_SUCCESS;
115 
116  s = sz_stack[sz_ptr];
117  ndim = s->ndim;
118 
119  /*
120  Free the MPI data types
121  */
122  if( s->compiled == AL_TRUE ){
123  for( i=0; i<ndim; i++){
124  MPI_Type_free(&(s->strided[i]));
125  MPI_Type_free(&(s->type_lr[i]));
126  MPI_Type_free(&(s->type_rl[i]));
127  MPI_Comm_free(&(s->oned_comm[i]));
128  }
129 
130  MPI_Type_free(&(s->gsubarr));
131  MPI_Type_free(&(s->lsubarr));
132  MPI_Comm_free(&(s->cart_comm));
133  }
134 
135  if( (s->begs != NULL) ) free(s->begs);
136 
137  /*
138  Begin by dellocating the SZ structure
139  */
140  free(sz_stack[sz_ptr]);
141 
142  /*
143  Now we can reset the integer pointer
144  */
145  stack_ptr[sz_ptr] = AL_STACK_FREE;
146  return (int) AL_SUCCESS;
147 }
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_TRUE
Definition: al_codes.h:28
#define AL_STACK_FREE
Definition: al_codes.h:24
int * begs
Definition: al_hidden.h:63
MPI_Datatype gsubarr
Definition: al_hidden.h:110
MPI_Comm oned_comm[AL_MAX_DIM]
Definition: al_hidden.h:50
MPI_Comm cart_comm
Definition: al_hidden.h:49
MPI_Datatype type_rl[AL_MAX_DIM]
Definition: al_hidden.h:106
MPI_Datatype lsubarr
Definition: al_hidden.h:111
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
MPI_Datatype strided[AL_MAX_DIM]
Definition: al_hidden.h:97
#define s
int i
Definition: analysis.c:2
Definition: al_hidden.h:38
int compiled
Definition: al_hidden.h:39
MPI_Datatype type_lr[AL_MAX_DIM]
Definition: al_hidden.h:104
int ndim
Definition: al_hidden.h:43

Here is the caller graph for this function:

int AL_Decompose ( int  sz_ptr,
int *  procs,
int  mode 
)

Create a distributed array descriptor and compile it

Parameters
[in]sz_ptrinteger pointer to the distributed array descriptor
[in,out]procsarray with the processor decomposition
[in]modeavailable value:
AL_AUTO_DECOMP (internal decomposition) [Only for powers of two dimensions];
AL_MPI_DECOMP (MPI decomposition);
AL_USER_DECOMP (user defined)

Bug fixed on Aug 26, 2012: When the global view of the file for a staggered variable is computed, the gdims[nd] has to be computed and then passed to the function AL_Type_create_subarray.

Since gdims[nd]=s->arrdim[nd], and in s->arrdim[nd] is taking into account that the variable is staggered, we comment the line "gdims[istag]++;" because it is no more needed.

Bug fixed on Aug 26, 2012: The following "if" has been modified:

OLD version: if (s->beg[istag] == s->bg[istag]) { ldims[istag]++; }else{ starts[istag]++; }

NEW version: if (s->beg[istag] == s->bg[istag]) { }else{ starts[istag]++; ldims[istag]–; }

  • "ldims[istag]++;" has been cancelled from the firt part of the if for the same motivation explaned before;
  • "ldims[istag]--;" has been added in the second part of the if because in PLUTO the index of the staggered variables start from -1, while in the ArrayLib they start from 0.

Bugs fixed on Aug 26, 2012: The following if has been modified:

OLD version: if (s->beg[istag] == s->bg[istag]) { ldims[istag]++; starts[istag]–; }

NEW version: if (s->beg[istag] == s->bg[istag]) { }else{ ldims[istag]–; starts[istag]++; }

  • "ldims[istag]++;" has been cancelled from the firt part of the if because when the local subarray for the MPI_Set_write_all for a staggered variable is computed, the ldims[nd] has to be computed and then passed to the function AL_Type_create_subarray. Since ldims[nd]=s->larrdim[nd], and in s->larrdim[nd] is taking into account that the variable is staggered, we remouved the line "ldims[istag]++;" because it is no more needed;
  • "starts[istag]--;" has been cancelled from the first part of the if because in PLUTO the indexes for the staggered variables start locally from -1, while in ArrayLib they start from 0;
  • the "else" has been added to take into account the motivation explaned in the point before.

Definition at line 39 of file al_decompose.c.

49 {
50  int myrank, nproc;
51  SZ *s;
52 
53  register int i, nd, nb;
54 
55  int ndim;
56  int gdims[AL_MAX_DIM], ldims[AL_MAX_DIM], starts[AL_MAX_DIM];
57  int count[AL_MAX_DIM], blocklen[AL_MAX_DIM], periods[AL_MAX_DIM];
58  int reord, nleft, nright;
59  int stagpt, type_size;
60  int stride[AL_MAX_DIM];
61  MPI_Datatype itype,ivector;
62  MPI_Comm comm, cart_comm;
63 
64  int Find_decomp = AL_TRUE; /* By default, we use the internal decomposition */
65 
66  /* DIAGNOSTICS
67  Check that sz_ptr points to an allocated SZ
68  */
69  if( stack_ptr[sz_ptr] == AL_STACK_FREE){
70  printf("AL_Decompose: wrong SZ pointer\n");
71  }
72 
73  s = sz_stack[sz_ptr];
74 
75  myrank = s->rank;
76  nproc = s->size;
77 
78  ndim = s->ndim;
79 
80  /*
81  Find if we need to find our own decomposition
82  */
83 
84  if( mode != AL_USER_DECOMP && mode != AL_AUTO_DECOMP && mode != AL_MPI_DECOMP ) {
85  printf("AL_Decompose: wrong mode. Using default AL_MPI_DECOMP\n");
86  mode = AL_MPI_DECOMP;
87  }
88 
89  if( mode == AL_USER_DECOMP ){ Find_decomp = AL_FALSE; }
90 
91  /*
92  If Find_decomp is TRUE, we proceed to find the
93  decomposition
94  */
95  if( Find_decomp == AL_TRUE ){
96  AL_Find_decomp_(sz_ptr, mode, procs);
97  } else {
98  for(nd=0;nd<ndim;nd++){
99  s->lsize[nd] = procs[nd];
100  }
101  }
102 
103  /*
104  Now we set the MPI datatypes for the ghost points
105  exchange
106  */
107  comm = s->comm;
108 
109  /*
110  First, we create the cartesian communicator
111  */
112  reord = 0;
113 
114  for(nd=0;nd<ndim;nd++){
115  ldims[nd] = s->lsize[nd];
116  periods[nd] = s->isperiodic[nd];
117  }
118 
119 #ifdef DEBUG
120  if(sz_ptr==4){ /* -- print SZ_stagx -- */
121  if(myrank==0) printf("MPI_Cart_create: %d x %d x %d, ndim %d, nproc %d\n",ldims[0],ldims[1],ldims[2],ndim,nproc);
122  }
123 #endif
124 
125  MPI_Cart_create(comm, ndim, ldims, periods, reord, &cart_comm);
126 
127  s->cart_comm = cart_comm;
128 
129  /*
130  Then we generate the cartesian topology
131  of the processors
132  */
133  MPI_Cart_coords(cart_comm, myrank, ndim, ldims);
134 
135  for(nd=0;nd<ndim;nd++){
136  s->lrank[nd] = ldims[nd];
137  MPI_Cart_shift(cart_comm, nd, 1, &nleft, &nright);
138  s->left[nd] = nleft;
139  s->right[nd] = nright;
140  }
141 
142  /*
143  We also generate a set of 1D communicators that
144  are useful when doing global operations along
145  single directions.
146  */
147  {
148  int rem_dims[AL_MAX_DIM];
149 
150  for(nd=0;nd<ndim;nd++){
151  for(i=0;i<ndim;i++){
152  if( i==nd ) {
153  rem_dims[i] = AL_TRUE;
154  } else {
155  rem_dims[i] = AL_FALSE;
156  }
157  }
158  MPI_Cart_sub(cart_comm, rem_dims, &(s->oned_comm[nd]));
159  }
160  }
161 
162  /*
163  Now we create the strided data types for the
164  exchange of the ghost points
165  */
166 
167  /*
168  Get the local array indexes (in global addressing mode)
169  */
170  AL_Global_to_local_(sz_ptr);
171 
172  /*
173  Set the size of the local buffer for this array
174  */
175  s->buffsize = 1;
176  for(nd=0;nd<ndim;nd++){
177  s->buffsize *=
178  ( (s->end[nd]) - (s->beg[nd]) + (s->bg[nd]) + (s->eg[nd]) + 1 );
179  }
180 
181 #ifdef AL_GATHER_INDEXES
182 
183  /*
184  If it has not yet been done, allocate the array of begs
185  and ends, then gather the indexing information from all
186  processors. We check that the array of 'begs' was not
187  previously allocated, in order to avoid problems with
188  duplicated descriptors.
189  */
190 /* if( (s->begs) ) free(s->begs); */
191 
192  if( !(s->begs = (int *)malloc(sizeof(int)*nproc*ndim*2))) {
193  printf("AL_Decompose: could not allocate s->begs\n");
194  }
195  s->ends = s->begs + nproc*ndim;
196 
197  /* Gather the 'beg' index information from all nodes */
198  MPI_Allgather(s->beg, ndim, MPI_INT, s->begs, ndim, MPI_INT, comm);
199 
200  /* Gather the 'end' index information from all nodes */
201  MPI_Allgather(s->end, ndim, MPI_INT, s->ends, ndim, MPI_INT, comm);
202 
203 #endif /* End #ifdef AL_GATHER_INDEXES */
204 
205  /* -- Set the type size of the array -- */
206 
207  MPI_Type_size(s->type, &(s->type_size));
208 
209  /*--------------------- Begin creation of strided data types -------*/
210 
211  /*
212  Set the count, blocklen, and stride elements for
213  the MPI_Vector calls.
214  This is VERY, VERY cryptic, but it seems to work ..
215  */
216 
217  /* Right->Left exchange */
218  for(nd=0;nd<ndim;nd++){
219  count[nd] = 1;
220  blocklen[nd] = s->bg[nd];
221  stride[nd] = s->larrdim_gp[nd];
222  for(nb=0;nb<nd;nb++){
223  blocklen[nd] *= s->larrdim_gp[nb];
224  stride[nd] *= s->larrdim_gp[nb];
225  }
226  for(nb=nd+1;nb<ndim;nb++){
227  count[nd] *= s->larrdim_gp[nb];
228  }
229  }
230 
231  itype = s->type;
232  for(nd=0;nd<ndim;nd++){
233  if( count[nd] != 1){
234  MPI_Type_vector( count[nd], blocklen[nd], stride[nd], itype, &ivector);
235  } else {
236  MPI_Type_contiguous( blocklen[nd], itype, &ivector);
237  }
238  MPI_Type_commit(&ivector);
239  s->strided[nd] = ivector;
240  s->type_rl[nd] = ivector;
241  }
242 
243  /* Left->Right exchange */
244  for(nd=0;nd<ndim;nd++){
245  count[nd] = 1;
246  /* If the array is staggered, we give
247  unique ownership of the overlapping
248  point to the LEFT processor */
249  if( s->isstaggered[nd] == AL_TRUE ){
250  blocklen[nd] = s->bg[nd]+1;
251  } else {
252  blocklen[nd] = s->bg[nd];
253  }
254 
255  stride[nd] = s->larrdim_gp[nd];
256  for(nb=0;nb<nd;nb++){
257  blocklen[nd] *= s->larrdim_gp[nb];
258  stride[nd] *= s->larrdim_gp[nb];
259  }
260  for(nb=nd+1;nb<ndim;nb++){
261  count[nd] *= s->larrdim_gp[nb];
262  }
263  }
264 
265  itype = s->type;
266  for(nd=0;nd<ndim;nd++){
267  if( count[nd] != 1){
268  MPI_Type_vector( count[nd], blocklen[nd], stride[nd], itype, &ivector);
269  } else {
270  MPI_Type_contiguous( blocklen[nd], itype, &ivector);
271  }
272  MPI_Type_commit(&ivector);
273  s->type_lr[nd] = ivector;
274  }
275 
276  /*--------------------- End creation of strided data types -------*/
277 
278  /*
279  Create the buffer pointers to the array
280  for the SendRecv operations
281  */
282  MPI_Type_size( s->type, &type_size);
283 
284  for(nd=0;nd<ndim;nd++){
285  s->tag1[nd] = nd*100;
286  s->tag2[nd] = nd*100+1;
287  /* This is the correction for a staggered mesh */
288  if( s->isstaggered[nd] == AL_TRUE ){
289  stagpt = 1;
290  } else {
291  stagpt = 0;
292  }
293 
294  s->sendb1[nd] = s->lbeg[nd] + stagpt;
295  s->recvb1[nd] = s->lend[nd]+1;
296  s->sendb2[nd] = s->lend[nd]-s->bg[nd]+1 - stagpt;
297  s->recvb2[nd] = 0;
298 
299  for(nb=0;nb<nd;nb++){
300  s->sendb1[nd] *= s->larrdim_gp[nb];
301  s->recvb1[nd] *= s->larrdim_gp[nb];
302  s->sendb2[nd] *= s->larrdim_gp[nb];
303  s->recvb2[nd] *= s->larrdim_gp[nb];
304  }
305 
306  /*
307  The buffer will be considered as char, so we
308  have to correct for the size of the data type
309  */
310  s->sendb1[nd] *= type_size;
311  s->recvb1[nd] *= type_size;
312  s->sendb2[nd] *= type_size;
313  s->recvb2[nd] *= type_size;
314 
315  }
316 
317  /*
318  Create the subarrays for the MPI-IO operations
319  */
320  {
321  int istag;
322 
323  MPI_Datatype igsubarr, ilsubarr;
324 
325  MPI_Datatype igsubarr_stag[AL_MAX_DIM];
326  MPI_Datatype ilsubarr_stag[AL_MAX_DIM];
327 
328  /* -----------------------------------------------------
329  Create the global subarray for MPI_Set_file_view
330  ----------------------------------------------------- */
331 
332  for(nd=0;nd<ndim;nd++){
333  gdims[nd] = s->arrdim[nd];
334  ldims[nd] = s->larrdim[nd];
335  starts[nd] = s->beg[nd] - s->bg[nd];
336  }
337 #ifdef DEBUG
338  if(sz_ptr==4){ /* print SZ_stagx */
339 printf("%d, gsubarr: gdims[0:2] %d,%d,%d, ldims[0:2] %d,%d,%d, starts[0:2] %d,%d,%d\n", s->rank, gdims[0], gdims[1], gdims[2], ldims[0], ldims[1], ldims[2], starts[0], starts[1], starts[2]);
340 }
341 #endif
342 
343  AL_Type_create_subarray(ndim, gdims, ldims, starts,
344  AL_ORDER_FORTRAN, s->type, &igsubarr);
345  MPI_Type_commit(&igsubarr);
346  s->gsubarr = igsubarr;
347 
348  /* -----------------------------------------------------
349  Create the global staggered subarrays
350  for MPI_Set_file_view
351  ----------------------------------------------------- */
352 
353  for (istag = 0; istag < ndim; istag++){
354  for(nd = 0; nd < ndim; nd++){
355  gdims[nd] = s->arrdim[nd];
356  ldims[nd] = s->larrdim[nd];
357  starts[nd] = s->beg[nd] - s->bg[nd];
358  }
359  /* --------------------------------------------------------------- */
360  /*! <b> Bug fixed on Aug 26, 2012:</b>
361  When the global view of the file for a staggered
362  variable is computed, the gdims[nd] has to be computed
363  and then passed to the function AL_Type_create_subarray.
364 
365  Since gdims[nd]=s->arrdim[nd], and in s->arrdim[nd] is
366  taking into account that the variable is staggered, we comment
367  the line "gdims[istag]++;" because it is no more needed. */
368  /* --------------------------------------------------------------- */
369  /* gdims[istag]++; */
370 
371 
372 
373  /* --------------------------------------------------------------- */
374  /*! <b> Bug fixed on Aug 26, 2012: </b>
375  The following "if" has been modified:
376 
377  OLD version:
378  if (s->beg[istag] == s->bg[istag]) {
379  ldims[istag]++;
380  }else{
381  starts[istag]++;
382  }
383 
384  NEW version:
385  if (s->beg[istag] == s->bg[istag]) {
386  }else{
387  starts[istag]++;
388  ldims[istag]--;
389  }
390 
391  \li "ldims[istag]++;" has been cancelled from the firt part
392  of the if for the same motivation explaned before;
393  \li "ldims[istag]--;" has been added in the second part of the
394  if because in PLUTO the index of the staggered variables
395  start from -1, while in the ArrayLib they start from 0. */
396  /* --------------------------------------------------------------- */
397 
398  if (s->beg[istag] == s->bg[istag]) {
399  }else{
400  starts[istag]++;
401  ldims[istag]--;
402  }
403 
404 #ifdef DEBUG
405  if(sz_ptr==4){ /* print SZ_stagx */
406  printf("%d, gsubarr_stag[%d]: gdims[0:2] %d,%d,%d, ldims[0:2] %d,%d,%d, starts[0:2] %d,%d,%d \n", s->rank, istag, gdims[0], gdims[1], gdims[2], ldims[0], ldims[1], ldims[2], starts[0], starts[1], starts[2]);
407 }
408 #endif
409 
410  AL_Type_create_subarray(ndim, gdims, ldims, starts,
411  AL_ORDER_FORTRAN, s->type, igsubarr_stag + istag);
412  MPI_Type_commit(igsubarr_stag + istag);
413  s->gsubarr_stag[istag] = igsubarr_stag[istag];
414  }
415 
416  /* -----------------------------------------------------
417  Create the local subarray for the MPI_Set_write_all
418  ----------------------------------------------------- */
419 
420  for(nd=0;nd<ndim;nd++){
421  gdims[nd] = s->larrdim_gp[nd];
422  ldims[nd] = s->larrdim[nd];
423  starts[nd] = s->lbeg[nd];
424  }
425 #ifdef DEBUG
426  if(sz_ptr==4){ /* print SZ_stagx */
427 printf("%d, lsubarr: gdims[0:2] %d,%d,%d, ldims[0:2] %d,%d,%d, starts[0:2] %d,%d,%d \n", s->rank, gdims[0], gdims[1], gdims[2], ldims[0], ldims[1], ldims[2], starts[0], starts[1], starts[2]);
428 }
429 #endif
430  AL_Type_create_subarray(ndim, gdims, ldims, starts,
431  AL_ORDER_FORTRAN, s->type, &ilsubarr);
432  MPI_Type_commit(&ilsubarr);
433  s->lsubarr = ilsubarr;
434 
435 
436  /* -----------------------------------------------------
437  Create the local staggered subarrays
438  for MPI_Set_write_all
439  ----------------------------------------------------- */
440 
441  for (istag = 0; istag < ndim; istag++){
442  for(nd = 0; nd < ndim; nd++){
443  gdims[nd] = s->larrdim_gp[nd];
444  ldims[nd] = s->larrdim[nd];
445  starts[nd] = s->lbeg[nd];
446  }
447  /* --------------------------------------------------------------- */
448  /*! <b> Bugs fixed on Aug 26, 2012: </b>
449  The following if has been modified:
450 
451  OLD version:
452  if (s->beg[istag] == s->bg[istag]) {
453  ldims[istag]++;
454  starts[istag]--;
455  }
456 
457  NEW version:
458  if (s->beg[istag] == s->bg[istag]) {
459  }else{
460  ldims[istag]--;
461  starts[istag]++;
462  }
463 
464  \li "ldims[istag]++;" has been cancelled from the firt part of
465  the if because when the local subarray for the MPI_Set_write_all
466  for a staggered variable is computed, the ldims[nd] has to be
467  computed and then passed to the function AL_Type_create_subarray.
468  Since ldims[nd]=s->larrdim[nd], and in s->larrdim[nd] is
469  taking into account that the variable is staggered, we remouved
470  the line "ldims[istag]++;" because it is no more needed;
471  \li "starts[istag]--;" has been cancelled from the first part of
472  the if because in PLUTO the indexes for the staggered variables
473  start locally from -1, while in ArrayLib they start from 0;
474  \li the "else" has been added to take into account the motivation
475  explaned in the point before. */
476  /* --------------------------------------------------------------- */
477 
478  if (s->beg[istag] == s->bg[istag]) {
479  }else{
480  ldims[istag]--;
481  starts[istag]++;
482  }
483 
484 
485 #ifdef DEBUG
486  if(sz_ptr==4){ /* print SZ_stagx */
487 printf("%d, lsubarr_stag[%d]: gdims[0:2] %d,%d,%d, ldims[0:2] %d,%d,%d, starts[0:2] %d,%d,%d \n", s->rank, istag, gdims[0], gdims[1], gdims[2], ldims[0], ldims[1], ldims[2], starts[0], starts[1], starts[2]);
488 }
489 #endif
490  AL_Type_create_subarray(ndim, gdims, ldims, starts,
491  AL_ORDER_FORTRAN, s->type, ilsubarr_stag + istag);
492  MPI_Type_commit(ilsubarr_stag + istag);
493  s->lsubarr_stag[istag] = ilsubarr_stag[istag];
494  }
495 
496  }
497 
498  /*
499  Now array is officially compiled
500  */
501  s->compiled = AL_TRUE;
502 
503  /* DIAGNOSTICS */
504 #ifdef DEBUG
505  if(sz_ptr==4){ /* print SZ_stagx */
506  printf("AL_Decompose: Decomposition successful\n");
507  for(nd=0;nd<ndim;nd++){
508  printf("AL_Decompose: %d %d\n", s->beg[nd], s->end[nd]);
509  }
510 }
511 #endif
512 
513  return (int) AL_SUCCESS;
514 }
int recvb2[AL_MAX_DIM]
Definition: al_hidden.h:95
int isperiodic[AL_MAX_DIM]
Definition: al_hidden.h:80
#define AL_MAX_DIM
Definition: al_codes.h:18
#define AL_MPI_DECOMP
Definition: al_codes.h:41
int end[AL_MAX_DIM]
Definition: al_hidden.h:60
int lrank[AL_MAX_DIM]
Definition: al_hidden.h:88
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_TRUE
Definition: al_codes.h:28
#define AL_STACK_FREE
Definition: al_codes.h:24
int AL_Find_decomp_(int sz_ptr, int mode, int *procs)
Definition: al_decompose.c:519
int type_size
Definition: al_hidden.h:42
int * begs
Definition: al_hidden.h:63
MPI_Datatype gsubarr
Definition: al_hidden.h:110
MPI_Comm oned_comm[AL_MAX_DIM]
Definition: al_hidden.h:50
int recvb1[AL_MAX_DIM]
Definition: al_hidden.h:94
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
#define AL_USER_DECOMP
Definition: al_codes.h:42
MPI_Comm cart_comm
Definition: al_hidden.h:49
int arrdim[AL_MAX_DIM]
Definition: al_hidden.h:53
int larrdim_gp[AL_MAX_DIM]
Definition: al_hidden.h:56
int sendb1[AL_MAX_DIM]
Definition: al_hidden.h:92
int larrdim[AL_MAX_DIM]
Definition: al_hidden.h:54
int AL_Type_create_subarray(int, int *, int *, int *, int, MPI_Datatype, MPI_Datatype *)
Definition: al_subarray_.c:18
int rank
Definition: al_hidden.h:44
#define AL_AUTO_DECOMP
Definition: al_codes.h:40
int * ends
Definition: al_hidden.h:63
MPI_Datatype type_rl[AL_MAX_DIM]
Definition: al_hidden.h:106
int bg[AL_MAX_DIM]
Definition: al_hidden.h:72
MPI_Datatype lsubarr
Definition: al_hidden.h:111
int left[AL_MAX_DIM]
Definition: al_hidden.h:84
int beg[AL_MAX_DIM]
Definition: al_hidden.h:58
MPI_Datatype lsubarr_stag[AL_MAX_DIM]
Definition: al_hidden.h:114
#define AL_ORDER_FORTRAN
Definition: al_codes.h:46
int sendb2[AL_MAX_DIM]
Definition: al_hidden.h:93
int isstaggered[AL_MAX_DIM]
Definition: al_hidden.h:82
int eg[AL_MAX_DIM]
Definition: al_hidden.h:74
int lsize[AL_MAX_DIM]
Definition: al_hidden.h:90
MPI_Datatype strided[AL_MAX_DIM]
Definition: al_hidden.h:97
#define s
int size
Definition: al_hidden.h:45
int tag2[AL_MAX_DIM]
Definition: al_hidden.h:96
int i
Definition: analysis.c:2
MPI_Comm comm
Definition: al_hidden.h:47
long int buffsize
Definition: al_hidden.h:52
MPI_Datatype type
Definition: al_hidden.h:40
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
#define AL_FALSE
Definition: al_codes.h:29
int AL_Global_to_local_(int sz_ptr)
Definition: al_decompose.c:682
Definition: al_hidden.h:38
MPI_Datatype gsubarr_stag[AL_MAX_DIM]
Definition: al_hidden.h:113
int compiled
Definition: al_hidden.h:39
int lbeg[AL_MAX_DIM]
Definition: al_hidden.h:68
MPI_Datatype type_lr[AL_MAX_DIM]
Definition: al_hidden.h:104
int lend[AL_MAX_DIM]
Definition: al_hidden.h:70
int ndim
Definition: al_hidden.h:43

Here is the call graph for this function:

Here is the caller graph for this function:

int AL_Exchange ( void *  vbuf,
int  sz_ptr 
)

Fill the ghost boundaries

Parameters
[in]vbufpointer to buffer
[in]sz_ptrinteger pointer to the distributed array descriptor

Definition at line 26 of file al_exchange.c.

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 }
int recvb2[AL_MAX_DIM]
Definition: al_hidden.h:95
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 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
MPI_Datatype type_lr[AL_MAX_DIM]
Definition: al_hidden.h:104
int ndim
Definition: al_hidden.h:43
int AL_Exchange_dim ( char *  buf,
int *  dims,
int  sz_ptr 
)

Fill the ghost boundaries along selected dimensions

Parameters
[in]bufpointer to buffer
[in]dimsif dims[i]=0, do not perform the exchange in this dimension (array if int)
[in]sz_ptrinteger pointer to the distributed array descriptor

Definition at line 24 of file al_exchange_dim.c.

33 {
34  register int nd;
35  int myrank, nproc;
36  int ndim, gp, nleft, nright, tag1, tag2;
37  int sendb, recvb;
38  MPI_Datatype itype;
39  MPI_Comm comm;
40  MPI_Status status;
41  SZ *s;
42 
43  /* DIAGNOSTICS
44  Check that sz_ptr points to an allocated SZ
45  */
46  if( stack_ptr[sz_ptr] == AL_STACK_FREE){
47  printf("AL_Decompose: wrong SZ pointer\n");
48  }
49 
50  s = sz_stack[sz_ptr];
51 
52  myrank = s->rank;
53  nproc = s->size;
54  comm = s->comm;
55  ndim = s->ndim;
56 
57  for(nd=0;nd<ndim;nd++){
58  gp = s->bg[nd];
59 
60  /* If gp=0, do nothing */
61  if( gp > 0 && dims[nd] != 0 ){
62  nleft = s->left[nd];
63  nright = s->right[nd];
64  itype = s->type_rl[nd];
65  tag1 = s->tag1[nd];
66 
67  sendb = s->sendb1[nd];
68  recvb = s->recvb1[nd];
69 
70  MPI_Sendrecv(&buf[sendb], 1, itype, nleft, tag1,
71  &buf[recvb], 1, itype, nright,tag1,
72  comm, &status);
73 
74  nleft = s->left[nd];
75  nright = s->right[nd];
76  itype = s->type_lr[nd];
77  tag2 = s->tag2[nd];
78 
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 }
int recvb2[AL_MAX_DIM]
Definition: al_hidden.h:95
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
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
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 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
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
#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
Definition: al_hidden.h:38
MPI_Datatype type_lr[AL_MAX_DIM]
Definition: al_hidden.h:104
int ndim
Definition: al_hidden.h:43

Here is the caller graph for this function:

int AL_Exchange_periods ( void *  vbuf,
int *  periods,
int  sz_ptr 
)

Same as AL_Exchange, but exchanges periodic boundaries at physical domain in the dim direction only if periods[dim] = 1. If a dimension is not periodic and periods[dim] = 1 nothing changes.

Parameters
[in]vbufpointer to buffer
[in]periods
[in]sz_ptrinteger pointer to the distributed array descriptor

Definition at line 98 of file al_exchange.c.

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 }
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
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 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
MPI_Datatype type_lr[AL_MAX_DIM]
Definition: al_hidden.h:104
int ndim
Definition: al_hidden.h:43

Here is the call graph for this function:

int AL_File_close ( int  sz_ptr)

Close a file associate with and array

Parameters
[in]sz_ptrinteger pointer to array descriptor

Definition at line 83 of file al_io.c.

89 {
90  int myrank;
91  int errcode;
92  SZ *s;
93 
94  s = sz_stack[sz_ptr];
95  myrank = s->rank;
96 
97  errcode = MPI_File_close(&(s->ifp));
98 
99  /* DIAGNOSTICS */
100 #ifdef DEBUG
101  int myid, len;
102  char es[128];
103  if( errcode ){
104  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
105  MPI_Error_string(errcode, es, &len);
106  printf("Errcode from MPI_File_close: %d | %s\n", errcode, es);
107  }
108  printf("myid %d, AL_File_close: Closed file, errcode %d \n",myid, errcode);
109 #endif
110 
111  return (int) AL_SUCCESS;
112 }
#define AL_SUCCESS
Definition: al_codes.h:32
int rank
Definition: al_hidden.h:44
MPI_File ifp
Definition: al_hidden.h:117
#define s
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
Definition: al_hidden.h:38

Here is the caller graph for this function:

int AL_File_open ( char *  filename,
int  sz_ptr 
)

Open a file associated with a distributed array

Parameters
[in]filenamename of the file
[in]sz_ptrinteger pointer to the distributed array descriptor
Returns
ifp pointer to integer file pointer

Definition at line 27 of file al_io.c.

36 {
37  int myrank, nproc;
38  int errcode;
39  MPI_File ifp;
40  MPI_Comm comm;
41  SZ *s;
42 
43  /* DIAGNOSTICS
44  Check that sz_ptr points to an allocated SZ
45  */
46  if( stack_ptr[sz_ptr] == AL_STACK_FREE){
47  printf("AL_Decompose: wrong SZ pointer\n");
48  }
49 
50  s = sz_stack[sz_ptr];
51 
52  myrank = s->rank;
53  nproc = s->size;
54 
55  comm = s->comm;
56 
57  MPI_Barrier(comm);
58 
59  errcode = MPI_File_open(comm, filename,
60  MPI_MODE_CREATE | MPI_MODE_RDWR | MPI_MODE_UNIQUE_OPEN,
61  MPI_INFO_NULL, &ifp);
62 
63  s->io_offset = 0;
64  s->ifp = ifp;
65 
66  /* DIAGNOSTICS */
67 #ifdef DEBUG
68  int myid, len;
69  char es[128];
70  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
71  if( errcode ){
72  MPI_Error_string(errcode, es, &len);
73  printf("Errcode from MPI_File_open: %d | %s\n", errcode, es);
74  }
75  printf("myid %d, AL_File_open: Opened file %s\n",myid,filename);
76 #endif
77 
78  return (int) AL_SUCCESS;
79 }
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_STACK_FREE
Definition: al_codes.h:24
int rank
Definition: al_hidden.h:44
MPI_File ifp
Definition: al_hidden.h:117
#define s
int size
Definition: al_hidden.h:45
MPI_Comm comm
Definition: al_hidden.h:47
MPI_Offset io_offset
Definition: al_hidden.h:116
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
Definition: al_hidden.h:38
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26

Here is the caller graph for this function:

int AL_Finalize ( )

Finalize the AL Tool. It contains a call to MPI_Finalize()

Definition at line 15 of file al_finalize.c.

20 {
21  int myrank, nproc, errcode;
22 
23  MPI_Comm_rank(MPI_COMM_WORLD, &myrank);
24  MPI_Comm_size(MPI_COMM_WORLD, &nproc);
25 
26  /* Synchronize just in case */
27  MPI_Barrier(MPI_COMM_WORLD);
28 
29  errcode = MPI_Finalize();
30 
31 #ifdef DEBUG
32  printf("AL_Finalize: Called MPI_Finalize: %d\n",errcode);
33 #endif
34 
35  return errcode;
36 }

Here is the caller graph for this function:

int AL_Free ( int  )
int AL_Get_bounds ( int  isz,
int *  beg,
int *  end,
int *  gp,
int  style 
)

Get the global indexes for the local portion of a distributed array

Parameters
[in]iszInteger pointer to the input array descriptor
[out]begArray of ndim integers containing the start points
[out]endArray of ndim integers containing the end points
[out]gpArray of ndim integers containing the ghost points
[in]styleIndex style: AL_C_INDEXES or AL_FORTRAN_INDEXES

Definition at line 707 of file al_sz_get.c.

717 {
718  register int i;
719  int ndim;
720  SZ *s;
721 
722  /*
723  Check that isz points to an allocated SZ
724  */
725  if( stack_ptr[isz] == AL_STACK_FREE ){
726  printf("AL_Get_gbounds: wrong SZ pointer\n");
727  return (int) AL_FAILURE;
728  }
729 
730  /*
731  Get the SZ structure isz is pointing at
732  */
733  s = sz_stack[isz];
734  ndim = s->ndim;
735  for(i=0;i<ndim;i++){
736  beg[i] = s->beg[i];
737  end[i] = s->end[i];
738  gp[i] = s->bg[i];
739  }
740 
741  return (int) AL_SUCCESS;
742 }
int end[AL_MAX_DIM]
Definition: al_hidden.h:60
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_STACK_FREE
Definition: al_codes.h:24
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
int bg[AL_MAX_DIM]
Definition: al_hidden.h:72
int beg[AL_MAX_DIM]
Definition: al_hidden.h:58
#define AL_FAILURE
Definition: al_codes.h:33
#define s
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
int i
Definition: analysis.c:2
Definition: al_hidden.h:38
int ndim
Definition: al_hidden.h:43

Here is the caller graph for this function:

int AL_Get_buffsize ( int  sz_ptr,
int *  buffsize 
)

Get the local buffer size for a distributed array, including the number of elements of the ghsot regions.

Parameters
[in]sz_ptrInteger pointer to the input array descriptor
[out]buffsizeBuffer size [in number of elements]

Definition at line 242 of file al_sz_get.c.

250 {
251  int isz;
252  SZ *s;
253 
254  isz = sz_ptr;
255 
256  /*
257  Check that isz points to an allocated SZ
258  */
259  if( stack_ptr[isz] == AL_STACK_FREE ){
260  printf("AL_Get_buffsize: wrong SZ pointer\n");
261  return (int) AL_FAILURE;
262  }
263 
264  /*
265  Get the SZ structure isz is pointing at
266  */
267  s = sz_stack[isz];
268 
269  *buffsize = s->buffsize;
270 
271  return (int) AL_SUCCESS;
272 }
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_STACK_FREE
Definition: al_codes.h:24
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
#define AL_FAILURE
Definition: al_codes.h:33
#define s
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
long int buffsize
Definition: al_hidden.h:52
Definition: al_hidden.h:38
int AL_Get_cart_comm ( int  isz,
MPI_Comm *  cart_comm 
)

Get the cartesian communicator for a distributed array

Parameters
[in]iszInteger pointer to the input array descriptor
[out]cart_commPointer to cartesian communicator

Definition at line 117 of file al_sz_get.c.

124 {
125  SZ *s;
126 
127  /*
128  Check that isz points to an allocated SZ
129  */
130  if( stack_ptr[isz] == AL_STACK_FREE ){
131  printf("AL_Get_comm: wrong SZ pointer\n");
132  return (int) AL_FAILURE;
133  }
134 
135  /*
136  Get the SZ structure isz is pointing at
137  */
138  s = sz_stack[isz];
139 
140  *cart_comm = s->cart_comm;
141 
142  return (int) AL_SUCCESS;
143 }
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_STACK_FREE
Definition: al_codes.h:24
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
MPI_Comm cart_comm
Definition: al_hidden.h:49
#define AL_FAILURE
Definition: al_codes.h:33
#define s
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
Definition: al_hidden.h:38

Here is the caller graph for this function:

int AL_Get_comm ( int  isz,
MPI_Comm *  comm 
)

Get the communicator for a distributed array

Parameters
[in]iszInteger pointer to the input array descriptor
[out]commPointer to communicator

Definition at line 87 of file al_sz_get.c.

94 {
95  SZ *s;
96 
97  /*
98  Check that isz points to an allocated SZ
99  */
100  if( stack_ptr[isz] == AL_STACK_FREE ){
101  printf("AL_Get_comm: wrong SZ pointer\n");
102  return (int) AL_FAILURE;
103  }
104 
105  /*
106  Get the SZ structure isz is pointing at
107  */
108  s = sz_stack[isz];
109 
110  *comm = s->comm;
111 
112  return (int) AL_SUCCESS;
113 }
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_STACK_FREE
Definition: al_codes.h:24
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
#define AL_FAILURE
Definition: al_codes.h:33
#define s
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
MPI_Comm comm
Definition: al_hidden.h:47
Definition: al_hidden.h:38
int AL_Get_dimensions ( int  isz,
int *  ndim 
)

Get the dimensions of the distributed array

Parameters
[in]iszInteger pointer to the input array descriptor
[out]ndimPointer to integer number of dimensions

Definition at line 148 of file al_sz_get.c.

155 {
156  SZ *s;
157 
158  /*
159  Check that isz points to an allocated SZ
160  */
161  if( stack_ptr[isz] == AL_STACK_FREE ){
162  printf("AL_Get_dimensions: wrong SZ pointer\n");
163  return (int) AL_FAILURE;
164  }
165 
166  /*
167  Get the SZ structure isz is pointing at
168  */
169  s = sz_stack[isz];
170 
171  *ndim = s->ndim;
172 
173  return (int) AL_SUCCESS;
174 }
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_STACK_FREE
Definition: al_codes.h:24
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
#define AL_FAILURE
Definition: al_codes.h:33
#define s
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
Definition: al_hidden.h:38
int ndim
Definition: al_hidden.h:43
int AL_Get_gbounds ( int  isz,
int *  gbeg,
int *  gend,
int *  gp,
int  style 
)

Get the global bounds of the global array

Parameters
[in]iszInteger pointer to the input array descriptor
[out]gbegArray of ndim integers containing the start points
[out]gendArray of ndim integers containing the end points
[out]gpArray of ndim integers containing the ghost points
[in]styleIndex style: AL_C_INDEXES or AL_FORTRAN_INDEXES

Definition at line 746 of file al_sz_get.c.

756 {
757  register int i;
758  int ndim;
759  SZ *s;
760 
761  /*
762  Check that isz points to an allocated SZ
763  */
764  if( stack_ptr[isz] == AL_STACK_FREE ){
765  printf("AL_Get_gbounds: wrong SZ pointer\n");
766  return (int) AL_FAILURE;
767  }
768 
769  /*
770  Get the SZ structure isz is pointing at
771  */
772  s = sz_stack[isz];
773  ndim = s->ndim;
774  for(i=0;i<ndim;i++){
775  gbeg[i] = s->bg[i];
776  gend[i] = s->arrdim[i]+s->bg[i]-1;
777  gp[i] = s->bg[i];
778  }
779 
780  return (int) AL_SUCCESS;
781 }
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_STACK_FREE
Definition: al_codes.h:24
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
int arrdim[AL_MAX_DIM]
Definition: al_hidden.h:53
int bg[AL_MAX_DIM]
Definition: al_hidden.h:72
#define AL_FAILURE
Definition: al_codes.h:33
#define s
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
int i
Definition: analysis.c:2
Definition: al_hidden.h:38
int ndim
Definition: al_hidden.h:43

Here is the caller graph for this function:

int AL_Get_ghosts ( int  isz,
int *  ghosts 
)

Get the ghost points of the distributed array

Parameters
[in]iszInteger pointer to the input array descriptor
[out]ghostArray of integers with size of ghost points in each dimension

Definition at line 633 of file al_sz_get.c.

641 {
642  register int i;
643  int ndim;
644  SZ *s;
645 
646  /*
647  Check that isz points to an allocated SZ
648  */
649  if( stack_ptr[isz] == AL_STACK_FREE ){
650  printf("AL_Get_ghosts: wrong SZ pointer\n");
651  return (int) AL_FAILURE;
652  }
653 
654  /*
655  Get the SZ structure isz is pointing at
656  */
657  s = sz_stack[isz];
658  ndim = s->ndim;
659  for(i=0;i<ndim;i++){
660  ghosts[i] = s->bg[i];
661  }
662 
663  return (int) AL_SUCCESS;
664 }
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_STACK_FREE
Definition: al_codes.h:24
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
int bg[AL_MAX_DIM]
Definition: al_hidden.h:72
#define AL_FAILURE
Definition: al_codes.h:33
#define s
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
int i
Definition: analysis.c:2
Definition: al_hidden.h:38
int ndim
Definition: al_hidden.h:43
int AL_Get_global_dim ( int  isz,
int *  gdims 
)

Get the global dimensions of the distributed array

Parameters
[in]iszInteger pointer to the input array descriptor
[out]gdimsArray of integers with global dimensions

Definition at line 276 of file al_sz_get.c.

283 {
284  register int i;
285  int ndim;
286  SZ *s;
287 
288  /*
289  Check that isz points to an allocated SZ
290  */
291  if( stack_ptr[isz] == AL_STACK_FREE ){
292  printf("AL_Get_global_dim: wrong SZ pointer\n");
293  return (int) AL_FAILURE;
294  }
295 
296  /*
297  Get the SZ structure isz is pointing at
298  */
299  s = sz_stack[isz];
300 
301  ndim = s->ndim;
302  for(i=0;i<ndim;i++){ gdims[i] = s->arrdim[i];}
303 
304  return (int) AL_SUCCESS;
305 }
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_STACK_FREE
Definition: al_codes.h:24
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
int arrdim[AL_MAX_DIM]
Definition: al_hidden.h:53
#define AL_FAILURE
Definition: al_codes.h:33
#define s
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
int i
Definition: analysis.c:2
Definition: al_hidden.h:38
int ndim
Definition: al_hidden.h:43
int AL_Get_lbounds ( int  isz,
int *  lbeg,
int *  lend,
int *  gp,
int  style 
)

Get the local indexes for the local portion of a distributed array

Parameters
[in]iszInteger pointer to the input array descriptor
[out]lbegArray of ndim integers containing the start points
[out]lendArray of ndim integers containing the end points
[out]gpArray of ndim integers containing the ghost points
[in]styleIndex style: AL_C_INDEXES or AL_FORTRAN_INDEXES

Definition at line 668 of file al_sz_get.c.

678 {
679  register int i;
680  int ndim;
681  SZ *s;
682 
683  /*
684  Check that isz points to an allocated SZ
685  */
686  if( stack_ptr[isz] == AL_STACK_FREE ){
687  printf("AL_Get_lbounds: wrong SZ pointer\n");
688  return (int) AL_FAILURE;
689  }
690 
691  /*
692  Get the SZ structure isz is pointing at
693  */
694  s = sz_stack[isz];
695  ndim = s->ndim;
696  for(i=0;i<ndim;i++){
697  lbeg[i] = s->lbeg[i];
698  lend[i] = s->lend[i];
699  gp[i] = s->bg[i];
700  }
701 
702  return (int) AL_SUCCESS;
703 }
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_STACK_FREE
Definition: al_codes.h:24
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
int bg[AL_MAX_DIM]
Definition: al_hidden.h:72
#define AL_FAILURE
Definition: al_codes.h:33
#define s
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
int i
Definition: analysis.c:2
Definition: al_hidden.h:38
int lbeg[AL_MAX_DIM]
Definition: al_hidden.h:68
int lend[AL_MAX_DIM]
Definition: al_hidden.h:70
int ndim
Definition: al_hidden.h:43

Here is the caller graph for this function:

int AL_Get_local_dim ( int  isz,
int *  ldims 
)

Get the local dimensions of a distributed array WITHOUT GHOST POINTS

Parameters
[in]iszInteger pointer to the input array descriptor
[out]ldimsArray of integers with local dimensions

Definition at line 309 of file al_sz_get.c.

316 {
317  register int i;
318  int ndim;
319  SZ *s;
320 
321  /*
322  Check that isz points to an allocated SZ
323  */
324  if( stack_ptr[isz] == AL_STACK_FREE ){
325  printf("AL_Get_local_dim: wrong SZ pointer\n");
326  return (int) AL_FAILURE;
327  }
328 
329  /*
330  Get the SZ structure isz is pointing at
331  */
332  s = sz_stack[isz];
333 
334  ndim = s->ndim;
335  for(i=0;i<ndim;i++){ ldims[i] = s->larrdim[i];}
336 
337  return (int) AL_SUCCESS;
338 }
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_STACK_FREE
Definition: al_codes.h:24
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
int larrdim[AL_MAX_DIM]
Definition: al_hidden.h:54
#define AL_FAILURE
Definition: al_codes.h:33
#define s
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
int i
Definition: analysis.c:2
Definition: al_hidden.h:38
int ndim
Definition: al_hidden.h:43

Here is the caller graph for this function:

int AL_Get_local_dim_gp ( int  sz_ptr,
int *  ldims_gp 
)

Get the local dimensions of a distributed array WITH GHOST POINTS

Parameters
[in]sz_ptrInteger pointer to the input array descriptor
[out]ldims_gpArray of integers with local dimensions (including ghost points)

Definition at line 341 of file al_sz_get.c.

349 {
350  register int i;
351  int isz;
352  int ndim;
353  SZ *s;
354 
355  isz = sz_ptr;
356 
357  /*
358  Check that isz points to an allocated SZ
359  */
360  if( stack_ptr[isz] == AL_STACK_FREE ){
361  printf("AL_Get_local_dim_gp: wrong SZ pointer\n");
362  return (int) AL_FAILURE;
363  }
364 
365  /*
366  Get the SZ structure isz is pointing at
367  */
368  s = sz_stack[isz];
369 
370  ndim = s->ndim;
371  for(i=0;i<ndim;i++){ ldims_gp[i] = s->larrdim_gp[i];}
372 
373  return (int) AL_SUCCESS;
374 }
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_STACK_FREE
Definition: al_codes.h:24
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
int larrdim_gp[AL_MAX_DIM]
Definition: al_hidden.h:56
#define AL_FAILURE
Definition: al_codes.h:33
#define s
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
int i
Definition: analysis.c:2
Definition: al_hidden.h:38
int ndim
Definition: al_hidden.h:43
long long AL_Get_offset ( int  )

Definition at line 327 of file al_io.c.

332 {
333  SZ *s;
334 
335  if( stack_ptr[sz_ptr] == AL_STACK_FREE){
336  printf("AL_Get_offset: wrong SZ pointer\n");
337  }
338  s = sz_stack[sz_ptr];
339  return s->io_offset;
340 }
#define AL_STACK_FREE
Definition: al_codes.h:24
#define s
MPI_Offset io_offset
Definition: al_hidden.h:116
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
Definition: al_hidden.h:38
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26

Here is the caller graph for this function:

int AL_Get_offsets ( int  sz_ptr,
int *  offset 
)

Get the offsets for the for loop index calculation

NOTE: The typical multidimensional C for loop will look like, for example:

... ;

AL_Get_offsets(sz_ptr, offset);

for(k=kbeg; k<=kend; k++){ koff = k*offset[2]; for(j=jbeg; j<=jend; j++){ joff = j*offset[1]; for(i=ibeg; i<=iend; i++){ ioff = i*offset[0]+joff+koff; a[ioff] = ... ; } } } ... ;

Where a[ioff] is a generic 3D array.

Parameters
[in]sz_ptrInteger pointer to the input array descriptor
[out]offsetArray of integers with offsets

Definition at line 378 of file al_sz_get.c.

406 {
407  register int i;
408  int isz;
409  int ndim;
410  SZ *s;
411 
412  isz = sz_ptr;
413 
414  /*
415  Check that isz points to an allocated SZ
416  */
417  if( stack_ptr[isz] == AL_STACK_FREE ){
418  printf("AL_Get_offsets: wrong SZ pointer\n");
419  return (int) AL_FAILURE;
420  }
421 
422  /*
423  Get the SZ structure isz is pointing at
424  */
425  s = sz_stack[isz];
426  ndim = s->ndim;
427  for(i=0;i<ndim;i++){
428  offset[i] = s->offset[i];
429  }
430 
431  return (int) AL_SUCCESS;
432 }
int offset[AL_MAX_DIM]
Definition: al_hidden.h:76
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_STACK_FREE
Definition: al_codes.h:24
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
#define AL_FAILURE
Definition: al_codes.h:33
#define s
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
int i
Definition: analysis.c:2
Definition: al_hidden.h:38
int ndim
Definition: al_hidden.h:43
int AL_Get_parallel_dim ( int  isz,
int *  pardims 
)

Get the parallel dimensions of a distributed array

Parameters
[in]iszInteger pointer to the input array descriptor
[out]pardimsArray of integers with parallel dimensions [AL_TRUE|AL_FALSE]

Definition at line 528 of file al_sz_get.c.

536 {
537  register int i;
538  int ndim;
539  SZ *s;
540 
541  /*
542  Check that isz points to an allocated SZ
543  */
544  if( stack_ptr[isz] == AL_STACK_FREE ){
545  printf("AL_Get_parallel_dim: wrong SZ pointer\n");
546  return (int) AL_FAILURE;
547  }
548 
549  /*
550  Get the SZ structure isz is pointing at
551  */
552  s = sz_stack[isz];
553  ndim = s->ndim;
554  for(i=0;i<ndim;i++){
555  pardims[i] = s->isparallel[i];
556  }
557 
558  return (int) AL_SUCCESS;
559 }
int isparallel[AL_MAX_DIM]
Definition: al_hidden.h:78
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_STACK_FREE
Definition: al_codes.h:24
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
#define AL_FAILURE
Definition: al_codes.h:33
#define s
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
int i
Definition: analysis.c:2
Definition: al_hidden.h:38
int ndim
Definition: al_hidden.h:43
int AL_Get_periodic_dim ( int  isz,
int *  periods 
)

Get the periodic dimensions of the distributed arrays

Parameters
[in]iszInteger pointer to the input array descriptor
[out]periodsArray of integers with periodic dimensions [AL_TRUE|AL_FALSE]

Definition at line 563 of file al_sz_get.c.

571 {
572  register int i;
573  int ndim;
574  SZ *s;
575 
576  /*
577  Check that isz points to an allocated SZ
578  */
579  if( stack_ptr[isz] == AL_STACK_FREE ){
580  printf("AL_Get_periodic_dim: wrong SZ pointer\n");
581  return (int) AL_FAILURE;
582  }
583 
584  /*
585  Get the SZ structure isz is pointing at
586  */
587  s = sz_stack[isz];
588  ndim = s->ndim;
589  for(i=0;i<ndim;i++){
590  periods[i] = s->isperiodic[i];
591  }
592 
593  return (int) AL_SUCCESS;
594 }
int isperiodic[AL_MAX_DIM]
Definition: al_hidden.h:80
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_STACK_FREE
Definition: al_codes.h:24
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
#define AL_FAILURE
Definition: al_codes.h:33
#define s
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
int i
Definition: analysis.c:2
Definition: al_hidden.h:38
int ndim
Definition: al_hidden.h:43
int AL_Get_size ( int  isz,
int *  size 
)

Get the size of the communicator associated with a distributed array

Parameters
[in]iszInteger pointer to the input array descriptor (input)
[in]sizeInteger pointer to size (output)

Definition at line 23 of file al_sz_get.c.

30 {
31  int nproc;
32  MPI_Comm comm;
33  SZ *s;
34 
35  /*
36  Check that isz points to an allocated SZ
37  */
38  if( stack_ptr[isz] == AL_STACK_FREE ){
39  printf("AL_Get_comm: wrong SZ pointer\n");
40  return (int) AL_FAILURE;
41  }
42 
43  /*
44  Get the SZ structure isz is pointing at
45  */
46  s = sz_stack[isz];
47 
48  comm = s->comm;
49 
50  MPI_Comm_size(comm, &nproc);
51 
52  *size = nproc;
53 
54  return (int) AL_SUCCESS;
55 }
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_STACK_FREE
Definition: al_codes.h:24
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
#define AL_FAILURE
Definition: al_codes.h:33
#define s
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
MPI_Comm comm
Definition: al_hidden.h:47
Definition: al_hidden.h:38
int AL_Get_staggered_dim ( int  isz,
int *  stagger 
)

Get the staggered dimensions of the distributed array

Parameters
[in]iszInteger pointer to the input array descriptor
[out]staggerArray of integers with staggered dimensions [AL_TRUE|AL_FALSE]

Definition at line 598 of file al_sz_get.c.

606 {
607  register int i;
608  int ndim;
609  SZ *s;
610 
611  /*
612  Check that isz points to an allocated SZ
613  */
614  if( stack_ptr[isz] == AL_STACK_FREE ){
615  printf("AL_Get_staggered_dim: wrong SZ pointer\n");
616  return (int) AL_FAILURE;
617  }
618 
619  /*
620  Get the SZ structure isz is pointing at
621  */
622  s = sz_stack[isz];
623  ndim = s->ndim;
624  for(i=0;i<ndim;i++){
625  stagger[i] = s->isstaggered[i];
626  }
627 
628  return (int) AL_SUCCESS;
629 }
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_STACK_FREE
Definition: al_codes.h:24
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
int isstaggered[AL_MAX_DIM]
Definition: al_hidden.h:82
#define AL_FAILURE
Definition: al_codes.h:33
#define s
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
int i
Definition: analysis.c:2
Definition: al_hidden.h:38
int ndim
Definition: al_hidden.h:43
int AL_Get_stride ( int  sz_ptr,
int *  stride 
)

Get the strides for the for loop index calculation NOTE: The typical multidimensional C for loop will look like, for example:

... ;

AL_Get_stride(sz_ptr, stride);

for(k=kbeg; k<=kend; k++){ koff = k*stride[2]; for(j=jbeg; j<=jend; j++){ joff = j*stride[1]; for(i=ibeg; i<=iend; i++){ ioff = i*stride[0]+joff+koff; a[ioff] = ... ; } } } ... ;

Where a[ioff] is a generic 3D array.

Parameters
[in]sz_ptrInteger pointer to the input array descriptor
[out]strideArray of strides

Definition at line 437 of file al_sz_get.c.

464 {
465  register int i;
466  int isz;
467  int ndim;
468  SZ *s;
469 
470  isz = sz_ptr;
471 
472  /*
473  Check that isz points to an allocated SZ
474  */
475  if( stack_ptr[isz] == AL_STACK_FREE ){
476  printf("AL_Get_offsets: wrong SZ pointer\n");
477  return (int) AL_FAILURE;
478  }
479 
480  /*
481  Get the SZ structure isz is pointing at
482  */
483  s = sz_stack[isz];
484  ndim = s->ndim;
485  for(i=0;i<ndim;i++){
486  stride[i] = s->stride[i];
487  }
488 
489  return (int) AL_SUCCESS;
490 }
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_STACK_FREE
Definition: al_codes.h:24
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
int stride[AL_MAX_DIM]
Definition: al_hidden.h:77
#define AL_FAILURE
Definition: al_codes.h:33
#define s
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
int i
Definition: analysis.c:2
Definition: al_hidden.h:38
int ndim
Definition: al_hidden.h:43
int AL_Get_type ( int  isz,
AL_Datatype type 
)

Get the basic data type of a distributed array

Parameters
[in]iszInteger pointer to the input array descriptor
[out]typePointer to datatype

Definition at line 177 of file al_sz_get.c.

184 {
185  SZ *s;
186 
187  /*
188  Check that isz points to an allocated SZ
189  */
190  if( stack_ptr[isz] == AL_STACK_FREE ){
191  printf("AL_Get_type: wrong SZ pointer\n");
192  return (int) AL_FAILURE;
193  }
194 
195  /*
196  Get the SZ structure isz is pointing at
197  */
198  s = sz_stack[isz];
199 
200  *type = s->type;
201 
202  return (int) AL_SUCCESS;
203 }
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_STACK_FREE
Definition: al_codes.h:24
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
#define AL_FAILURE
Definition: al_codes.h:33
#define s
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
MPI_Datatype type
Definition: al_hidden.h:40
Definition: al_hidden.h:38
int AL_Init ( int *  argc,
char ***  argv 
)

Initialize the AL Tool. It contains a call to MPI_Init().

Parameters
[in]argcinteger pointer to number of arguments
[in]argvpointer to argv list

Definition at line 18 of file al_init.c.

25 {
26  int myrank, nproc, errcode;
27  int flag;
28 
29  errcode = MPI_Initialized(&flag);
30 
31  if( !flag )
32  errcode = MPI_Init(argc, argv);
33 
34  MPI_Comm_rank(MPI_COMM_WORLD, &myrank);
35  MPI_Comm_size(MPI_COMM_WORLD, &nproc);
36 
37 #ifdef DEBUG
38  printf("AL_Init: Called MPI_init from C: %d\n",errcode);
39 #endif
40 
41  /* Initialize the SZ pointers stack */
43 
44  /* Synchronize just in case */
45  MPI_Barrier(MPI_COMM_WORLD);
46 
48 
49  return errcode;
50 }
#define AL_TRUE
Definition: al_codes.h:28
static int AL_initialized
Definition: al_init.c:15
int AL_Init_stack_()
Definition: al_szptr_.c:33

Here is the call graph for this function:

Here is the caller graph for this function:

int AL_Init_stack_ ( )

Initialize the stack of SZ pointers This routine is called internally by AL_Init.

Definition at line 33 of file al_szptr_.c.

38 {
39  register int i;
40  int myrank;
41 
42  stack_top = 0;
43  stack_used = 0;
44 
45  for( i=0; i<AL_MAX_ARRAYS;i++){ stack_ptr[i] = AL_STACK_FREE ;}
46 
47  MPI_Comm_rank(MPI_COMM_WORLD, &myrank);
48 
49 #ifdef DEBUG
50  printf("AL_Init_stack_: SZ stack initialized\n");
51 #endif
52 
53  return 0;
54 }
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
#define AL_STACK_FREE
Definition: al_codes.h:24
static int stack_used
Definition: al_szptr_.c:30
static int stack_top
Definition: al_szptr_.c:29
int i
Definition: analysis.c:2
#define AL_MAX_ARRAYS
Definition: al_codes.h:21

Here is the caller graph for this function:

int AL_Initialized ( )

Test whether or not AL was initialized

Returns
This function returns AL_TRUE if AL is initialized, AL_FALSE otherwise.

Definition at line 53 of file al_init.c.

58 {
59  return AL_initialized;
60 }
static int AL_initialized
Definition: al_init.c:15
int AL_Is_boundary ( int  sz_ptr,
int *  is_gbeg,
int *  is_gend 
)

This routine is useful when implementing the physical boundary conditions on a problem. It returns two arrays that are set to AL_TRUE or AL_FALSE depending on wether or not the local beginning (ending) address for the array in each direction is actually a global one.

Parameters
[in]sz_ptrpointer to a distributed array descriptor (integer)
[out]is_gbegint array set to AL_TRUE or AL_FALSE for global beginning
[out]is_gendint array set to AL_TRUE or AL_FALSE for global ending

Definition at line 28 of file al_boundary.c.

42 {
43  register int i;
44  int myrank, nproc,ndims;
45  MPI_Comm comm;
46  SZ *s;
47 
48  /* DIAGNOSTICS
49  Check that sz_ptr points to an allocated SZ
50  */
51  if( stack_ptr[sz_ptr] == AL_STACK_FREE){
52  printf("AL_Is_boundary: wrong SZ pointer\n");
53  }
54 
55  s = sz_stack[sz_ptr];
56 
57  myrank = s->rank;
58  nproc = s->size;
59  comm = s->comm;
60  ndims = s->ndim;
61 
62  for(i=0;i<ndims;i++){
63  if( s->beg[i] == s->bg[i] ){
64  is_gbeg[i] = AL_TRUE;
65  } else {
66  is_gbeg[i] = AL_FALSE;
67  }
68  if( s->end[i] == s->arrdim[i]+s->bg[i]-1 ){
69  is_gend[i] = AL_TRUE;
70  } else {
71  is_gend[i] = AL_FALSE;
72  }
73  }
74 
75  /* DIAGNOSTICS */
76 #ifdef DEBUG
77  printf("AL_Is_boundary: completed\n");
78 #endif
79 
80  return (int) AL_SUCCESS;
81 }
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
int end[AL_MAX_DIM]
Definition: al_hidden.h:60
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_TRUE
Definition: al_codes.h:28
#define AL_STACK_FREE
Definition: al_codes.h:24
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
int arrdim[AL_MAX_DIM]
Definition: al_hidden.h:53
int rank
Definition: al_hidden.h:44
int bg[AL_MAX_DIM]
Definition: al_hidden.h:72
int beg[AL_MAX_DIM]
Definition: al_hidden.h:58
#define s
int size
Definition: al_hidden.h:45
int i
Definition: analysis.c:2
MPI_Comm comm
Definition: al_hidden.h:47
#define AL_FALSE
Definition: al_codes.h:29
Definition: al_hidden.h:38
int ndim
Definition: al_hidden.h:43

Here is the caller graph for this function:

int AL_Read_array ( void *  va,
int  sz_ptr,
int  istag 
)

Read a distributed array to a file in parallel using synchronous and collective IO operations

Parameters
[in]bufferpointer to the buffer to write
[in]sz_ptrinteger pointer to a distributed array descriptor
[in]istagset it to -1 for centred variables, set to 0,1,2 for staggered field in the x,y,z direction

Bugs fixed on Aug 26, 2012:

  • nelem has been declared long long in order to manage dbl and flt output 'single_file' in which each PLUTO variable is >= 4GB
  • the updating of nelem after the writing of a variable has been modified from ''nelem *= (long long)(s->arrdim[i] +(istag == i));'' to ''nelem *= (long long)(s->arrdim[i]);'' because when a staggered variable is written, in the sz_ptr descriptor there is the information about the really number of element of this variable. This appens because the function AL_Set_staggered_dim is called for a staggered array.

Definition at line 209 of file al_io.c.

220 {
221  char *a;
222  register int i;
223  int myrank;
224  long long nelem;
225  SZ *s;
226  MPI_Status status;
227  AL_Datatype gsub_arr, lsub_arr;
228 
229 
230  MPI_File ifp;
231  MPI_Offset offset;
232  int size;
233 
234  a = (void *)va;
235 
236  s = sz_stack[sz_ptr];
237  myrank = s->rank;
238 
239  ifp = s->ifp;
240 
241  offset = s->io_offset;
242  if (istag == -1){
243  gsub_arr = s->gsubarr;
244  lsub_arr = s->lsubarr;
245  }else{
246  gsub_arr = s->gsubarr_stag[istag];
247  lsub_arr = s->lsubarr_stag[istag];
248  }
249 
250  MPI_Barrier(s->comm);
251 
252  MPI_File_set_view(ifp, offset, MPI_BYTE, gsub_arr,
253  "native", MPI_INFO_NULL);
254  MPI_File_read_all(ifp, a, 1, lsub_arr, &status);
255 
256  MPI_Type_size( s->type, &size);
257 
258  nelem = 1;
259  for(i = 0; i < s->ndim; i++) {
260 /* -------------------------------------------------------------------- */
261 /*! Bugs fixed on Aug 26, 2012:
262  \li nelem has been declared long long in order to manage dbl and flt
263  output 'single_file' in which each PLUTO variable is >= 4GB
264  \li the updating of nelem after the writing of a variable has been
265  modified from ''nelem *= (long long)(s->arrdim[i] +(istag == i));''
266  to ''nelem *= (long long)(s->arrdim[i]);'' because when a staggered
267  variable is written, in the sz_ptr descriptor there is the
268  information about the really number of element of this variable.
269  This appens because the function AL_Set_staggered_dim is called
270  for a staggered array. */
271 /* -------------------------------------------------------------------- */
272 /* nelem *= (long long)(s->arrdim[i] + (istag == i)); */
273  nelem *= (long long)(s->arrdim[i]);
274  }
275  s->io_offset += (long long)(size)*nelem;
276 
277  return (int) AL_SUCCESS;
278 }
static double a
Definition: init.c:135
#define AL_Datatype
Definition: al_defs.h:22
#define AL_SUCCESS
Definition: al_codes.h:32
MPI_Datatype gsubarr
Definition: al_hidden.h:110
int arrdim[AL_MAX_DIM]
Definition: al_hidden.h:53
int rank
Definition: al_hidden.h:44
MPI_File ifp
Definition: al_hidden.h:117
MPI_Datatype lsubarr
Definition: al_hidden.h:111
MPI_Datatype lsubarr_stag[AL_MAX_DIM]
Definition: al_hidden.h:114
#define s
int i
Definition: analysis.c:2
MPI_Comm comm
Definition: al_hidden.h:47
MPI_Datatype type
Definition: al_hidden.h:40
MPI_Offset io_offset
Definition: al_hidden.h:116
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
Definition: al_hidden.h:38
MPI_Datatype gsubarr_stag[AL_MAX_DIM]
Definition: al_hidden.h:113
int ndim
Definition: al_hidden.h:43

Here is the caller graph for this function:

int AL_Read_common ( void *  ,
int  ,
AL_Datatype  ,
int   
)
int AL_Set_comm ( MPI_Comm  comm,
int  isz 
)

Set the communicator for a distributed array

Parameters
[in]commMPI communicator the array is associated with
[out]iszInteger pointer to the input array descriptor
Returns
AL_SUCCESS if the communicator is set correctly, AL_FAILURE otherwise.

Definition at line 24 of file al_sz_set.c.

34 {
35  int myrank, nproc;
36  SZ *s;
37 
38  /*
39  Check that isz points to an allocated SZ
40  */
41  if( stack_ptr[isz] == AL_STACK_FREE ){
42  printf("AL_Set_comm: wrong SZ pointer\n");
43  return (int) AL_FAILURE;
44  }
45 
46  /*
47  Get the SZ structure isz is pointing at
48  */
49  s = sz_stack[isz];
50 
51  MPI_Comm_rank(comm, &myrank);
52  MPI_Comm_size(comm, &nproc);
53 
54  s->comm = comm;
55  s->size = nproc;
56  /*
57  Patch for the case in which this node does not
58  belong to this communicator
59  */
60  if( nproc > 0 ) {
61  s->rank = myrank;
62  } else { s->rank = MPI_UNDEFINED; }
63 
64  return (int) AL_SUCCESS;
65 }
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_STACK_FREE
Definition: al_codes.h:24
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
int rank
Definition: al_hidden.h:44
#define AL_FAILURE
Definition: al_codes.h:33
#define s
int size
Definition: al_hidden.h:45
MPI_Comm comm
Definition: al_hidden.h:47
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
Definition: al_hidden.h:38
int AL_Set_dimensions ( int  ndim,
int  isz 
)

Set the dimensions of a distributed array

Parameters
[in]ndimNumber of dimensions (1 to 5)
[out]iszInteger pointer to the input array descriptor
Returns
AL_SUCCESS if the dimensions are set correctly, AL_FAILURE otherwise.

Definition at line 69 of file al_sz_set.c.

79 {
80  SZ *s;
81 
82  /*
83  Check that isz points to an allocated SZ
84  */
85  if( stack_ptr[isz] == AL_STACK_FREE ){
86  printf("AL_Set_dimensions: wrong SZ pointer\n");
87  return (int) AL_FAILURE;
88  }
89 
90  /*
91  Get the SZ structure isz is pointing at
92  */
93  s = sz_stack[isz];
94 
95  s->ndim = ndim;
96 
97  return (int) AL_SUCCESS;
98 }
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_STACK_FREE
Definition: al_codes.h:24
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
#define AL_FAILURE
Definition: al_codes.h:33
#define s
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
Definition: al_hidden.h:38
int ndim
Definition: al_hidden.h:43

Here is the caller graph for this function:

int AL_Set_ghosts ( int *  ghosts,
int  isz 
)

Set the ghost points of a distributed array

Parameters
[in]ghostArray of integers with size of ghost points in each dimension
[out]iszInteger pointer to the input array descriptor
Returns
AL_SUCCESS if the ghost points are set correctly, AL_FAILURE otherwise.

Definition at line 333 of file al_sz_set.c.

344 {
345  register int i;
346  int ndim;
347  SZ *s;
348 
349  /*
350  Check that isz points to an allocated SZ
351  */
352  if( stack_ptr[isz] == AL_STACK_FREE ){
353  printf("AL_Set_ghosts: wrong SZ pointer\n");
354  return (int) AL_FAILURE;
355  }
356 
357  /*
358  Get the SZ structure isz is pointing at
359  */
360  s = sz_stack[isz];
361 
362  ndim = s->ndim;
363  for(i=0;i<ndim;i++){
364  s->bg[i] = ghosts[i];
365  s->eg[i] = ghosts[i];
366  }
367 
368  return (int) AL_SUCCESS;
369 }
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_STACK_FREE
Definition: al_codes.h:24
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
int bg[AL_MAX_DIM]
Definition: al_hidden.h:72
int eg[AL_MAX_DIM]
Definition: al_hidden.h:74
#define AL_FAILURE
Definition: al_codes.h:33
#define s
int i
Definition: analysis.c:2
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
Definition: al_hidden.h:38
int ndim
Definition: al_hidden.h:43

Here is the caller graph for this function:

int AL_Set_global_dim ( int *  gdims,
int  isz 
)

Set the global dimensions of a distributed array

Parameters
[in]gdimsArray of integers with global dimensions
[out]iszInteger pointer to the input array descriptor
Returns
AL_SUCCESS if the global dimension are set correctly, AL_FAILURE otherwise.

Definition at line 141 of file al_sz_set.c.

151 {
152  register int i;
153  int ndim;
154 
155  SZ *s;
156 
157  /*
158  Check that isz points to an allocated SZ
159  */
160  if( stack_ptr[isz] == AL_STACK_FREE ){
161  printf("AL_Set_global_dim: wrong SZ pointer\n");
162  return (int) AL_FAILURE;
163  }
164 
165  /*
166  Get the SZ structure isz is pointing at
167  */
168  s = sz_stack[isz];
169 
170  ndim = s->ndim;
171  for(i=0;i<ndim;i++){ s->arrdim[i] = gdims[i] ;}
172 
173  return (int) AL_SUCCESS;
174 }
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_STACK_FREE
Definition: al_codes.h:24
int arrdim[AL_MAX_DIM]
Definition: al_hidden.h:53
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
#define AL_FAILURE
Definition: al_codes.h:33
#define s
int i
Definition: analysis.c:2
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
Definition: al_hidden.h:38
int ndim
Definition: al_hidden.h:43

Here is the caller graph for this function:

int AL_Set_local_dim ( int *  ldims,
int  isz 
)

Set the local dimensions of a distributed array

Parameters
[in]ldimsArray of integers with local dimensions
[out]iszInteger pointer to the input array descriptor
Returns
AL_SUCCESS if the local dimension are set correctly, AL_FAILURE otherwise.

Definition at line 177 of file al_sz_set.c.

187 {
188  register int i;
189  int ndim;
190  SZ *s;
191 
192  /*
193  Check that isz points to an allocated SZ
194  */
195  if( stack_ptr[isz] == AL_STACK_FREE ){
196  printf("AL_Set_local_dim: wrong SZ pointer\n");
197  return (int) AL_FAILURE;
198  }
199 
200  /*
201  Get the SZ structure isz is pointing at
202  */
203  s = sz_stack[isz];
204 
205  ndim = s->ndim;
206  for(i=0;i<ndim;i++){
207  s->larrdim[i] = ldims[i];
208  s->lbeg[i] = 0;
209  s->lend[i] = ldims[i]-1;
210  }
211 
212  return (int) AL_SUCCESS;
213 }
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_STACK_FREE
Definition: al_codes.h:24
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
int larrdim[AL_MAX_DIM]
Definition: al_hidden.h:54
#define AL_FAILURE
Definition: al_codes.h:33
#define s
int i
Definition: analysis.c:2
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
Definition: al_hidden.h:38
int lbeg[AL_MAX_DIM]
Definition: al_hidden.h:68
int lend[AL_MAX_DIM]
Definition: al_hidden.h:70
int ndim
Definition: al_hidden.h:43
int AL_Set_offset ( int  ,
long  long 
)

Definition at line 342 of file al_io.c.

347 {
348  SZ *s;
349 
350  if( stack_ptr[sz_ptr] == AL_STACK_FREE){
351  printf("AL_Get_offset: wrong SZ pointer\n");
352  }
353  s = sz_stack[sz_ptr];
354  s->io_offset = offset;
355  return (int)AL_SUCCESS;
356 }
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_STACK_FREE
Definition: al_codes.h:24
#define s
MPI_Offset io_offset
Definition: al_hidden.h:116
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
Definition: al_hidden.h:38
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26

Here is the caller graph for this function:

int AL_Set_parallel_dim ( int *  pardims,
int  isz 
)

Set the parallel dimensions of a distributed array

Parameters
[in]pardimsArray of integers with parallel dimensions [AL_TRUE|AL_FALSE]
[out]iszInteger pointer to the input array descriptor
Returns
AL_SUCCESS if the parallel dimension are set correctly, AL_FAILURE otherwise.

Definition at line 216 of file al_sz_set.c.

226 {
227  register int i;
228  int ndim;
229  SZ *s;
230 
231  /*
232  Check that isz points to an allocated SZ
233  */
234  if( stack_ptr[isz] == AL_STACK_FREE ){
235  printf("AL_Set_parallel_dim: wrong SZ pointer\n");
236  return (int) AL_FAILURE;
237  }
238 
239  /*
240  Get the SZ structure isz is pointing at
241  */
242  s = sz_stack[isz];
243 
244  ndim = s->ndim;
245  for(i=0;i<ndim;i++){
246  s->isparallel[i] = pardims[i];
247  }
248 
249  return (int) AL_SUCCESS;
250 }
int isparallel[AL_MAX_DIM]
Definition: al_hidden.h:78
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_STACK_FREE
Definition: al_codes.h:24
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
#define AL_FAILURE
Definition: al_codes.h:33
#define s
int i
Definition: analysis.c:2
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
Definition: al_hidden.h:38
int ndim
Definition: al_hidden.h:43

Here is the caller graph for this function:

int AL_Set_periodic_dim ( int *  periods,
int  isz 
)

Set the periodic dimensions of a distributed array

Parameters
[in]periodsArray of integers with periodic dimensions [AL_TRUE|AL_FALSE]
[out]iszInteger pointer to the input array descriptor
Returns
AL_SUCCESS if the periodic dimension are set correctly, AL_FAILURE otherwise.

Definition at line 253 of file al_sz_set.c.

263 {
264  register int i;
265  int ndim;
266  SZ *s;
267 
268  /*
269  Check that isz points to an allocated SZ
270  */
271  if( stack_ptr[isz] == AL_STACK_FREE ){
272  printf("AL_Set_periodic_dim: wrong SZ pointer\n");
273  return (int) AL_FAILURE;
274  }
275 
276  /*
277  Get the SZ structure isz is pointing at
278  */
279  s = sz_stack[isz];
280 
281  ndim = s->ndim;
282  for(i=0;i<ndim;i++){
283  if( periods[i] != AL_TRUE && periods[i] != AL_FALSE ){
284  printf("Warning: periods has illegal values in AL_Set_periodic_dim\n");
285  }
286  s->isperiodic[i] = periods[i];
287  }
288 
289  return (int) AL_SUCCESS;
290 }
int isperiodic[AL_MAX_DIM]
Definition: al_hidden.h:80
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_TRUE
Definition: al_codes.h:28
#define AL_STACK_FREE
Definition: al_codes.h:24
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
#define AL_FAILURE
Definition: al_codes.h:33
#define s
int i
Definition: analysis.c:2
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
#define AL_FALSE
Definition: al_codes.h:29
Definition: al_hidden.h:38
int ndim
Definition: al_hidden.h:43

Here is the caller graph for this function:

int AL_Set_staggered_dim ( int *  stagger,
int  isz 
)

Set the staggered dimensions of a distributed array

Parameters
[in]staggerArray of integers with staggered dimensions [AL_TRUE|AL_FALSE]
[out]iszInteger pointer to the input array descriptor
Returns
AL_SUCCESS if the staggered dimension are set correctly, AL_FAILURE otherwise.

Definition at line 293 of file al_sz_set.c.

303 {
304  register int i;
305  int ndim;
306  SZ *s;
307 
308  /*
309  Check that isz points to an allocated SZ
310  */
311  if( stack_ptr[isz] == AL_STACK_FREE ){
312  printf("AL_Set_staggered_dim: wrong SZ pointer\n");
313  return (int) AL_FAILURE;
314  }
315 
316  /*
317  Get the SZ structure isz is pointing at
318  */
319  s = sz_stack[isz];
320 
321  ndim = s->ndim;
322  for(i=0;i<ndim;i++){
323  if( stagger[i] != AL_TRUE && stagger[i] != AL_FALSE ){
324  printf("Warning: stagger has illegal values in AL_Set_staggered_dim\n");
325  }
326  s->isstaggered[i] = stagger[i];
327  }
328 
329  return (int) AL_SUCCESS;
330 }
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_TRUE
Definition: al_codes.h:28
#define AL_STACK_FREE
Definition: al_codes.h:24
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
int isstaggered[AL_MAX_DIM]
Definition: al_hidden.h:82
#define AL_FAILURE
Definition: al_codes.h:33
#define s
int i
Definition: analysis.c:2
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
#define AL_FALSE
Definition: al_codes.h:29
Definition: al_hidden.h:38
int ndim
Definition: al_hidden.h:43

Here is the caller graph for this function:

int AL_Set_type ( AL_Datatype  type,
int  nelem,
int  isz 
)

Set the basic data type of a distributed array. The data types are identical to the MPI datatypes, and they can be defined as AL_<type> or MPI_<type> (e.g. AL_FLOAT or MPI_FLOAT).

Parameters
[in]typeDatatype (AL_Datatype or MPI_Datatype)
[in]nelemNumber of type 'type' elements in array elements
[out]iszInteger pointer to the input array descriptor
Returns
AL_SUCCESS if the type is set correctly, AL_FAILURE otherwise.

Definition at line 101 of file al_sz_set.c.

112 {
113  SZ *s;
114  AL_Datatype ivector;
115 
116  /*
117  Check that isz points to an allocated SZ
118  */
119  if( stack_ptr[isz] == AL_STACK_FREE ){
120  printf("AL_Set_type: wrong SZ pointer\n");
121  return (int) AL_FAILURE;
122  }
123 
124  /*
125  Get the SZ structure isz is pointing at
126  */
127  s = sz_stack[isz];
128 
129  if( nelem <= 1 ){
130  s->type = type;
131  } else {
132  MPI_Type_contiguous( nelem, type, &ivector);
133  MPI_Type_commit(&ivector);
134  s->type = ivector;
135  }
136 
137  return (int) AL_SUCCESS;
138 }
#define AL_Datatype
Definition: al_defs.h:22
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_STACK_FREE
Definition: al_codes.h:24
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
#define AL_FAILURE
Definition: al_codes.h:33
#define s
MPI_Datatype type
Definition: al_hidden.h:40
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
Definition: al_hidden.h:38

Here is the caller graph for this function:

int AL_Sort_ ( int  n,
int *  in,
int *  ind 
)

Sort and array of integers

This is a really simple implementation, since we do not really use this for large arrays.

Parameters
[in]nsize of input array (integer)
[in]ininput array
[in]indarray of the sorted index arrays (max to min)

Definition at line 15 of file al_sort_.c.

26 {
27  register int i, temp;
28  int m, r, s;
29 
30  /* Start with the normal index ordering */
31  for(i=0;i<n;i++) ind[i]=i;
32 
33  m=0;
34  while(m<n) {
35  r = m;
36  s = in[m];
37  for(i=m+1;i<n;i++) if( in[ind[i]] > s ) { r=i; s=in[ind[i]]; }
38  swapi_(ind[m],ind[r]);
39  m++;
40  }
41  return 0;
42 }
static int n
Definition: analysis.c:3
#define swapi_(a, b)
Definition: al_sort_.c:12
#define s
int i
Definition: analysis.c:2

Here is the caller graph for this function:

int AL_Sz_free ( int  sz_ptr)

Dellocate a distributed array descriptor

Parameters
[in]sz_ptrInteger pointer to the array descriptor

Definition at line 23 of file al_sz_free.c.

29 {
30  /*
31  Deallocate the SZ structure
32  */
33  return AL_Deallocate_sz_(sz_ptr);
34 
35 }
int AL_Deallocate_sz_(int)
Definition: al_szptr_.c:102

Here is the call graph for this function:

int AL_Sz_init ( MPI_Comm  comm,
int *  sz_ptr 
)

Allocate and initialize a distributed array descriptor

Parameters
[in]commMPI communicator the array is associated with
[out]sz_ptrInteger pointer to the array descriptor
Returns
AL_SUCCESS if the array descriptor is correct initialized, AL_FAILURE otherwise.

Definition at line 22 of file al_sz_init.c.

32 {
33  int myrank, nproc;
34  register int i;
35 
36  /*
37  If this process does not belong to the communicator,
38  or if the communicator is MPI_COMM_NULL, return a
39  AL_UNDEFINED pointer
40  */
41  if( comm == MPI_COMM_NULL ) {
42  *sz_ptr = AL_UNDEFINED;
43  return AL_FAILURE;
44  }
45 
46  /*
47  Allocate the SZ structure
48  */
49  *sz_ptr = AL_Allocate_sz_();
50 
51  MPI_Comm_rank(comm, &myrank);
52  MPI_Comm_size(comm, &nproc);
53 
54  /*
55  Set the initial entries
56  */
57  sz_stack[*sz_ptr]->compiled = AL_FALSE; /* The array is not compiled yet */
58  sz_stack[*sz_ptr]->comm = comm;
59  sz_stack[*sz_ptr]->rank = myrank;
60  sz_stack[*sz_ptr]->size = nproc;
61 
62  /*
63  Set some defaults
64  */
65  sz_stack[*sz_ptr]->type = MPI_DOUBLE;
66  for(i=0;i<AL_MAX_DIM;i++){
67  sz_stack[*sz_ptr]->isparallel[i] = AL_TRUE;
68  sz_stack[*sz_ptr]->isperiodic[i] = AL_TRUE;
69  sz_stack[*sz_ptr]->isstaggered[i] = AL_FALSE;
70  sz_stack[*sz_ptr]->larrdim_gp[i]=1;
71  sz_stack[*sz_ptr]->larrdim[i]=1;
72  sz_stack[*sz_ptr]->arrdim[i]=1;
73  sz_stack[*sz_ptr]->beg[i]=0;
74  sz_stack[*sz_ptr]->end[i]=0;
75  sz_stack[*sz_ptr]->lbeg[i]=0;
76  sz_stack[*sz_ptr]->lend[i]=0;
77  sz_stack[*sz_ptr]->bg[i]=0;
78  sz_stack[*sz_ptr]->eg[i]=0;
79  sz_stack[*sz_ptr]->offset[i]=1;
80  sz_stack[*sz_ptr]->stride[i]=1;
81  }
82 
83  sz_stack[*sz_ptr]->begs = NULL;
84  sz_stack[*sz_ptr]->ends = NULL;
85 
86  return (int) AL_SUCCESS;
87 }
int isperiodic[AL_MAX_DIM]
Definition: al_hidden.h:80
#define AL_MAX_DIM
Definition: al_codes.h:18
int end[AL_MAX_DIM]
Definition: al_hidden.h:60
int isparallel[AL_MAX_DIM]
Definition: al_hidden.h:78
int AL_Allocate_sz_()
Definition: al_szptr_.c:58
int offset[AL_MAX_DIM]
Definition: al_hidden.h:76
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_TRUE
Definition: al_codes.h:28
int * begs
Definition: al_hidden.h:63
int arrdim[AL_MAX_DIM]
Definition: al_hidden.h:53
int larrdim_gp[AL_MAX_DIM]
Definition: al_hidden.h:56
int larrdim[AL_MAX_DIM]
Definition: al_hidden.h:54
int rank
Definition: al_hidden.h:44
int * ends
Definition: al_hidden.h:63
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
int stride[AL_MAX_DIM]
Definition: al_hidden.h:77
int bg[AL_MAX_DIM]
Definition: al_hidden.h:72
int beg[AL_MAX_DIM]
Definition: al_hidden.h:58
int isstaggered[AL_MAX_DIM]
Definition: al_hidden.h:82
int eg[AL_MAX_DIM]
Definition: al_hidden.h:74
#define AL_FAILURE
Definition: al_codes.h:33
int size
Definition: al_hidden.h:45
int i
Definition: analysis.c:2
MPI_Comm comm
Definition: al_hidden.h:47
MPI_Datatype type
Definition: al_hidden.h:40
#define AL_FALSE
Definition: al_codes.h:29
#define AL_UNDEFINED
Definition: al_codes.h:127
int compiled
Definition: al_hidden.h:39
int lbeg[AL_MAX_DIM]
Definition: al_hidden.h:68
int lend[AL_MAX_DIM]
Definition: al_hidden.h:70

Here is the call graph for this function:

Here is the caller graph for this function:

int AL_Type_create_subarray ( int  ndims,
int *  array_of_sizes,
int *  array_of_subsizes,
int *  array_of_starts,
int  order,
MPI_Datatype  oldtype,
MPI_Datatype *  newtype 
)

Creates a datatype describing a subarray of a multidimensional array NOTE: This routine has been modified from R. Thakur's ROMIO implementation of MPI_Type_create_subarry. The reason for keeping a local copy is mainly to reduce this package's dependence on ROMIO. It only affects the non I/O components of the library (e.g. AL_Allgather).

Parameters
[in]ndimsnumber of array dimensions
[in]array_of_sizesnumber of elements of type oldtype in each dimension of the full array
[in]array_of_subsizesnumber of elements of type oldtype in each dimension of the subarray
[in]array_of_startsstarting coordinates of the subarray in each dimension
[in]orderarray storage order flag
[in]oldtypeold datatype (handle)
[out]newtypenew datatype (handle)

Definition at line 18 of file al_subarray_.c.

43 {
44  MPI_Aint extent, disps[AL_MAX_DIM], size, size_with_aint;
45  int i, blklens[AL_MAX_DIM];
46  MPI_Datatype tmp1, tmp2, types[AL_MAX_DIM];
47  MPI_Offset size_with_offset;
48 
49  if (ndims <= 0) {
50  printf("MPI_Type_create_subarray: Invalid ndims argument\n");
51  MPI_Abort(MPI_COMM_WORLD, 1);
52  }
53  if (array_of_sizes <= (int *) 0) {
54  printf("MPI_Type_create_subarray: array_of_sizes is an invalid address\n");
55  MPI_Abort(MPI_COMM_WORLD, 1);
56  }
57  if (array_of_subsizes <= (int *) 0) {
58  printf("MPI_Type_create_subarray: array_of_subsizes is an invalid address\n");
59  MPI_Abort(MPI_COMM_WORLD, 1);
60  }
61  if (array_of_starts <= (int *) 0) {
62  printf("MPI_Type_create_subarray: array_of_starts is an invalid address\n");
63  MPI_Abort(MPI_COMM_WORLD, 1);
64  }
65 
66  for (i=0; i<ndims; i++) {
67  if (array_of_sizes[i] <= 0) {
68  printf("MPI_Type_create_subarray: Invalid value in array_of_sizes\n");
69  MPI_Abort(MPI_COMM_WORLD, 1);
70  }
71  if (array_of_subsizes[i] <= 0) {
72  printf("MPI_Type_create_subarray: Invalid value in array_of_subsizes\n");
73  MPI_Abort(MPI_COMM_WORLD, 1);
74  }
75  if (array_of_starts[i] < 0) {
76  printf("MPI_Type_create_subarray: Invalid value in array_of_starts\n");
77  MPI_Abort(MPI_COMM_WORLD, 1);
78  }
79  }
80 
81  /* order argument checked below */
82 
83  if (oldtype == MPI_DATATYPE_NULL) {
84  printf("MPI_Type_create_subarray: oldtype is an invalid datatype\n");
85  MPI_Abort(MPI_COMM_WORLD, 1);
86  }
87 
88  MPI_Type_extent(oldtype, (MPI_Aint *) &extent);
89 
90 /* check if MPI_Aint is large enough for size of global array.
91  if not, complain. */
92 
93  size_with_aint = extent;
94  for (i=0; i<ndims; i++) size_with_aint *= array_of_sizes[i];
95  size_with_offset = extent;
96  for (i=0; i<ndims; i++) size_with_offset *= array_of_sizes[i];
97  if (size_with_aint != size_with_offset) {
98  printf("MPI_Type_create_subarray: Can't use an array of this size unless the MPI implementation defines a 64-bit MPI_Aint\n");
99  MPI_Abort(MPI_COMM_WORLD, 1);
100  }
101 
102  if (order == AL_ORDER_FORTRAN) {
103  /* dimension 0 changes fastest */
104  if (ndims == 1)
105  MPI_Type_contiguous(array_of_subsizes[0], oldtype, &tmp1);
106  else {
107  MPI_Type_vector(array_of_subsizes[1], array_of_subsizes[0],
108  array_of_sizes[0], oldtype, &tmp1);
109 
110  size = array_of_sizes[0]*extent;
111  for (i=2; i<ndims; i++) {
112  size *= array_of_sizes[i-1];
113  MPI_Type_hvector(array_of_subsizes[i], 1, size, tmp1, &tmp2);
114  MPI_Type_free(&tmp1);
115  tmp1 = tmp2;
116  }
117  }
118 
119  /* add displacement and UB */
120 
121  disps[1] = array_of_starts[0];
122  size = 1;
123  for (i=1; i<ndims; i++) {
124  size *= array_of_sizes[i-1];
125  disps[1] += size*array_of_starts[i];
126  }
127  /* rest done below for both Fortran and C order */
128  }
129 
130  else if (order == AL_ORDER_C) {
131  /* dimension ndims-1 changes fastest */
132  if (ndims == 1)
133  MPI_Type_contiguous(array_of_subsizes[0], oldtype, &tmp1);
134  else {
135  MPI_Type_vector(array_of_subsizes[ndims-2],
136  array_of_subsizes[ndims-1],
137  array_of_sizes[ndims-1], oldtype, &tmp1);
138 
139  size = array_of_sizes[ndims-1]*extent;
140  for (i=ndims-3; i>=0; i--) {
141  size *= array_of_sizes[i+1];
142  MPI_Type_hvector(array_of_subsizes[i], 1, size, tmp1, &tmp2);
143  MPI_Type_free(&tmp1);
144  tmp1 = tmp2;
145  }
146  }
147 
148  /* add displacement and UB */
149 
150  disps[1] = array_of_starts[ndims-1];
151  size = 1;
152  for (i=ndims-2; i>=0; i--) {
153  size *= array_of_sizes[i+1];
154  disps[1] += size*array_of_starts[i];
155  }
156  }
157  else {
158  printf("MPI_Type_create_subarray: Invalid order argument\n");
159  MPI_Abort(MPI_COMM_WORLD, 1);
160  }
161 
162  disps[1] *= extent;
163 
164  disps[2] = extent;
165  for (i=0; i<ndims; i++) disps[2] *= array_of_sizes[i];
166 
167  disps[0] = 0;
168  blklens[0] = blklens[1] = blklens[2] = 1;
169  types[0] = MPI_LB;
170  types[1] = tmp1;
171  types[2] = MPI_UB;
172 
173  MPI_Type_struct(3, blklens, disps, types, newtype);
174 
175  MPI_Type_free(&tmp1);
176 
177  return MPI_SUCCESS;
178 }
#define AL_MAX_DIM
Definition: al_codes.h:18
#define AL_ORDER_C
Definition: al_codes.h:45
#define AL_ORDER_FORTRAN
Definition: al_codes.h:46
int i
Definition: analysis.c:2

Here is the caller graph for this function:

int AL_Valid_ptr ( int  sz_ptr)

Return AL_TRUE if the input pointer points to an allocated distributed array. Return AL_FALSE otherwise.

Parameters
[in]sz_ptrInteger pointer to an entry in the SZ stack

Definition at line 150 of file al_szptr_.c.

158 {
159  if(stack_ptr[sz_ptr]!=AL_STACK_USED ) return AL_FALSE;
160  return sz_stack[sz_ptr]->compiled;
161 }
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
#define AL_STACK_USED
Definition: al_codes.h:25
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
#define AL_FALSE
Definition: al_codes.h:29
int compiled
Definition: al_hidden.h:39
int AL_Write_array ( void *  va,
int  sz_ptr,
int  istag 
)

Write a distributed array to a file in parallel using synchronous and collective IO operations

Parameters
[in]bufferpointer to the buffer to write
[in]sz_ptrinteger pointer to the distributed array descriptor
[in]istagset it to -1 for centred variables, set to 0,1,2 for staggered field in the x,y,z direction

Bugs fixed on Aug 26, 2012:

  • nelem has been declared long long in order to manage dbl and flt output 'single_file' in which each PLUTO variable is >= 4GB
  • the updating of nelem after the writing of a variable has been modified from ''nelem *= (long long)(s->arrdim[i] +(istag == i));'' to ''nelem *= (long long)(s->arrdim[i]);'' because when a staggered variable is written, in the sz_ptr descriptor there is the information about the really number of element of this variable. This appens because the function AL_Set_staggered_dim is called for a staggered array.

Definition at line 115 of file al_io.c.

126 {
127  char *a;
128  register int i;
129  int myrank;
130  long long nelem;
131  int errcode;
132  SZ *s;
133 
134  MPI_Status status;
135  AL_Datatype gsub_arr, lsub_arr;
136 
137  MPI_File ifp;
138  MPI_Offset offset;
139  int size;
140 
141 
142  a = (char *) va;
143 
144  s = sz_stack[sz_ptr];
145  myrank = s->rank;
146  ifp = s->ifp;
147 
148  offset = s->io_offset;
149  if (istag == -1){
150  gsub_arr = s->gsubarr;
151  lsub_arr = s->lsubarr;
152  }else{
153  gsub_arr = s->gsubarr_stag[istag];
154  lsub_arr = s->lsubarr_stag[istag];
155  }
156 
157  MPI_Barrier(s->comm);
158 
159  errcode = MPI_File_set_view(ifp, offset, MPI_BYTE, gsub_arr,
160  "native", MPI_INFO_NULL);
161 
162 #ifdef DEBUG
163  int myid, len;
164  char es[256];
165  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
166  if( errcode ){
167  MPI_Error_string(errcode, es, &len);
168  printf("Errcode from MPI_File_set_view: %d | %s\n", errcode, es);
169  }
170  printf("myid %d, offset file_set_view %lld\n", myid,offset);
171 #endif
172 
173  errcode = MPI_File_write_all(ifp, a, 1, lsub_arr, &status);
174 
175 #ifdef DEBUG
176  if( errcode ){
177  MPI_Error_string(errcode, es, &len);
178  printf("Errcode from MPI_File_write_all: %d | %s\n", errcode, es);
179  }
180 #endif
181 
182  MPI_Type_size( s->type, &size);
183 
184  nelem = 1;
185  for(i = 0; i < s->ndim; i++) {
186 /* -------------------------------------------------------------------- */
187 /*! Bugs fixed on Aug 26, 2012:
188  \li nelem has been declared long long in order to manage dbl and flt
189  output 'single_file' in which each PLUTO variable is >= 4GB
190  \li the updating of nelem after the writing of a variable has been
191  modified from ''nelem *= (long long)(s->arrdim[i] +(istag == i));''
192  to ''nelem *= (long long)(s->arrdim[i]);'' because when a staggered
193  variable is written, in the sz_ptr descriptor there is the
194  information about the really number of element of this variable.
195  This appens because the function AL_Set_staggered_dim is called
196  for a staggered array. */
197 /* -------------------------------------------------------------------- */
198 /* nelem *= (long long)(s->arrdim[i] +(istag == i)); */
199  nelem *= (long long)(s->arrdim[i]);
200  }
201 
202  s->io_offset += (long long)(size)*nelem;
203 
204  return (int) AL_SUCCESS;
205 }
static double a
Definition: init.c:135
#define AL_Datatype
Definition: al_defs.h:22
#define AL_SUCCESS
Definition: al_codes.h:32
MPI_Datatype gsubarr
Definition: al_hidden.h:110
int arrdim[AL_MAX_DIM]
Definition: al_hidden.h:53
int rank
Definition: al_hidden.h:44
MPI_File ifp
Definition: al_hidden.h:117
MPI_Datatype lsubarr
Definition: al_hidden.h:111
MPI_Datatype lsubarr_stag[AL_MAX_DIM]
Definition: al_hidden.h:114
#define s
int i
Definition: analysis.c:2
MPI_Comm comm
Definition: al_hidden.h:47
MPI_Datatype type
Definition: al_hidden.h:40
MPI_Offset io_offset
Definition: al_hidden.h:116
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
Definition: al_hidden.h:38
MPI_Datatype gsubarr_stag[AL_MAX_DIM]
Definition: al_hidden.h:113
int ndim
Definition: al_hidden.h:43

Here is the caller graph for this function:

int AL_Write_array_begin ( void *  va,
int  sz_ptr,
int *  output_stag,
int *  output_dump,
int  output_nvar 
)

Write a distributed array to a parallel file by using asynchronous MPI-IO

Parameters
[in]bufferpointer to the buffer to write
[in]sz_ptrinteger pointer to the distributed array descriptor
[in]output_stagvector sets to -1 for centred variables, and sets to 0,1,2 for staggered field in the x,y,z direction
[in]output_dumpvector sets to 1 if the variable has to be dumped, 0 in the contrary case
[in]output_nvartotal number of variables in PLUTO simulation

Definition at line 26 of file al_write_array_async.c.

41 {
42  char *a;
43  register int i;
44  SZ *s;
45 
46  MPI_File ifp;
47 
48  int errcode;
49 
50  int myrank;
51  MPI_Comm_rank(MPI_COMM_WORLD, &myrank);
52 
53  a = (char *) va;
54 
55  s = sz_stack[sz_ptr];
56  ifp = s->ifp;
57 
58  /* Definition of a "global" filetype given by nvar times gsubarr for the MPI_File_set_view */
59  MPI_Datatype filetype_gsub_arr;
60  int count; /* number of variable to be dumped */
61  count=0;
62  for (i=0; i<output_nvar; i++){
63  if (!(output_dump[i]== 0)) count = count+1;
64  }
65  errcode = MPI_Type_contiguous(count, s->gsubarr, &filetype_gsub_arr);
66  errcode = MPI_Type_commit(&filetype_gsub_arr);
67 
68  /* Definition of a "global" datatype for the MPI_File_write_all_begin */
69  MPI_Datatype dtype_lsub_arr;
70 
71 #ifdef DEBUG
72  if(myrank==0) printf("count=%d, output_nvar=%d\n",count,output_nvar);
73 #endif
74 
75  int *block_leng, *displ;
76  int k=0;
77 
78  block_leng= (int*) malloc(count*sizeof(int));
79  displ=(int*) malloc(count*sizeof(int));
80 
81  for (i=0; i<output_nvar; i++){
82  if (!(output_dump[i]==0)){ /* cioè devo scrivere la variabile */
83  block_leng[k] = 1;
84  displ[k] = i;
85  k=k+1;
86  }
87  }
88 
89  errcode = MPI_Type_indexed(count,block_leng, displ, s->lsubarr, &dtype_lsub_arr);
90  errcode = MPI_Type_commit(&dtype_lsub_arr);
91 /*
92  MPI_Type_get_extent(dtype_lsub_arr,lb,extent);
93  printf("myid, %d, lb %d, extent %d of dtype_lsub_arr\n",myrank,lb,extent);
94 */
95 /* DIAGNOSTICS */
96 #ifdef DEBUG
97  if (myrank==0) {
98  for (i=0; i<output_nvar; i++) printf("myid, %d, count %d, output_nvar %d, i %d, output_dump[i] %d, block_leng[i] %d, displ[i] %d\n", myrank, count, output_nvar, i, output_dump[i], block_leng[i], displ[i]);
99  }
100 #endif
101 
102  free(block_leng);
103  free(displ);
104  MPI_Barrier(s->comm);
105 
106  /* Setting of the new view. Each procs see all the file at the beginning of the writing */
107  errcode = MPI_File_set_view(ifp, 0, MPI_BYTE, filetype_gsub_arr,
108  "native", MPI_INFO_NULL);
109  MPI_Type_free(&filetype_gsub_arr);
110 
111 /* DIAGNOSTICS */
112 #ifdef DEBUG
113  int len;
114  char es[MPI_MAX_ERROR_STRING];
115  MPI_Error_string(errcode, es, &len);
116  printf("myid %d, Errcode from MPI_File_set_view: %d | %s\n", myrank, errcode, es);
117 #endif
118 
119  errcode = MPI_File_write_all_begin(ifp, va, 1,dtype_lsub_arr);
120  MPI_Type_free(&dtype_lsub_arr);
121 
122 /* DIAGNOSTICS */
123 #ifdef DEBUG
124  printf("myid %d\n", myrank);
125  MPI_Error_string(errcode, es, &len);
126  printf("myid %d, Errcode from MPI_File_write_all_begin: %d | %s\n", myrank, errcode, es);
127 #endif
128 
129  return (int) AL_SUCCESS;
130 }
static double a
Definition: init.c:135
#define AL_SUCCESS
Definition: al_codes.h:32
MPI_Datatype gsubarr
Definition: al_hidden.h:110
MPI_File ifp
Definition: al_hidden.h:117
MPI_Datatype lsubarr
Definition: al_hidden.h:111
int k
Definition: analysis.c:2
#define s
int i
Definition: analysis.c:2
MPI_Comm comm
Definition: al_hidden.h:47
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
Definition: al_hidden.h:38
int AL_Write_array_end ( void *  va,
int  sz_ptr 
)

Completition of writing of a distributed array to a parallel file by using asynchronous MPI-IO

Parameters
[in]bufferpointer to the buffer to write
[in]sz_ptrinteger pointer to the distributed array descriptor

Definition at line 133 of file al_write_array_async.c.

141 {
142  char *a;
143  SZ *s;
144 
145  MPI_Status status;
146 
147  MPI_File ifp;
148 
149  int errcode;
150 
151  int myrank;
152  MPI_Comm_rank(MPI_COMM_WORLD, &myrank);
153 
154  a = (char *) va;
155 
156  s = sz_stack[sz_ptr];
157  ifp = s->ifp;
158 
159  errcode=MPI_File_write_all_end(ifp, a, &status);
160 
161  /* DIAGNOSTICS */
162 #ifdef DEBUG
163  int len;
164  char es[128];
165  if( errcode ){
166  MPI_Error_string(errcode, es, &len);
167  printf("Errcode from MPI_File_write_all_end: %d | %s\n", errcode, es);
168  }
169 #endif
170 
171  MPI_File_sync(ifp);
172 
173  return (int) AL_SUCCESS;
174 }
static double a
Definition: init.c:135
#define AL_SUCCESS
Definition: al_codes.h:32
MPI_File ifp
Definition: al_hidden.h:117
#define s
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
Definition: al_hidden.h:38
int AL_Write_common ( void *  ,
int  ,
AL_Datatype  ,
int   
)
int AL_Write_header ( void *  ,
int  ,
AL_Datatype  ,
int   
)

Definition at line 285 of file al_io.c.

293 {
294  char *buffer;
295  int myrank;
296  int nbytes, size;
297  SZ *s;
298  MPI_Status status;
299  MPI_File ifp;
300 
301 
302  MPI_Offset offset;
303  MPI_Type_size(type, &size);
304  buffer = (char *) vbuffer;
305 
306  nbytes = nelem*size;
307  s = sz_stack[sz_ptr];
308  myrank = s->rank;
309 
310  ifp = s->ifp;
311  offset = s->io_offset;
312  MPI_Barrier(s->comm);
313  MPI_File_set_view(ifp, offset, MPI_BYTE, MPI_CHAR,
314  "native", MPI_INFO_NULL);
315  if( myrank == 0 ){
316 /* MPI_File_write(ifp, buffer, strlen(buffer), MPI_CHAR, &status); */
317  MPI_File_write(ifp, buffer, nelem, type, &status);
318  }
319 
320  s->io_offset += nbytes;
321 
322 
323  return (int) AL_SUCCESS;
324 }
#define AL_SUCCESS
Definition: al_codes.h:32
int rank
Definition: al_hidden.h:44
MPI_File ifp
Definition: al_hidden.h:117
#define s
MPI_Comm comm
Definition: al_hidden.h:47
MPI_Offset io_offset
Definition: al_hidden.h:116
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
Definition: al_hidden.h:38

Here is the caller graph for this function: