HyPar  1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
MPIPartitionArraynD.c File Reference

Partition a global n-dimensional array. More...

#include <stdio.h>
#include <stdlib.h>
#include <basic.h>
#include <arrayfunctions.h>
#include <mpivars.h>

Go to the source code of this file.

Functions

int MPIPartitionArraynD (int ndims, void *m, double *xg, double *x, int *dim_global, int *dim_local, int ghosts, int nvars)
 
int MPIPartitionArraynDwGhosts (int ndims, void *m, double *xg, double *x, int *dim_global, int *dim_local, int ghosts, int nvars)
 

Detailed Description

Partition a global n-dimensional array.

Author
Debojyoti Ghosh

Definition in file MPIPartitionArraynD.c.

Function Documentation

◆ MPIPartitionArraynD()

int MPIPartitionArraynD ( int  ndims,
void *  m,
double *  xg,
double *  x,
int *  dim_global,
int *  dim_local,
int  ghosts,
int  nvars 
)

Partitions the contents of a global n-dimensional array on rank 0 (root) to local n-dimensional arrays on all the MPI ranks. See documentation of MPIExchangeBoundariesnD() for how the n-dimensional array is stored in the memory as a single-index array.

Notes:

  • The global array must have no ghost points.
  • The global array must be allocated only on rank 0. On other ranks, it must be NULL.
Parameters
ndimsNumber of spatial dimensions
mMPI object of type MPIVariables
xgGlobal array (preallocated) without ghost points
xLocal array
dim_globalInteger array with elements as global size along each spatial dimension
dim_localInteger array with elements as local size along each spatial dimension
ghostsNumber of ghost points
nvarsNumber of variables (vector components)

Definition at line 21 of file MPIPartitionArraynD.c.

31 {
32  MPIVariables *mpi = (MPIVariables*) m;
34 
35  int is[ndims], ie[ndims], index[ndims], bounds[ndims];
36 
37  /* xg should be non-null only on root */
38  if (mpi->rank && xg) {
39  fprintf(stderr,"Error in MPIPartitionArraynD(): global array exists on non-root processors (rank %d).\n",
40  mpi->rank);
41  return(1);
42  }
43  if ((!mpi->rank) && (!xg)) {
44  fprintf(stderr,"Error in MPIPartitionArraynD(): global array is not allocated on root processor.\n");
45  return(1);
46  }
47 
48  if (!mpi->rank) {
49  int proc;
50  for (proc = 0; proc < mpi->nproc; proc++) {
51  int d,done,size;
52  /* Find out the domain limits for each process */
53  IERR MPILocalDomainLimits(ndims,proc,mpi,dim_global,is,ie); CHECKERR(ierr);
54  size = 1;
55  for (d=0; d<ndims; d++) {
56  size *= (ie[d]-is[d]);
57  bounds[d] = ie[d] - is[d];
58  }
59  double *buffer = (double*) calloc (size*nvars, sizeof(double));
60  done = 0; _ArraySetValue_(index,ndims,0);
61  while (!done) {
62  int p1; _ArrayIndex1DWO_(ndims,dim_global,index,is,0,p1);
63  int p2; _ArrayIndex1D_(ndims,bounds,index,0,p2);
64  _ArrayCopy1D_((xg+nvars*p1),(buffer+nvars*p2),nvars);
65  _ArrayIncrementIndex_(ndims,bounds,index,done);
66  }
67  if (proc) {
68 #ifndef serial
69  MPI_Send(buffer,size*nvars,MPI_DOUBLE,proc,1538,mpi->world);
70 #endif
71  } else {
72  done = 0; _ArraySetValue_(index,ndims,0);
73  while (!done) {
74  int p1; _ArrayIndex1D_(ndims,dim_local,index,ghosts,p1);
75  int p2; _ArrayIndex1D_(ndims,dim_local,index,0,p2);
76  _ArrayCopy1D_((buffer+nvars*p2),(x+nvars*p1),nvars);
77  _ArrayIncrementIndex_(ndims,dim_local,index,done);
78  }
79  }
80  free(buffer);
81  }
82 
83  } else {
84 #ifndef serial
85  /* Meanwhile, on other processes */
86  MPI_Status status;
87  int d, done, size;
88  size = 1; for (d=0; d<ndims; d++) size *= dim_local[d];
89  double *buffer = (double*) calloc (size*nvars, sizeof(double));
90  MPI_Recv(buffer,size*nvars,MPI_DOUBLE,0,1538,mpi->world,&status);
91  done = 0; _ArraySetValue_(index,ndims,0);
92  while (!done) {
93  int p1; _ArrayIndex1D_(ndims,dim_local,index,ghosts,p1);
94  int p2; _ArrayIndex1D_(ndims,dim_local,index,0,p2);
95  _ArrayCopy1D_((buffer+nvars*p2),(x+nvars*p1),nvars);
96  _ArrayIncrementIndex_(ndims,dim_local,index,done);
97  }
98  free(buffer);
99 #endif
100  }
101  return(0);
102 }
#define IERR
Definition: basic.h:16
#define CHECKERR(ierr)
Definition: basic.h:18
MPI_Comm world
#define _ArrayIndex1D_(N, imax, i, ghost, index)
#define _ArrayIndex1DWO_(N, imax, i, offset, ghost, index)
#define _ArraySetValue_(x, size, value)
#define _ArrayIncrementIndex_(N, imax, i, done)
int MPILocalDomainLimits(int, int, void *, int *, int *, int *)
Structure of MPI-related variables.
#define _DECLARE_IERR_
Definition: basic.h:17
#define _ArrayCopy1D_(x, y, size)

