source: CIVL/examples/messagePassing/diffusion1d.cvl@ f8b6015

1.23 2.0 main test-branch
Last change on this file since f8b6015 was 04e9b23, checked in by Manchun Zheng <zmanchun@…>, 12 years ago

implemented CIVL pragma transformer, added test accordingly; minor correction of examples.

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

  • Property mode set to 100644
File size: 8.1 KB
Line 
1/* Composite model of sequential and parallel 1d-diffusion solvers.
2 * Compares the results of the sequential and parallel computation to
3 * determine functional equivalence of the two versions.
4 *
5 * The initial values are taken as inputs. The two global boundary
6 * points are fixed.
7 *
8 * Command line example:
9 *
10 * civl verify -inputNPROCSB=3 -inputNSTEPSB=3 \
11 -inputNXB=6 diffusion1d_good.cvl
12 */
13#include <civlc.h>
14#include <stdio.h>
15#include <stdlib.h>
16#include <assert.h>
17#include "mpi.cvl"
18
19// Definitions from programs:
20#define OWNER(index) ((nprocs*(index+1)-1)/nx)
21
22// inputs:
23$input int NPROCS; // number of processes for parallel version
24$input int NPROCSB; // upper bound on nprocs
25$input double K; // k = alpha^2 * dt/(dx^2)
26$input int NSTEPS; // number of time steps
27$input int WSTEP; // write every this many time steps
28$input int NSTEPSB; // upperbound on nsteps
29$input int NX; // global number of discrete spatial points
30$input int NXB; // upper bound on nx
31$input double u_init[NX]; // initial values for temperature (u)
32
33// assumptions:
34$assume 2 <= NX && NX <= NXB;
35$assume K>0 && K<.5;
36$assume WSTEP >= 1 && WSTEP <= NSTEPS;
37// $assume NX >= NPROCS;
38$assume 1 <= NPROCS && NPROCS <= NPROCSB;
39$assume 1 <= NSTEPS && NSTEPS <= NSTEPSB;
40
41// global variables:
42int output_seq[NX]; // final (color) result of sequential computation
43int output_par[NX]; // final (color) result of parallel computation
44CMPI_Gcomm MPI_GCOMM_WORLD; // the global communicator object
45
46/* Abstract function representing conversion from temperature to color */
47$abstract int colorOf(double x);
48
49/* Sequential algorithm: performs simulation storing result of final
50 * time step in output_seq */
51void run_seq() {
52 $scope seq_scope = $here;
53 int nx, nsteps, wstep;
54 double k, *u;
55
56 void init_seq() {
57 int i;
58 int pos = 0; // file position
59
60 // in place of reading these from file:
61 nx = NX;
62 nsteps = NSTEPS;
63 wstep = WSTEP;
64 k = K;
65 printf("Diffusion1d (seq) with nx=%d, k=%f, nsteps=%d, wstep=%d\n",
66 nx, k, nsteps, wstep);
67 assert(nx>=2);
68 assert(k>0 && k<.5);
69 assert(nsteps >= 1);
70 assert(wstep >= 1 && wstep <=nsteps);
71 u = (double*)$malloc(seq_scope, nx*sizeof(double));
72 assert(u);
73 for (i = 0; i < nx; i++) u[i] = u_init[pos++];
74 $free(u);
75 }
76
77 void write_frame(int time) {
78 for (int i=0; i<nx; i++)
79 output_seq[i] = colorOf(u[i]);
80 }
81
82 void update() {
83 int i;
84 double u_new[nx];
85
86 for (i = 1; i < nx-1; i++)
87 u_new[i] = u[i] + k*(u[i+1] + u[i-1] -2*u[i]);
88 for (i = 1; i < nx-1; i++)
89 u[i] = u_new[i];
90 }
91
92 void _main() {
93 int iter;
94
95 init_seq();
96 write_frame(0);
97 for (iter = 1; iter <= nsteps; iter++) {
98 update();
99 if (iter%wstep==0) write_frame(iter);
100 }
101 free(u);
102 }
103
104 _main();
105}
106
107/* MPI Process in parallel algorithm */
108void MPI_Process (int _rank) {
109 $scope proc_scope = $here; /* this scope */
110 MPI_Comm MPI_COMM_WORLD;
111
112 int nx = -1; /* number of discrete points including endpoints */
113 double k = -1; /* D*dt/(dx*dx) */
114 int nsteps = -1; /* number of time steps */
115 int wstep = -1; /* write frame every this many time steps */
116 double *u; /* temperature function */
117 double *u_new; /* temp. used to update u */
118
119 int nprocs; /* number of processes */
120 int rank; /* the rank of this process */
121 int left; /* rank of left neighbor */
122 int right; /* rank of right neighbor on torus */
123 int nxl; /* horizontal extent of one process */
124 int first; /* global index for local index 0 */
125 int start; /* first local index to update */
126 int stop; /* last local index to update */
127 double *buf; /* temp. buffer used on proc 0 only */
128
129 int firstForProc(int rank) {
130 return (rank*nx)/nprocs;
131 }
132
133 int countForProc(int rank) {
134 int a;
135 int b;
136
137 a = firstForProc(rank+1);
138 b = firstForProc(rank);
139 return a-b;
140 }
141
142 /* init: initializes global variables. */
143 void init() {
144 int i, j;
145 int pos = 0; // position in input file
146
147 MPI_Comm_size(MPI_COMM_WORLD, &nprocs);
148 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
149 if (rank == 0) {
150 // in place of reading these from file:
151 nx = NX;
152 k = K;
153 nsteps = NSTEPS;
154 wstep = WSTEP;
155 printf("Diffusion1d (par) with nx=%d, k=%f, nsteps=%d, wstep=%d nprocs=%d\n",
156 nx, k, nsteps, wstep, nprocs);
157 }
158 MPI_Bcast(&nx, 1, MPI_INT, 0, MPI_COMM_WORLD);
159 MPI_Bcast(&k, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD);
160 MPI_Bcast(&nsteps, 1, MPI_INT, 0, MPI_COMM_WORLD);
161 MPI_Bcast(&wstep, 1, MPI_INT, 0, MPI_COMM_WORLD);
162 // assert(nx>=nprocs);
163 assert(k>0 && k<.5);
164 assert(nx>=2);
165 assert(nsteps>=1);
166 first = firstForProc(rank);
167 nxl = countForProc(rank);
168 if (first == 0 || nxl == 0)
169 left = MPI_PROC_NULL;
170 else
171 left = OWNER(first-1);
172 if (first+nxl >= nx || nxl == 0)
173 right = MPI_PROC_NULL;
174 else
175 right = OWNER(first+nxl);
176 u = (double*)$malloc(proc_scope, (nxl+2)*sizeof(double));
177 assert(u != NULL);
178 u_new = (double*)$malloc(proc_scope, (nxl+2)*sizeof(double));
179 assert(u_new != NULL);
180 if (rank == 0) {
181 buf = (double*)$malloc(proc_scope, (1+nx/nprocs)*sizeof(double));
182 for (i=1; i <= nxl; i++)
183 u[i] = u_init[pos++]; // instead of reading from file
184 for (i=1; i < nprocs; i++) {
185 int count_i = countForProc(i);
186
187 for (j=0; j<count_i; j++)
188 buf[j] = u_init[pos++];
189 MPI_Send(buf, count_i, MPI_DOUBLE, i, 0, MPI_COMM_WORLD);
190 }
191 $free(buf);
192 }
193 else {
194 buf = NULL;
195 MPI_Recv(u+1, nxl, MPI_DOUBLE, 0, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
196 }
197 if (rank == OWNER(0)) {
198 start = 2;
199 }
200 else {
201 start = 1;
202 }
203 if (rank == OWNER(nx-1)) {
204 stop = nxl - 1;
205 }
206 else {
207 stop = nxl;
208 }
209 }
210
211 void write_frame(int time) {
212 if (rank != 0) {
213 MPI_Send(u+1, nxl, MPI_DOUBLE, 0, 0, MPI_COMM_WORLD);
214 } else {
215 int source, i, j, count, global_index;
216 MPI_Status status;
217
218 global_index = 0;
219 for (source = 0; source < nprocs; source++) {
220 if (source != 0) {
221 MPI_Recv(buf, 1+nx/nprocs, MPI_DOUBLE, source, 0,
222 MPI_COMM_WORLD, &status);
223 MPI_Get_count(&status, MPI_DOUBLE, &count);
224 } else {
225 for (i = 1; i <= nxl; i++) buf[i-1] = u[i];
226 count = nxl;
227 }
228 for (i = 0; i < count; i++) {
229 output_par[global_index] = colorOf(buf[i]);
230 global_index++;
231 }
232 }
233 }
234 }
235
236 /* exchange_ghost_cells: updates ghost cells using MPI communication */
237 void exchange_ghost_cells() {
238 MPI_Sendrecv(&u[1], 1, MPI_DOUBLE, left, 0,
239 &u[nxl+1], 1, MPI_DOUBLE, right, 0,
240 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
241 MPI_Sendrecv(&u[nxl], 1, MPI_DOUBLE, right, 0,
242 &u[0], 1, MPI_DOUBLE, left, 0,
243 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
244 }
245
246 /* update: updates u. Uses ghost cells. Purely local operation. */
247 void update() {
248 int i;
249
250 for (i = start; i <= stop; i++)
251 u_new[i] = u[i] + k*(u[i+1] + u[i-1] - 2*u[i]);
252 for (i = start; i <= stop; i++)
253 u[i] = u_new[i];
254 }
255
256 void _main() {
257 int iter;
258
259 MPI_Init();
260 MPI_COMM_WORLD = MPI_Comm_create(proc_scope, MPI_GCOMM_WORLD, _rank);
261 init();
262 write_frame(0);
263 for (iter = 1; iter <= nsteps; iter++) {
264 exchange_ghost_cells();
265 update();
266 if (iter%wstep==0) write_frame(iter);
267 }
268 MPI_Finalize();
269 MPI_Comm_destroy(MPI_COMM_WORLD);
270 free(u);
271 free(u_new);
272 if (rank == 0)
273 free(buf);
274 }
275
276 _main();
277}
278
279void run_par() {
280 $proc procs[NPROCS]; // the MPI processes
281
282 MPI_GCOMM_WORLD = CMPI_Gcomm_create($root, NPROCS);
283 for (int i=0; i<NPROCS; i++) procs[i] = $spawn MPI_Process(i);
284 for (int i=0; i<NPROCS; i++) $wait(procs[i]);
285 CMPI_Gcomm_destroy(MPI_GCOMM_WORLD);
286}
287
288void main() {
289 for (int i=0; i<NX; i++) ;
290 for (int i=0; i<NSTEPS; i++) ;
291 for (int i=0; i<WSTEP; i++) ;
292 for (int i=0; i<NPROCS; i++) ;
293
294 run_seq();
295 run_par();
296 for (int i=0; i<NX; i++) {
297 printf("output[%d] = %d\n", i, output_seq[i]);
298 $assert(output_seq[i]==output_par[i]);
299 }
300}
Note: See TracBrowser for help on using the repository browser.