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

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

add a broken test for compare diffusion1d

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

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