◆ MPIPartitionArraynDwGhosts()

int MPIPartitionArraynDwGhosts ( int  ndims,
void *  m,
double *  xg,
double *  x,
int *  dim_global,
int *  dim_local,
int  ghosts,
int  nvars 
)

Partitions the contents of a global n-dimensional array on rank 0 (root) to local n-dimensional arrays on all the MPI ranks. See documentation of MPIExchangeBoundariesnD() for how the n-dimensional array is stored in the memory as a single-index array. This is same as MPIPartitionArraynD() but where the global array has ghost points.

Notes:

  • The global array must have ghost points (the same number as the local arrays).
  • The global array must be allocated only on rank 0. On other ranks, it must be NULL.
Parameters
ndimsNumber of spatial dimensions
mMPI object of type MPIVariables
xgGlobal array (preallocated) without ghost points
xLocal array
dim_globalInteger array with elements as global size along each spatial dimension
dim_localInteger array with elements as local size along each spatial dimension
ghostsNumber of ghost points
nvarsNumber of variables (vector components)

Definition at line 114 of file MPIPartitionArraynD.c.

124 {
125  MPIVariables *mpi = (MPIVariables*) m;
126  int d;
128 
129  int is[ndims], ie[ndims], index[ndims], bounds[ndims];
130  int dim_global_wghosts[ndims];
131  for (d = 0; d < ndims; d++) dim_global_wghosts[d] = dim_global[d] + 2*ghosts;
132 
133  /* xg should be non-null only on root */
134  if (mpi->rank && xg) {
135  fprintf(stderr,"Error in MPIPartitionArraynD(): global array exists on non-root processors (rank %d).\n",
136  mpi->rank);
137  return(1);
138  }
139  if ((!mpi->rank) && (!xg)) {
140  fprintf(stderr,"Error in MPIPartitionArraynD(): global array is not allocated on root processor.\n");
141  return(1);
142  }
143 
144  if (!mpi->rank) {
145  int proc;
146  for (proc = 0; proc < mpi->nproc; proc++) {
147  int done,size;
148  /* Find out the domain limits for each process */
149  IERR MPILocalDomainLimits(ndims,proc,mpi,dim_global,is,ie); CHECKERR(ierr);
150  size = 1;
151  for (d=0; d<ndims; d++) {
152  size *= (ie[d]-is[d]+2*ghosts);
153  bounds[d] = ie[d] - is[d] + 2*ghosts;
154  }
155  double *buffer = (double*) calloc (size*nvars, sizeof(double));
156  done = 0; _ArraySetValue_(index,ndims,0);
157  while (!done) {
158  int p1; _ArrayIndex1DWO_(ndims,dim_global_wghosts,index,is,0,p1);
159  int p2; _ArrayIndex1D_(ndims,bounds,index,0,p2);
160  _ArrayCopy1D_((xg+nvars*p1),(buffer+nvars*p2),nvars);
161  _ArrayIncrementIndex_(ndims,bounds,index,done);
162  }
163  if (proc) {
164 #ifndef serial
165  MPI_Send(buffer,size*nvars,MPI_DOUBLE,proc,1538,mpi->world);
166 #endif
167  } else {
168  done = 0; _ArraySetValue_(index,ndims,0);
169  _ArrayCopy1D_(buffer, x, size*nvars);
170  }
171  free(buffer);
172  }
173 
174  } else {
175 #ifndef serial
176  /* Meanwhile, on other processes */
177  MPI_Status status;
178  int done, size;
179  size = 1; for (d=0; d<ndims; d++) size *= (dim_local[d]+2*ghosts);
180  double *buffer = (double*) calloc (size*nvars, sizeof(double));
181  MPI_Recv(buffer,size*nvars,MPI_DOUBLE,0,1538,mpi->world,&status);
182  _ArrayCopy1D_(buffer, x, size*nvars);
183  free(buffer);
184 #endif
185  }
186  return(0);
187 }
#define IERR
Definition: basic.h:16
#define CHECKERR(ierr)
Definition: basic.h:18
MPI_Comm world
#define _ArrayIndex1D_(N, imax, i, ghost, index)
#define _ArrayIndex1DWO_(N, imax, i, offset, ghost, index)
#define _ArraySetValue_(x, size, value)
#define _ArrayIncrementIndex_(N, imax, i, done)
int MPILocalDomainLimits(int, int, void *, int *, int *, int *)
Structure of MPI-related variables.
#define _DECLARE_IERR_
Definition: basic.h:17
#define _ArrayCopy1D_(x, y, size)