source: CIVL/examples/mpi/wave1d.c@ 887a26b4

1.23 2.0 main test-branch
Last change on this file since 887a26b4 was 2321281, checked in by Manchun Zheng <zmanchun@…>, 11 years ago

fixed $elaborate for examples.

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

  • Property mode set to 100644
File size: 6.7 KB
Line 
1/* wave1d.c: parallel 1d-wave equation solver with constant boundary.
2 * To execute: mpicc wave1d.c ; mpiexec -n 4 ./a.out
3 * Or replace "4" with however many procs you want to use.
4 * To verify: civl verify wave1d.c
5 */
6#include <stdio.h>
7#include <stdlib.h>
8#include <string.h>
9#include <assert.h>
10#include <math.h>
11#include <mpi.h>
12
13#define SQR(x) ((x)*(x))
14#define OWNER(index) ((nprocs*(index+1)-1)/nx)
15
16/* MPI message tags */
17#define FROMLEFT 1
18#define FROMRIGHT 2
19#define DATAPASS 3
20
21/* Input parameters */
22#ifdef _CIVL
23#include <civlc.cvh>
24
25$input int NXB = 5; /* upper bound on nx */
26$input int nx; /* number of discrete points including endpoints */
27$assume(2 <= nx && nx <= NXB); /* setting bounds */
28$input double c; /* physical constant to do with string */
29$assume(c > 0.0);
30$input int NSTEPSB = 5;
31$input int nsteps; /* number of time steps */
32$assume(0 < nsteps && nsteps <= NSTEPSB);
33$input int wstep = 1; /* number of time steps between printing */
34$input int _NPROCS_LOWER_BOUND = 1;
35$input int _NPROCS_UPPER_BOUND = 4;
36double oracle[nsteps + 1][nx + 2];/* result of sequential run at every time step */
37$input double u_init[nx]; /* arbitrary initial position of string */
38
39#else
40
41int nx, height_init, width_init;
42int nsteps, wstep;
43double c;
44
45#endif
46
47/* Global varibales */
48double *u_prev, *u_curr, *u_next;
49double k;
50int nprocs, nxl, rank;
51int first; /* global index of first cell owned by this process */
52int left, right; /* ranks of left neighbor and right neighbor */
53
54/* Returns the global index of the first cell owned
55 * by the process with given rank */
56int firstForProc(int rank) {
57 return (rank*nx)/nprocs;
58}
59
60/* Returns the number of cells the given process owns */
61int countForProc(int rank) {
62 int a = firstForProc(rank);
63 int b = firstForProc(rank + 1);
64
65 return b - a;
66}
67
68#ifndef _CIVL
69
70/* Some particular initial condition for testing only */
71void init(int first, int nxl) {
72 int i;
73 double e = exp(1.0);
74
75 for(i = 1; i < nxl+1; i++) {
76 int idx = first + i - 1;
77
78 if(idx == 1 || idx >= width_init)
79 u_prev[i] = 0.0;
80 else
81 u_prev[i] = height_init * e *
82 exp(-1.0/(1-SQR(2.0*(idx-width_init/2.0)/width_init)));
83 }
84}
85
86#endif
87
88/* Update cells owned by processes */
89void update() {
90 int i;
91 double *tmp;
92
93 for (i = 1; i < nxl + 1; i++){
94 u_next[i] = 2.0*u_curr[i] - u_prev[i] +
95 k*(u_curr[i+1] + u_curr[i-1] -2.0*u_curr[i]);
96 }
97 //cycle pointers:
98 tmp = u_prev;
99 u_prev = u_curr;
100 u_curr = u_next;
101 u_next = tmp;
102}
103
104/* Initialization function, initializes all parameters and data array.
105 * In CIVL mode, process 0 runs the complete problem sequentially to
106 * create the oracle which will be compared with the results of the
107 * parallel run.
108 */
109void initialization() {
110 int i, j;
111
112#ifndef _CIVL
113
114 nx = 50;
115 c = 0.3;
116 height_init = 10;
117 width_init = 10;
118 nsteps = 500;
119 wstep = 5;
120
121#endif
122
123 assert(nx >= 2);
124 assert(c > 0);
125 assert(nsteps >= 1);
126 assert(wstep >= 1 && wstep <= nsteps);
127 k = c * c;
128 printf("Wave1d with nx=%d, c=%f, nsteps=%d, wstep=%d\n", nx, c, nsteps, wstep);
129
130#ifdef _CIVL
131
132 // If in CIVL verification mode and rank is 0,
133 // do a sequential run and store result in "oracle"
134 // for comparison later
135 if(rank == 0) {
136 // sets constant boundaries for first string and
137 // it's virtual previous string
138 oracle[0][0] = 0;
139 oracle[0][nx + 1] = 0;
140 // initialization for first string and
141 // it's virtual previous string
142 for(i = 1; i < nx + 1; i++) {
143 oracle[0][i] = u_init[i - 1];
144 oracle[1][i] = u_init[i - 1];
145 }
146 for(i = 1; i < nsteps; i++){
147 oracle[i][0] = 0;
148 oracle[i][nx + 1] = 0;
149 // update
150 for (j = 1; j < nx + 1; j++)
151 oracle[i+1][j] = 2.0*oracle[i][j] - oracle[i-1][j] +
152 k*(oracle[i][j+1] + oracle[i][j-1] -2.0*oracle[i][j]);
153 }
154 }
155
156#endif
157
158 nxl = countForProc(rank);
159 first = firstForProc(rank);
160 u_prev = (double *)malloc((nxl + 2) * sizeof(double));
161 assert(u_prev);
162 u_curr = (double *)malloc((nxl + 2) * sizeof(double));
163 assert(u_curr);
164 u_next = (double *)malloc((nxl + 2) * sizeof(double));
165 assert(u_next);
166 // sets constant boundaries
167 u_prev[0] = 0;
168 u_curr[0] = 0;
169 u_next[0] = 0;
170 u_prev[nxl + 1] = 0;
171 u_curr[nxl + 1] = 0;
172 u_next[nxl + 1] = 0;
173 // skip processes which are allocated no cells
174 if(first != 0 && nxl != 0)
175 left = OWNER(first - 1);
176 else
177 left = MPI_PROC_NULL;
178 if(first + nxl < nx && nxl != 0)
179 right = OWNER(first + nxl);
180 else
181 right = MPI_PROC_NULL;
182}
183
184/* Print out the value of data cells;
185 Do comparison in CIVL mode */
186void printData (int time, int first, int length, double * buf) {
187 int i;
188
189 for(i = 0; i < length; i++){
190 printf("u_curr[%d]=%8.8f ", first + i, buf[i]);
191#ifdef _CIVL
192
193 $assert((oracle[time + 1][first + i + 1] == buf[i]), \
194 "Error: disagreement at time %d position %d: saw %lf, expected %lf", \
195 time, first + i, buf[i], oracle[time + 1][first + i + 1]);
196
197#endif
198 printf("\n");
199 }
200}
201
202/* receives data from other processes and wirte frames */
203void write_frame (int time) {
204 if(rank == 0) {
205 double buf[nx + 2];
206 MPI_Status status;
207
208 printf("\n======= Time %d =======\n", time);
209 printData(time, first, nxl, &u_curr[1]);
210 for(int i=1; i < nprocs; i++) {
211 int displ = firstForProc(i);
212 int count;
213
214 MPI_Recv(buf, nx, MPI_DOUBLE, i, DATAPASS, MPI_COMM_WORLD, &status);
215 MPI_Get_count(&status, MPI_DOUBLE, &count);
216 printData(time, displ, count, buf);
217 }
218 printf("\n");
219 } else
220 MPI_Send(&u_curr[1], nxl, MPI_DOUBLE, 0, DATAPASS, MPI_COMM_WORLD);
221}
222
223/* Exchanging ghost cells */
224void exchange(){
225 MPI_Sendrecv(&u_curr[1], 1, MPI_DOUBLE, left, FROMRIGHT, &u_curr[nxl+1], 1, MPI_DOUBLE,
226 right, FROMRIGHT, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
227 MPI_Sendrecv(&u_curr[nxl], 1, MPI_DOUBLE, right, FROMLEFT, &u_curr[0], 1,
228 MPI_DOUBLE, left, FROMLEFT, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
229}
230
231int main(int argc, char * argv[]) {
232 int iter;
233
234
235#ifdef _CIVL
236
237 // elaborate nx to concrete value...
238 $elaborate(nx);
239
240#endif
241 MPI_Init(&argc, &argv);
242 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
243 MPI_Comm_size(MPI_COMM_WORLD, &nprocs);
244 initialization();
245#ifdef _CIVL
246
247 if(nxl != 0) {
248 memcpy(&u_curr[1], &u_init[first], sizeof(double) * nxl);
249 memcpy(&u_prev[1], &u_curr[1], sizeof(double) * nxl);
250 }
251
252#else
253
254 if(nxl != 0) {
255 init(first, nxl);
256 memcpy(&u_curr[1], &u_prev[1], nxl * sizeof(double));
257 }
258
259#endif
260 for(iter = 0; iter < nsteps; iter++) {
261 if(iter % wstep == 0)
262 write_frame(iter);
263 if(nxl != 0){
264 exchange();
265 update();
266 }
267 }
268 free(u_curr);
269 free(u_prev);
270 free(u_next);
271 MPI_Finalize();
272 return 0;
273}
Note: See TracBrowser for help on using the repository browser.