source: CIVL/examples/compare/diffusion1d/diffusion1d_par.c@ 631a953

1.23 2.0 main test-branch
Last change on this file since 631a953 was 3ad5f25, checked in by Manchun Zheng <zmanchun@…>, 11 years ago

cleaned up elaborate calls

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

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