source: CIVL/examples/mpi/diffusion1d_ab.c@ 397ae5f

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