| [276a7df] | 1 | /* diffusion2d.c: parallel 2d-diffusion equation solver with constant boundaries
|
|---|
| 2 | * slicing matrix as a checker board.
|
|---|
| 3 | * To execute: mpicc diffusion2d.c ; mpiexec -n 4 ./a.out 2 2
|
|---|
| 4 | * To verify: civl verify diffusion2d.c
|
|---|
| 5 | */
|
|---|
| 6 | #include<stdio.h>
|
|---|
| 7 | #include<stdlib.h>
|
|---|
| 8 | #include<assert.h>
|
|---|
| 9 | #include<string.h>
|
|---|
| 10 | #include<mpi.h>
|
|---|
| 11 |
|
|---|
| 12 | /* Message tags */
|
|---|
| 13 | #define FROMLEFT 0
|
|---|
| 14 | #define FROMRIGHT 1
|
|---|
| 15 | #define FROMTOP 2
|
|---|
| 16 | #define FROMBOTTOM 3
|
|---|
| 17 | #define DATAPASS 4
|
|---|
| 18 | #define comm MPI_COMM_WORLD
|
|---|
| 19 |
|
|---|
| 20 | #ifdef _CIVL
|
|---|
| [6f8bc6ea] | 21 | #include <civlc.cvh>
|
|---|
| 22 |
|
|---|
| [276a7df] | 23 | $input long NXB = 5; // nx upper bound
|
|---|
| 24 | $input long nx; // global number of columns in matrix
|
|---|
| [3ff27cf] | 25 | $assume(1 <= nx && nx <= NXB);
|
|---|
| [276a7df] | 26 | $input long NYB = 5; // ny upper bound
|
|---|
| 27 | $input long ny; // global number of rows of matrix
|
|---|
| [3ff27cf] | 28 | $assume(1 <= ny && ny <= NYB);
|
|---|
| [276a7df] | 29 | $input double u_init[ny+2][nx+2]; // initial value of temperatures, including boundaries
|
|---|
| 30 | $input double k; // constant coefficient
|
|---|
| [3ff27cf] | 31 | $assume(k > 0.0 && k < 0.5);
|
|---|
| [276a7df] | 32 | $input int NSTEPSB = 3; // upper bound for nsteps
|
|---|
| 33 | $input int nsteps; // number of steps
|
|---|
| [3ff27cf] | 34 | $assume(1<=nsteps && nsteps<=NSTEPSB);
|
|---|
| [276a7df] | 35 | $input int wstep = 1; // write frame every this many time steps
|
|---|
| 36 | double oracle[nsteps][ny+2][nx+2]; // solution computed sequentially, done by proc 0 only
|
|---|
| 37 | $input int NPROCSXB = 2; // upper bound for NPROCSX
|
|---|
| 38 | $input int NPROCSX = 2; // number of procs in x direction
|
|---|
| [3ff27cf] | 39 | $assume(NPROCSX >= 1 && NPROCSX <= NPROCSXB);
|
|---|
| [276a7df] | 40 | $input int NPROCSYB = 2; // upper bound for NPROCSY
|
|---|
| 41 | $input int NPROCSY = 2; // number of procs in y direction
|
|---|
| [3ff27cf] | 42 | $assume(NPROCSY >= 1 && NPROCSY <= NPROCSYB);
|
|---|
| [024a9eb] | 43 | $input int _mpi_nprocs = NPROCSX * NPROCSY;
|
|---|
| 44 | $assume(_mpi_nprocs == NPROCSX * NPROCSY);
|
|---|
| [276a7df] | 45 | #else
|
|---|
| 46 | long nx, ny;
|
|---|
| 47 | int nsteps, wstep;
|
|---|
| 48 | int NPROCSX, NPROCSY;
|
|---|
| 49 | double constTemp; // value of constant boundaries for test
|
|---|
| 50 | double initTemp; // value of initial temperature for test
|
|---|
| 51 | double k;
|
|---|
| 52 | #endif
|
|---|
| 53 |
|
|---|
| 54 | /* Global variables */
|
|---|
| 55 | int nprocs, rank;
|
|---|
| 56 | int left, right, top, bottom;// ranks of neighbors
|
|---|
| 57 | int nxl; // local dimension of x
|
|---|
| 58 | int nyl; // local dimension of y
|
|---|
| 59 | long firstCol; // index of the first column in global array
|
|---|
| 60 | long firstRow; // index of the first row in global array
|
|---|
| 61 | /* dynamically-allocated 2d array of doubles specifying state of local
|
|---|
| 62 | * grid in current time step. Dimensions are u_curr[nyl+2][nxl+2].
|
|---|
| 63 | */
|
|---|
| 64 | double ** u_curr;
|
|---|
| 65 | /* dynamically-allocated 2d array of doubles specifying state of local
|
|---|
| 66 | * grid in next time step. Dimensions are u_next[nyl+2][nxl+2].
|
|---|
| 67 | */
|
|---|
| 68 | double ** u_next;
|
|---|
| 69 | /* dynamically-allocated 1d temporary buffer. Length is recvbuf[nxl+1]
|
|---|
| 70 | */
|
|---|
| 71 | double * recvbuf;
|
|---|
| 72 |
|
|---|
| 73 | /* The processes are arranged geometrically as follows for the case
|
|---|
| 74 | * NPROCSX = 3:
|
|---|
| 75 | * row 0: 0 1 2
|
|---|
| 76 | * row 1: 3 4 5
|
|---|
| 77 | * ...
|
|---|
| 78 | * This function computes the global index of the first column
|
|---|
| 79 | * owned by the process of the given rank. rank must be in the
|
|---|
| [024a9eb] | 80 | * range [0, _mpi_nprocs-1]. The result will in the range [0, nx-1].
|
|---|
| [276a7df] | 81 | */
|
|---|
| 82 | long firstColForProc(int rank) {
|
|---|
| 83 | long tmp = (rank - (rank / NPROCSX)*NPROCSX)*nx;
|
|---|
| 84 |
|
|---|
| 85 | return tmp/NPROCSX;
|
|---|
| 86 | }
|
|---|
| 87 |
|
|---|
| 88 | /* This function computes the global index of the first row owned by
|
|---|
| 89 | the process of the given rank. rank must be in the range[0,
|
|---|
| [024a9eb] | 90 | _mpi_nprocs-1]. The result will in the range [0, ny-1]. */
|
|---|
| [276a7df] | 91 | long firstRowForProc(int rank) {
|
|---|
| 92 | long tmp = ((rank / NPROCSX)*ny);
|
|---|
| 93 |
|
|---|
| 94 | return tmp/NPROCSY;
|
|---|
| 95 | }
|
|---|
| 96 |
|
|---|
| 97 | /* Computes the number of columns owned by the process. The result
|
|---|
| 98 | will be in the range [0, nx]. */
|
|---|
| 99 | int countColForProc(int rank) {
|
|---|
| 100 | long a = firstColForProc(rank);
|
|---|
| 101 | long b;
|
|---|
| 102 |
|
|---|
| 103 | if ((rank / NPROCSX) == ((rank+1) / NPROCSX))
|
|---|
| 104 | b = firstColForProc(rank+1);
|
|---|
| 105 | else
|
|---|
| 106 | b = nx;
|
|---|
| 107 | return b - a;
|
|---|
| 108 | }
|
|---|
| 109 |
|
|---|
| 110 | /* Computes the number of rows owned by the process. The result
|
|---|
| 111 | will be in the range [0, ny]. */
|
|---|
| 112 | int countRowForProc(int rank) {
|
|---|
| 113 | long a = firstRowForProc(rank);
|
|---|
| 114 | long b = firstRowForProc(rank+NPROCSX);
|
|---|
| 115 |
|
|---|
| 116 | return b - a;
|
|---|
| 117 | }
|
|---|
| 118 |
|
|---|
| 119 | /* Get the owner process of the given cell */
|
|---|
| 120 | int OWNER(long col, long row) {
|
|---|
| 121 | long tmp = ((NPROCSY * (row+1))-1);
|
|---|
| 122 | int procRow = tmp / ny;
|
|---|
| 123 | int procCol;
|
|---|
| 124 |
|
|---|
| 125 | tmp = ((col + 1) * NPROCSX - 1);
|
|---|
| 126 | procCol = tmp / nx;
|
|---|
| 127 | tmp = procRow * NPROCSX;
|
|---|
| 128 | return tmp + procCol;
|
|---|
| 129 | }
|
|---|
| 130 |
|
|---|
| 131 |
|
|---|
| 132 | /* initialize all data values owned by each process */
|
|---|
| 133 | void initData() {
|
|---|
| 134 | #ifdef _CIVL
|
|---|
| 135 | // Data is initialized with totally unconstrained values
|
|---|
| 136 | // for verification.
|
|---|
| 137 | // set the vertical boundary cells:
|
|---|
| 138 | if (left == MPI_PROC_NULL)
|
|---|
| 139 | for (int i=0; i<nyl+2; i++)
|
|---|
| 140 | u_next[i][0] = u_curr[i][0] = u_init[i + firstRow][0];
|
|---|
| 141 | if (right == MPI_PROC_NULL)
|
|---|
| 142 | for (int i=0; i<nyl+2; i++)
|
|---|
| 143 | u_next[i][nxl+1] = u_curr[i][nxl+1] = u_init[i + firstRow][nx+1];
|
|---|
| 144 | // set the horizontal boundary cells:
|
|---|
| 145 | if (top == MPI_PROC_NULL)
|
|---|
| 146 | for (int i=0; i<nxl+2; i++)
|
|---|
| 147 | u_next[0][i] = u_curr[0][i] = u_init[0][i + firstCol];
|
|---|
| 148 | if (bottom == MPI_PROC_NULL)
|
|---|
| 149 | for (int i=0; i<nxl+2; i++)
|
|---|
| 150 | u_next[nyl+1][i] = u_curr[nyl+1][i] = u_init[ny+1][i + firstCol];
|
|---|
| 151 | for (int i=1; i<nyl+1; i++)
|
|---|
| 152 | memcpy(&u_curr[i][1], &u_init[firstRow + i][firstCol + 1], nxl * sizeof(double));
|
|---|
| 153 | #else
|
|---|
| 154 | // All boundary cells are set to the same value constTemp.
|
|---|
| 155 | // All cells in the interior region are set to the same value
|
|---|
| 156 | // initTemp.
|
|---|
| 157 | for (int i=0; i < nyl+2; i++)
|
|---|
| 158 | for (int j=0; j < nxl+2; j++)
|
|---|
| 159 | if (i == 0 || j == 0 || i == nyl+1 || j == nxl+1)
|
|---|
| 160 | u_next[i][j] = u_curr[i][j] = constTemp;
|
|---|
| 161 | else
|
|---|
| 162 | u_curr[i][j] = initTemp;
|
|---|
| 163 | #endif
|
|---|
| 164 | }
|
|---|
| 165 |
|
|---|
| 166 | /* Initialize all global variables, allocate memory spaces for
|
|---|
| 167 | * pointers and proc 0 will do a sequential run. The results of the
|
|---|
| 168 | * sequential run will be used to compare with parallel run later. */
|
|---|
| 169 | void initialization(int argc, char * argv[]) {
|
|---|
| 170 | #ifndef _CIVL
|
|---|
| 171 | nsteps = 300;
|
|---|
| 172 | wstep = 5;
|
|---|
| 173 | nx = 15;
|
|---|
| 174 | ny = 15;
|
|---|
| 175 | if (argc != 3) {
|
|---|
| 176 | printf("Usage: mpiexec -n NPROCS diffusion2d NPROCSX NPROCSY\n"
|
|---|
| 177 | " NPROCSX: number of processes in x direction\n"
|
|---|
| 178 | " NPROCSY: number of processes in y direction\n"
|
|---|
| 179 | " NPROCS: the product of NPROCSX and NPROCSY\n");
|
|---|
| 180 | exit(1);
|
|---|
| 181 | }
|
|---|
| 182 | NPROCSX = atoi(argv[1]);
|
|---|
| 183 | NPROCSY = atoi(argv[2]);
|
|---|
| 184 | assert(NPROCSX * NPROCSY == nprocs);
|
|---|
| 185 | constTemp = 0.0;
|
|---|
| 186 | initTemp = 100.0;
|
|---|
| 187 | k = 0.13;
|
|---|
| 188 | #endif
|
|---|
| 189 | printf("Diffusion2d with k=%f, nx=%ld, ny=%ld, nsteps=%d, wstep=%d\n",
|
|---|
| 190 | k, nx, ny, nsteps, wstep);
|
|---|
| 191 | nxl = countColForProc(rank);
|
|---|
| 192 | nyl = countRowForProc(rank);
|
|---|
| 193 | if (rank == 0)
|
|---|
| 194 | recvbuf = (double *)malloc((nxl + 1) * sizeof(double));
|
|---|
| 195 | u_curr = (double **)malloc((nyl + 2) * sizeof(double *));
|
|---|
| 196 | assert(u_curr);
|
|---|
| 197 | u_next = (double **)malloc((nyl + 2) * sizeof(double *));
|
|---|
| 198 | assert(u_next);
|
|---|
| 199 | for (int i=0; i < nyl + 2; i++){
|
|---|
| 200 | u_curr[i] = (double *)malloc((nxl + 2) * sizeof(double));
|
|---|
| 201 | assert(u_curr[i]);
|
|---|
| 202 | u_next[i] = (double *)malloc((nxl + 2) * sizeof(double));
|
|---|
| 203 | assert(u_next[i]);
|
|---|
| 204 | }
|
|---|
| 205 | firstCol = firstColForProc(rank);
|
|---|
| 206 | firstRow = firstRowForProc(rank);
|
|---|
| 207 | // computes neighbors
|
|---|
| 208 | if (firstCol != 0)
|
|---|
| 209 | left = OWNER(firstCol - 1, firstRow);
|
|---|
| 210 | else
|
|---|
| 211 | left = MPI_PROC_NULL;
|
|---|
| 212 | if (firstRow != 0)
|
|---|
| 213 | top = OWNER(firstCol, firstRow - 1);
|
|---|
| 214 | else
|
|---|
| 215 | top = MPI_PROC_NULL;
|
|---|
| 216 | if (firstCol + nxl < nx)
|
|---|
| 217 | right = OWNER(firstCol + nxl, firstRow);
|
|---|
| 218 | else
|
|---|
| 219 | right = MPI_PROC_NULL;
|
|---|
| 220 | if (firstRow + nyl < ny)
|
|---|
| 221 | bottom = OWNER(firstCol, firstRow + nyl);
|
|---|
| 222 | else
|
|---|
| 223 | bottom = MPI_PROC_NULL;
|
|---|
| 224 | #ifdef _CIVL
|
|---|
| 225 | /* In CIVL mode process with rank 0 will be responsible for
|
|---|
| 226 | * computing the diffusion2d equation sequentially such that the
|
|---|
| 227 | * results can be used to compare with the ones of parallel run. */
|
|---|
| 228 | if (rank == 0) {
|
|---|
| 229 | for (long i = 0; i < ny + 2; i++)
|
|---|
| 230 | for (long j = 0; j < nx + 2; j++)
|
|---|
| 231 | oracle[0][i][j] = u_init[i][j];
|
|---|
| 232 | for (int t=1; t < nsteps; t++)
|
|---|
| 233 | for (long i = 0; i < ny + 2; i++)
|
|---|
| 234 | for (long j = 0; j < nx + 2; j++)
|
|---|
| 235 | if (i==0 || j==0 || i == ny + 1 || j == nx + 1)
|
|---|
| 236 | oracle[t][i][j] = oracle[t-1][i][j];
|
|---|
| 237 | else
|
|---|
| 238 | oracle[t][i][j] = oracle[t-1][i][j] +
|
|---|
| 239 | k*(oracle[t-1][i+1][j] + oracle[t-1][i-1][j] +
|
|---|
| 240 | oracle[t-1][i][j+1] + oracle[t-1][i][j-1] -
|
|---|
| 241 | 4*oracle[t-1][i][j]);
|
|---|
| 242 | }
|
|---|
| 243 | #endif
|
|---|
| 244 | }
|
|---|
| 245 |
|
|---|
| 246 | /* Update local cells owned by process */
|
|---|
| 247 | void update() {
|
|---|
| 248 | double **tmp;
|
|---|
| 249 |
|
|---|
| 250 | for (int i = 1; i < nyl + 1; i++)
|
|---|
| 251 | for (int j = 1; j < nxl + 1; j++) {
|
|---|
| 252 | u_next[i][j] = u_curr[i][j] +
|
|---|
| 253 | k*(u_curr[i+1][j] + u_curr[i-1][j] +
|
|---|
| 254 | u_curr[i][j+1] + u_curr[i][j-1] - 4*u_curr[i][j]);
|
|---|
| 255 | }
|
|---|
| 256 | // swap two pointers
|
|---|
| 257 | tmp = u_curr;
|
|---|
| 258 | u_curr = u_next;
|
|---|
| 259 | u_next = tmp;
|
|---|
| 260 | }
|
|---|
| 261 |
|
|---|
| 262 | /* Exchange ghost cells between proceeses */
|
|---|
| 263 | void exchange() {
|
|---|
| 264 | double sendbuf[nyl];
|
|---|
| 265 | double recvbuf[nyl];
|
|---|
| 266 |
|
|---|
| 267 | // sends top border row, receives into bottom ghost cell row
|
|---|
| 268 | MPI_Sendrecv(&u_curr[1][1], nxl, MPI_DOUBLE, top, FROMBOTTOM, &u_curr[nyl+1][1], nxl,
|
|---|
| 269 | MPI_DOUBLE, bottom, FROMBOTTOM, comm, MPI_STATUS_IGNORE);
|
|---|
| 270 | // sends bottom border row, receives into top ghost cell row
|
|---|
| 271 | MPI_Sendrecv(&u_curr[nyl][1], nxl, MPI_DOUBLE, bottom, FROMTOP, &u_curr[0][1], nxl,
|
|---|
| 272 | MPI_DOUBLE, top, FROMTOP, comm, MPI_STATUS_IGNORE);
|
|---|
| 273 | // sends left border column, receives into temporary buffer
|
|---|
| 274 | for (int i = 0; i < nyl; i++) sendbuf[i] = u_curr[i+1][1];
|
|---|
| 275 | MPI_Sendrecv(sendbuf, nyl, MPI_DOUBLE, left, FROMRIGHT, recvbuf, nyl,
|
|---|
| 276 | MPI_DOUBLE, right, FROMRIGHT, comm, MPI_STATUS_IGNORE);
|
|---|
| 277 | // copies temporary buffer into right ghost cell column
|
|---|
| 278 | if (right != MPI_PROC_NULL)
|
|---|
| 279 | for (int i = 0; i < nyl; i++) u_curr[i+1][nxl+1] = recvbuf[i];
|
|---|
| 280 | // sends right border column, receives into temporary buffer
|
|---|
| 281 | for (int i = 0; i < nyl; i++) sendbuf[i] = u_curr[i+1][nxl];
|
|---|
| 282 | MPI_Sendrecv(sendbuf, nyl, MPI_DOUBLE, right, FROMLEFT, recvbuf, nyl,
|
|---|
| 283 | MPI_DOUBLE, left, FROMLEFT, comm, MPI_STATUS_IGNORE);
|
|---|
| 284 | // copies temporary buffer into left ghost cell column
|
|---|
| 285 | if (left != MPI_PROC_NULL)
|
|---|
| 286 | for (int i = 0; i < nyl; i++) u_curr[i+1][0] = recvbuf[i];
|
|---|
| 287 | }
|
|---|
| 288 |
|
|---|
| 289 | /* Helper function for write_frame(int). In CIVL mode, it takes the
|
|---|
| 290 | index of the first column, the number of columns owned by the
|
|---|
| 291 | process and the index of current row to locate the corresponding
|
|---|
| 292 | cell in oracle and compare it with the given cell's value computed
|
|---|
| 293 | in parallel */
|
|---|
| 294 | void printData(int time, int firstCol, int nxl, int currRow, double * buf) {
|
|---|
| 295 | for (int i=0; i<nxl; i++) {
|
|---|
| 296 | printf("%6.2f", *(buf + i));
|
|---|
| 297 | #ifdef _CIVL
|
|---|
| [3ff27cf] | 298 | $assert((*(buf + i) == oracle[time][currRow + 1][firstCol + i + 1]), \
|
|---|
| [276a7df] | 299 | "Error: disagreement at time %d position [%d][%d]: saw %lf, expected %lf", \
|
|---|
| 300 | time, currRow, firstCol + i,
|
|---|
| [3ff27cf] | 301 | *(buf + i), oracle[time][currRow + 1][firstCol + i + 1]);
|
|---|
| [276a7df] | 302 | #endif
|
|---|
| 303 | }
|
|---|
| 304 | }
|
|---|
| 305 |
|
|---|
| 306 | /* Print the computed matrix at the given time step all processes
|
|---|
| 307 | * should send their local data to process rank 0 which is responsible
|
|---|
| 308 | * for printing */
|
|---|
| 309 | void write_frame(int time) {
|
|---|
| 310 | // sends data row by row
|
|---|
| 311 | if (rank != 0) {
|
|---|
| 312 | for (int i=0; i<nyl; i++)
|
|---|
| 313 | MPI_Send(&u_curr[i+1][1], nxl, MPI_DOUBLE, 0, DATAPASS, comm);
|
|---|
| 314 | } else {
|
|---|
| 315 | printf("\n-------------------- time step:%d --------------------\n", time);
|
|---|
| 316 | for (int i=0; i < NPROCSY; i++) {
|
|---|
| 317 | int numRows = countRowForProc(i*NPROCSX);
|
|---|
| 318 |
|
|---|
| 319 | for (int j=0; j < numRows; j++) {
|
|---|
| 320 | for (int k=0; k < NPROCSX; k++) {
|
|---|
| 321 | int curr_rank = i*NPROCSX + k;
|
|---|
| 322 |
|
|---|
| 323 | if (curr_rank!=0) {
|
|---|
| 324 | int senderx = firstColForProc(curr_rank);
|
|---|
| 325 | int sendery = firstRowForProc(curr_rank);
|
|---|
| 326 | int senderNxl = countColForProc(curr_rank);
|
|---|
| 327 |
|
|---|
| 328 | MPI_Recv(recvbuf, senderNxl, MPI_DOUBLE, curr_rank, DATAPASS, comm, MPI_STATUS_IGNORE);
|
|---|
| 329 | printData(time, senderx, senderNxl, sendery+j, recvbuf);
|
|---|
| 330 | } else {
|
|---|
| 331 | printData(time, firstCol, nxl, firstRow+j, &u_curr[j+1][1]);
|
|---|
| 332 | }
|
|---|
| 333 | }
|
|---|
| 334 | printf("\n");
|
|---|
| 335 | }
|
|---|
| 336 | }
|
|---|
| 337 | }
|
|---|
| 338 | }
|
|---|
| 339 |
|
|---|
| 340 | int main(int argc, char * argv[]) {
|
|---|
| 341 | int i,j;
|
|---|
| 342 |
|
|---|
| 343 | #ifdef _CIVL
|
|---|
| 344 | // elaborating nx, ny, NPROCSX and NPROCSY...
|
|---|
| [841d45a] | 345 | //$elaborate(NPROCSY);
|
|---|
| 346 | //$elaborate(nx);
|
|---|
| 347 | //$elaborate(ny);
|
|---|
| 348 | //$elaborate(NPROCSX);
|
|---|
| [276a7df] | 349 | #endif
|
|---|
| 350 | MPI_Init(&argc, &argv);
|
|---|
| 351 | MPI_Comm_rank(comm, &rank);
|
|---|
| 352 | MPI_Comm_size(comm, &nprocs);
|
|---|
| 353 | initialization(argc, argv);
|
|---|
| 354 | initData();
|
|---|
| 355 | for (i=0; i<nsteps; i++) {
|
|---|
| 356 | if (nxl != 0 && nyl != 0) {
|
|---|
| 357 | if (i%wstep == 0)
|
|---|
| 358 | write_frame(i);
|
|---|
| 359 | exchange();
|
|---|
| 360 | update();
|
|---|
| 361 | }
|
|---|
| 362 | }
|
|---|
| 363 | for (i=0; i<nyl+2; i++) {
|
|---|
| 364 | free(u_curr[i]);
|
|---|
| 365 | free(u_next[i]);
|
|---|
| 366 | }
|
|---|
| 367 | free(u_curr);
|
|---|
| 368 | free(u_next);
|
|---|
| 369 | if (rank == 0)
|
|---|
| 370 | free(recvbuf);
|
|---|
| 371 | return 0;
|
|---|
| 372 | }
|
|---|