PLUTO
al_sz_get.c
Go to the documentation of this file.
1 /* ///////////////////////////////////////////////////////////////////// */
2 /*!
3  \file
4  \brief Miscellaneous functions to obtain informations about the
5  properties of the distributed array
6 
7  \author A. Malagoli (University of Chicago)
8  \date Jul 17, 1999
9 */
10 /* ///////////////////////////////////////////////////////////////////// */
11 #include "al_hidden.h" /*I "al_hidden.h" I*/
12 
13 /*
14  The SZ structure stack is defined and maintained
15  in al_szptr_.c
16  Here we include an external reference to it in
17  order to be able to make internal references to it.
18 */
19 extern SZ *sz_stack[AL_MAX_ARRAYS];
20 extern int stack_ptr[AL_MAX_ARRAYS];
21 
22 /* ********************************************************************* */
23 int AL_Get_size(int isz, int *size)
24 /*!
25  * Get the size of the communicator associated with a distributed array
26  *
27  * \param [in] isz Integer pointer to the input array descriptor (input)
28  * \param [in] size Integer pointer to size (output)
29  *********************************************************************** */
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 }
56 
57 
58 /* ********************************************************************* */
59 int AL_Get_rank(int isz, int *rank)
60 /*!
61  * Get the rank of this node in the communicator associated with a
62  * distributed array
63  *
64  * \param [in] isz Integer pointer to the input array descriptor
65  * \param [out] rank Integer pointer to rank
66  *
67  * \return If the current node does not belong to the communicator,
68  * AL_UNDEFINED (same as MPI_UNDEFINED) is returned.
69  * Otherwise AL_SUCCESS is returned.
70  *********************************************************************** */
71 {
72 
73  /*
74  Check that isz points to an allocated SZ
75  */
76  if( stack_ptr[isz] == AL_STACK_FREE ){
77  printf("AL_Get_comm: wrong SZ pointer\n");
78  return (int) AL_FAILURE;
79  }
80 
81  *rank = sz_stack[isz]->rank;
82 
83  return (int) AL_SUCCESS;
84 }
85 
86 /* ********************************************************************* */
87 int AL_Get_comm(int isz, MPI_Comm *comm)
88 /*!
89  * Get the communicator for a distributed array
90  *
91  * \param [in] isz Integer pointer to the input array descriptor
92  * \param [out] comm Pointer to communicator
93  *********************************************************************** */
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 }
114 
115 
116 /* ********************************************************************* */
117 int AL_Get_cart_comm(int isz, MPI_Comm *cart_comm)
118 /*!
119  * Get the cartesian communicator for a distributed array
120  *
121  * \param [in] isz Integer pointer to the input array descriptor
122  * \param [out] cart_comm Pointer to cartesian communicator
123  *********************************************************************** */
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 }
144 
145 
146 
147 /* ********************************************************************* */
148 int AL_Get_dimensions(int isz, int *ndim)
149 /*!
150  * Get the dimensions of the distributed array
151  *
152  * \param [in] isz Integer pointer to the input array descriptor
153  * \param [out] ndim Pointer to integer number of dimensions
154  *********************************************************************** */
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 }
175 
176 /* ********************************************************************* */
177 int AL_Get_type(int isz, AL_Datatype *type)
178 /*!
179  * Get the basic data type of a distributed array
180  *
181  * \param [in] isz Integer pointer to the input array descriptor
182  * \param [out] type Pointer to datatype
183  *********************************************************************** */
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 }
204 
205 /* ********************************************************************* */
206 int AL_Get_type_size(int isz, int *tsize)
207 /*!
208  * Get size of the basic data type of a distributed array
209  *
210  * \param [in] isz Integer pointer to the input array descriptor
211  * \param [out] tsize Integer pointer to size of datatype.
212  *********************************************************************** */
213 {
214  int itsize;
215  AL_Datatype type;
216  SZ *s;
217 
218  /*
219  Check that isz points to an allocated SZ
220  */
221  if( stack_ptr[isz] == AL_STACK_FREE ){
222  printf("AL_Get_type_size: wrong SZ pointer\n");
223  return (int) AL_FAILURE;
224  }
225 
226  /*
227  Get the SZ structure isz is pointing at
228  */
229  s = sz_stack[isz];
230 
231  type = s->type;
232 
233  MPI_Type_size( type, &itsize);
234 
235  *tsize = itsize;
236 
237  return (int) AL_SUCCESS;
238 }
239 
240 
241 /* ********************************************************************* */
242 int AL_Get_buffsize(int sz_ptr, int *buffsize)
243 /*!
244  * Get the local buffer size for a distributed array, including
245  * the number of elements of the ghsot regions.
246  *
247  * \param [in] sz_ptr Integer pointer to the input array descriptor
248  * \param [out] buffsize Buffer size [in number of elements]
249  *********************************************************************** */
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 }
273 
274 
275 /* ********************************************************************* */
276 int AL_Get_global_dim(int isz, int *gdims)
277 /*!
278  * Get the global dimensions of the distributed array
279  *
280  * \param [in] isz Integer pointer to the input array descriptor
281  * \param [out] gdims Array of integers with global dimensions
282  *********************************************************************** */
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 }
306 
307 
308 /* ********************************************************************* */
309 int AL_Get_local_dim(int isz, int *ldims)
310 /*!
311  * Get the local dimensions of a distributed array WITHOUT GHOST POINTS
312  *
313  * \param [in] isz Integer pointer to the input array descriptor
314  * \param [out] ldims Array of integers with local dimensions
315  *********************************************************************** */
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 }
339 
340 /* ********************************************************************* */
341 int AL_Get_local_dim_gp(int sz_ptr, int *ldims_gp)
342 /*!
343  * Get the local dimensions of a distributed array WITH GHOST POINTS
344  *
345  * \param [in] sz_ptr Integer pointer to the input array descriptor
346  * \param [out] ldims_gp Array of integers with local dimensions
347  * (including ghost points)
348  *********************************************************************** */
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 }
375 
376 
377 /* ********************************************************************* */
378 int AL_Get_offsets(int sz_ptr, int *offset)
379 /*!
380  * Get the offsets for the for loop index calculation
381  *
382  * NOTE:
383  * The typical multidimensional C for loop will look like, for example:
384  *
385  * ... ;
386  *
387  * AL_Get_offsets(sz_ptr, offset);
388  *
389  * for(k=kbeg; k<=kend; k++){
390  * koff = k*offset[2];
391  * for(j=jbeg; j<=jend; j++){
392  * joff = j*offset[1];
393  * for(i=ibeg; i<=iend; i++){
394  * ioff = i*offset[0]+joff+koff;
395  * a[ioff] = ... ;
396  * }
397  * }
398  * }
399  * ... ;
400  *
401  * Where a[ioff] is a generic 3D array.
402  *
403  * \param [in] sz_ptr Integer pointer to the input array descriptor
404  * \param [out] offset Array of integers with offsets
405  *********************************************************************** */
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 }
433 
434 
435 
436 /* ********************************************************************* */
437 int AL_Get_stride(int sz_ptr, int *stride)
438 /*!
439  * Get the strides for the for loop index calculation
440  * NOTE:
441  * The typical multidimensional C for loop will look like, for example:
442  *
443  * ... ;
444  *
445  * AL_Get_stride(sz_ptr, stride);
446  *
447  * for(k=kbeg; k<=kend; k++){
448  * koff = k*stride[2];
449  * for(j=jbeg; j<=jend; j++){
450  * joff = j*stride[1];
451  * for(i=ibeg; i<=iend; i++){
452  * ioff = i*stride[0]+joff+koff;
453  * a[ioff] = ... ;
454  * }
455  * }
456  * }
457  * ... ;
458  *
459  * Where a[ioff] is a generic 3D array.
460  *
461  * \param [in] sz_ptr Integer pointer to the input array descriptor
462  * \param [out] stride Array of strides
463  *********************************************************************** */
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 }
491 
492 
493 /* ********************************************************************* */
494 int AL_Get_cart_sizes(int isz, int *cart_sizes)
495 /*!
496  * Get the sizes of the cartesian communicator
497  *
498  * \param [in] isz Integer pointer to the input array descriptor
499  * \param [out] cart_sizes Array of integers with cartesian sizes
500  *********************************************************************** */
501 {
502  register int i;
503  int ndim;
504  SZ *s;
505 
506  /*
507  Check that isz points to an allocated SZ
508  */
509  if( stack_ptr[isz] == AL_STACK_FREE ){
510  printf("AL_Get_cart_sizes: wrong SZ pointer\n");
511  return (int) AL_FAILURE;
512  }
513 
514  /*
515  Get the SZ structure isz is pointing at
516  */
517  s = sz_stack[isz];
518  ndim = s->ndim;
519  for(i=0;i<ndim;i++){
520  cart_sizes[i] = s->lsize[i];
521  }
522 
523  return (int) AL_SUCCESS;
524 }
525 
526 
527 /* ********************************************************************* */
528 int AL_Get_parallel_dim(int isz, int *pardims)
529 /*!
530  * Get the parallel dimensions of a distributed array
531  *
532  * \param [in] isz Integer pointer to the input array descriptor
533  * \param [out] pardims Array of integers with parallel dimensions
534  * [AL_TRUE|AL_FALSE]
535  *********************************************************************** */
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 }
560 
561 
562 /* ********************************************************************* */
563 int AL_Get_periodic_dim(int isz, int *periods)
564 /*!
565  * Get the periodic dimensions of the distributed arrays
566  *
567  * \param [in] isz Integer pointer to the input array descriptor
568  * \param [out] periods Array of integers with periodic dimensions
569  * [AL_TRUE|AL_FALSE]
570  *********************************************************************** */
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 }
595 
596 
597 /* ********************************************************************* */
598 int AL_Get_staggered_dim(int isz, int *stagger)
599 /*!
600  * Get the staggered dimensions of the distributed array
601  *
602  * \param [in] isz Integer pointer to the input array descriptor
603  * \param [out] stagger Array of integers with staggered dimensions
604  * [AL_TRUE|AL_FALSE]
605  *********************************************************************** */
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 }
630 
631 
632 /* ********************************************************************* */
633 int AL_Get_ghosts(int isz, int *ghosts)
634 /*!
635  * Get the ghost points of the distributed array
636  *
637  * \param [in] isz Integer pointer to the input array descriptor
638  * \param [out] ghost Array of integers with size of ghost points
639  * in each dimension
640  *********************************************************************** */
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 }
665 
666 
667 /* ********************************************************************* */
668 int AL_Get_lbounds(int isz, int *lbeg, int *lend, int *gp, int style)
669 /*!
670  * Get the local indexes for the local portion of a distributed array
671  *
672  * \param [in] isz Integer pointer to the input array descriptor
673  * \param [out] lbeg Array of ndim integers containing the start points
674  * \param [out] lend Array of ndim integers containing the end points
675  * \param [out] gp Array of ndim integers containing the ghost points
676  * \param [in] style Index style: AL_C_INDEXES or AL_FORTRAN_INDEXES
677  *********************************************************************** */
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 }
704 
705 
706 /* ********************************************************************* */
707 int AL_Get_bounds(int isz, int *beg, int *end, int *gp, int style)
708 /*!
709  * Get the global indexes for the local portion of a distributed array
710  *
711  * \param [in] isz Integer pointer to the input array descriptor
712  * \param [out] beg Array of ndim integers containing the start points
713  * \param [out] end Array of ndim integers containing the end points
714  * \param [out] gp Array of ndim integers containing the ghost points
715  * \param [in] style Index style: AL_C_INDEXES or AL_FORTRAN_INDEXES
716  *********************************************************************** */
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 }
743 
744 
745 /* ********************************************************************* */
746 int AL_Get_gbounds(int isz, int *gbeg, int *gend, int *gp, int style)
747 /*!
748  * Get the global bounds of the global array
749  *
750  * \param [in] isz Integer pointer to the input array descriptor
751  * \param [out] gbeg Array of ndim integers containing the start points
752  * \param [out] gend Array of ndim integers containing the end points
753  * \param [out] gp Array of ndim integers containing the ghost points
754  * \param [in] style Index style: AL_C_INDEXES or AL_FORTRAN_INDEXES
755  *********************************************************************** */
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 }
int isperiodic[AL_MAX_DIM]
Definition: al_hidden.h:80
int AL_Get_lbounds(int isz, int *lbeg, int *lend, int *gp, int style)
Definition: al_sz_get.c:668
int end[AL_MAX_DIM]
Definition: al_hidden.h:60
int AL_Get_periodic_dim(int isz, int *periods)
Definition: al_sz_get.c:563
int isparallel[AL_MAX_DIM]
Definition: al_hidden.h:78
#define AL_Datatype
Definition: al_defs.h:22
int AL_Get_offsets(int sz_ptr, int *offset)
Definition: al_sz_get.c:378
int offset[AL_MAX_DIM]
Definition: al_hidden.h:76
int AL_Get_size(int isz, int *size)
Definition: al_sz_get.c:23
int AL_Get_comm(int isz, MPI_Comm *comm)
Definition: al_sz_get.c:87
int AL_Get_type(int isz, AL_Datatype *type)
Definition: al_sz_get.c:177
#define AL_SUCCESS
Definition: al_codes.h:32
#define AL_STACK_FREE
Definition: al_codes.h:24
int AL_Get_cart_comm(int isz, MPI_Comm *cart_comm)
Definition: al_sz_get.c:117
int stack_ptr[AL_MAX_ARRAYS]
Definition: al_szptr_.c:26
MPI_Comm cart_comm
Definition: al_hidden.h:49
int arrdim[AL_MAX_DIM]
Definition: al_hidden.h:53
int AL_Get_local_dim(int isz, int *ldims)
Definition: al_sz_get.c:309
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 AL_Get_type_size(int isz, int *tsize)
Definition: al_sz_get.c:206
int stride[AL_MAX_DIM]
Definition: al_hidden.h:77
int AL_Get_stride(int sz_ptr, int *stride)
Definition: al_sz_get.c:437
int bg[AL_MAX_DIM]
Definition: al_hidden.h:72
int AL_Get_cart_sizes(int isz, int *cart_sizes)
Definition: al_sz_get.c:494
int beg[AL_MAX_DIM]
Definition: al_hidden.h:58
int isstaggered[AL_MAX_DIM]
Definition: al_hidden.h:82
int lsize[AL_MAX_DIM]
Definition: al_hidden.h:90
#define AL_FAILURE
Definition: al_codes.h:33
#define s
int AL_Get_gbounds(int isz, int *gbeg, int *gend, int *gp, int style)
Definition: al_sz_get.c:746
int AL_Get_dimensions(int isz, int *ndim)
Definition: al_sz_get.c:148
SZ * sz_stack[AL_MAX_ARRAYS]
Definition: al_szptr_.c:18
int AL_Get_parallel_dim(int isz, int *pardims)
Definition: al_sz_get.c:528
int i
Definition: analysis.c:2
MPI_Comm comm
Definition: al_hidden.h:47
long int buffsize
Definition: al_hidden.h:52
int AL_Get_bounds(int isz, int *beg, int *end, int *gp, int style)
Definition: al_sz_get.c:707
MPI_Datatype type
Definition: al_hidden.h:40
int AL_Get_global_dim(int isz, int *gdims)
Definition: al_sz_get.c:276
int AL_Get_ghosts(int isz, int *ghosts)
Definition: al_sz_get.c:633
int AL_Get_buffsize(int sz_ptr, int *buffsize)
Definition: al_sz_get.c:242
int AL_Get_local_dim_gp(int sz_ptr, int *ldims_gp)
Definition: al_sz_get.c:341
Definition: al_hidden.h:38
int AL_Get_rank(int isz, int *rank)
Definition: al_sz_get.c:59
Internal include file for the ArrayLib.
int lbeg[AL_MAX_DIM]
Definition: al_hidden.h:68
#define AL_MAX_ARRAYS
Definition: al_codes.h:21
int AL_Get_staggered_dim(int isz, int *stagger)
Definition: al_sz_get.c:598
int lend[AL_MAX_DIM]
Definition: al_hidden.h:70
int ndim
Definition: al_hidden.h:43