source: CIVL/examples/experimental/diff1d_par.c@ 7320499

main test-branch
Last change on this file since 7320499 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: 5.1 KB
Line 
1#ifdef _CIVL
2#include <civlc.cvh>
3#endif
4#include <stdlib.h>
5#include <stdio.h>
6#include <assert.h>
7#include <mpi.h>
8#define OWNER(index) ((nprocs*(index+1)-1)/nx)
9
10$input int NXB = 5; // upper bound on nx
11$input int nx; // global number of points excl. boundary
12$assume(1<=nx && nx<=NXB);
13$input double U_INIT[nx+2]; // initial values for temperature incl. boundary
14$input double k; // the constant D*dt/(dx*dx)
15$assume(k>0 && k<.5);
16$input int NSTEPS_BOUND=5; // upper bound on nsteps
17$input int nsteps; // number of time steps
18$assume(1<=nsteps && nsteps<=NSTEPS_BOUND);
19$input int wstep = 1; // write frame every this many time steps
20//$assume(1<=wstep && wstep<=nsteps);
21//$assume(nsteps%wstep == 0 && wstep != 0);
22$output double output[nsteps][nx+2]; // solution computed sequentially, proc 0 only
23
24int _NPROCS_LOWER_BOUND = 1;
25int _NPROCS_UPPER_BOUND = 3;
26
27/* Global variables */
28
29double lbound; /* left fixed boundary value */
30double rbound; /* right fixed boundary value */
31double *u; /* temperature function, local */
32double *u_new; /* second copy of temperature function, local */
33int nprocs; /* number of processes */
34int rank; /* the rank of this process */
35int left; /* rank of left neighbor */
36int right; /* rank of right neighbor on torus */
37int nxl; /* horizontal extent of one process */
38int first; /* global index for local index 0 */
39double *buf; /* temp. buffer used on proc 0 only */
40int print_pos; /* number of cells printed on current line */
41int time=0; /* current time step */
42
43/* Returns the global index of the first cell owned
44 * by the process with given rank */
45int firstForProc(int rank) {
46 return (rank*nx)/nprocs;
47}
48
49/* Returns the number of cells owned by the process
50 * of the given rank (excluding ghosts) */
51int countForProc(int rank) {
52 int a = firstForProc(rank+1);
53 int b = firstForProc(rank);
54
55 return a-b;
56}
57
58/* Initializes the global variables.
59 * Precondition: the configuration parameters have
60 * already been set. */
61void init_globals() {
62 MPI_Comm_size(MPI_COMM_WORLD, &nprocs);
63 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
64 // nxl: number actual points (incl. end-points)
65 // nxl+2: size of array (incl. ghost cells)
66 first = firstForProc(rank);
67 nxl = countForProc(rank);
68 left = first==0 || nxl==0 ? MPI_PROC_NULL : OWNER(first-1);
69 right = first+nxl >= nx || nxl == 0 ? MPI_PROC_NULL : OWNER(first+nxl);
70 u = (double*)malloc((nxl+2)*sizeof(double));
71 assert(u);
72 u_new = (double*)malloc((nxl+2)*sizeof(double));
73 assert(u_new);
74 if (rank == 0)
75 buf = (double*)malloc((1+nx/nprocs)*sizeof(double));
76}
77
78void initialize() {
79 // initialize globals and u...
80 init_globals();
81 lbound = U_INIT[0];
82 rbound = U_INIT[nx+1];
83 for (int i=1; i<=nxl; i++)
84 u[i] = U_INIT[first+i];
85 if (nx>=1 && rank == OWNER(0))
86 u[0] = u_new[0] = lbound;
87 if (nx>=1 && rank == OWNER(nx-1))
88 u[nxl+1] = u_new[nxl+1] = rbound;
89 //if (rank == 0)
90 //printf("nx=%d, k=%lf, nsteps=%d, wstep=%d, nprocs=%d\n",
91 // nx, k, nsteps, wstep, nprocs);
92}
93
94/* Prints header for time step. Called by proc 0 only */
95void print_time_header() {
96 //printf("======= Time %d =======\n", time);
97 print_pos = 0;
98}
99
100/* Prints one cell. Called by proc 0 only. */
101void print_cell(double value) {
102 printf("%7.2f\n", value);
103 output[time][print_pos] = value;
104 print_pos++;
105}
106
107/* Prints the current values of u. */
108void write_frame() {
109 if (rank != 0) {
110 MPI_Send(u+1, nxl, MPI_DOUBLE, 0, 0, MPI_COMM_WORLD);
111 } else {
112 print_time_header();
113 print_cell(lbound); // left boundary
114 for (int source = 0; source < nprocs; source++) {
115 int count;
116
117 if (source != 0) {
118 MPI_Status status;
119
120 MPI_Recv(buf, 1+nx/nprocs, MPI_DOUBLE, source, 0, MPI_COMM_WORLD,
121 &status);
122 MPI_Get_count(&status, MPI_DOUBLE, &count);
123 } else {
124 for (int i = 1; i <= nxl; i++)
125 buf[i-1] = u[i];
126 count = nxl;
127 }
128 for (int i = 0; i < count; i++)
129 print_cell(buf[i]);
130 }
131 print_cell(rbound); // right boundary
132 //printf("\n");
133 }
134}
135
136/* exchange_ghost_cells: updates ghost cells using MPI communication */
137void exchange_ghost_cells() {
138 MPI_Sendrecv(&u[1], 1, MPI_DOUBLE, left, 0,
139 &u[nxl+1], 1, MPI_DOUBLE, right, 0,
140 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
141 MPI_Sendrecv(&u[nxl], 1, MPI_DOUBLE, right, 0,
142 &u[0], 1, MPI_DOUBLE, left, 0,
143 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
144}
145
146/* Updates u_new using u, then swaps u and u_new.
147 * Reads the ghost cells in u. Purely local operation. */
148void update() {
149 for (int i = 1; i <= nxl; i++)
150 u_new[i] = u[i] + k*(u[i+1] + u[i-1] - 2*u[i]);
151 double * tmp = u_new; u_new=u; u=tmp;
152}
153
154/* Executes the simulation. */
155int main(int argc, char *argv[]) {
156 //elaborate nx to concrete value...
157 $elaborate(nx);
158 $elaborate(nsteps);
159 MPI_Init(&argc, &argv);
160 initialize();
161 write_frame();
162 printf("nx = %d", nx);
163 for (time=1; time < nsteps; time++) {
164 exchange_ghost_cells();
165 update();
166 if (time%wstep==0)
167 write_frame();
168 }
169 MPI_Finalize();
170 free(u);
171 free(u_new);
172 if (rank == 0)
173 free(buf);
174 return 0;
175}
176
177
Note: See TracBrowser for help on using the repository browser.