source: CIVL/examples/mpi/diffusion2d.c@ 486df51

1.23 2.0 main test-branch
Last change on this file since 486df51 was 23bf41d, checked in by Ziqing Luo <ziqing@…>, 11 years ago

fixed a tiny bug in diff2d

git-svn-id: svn://vsl.cis.udel.edu/civl/trunk@1899 fb995dde-84ed-4084-dfe6-e5aef3e2452c

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