source: CIVL/examples/mpi/diffusion1d.c@ 7d77e64

main test-branch
Last change on this file since 7d77e64 was ea777aa, checked in by Alex Wilton <awilton@…>, 3 years ago

Moved examples, include, build_default.properties, common.xml, and README out from dev.civl.com into the root of the repo.

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

  • Property mode set to 100644
File size: 6.3 KB
Line 
1/*******************************************************************
2 * diffusion1d.c: parallel 1d-diffusion solver with constant
3 * boundaries.
4 *
5 * This example contains a sequential 1d-diffusion 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 diffusion1d.c ; mpiexec -n 4 ./a.out Or replace
10 * "4" with however many procs you want to use.
11 *
12 * To verify: civl verify diffusion1d.
13 * Author: Ziqing Luo
14 ********************************************************************/
15#include <stdlib.h>
16#include <stdio.h>
17#include <assert.h>
18#include <mpi.h>
19
20#define OWNER(index) ((nprocs*(index+1)-1)/nx)
21
22#ifdef _CIVL
23#include <civlc.cvh>
24$input int NXB = 5; // upper bound on nx
25$input int nx; // global number of points excl. boundary
26$assume(1<=nx && nx<=NXB);
27$input double U_INIT[nx+2]; // initial values for temperature incl. boundary
28$input double k; // the constant D*dt/(dx*dx)
29$assume(k>0 && k<.5);
30$input int NSTEPS_BOUND=5; // upper bound on nsteps
31$input int nsteps; // number of time steps
32$assume(1<=nsteps && nsteps<=NSTEPS_BOUND);
33$input int wstep; // write frame every this many time steps
34$assume(1<=wstep && wstep<=nsteps);
35double oracle[nsteps][nx+2]; // solution computed sequentially, proc 0 only
36int _mpi_nprocs_lo = 1;
37int _mpi_nprocs_hi = 3;
38
39#else
40
41int nx;
42double k;
43int nsteps;
44int wstep;
45
46#endif
47
48/* Global variables */
49
50double lbound; /* left fixed boundary value */
51double rbound; /* right fixed boundary value */
52double *u; /* temperature function, local */
53double *u_new; /* second copy of temperature function, local */
54int nprocs; /* number of processes */
55int rank; /* the rank of this process */
56int left; /* rank of left neighbor */
57int right; /* rank of right neighbor on torus */
58int nxl; /* horizontal extent of one process */
59int first; /* global index for local index 0 */
60double *buf; /* temp. buffer used on proc 0 only */
61int print_pos; /* number of cells printed on current line */
62int time=0; /* current time step */
63
64
65/* Returns the global index of the first cell owned
66 * by the process with given rank */
67int firstForProc(int rank) {
68 return (rank*nx)/nprocs;
69}
70
71/* Returns the number of cells owned by the process
72 * of the given rank (excluding ghosts) */
73int countForProc(int rank) {
74 int a = firstForProc(rank+1);
75 int b = firstForProc(rank);
76
77 return a-b;
78}
79
80/* Initializes the global variables.
81 * Precondition: the configuration parameters have
82 * already been set. */
83void init_globals() {
84 MPI_Comm_size(MPI_COMM_WORLD, &nprocs);
85 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
86 // nxl: number actual points (incl. end-points)
87 // nxl+2: size of array (incl. ghost cells)
88 first = firstForProc(rank);
89 nxl = countForProc(rank);
90 left = first==0 || nxl==0 ? MPI_PROC_NULL : OWNER(first-1);
91 right = first+nxl >= nx || nxl == 0 ? MPI_PROC_NULL : OWNER(first+nxl);
92 u = (double*)malloc((nxl+2)*sizeof(double));
93 assert(u);
94 u_new = (double*)malloc((nxl+2)*sizeof(double));
95 assert(u_new);
96 if (rank == 0)
97 buf = (double*)malloc((1+nx/nprocs)*sizeof(double));
98}
99
100void initialize() {
101#ifdef _CIVL
102 // initialize globals and u...
103 init_globals();
104 lbound = U_INIT[0];
105 rbound = U_INIT[nx+1];
106 for (int i=1; i<=nxl; i++)
107 u[i] = U_INIT[first+i];
108 if (rank == 0) {
109 // compute the oracle...
110 for (int i=0; i<nx+2; i++)
111 oracle[0][i]=U_INIT[i];
112 for (int t=1; t<nsteps; t++) {
113 oracle[t][0] = lbound;
114 for (int i=1; i<=nx; i++)
115 oracle[t][i] = oracle[t-1][i] +
116 k*(oracle[t-1][i+1] + oracle[t-1][i-1]
117 - 2*oracle[t-1][i]);
118 oracle[t][nx+1] = rbound;
119 }
120 }
121#else
122 nx = 10;
123 k = 0.2;
124 nsteps = 20;
125 wstep = 2;
126 lbound = rbound = 0.0;
127 init_globals();
128 for (int i=1; i<=nxl; i++)
129 u[i]=100.0;
130#endif
131 if (nx>=1 && rank == OWNER(0))
132 u[0] = u_new[0] = lbound;
133 if (nx>=1 && rank == OWNER(nx-1))
134 u[nxl+1] = u_new[nxl+1] = rbound;
135 if (rank == 0)
136 printf("nx=%d, k=%lf, nsteps=%d, wstep=%d, nprocs=%d\n",
137 nx, k, nsteps, wstep, nprocs);
138}
139
140/* Prints header for time step. Called by proc 0 only */
141void print_time_header() {
142 printf("======= Time %d =======\n", time);
143 print_pos = 0;
144}
145
146/* Prints one cell. Called by proc 0 only. */
147void print_cell(double value) {
148 printf("%7.2f ", value);
149#pragma CIVL $assert(value == oracle[time][print_pos], \
150 "Error: disagreement at time %d position %d: saw %lf, expected %lf", \
151 time, print_pos, value, oracle[time][print_pos]);
152 print_pos++;
153}
154
155/* Prints the current values of u. */
156void write_frame() {
157 if (rank != 0) {
158 MPI_Send(u+1, nxl, MPI_DOUBLE, 0, 0, MPI_COMM_WORLD);
159 } else {
160 print_time_header();
161 print_cell(lbound); // left boundary
162 for (int source = 0; source < nprocs; source++) {
163 int count;
164
165 if (source != 0) {
166 MPI_Status status;
167
168 MPI_Recv(buf, 1+nx/nprocs, MPI_DOUBLE, source, 0, MPI_COMM_WORLD,
169 &status);
170 MPI_Get_count(&status, MPI_DOUBLE, &count);
171 } else {
172 for (int i = 1; i <= nxl; i++)
173 buf[i-1] = u[i];
174 count = nxl;
175 }
176 for (int i = 0; i < count; i++)
177 print_cell(buf[i]);
178 }
179 print_cell(rbound); // right boundary
180 printf("\n");
181 }
182}
183
184/* exchange_ghost_cells: updates ghost cells using MPI communication */
185void exchange_ghost_cells() {
186 MPI_Sendrecv(&u[1], 1, MPI_DOUBLE, left, 0,
187 &u[nxl+1], 1, MPI_DOUBLE, right, 0,
188 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
189 MPI_Sendrecv(&u[nxl], 1, MPI_DOUBLE, right, 0,
190 &u[0], 1, MPI_DOUBLE, left, 0,
191 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
192}
193
194/* Updates u_new using u, then swaps u and u_new.
195 * Reads the ghost cells in u. Purely local operation. */
196void update() {
197 for (int i = 1; i <= nxl; i++)
198 u_new[i] = u[i] + k*(u[i+1] + u[i-1] - 2*u[i]);
199 double * tmp = u_new; u_new=u; u=tmp;
200}
201
202/* Executes the simulation. */
203int main(int argc, char *argv[]) {
204 MPI_Init(&argc, &argv);
205 initialize();
206#ifdef _CIVL
207 $elaborate(nxl);
208#endif
209 write_frame();
210 for (time=1; time < nsteps; time++) {
211 exchange_ghost_cells();
212 update();
213 if (time%wstep==0)
214 write_frame();
215 }
216 MPI_Finalize();
217 free(u);
218 free(u_new);
219 if (rank == 0)
220 free(buf);
221 return 0;
222}
Note: See TracBrowser for help on using the repository browser.