source: CIVL/examples/mpi/wave1d.c@ 83af34d

1.23 2.0 main test-branch
Last change on this file since 83af34d was 276a7df, checked in by Ziqing Luo <ziqing@…>, 11 years ago

let gcomm_destroy check if there is any message remaining in buffer changed programs correspondingly

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

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