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