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

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

updates CIVL to use the new ABC FrontEnd. tests updated accordingly when necessary.

git-svn-id: svn://vsl.cis.udel.edu/civl/trunk@1380 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 <comm.cvh>
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++)
74 u[i] = u_init[pos++];
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
83 void update() {
84 int i;
85 double u_new[nx];
86
87 for (i = 1; i < nx-1; i++)
88 u_new[i] = u[i] + k*(u[i+1] + u[i-1] -2*u[i]);
89 for (i = 1; i < nx-1; i++)
90 u[i] = u_new[i];
91 }
92
93 void _main() {
94 int iter;
95
96 init_seq();
97 write_frame(0);
98 for (iter = 1; iter <= nsteps; iter++) {
99 update();
100 if (iter%wstep==0) write_frame(iter);
101 }
102 free(u);
103 }
104
105 _main();
106}
107
108/* MPI Process in parallel algorithm */
109void MPI_Process (int _rank) {
110 $scope proc_scope = $here; /* this scope */
111 MPI_Comm MPI_COMM_WORLD;
112
113 int nx = -1; /* number of discrete points including endpoints */
114 double k = -1; /* D*dt/(dx*dx) */
115 int nsteps = -1; /* number of time steps */
116 int wstep = -1; /* write frame every this many time steps */
117 double *u; /* temperature function */
118 double *u_new; /* temp. used to update u */
119
120 int nprocs; /* number of processes */
121 int rank; /* the rank of this process */
122 int left; /* rank of left neighbor */
123 int right; /* rank of right neighbor on torus */
124 int nxl; /* horizontal extent of one process */
125 int first; /* global index for local index 0 */
126 int start; /* first local index to update */
127 int stop; /* last local index to update */
128 double *buf; /* temp. buffer used on proc 0 only */
129
130 int firstForProc(int rank) {
131 return (rank*nx)/nprocs;
132 }
133
134 int countForProc(int rank) {
135 int a;
136 int b;
137
138 a = firstForProc(rank+1);
139 b = firstForProc(rank);
140 return a-b;
141 }
142
143 /* init: initializes global variables. */
144 void init() {
145 int i, j;
146 int pos = 0; // position in input file
147
148 MPI_Comm_size(MPI_COMM_WORLD, &nprocs);
149 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
150 if (rank == 0) {
151 // in place of reading these from file:
152 nx = NX;
153 k = K;
154 nsteps = NSTEPS;
155 wstep = WSTEP;
156 printf("Diffusion1d (par) with nx=%d, k=%f, nsteps=%d, wstep=%d nprocs=%d\n",
157 nx, k, nsteps, wstep, nprocs);
158 }
159 MPI_Bcast(&nx, 1, MPI_INT, 0, MPI_COMM_WORLD);
160 MPI_Bcast(&k, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD);
161 MPI_Bcast(&nsteps, 1, MPI_INT, 0, MPI_COMM_WORLD);
162 MPI_Bcast(&wstep, 1, MPI_INT, 0, MPI_COMM_WORLD);
163 // assert(nx>=nprocs);
164 assert(k>0 && k<.5);
165 assert(nx>=2);
166 assert(nsteps>=1);
167 first = firstForProc(rank);
168 nxl = countForProc(rank);
169 if (first == 0 || nxl == 0)
170 left = MPI_PROC_NULL;
171 else
172 left = OWNER(first-1);
173 if (first+nxl >= nx || nxl == 0)
174 right = MPI_PROC_NULL;
175 else
176 right = OWNER(first+nxl);
177 u = (double*)$malloc(proc_scope, (nxl+2)*sizeof(double));
178 assert(u != NULL);
179 u_new = (double*)$malloc(proc_scope, (nxl+2)*sizeof(double));
180 assert(u_new != NULL);
181 if (rank == 0) {
182 buf = (double*)$malloc(proc_scope, (1+nx/nprocs)*sizeof(double));
183 for (i=1; i <= nxl; i++)
184 u[i] = u_init[pos++]; // instead of reading from file
185 for (i=1; i < nprocs; i++) {
186 int count_i = countForProc(i);
187
188 for (j=0; j<count_i; j++)
189 buf[j] = u_init[pos++];
190 MPI_Send(buf, count_i, MPI_DOUBLE, i, 0, MPI_COMM_WORLD);
191 }
192 // $free(buf);
193 }
194 else {
195 buf = NULL;
196 MPI_Recv(u+1, nxl, MPI_DOUBLE, 0, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
197 }
198 if (rank == OWNER(0)) {
199 start = 2;
200 }
201 else {
202 start = 1;
203 }
204 if (rank == OWNER(nx-1)) {
205 stop = nxl - 1;
206 }
207 else {
208 stop = nxl;
209 }
210 }
211
212 void write_frame(int time) {
213 if (rank != 0) {
214 MPI_Send(u+1, nxl, MPI_DOUBLE, 0, 0, MPI_COMM_WORLD);
215 } else {
216 int source, i, j, count, global_index;
217 MPI_Status status;
218
219 global_index = 0;
220 for (source = 0; source < nprocs; source++) {
221 if (source != 0) {
222 MPI_Recv(buf, 1+nx/nprocs, MPI_DOUBLE, source, 0,
223 MPI_COMM_WORLD, &status);
224 MPI_Get_count(&status, MPI_DOUBLE, &count);
225 } else {
226 for (i = 1; i <= nxl; i++) buf[i-1] = u[i];
227 count = nxl;
228 }
229 for (i = 0; i < count; i++) {
230 output_par[global_index] = colorOf(buf[i]);
231 global_index++;
232 }
233 }
234 }
235 }
236
237 /* exchange_ghost_cells: updates ghost cells using MPI communication */
238 void exchange_ghost_cells() {
239 MPI_Sendrecv(&u[1], 1, MPI_DOUBLE, left, 0,
240 &u[nxl+1], 1, MPI_DOUBLE, right, 0,
241 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
242 MPI_Sendrecv(&u[nxl], 1, MPI_DOUBLE, right, 0,
243 &u[0], 1, MPI_DOUBLE, left, 0,
244 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
245 }
246
247 /* update: updates u. Uses ghost cells. Purely local operation. */
248 void update() {
249 int i;
250
251 for (i = start; i <= stop; i++)
252 u_new[i] = u[i] + k*(u[i+1] + u[i-1] - 2*u[i]);
253 for (i = start; i <= stop; i++)
254 u[i] = u_new[i];
255 }
256
257 void _main() {
258 int iter;
259
260 MPI_Init();
261 MPI_COMM_WORLD = MPI_Comm_create(proc_scope, MPI_GCOMM_WORLD, _rank);
262 init();
263 write_frame(0);
264 for (iter = 1; iter <= nsteps; iter++) {
265 exchange_ghost_cells();
266 update();
267 if (iter%wstep==0) write_frame(iter);
268 }
269 MPI_Finalize();
270 MPI_Comm_destroy(MPI_COMM_WORLD);
271 free(u);
272 free(u_new);
273 if (rank == 0)
274 free(buf);
275 }
276
277 _main();
278}
279
280void run_par() {
281 $proc procs[NPROCS]; // the MPI processes
282
283 MPI_GCOMM_WORLD = CMPI_Gcomm_create($root, NPROCS);
284 for (int i=0; i<NPROCS; i++) procs[i] = $spawn MPI_Process(i);
285 for (int i=0; i<NPROCS; i++) $wait(procs[i]);
286 CMPI_Gcomm_destroy(MPI_GCOMM_WORLD);
287}
288
289void main() {
290 for (int i=0; i<NX; i++) ;
291 for (int i=0; i<NSTEPS; i++) ;
292 for (int i=0; i<WSTEP; i++) ;
293 for (int i=0; i<NPROCS; i++) ;
294
295 run_seq();
296 run_par();
297 for (int i=0; i<NX; i++) {
298 printf("output[%d] = %d\n", i, output_seq[i]);
299 $assert(output_seq[i]==output_par[i]);
300 }
301}
Note: See TracBrowser for help on using the repository browser.