source: CIVL/examples/compare/diffusion1d/diffusion1d_par.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: 4.8 KB
RevLine 
[db66e7f]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);
[95b3837]63 $elaborate(nxl);
[db66e7f]64 left = first==0 || nxl==0 ? MPI_PROC_NULL : OWNER(first-1);
65 right = first+nxl >= nx || nxl == 0 ? MPI_PROC_NULL : OWNER(first+nxl);
66 u = (double*)malloc((nxl+2)*sizeof(double));
67 assert(u);
68 u_new = (double*)malloc((nxl+2)*sizeof(double));
69 assert(u_new);
70 if (rank == 0)
71 buf = (double*)malloc((1+nx/nprocs)*sizeof(double));
72}
73
74void initialize() {
75 // initialize globals and u...
76 init_globals();
77 lbound = U_INIT[0];
78 rbound = U_INIT[nx+1];
79 for (int i=1; i<=nxl; i++)
80 u[i] = U_INIT[first+i];
81 if (nx>=1 && rank == OWNER(0))
82 u[0] = u_new[0] = lbound;
83 if (nx>=1 && rank == OWNER(nx-1))
84 u[nxl+1] = u_new[nxl+1] = rbound;
85}
86
87/* Prints header for time step. Called by proc 0 only */
88void print_time_header() {
89 //printf("======= Time %d =======\n", time);
90 print_pos = 0;
91}
92
93/* Prints one cell. Called by proc 0 only. */
94void print_cell(double value) {
95 output[time][print_pos] = value;
96 print_pos++;
97}
98
99/* Prints the current values of u. */
100void write_frame() {
101 if (rank != 0) {
102 MPI_Send(u+1, nxl, MPI_DOUBLE, 0, 0, MPI_COMM_WORLD);
103 } else {
104 print_time_header();
105 print_cell(lbound); // left boundary
106 for (int source = 0; source < nprocs; source++) {
107 int count;
108
109 if (source != 0) {
110 MPI_Status status;
111
112 MPI_Recv(buf, 1+nx/nprocs, MPI_DOUBLE, source, 0, MPI_COMM_WORLD,
113 &status);
114 MPI_Get_count(&status, MPI_DOUBLE, &count);
115 } else {
116 for (int i = 1; i <= nxl; i++)
117 buf[i-1] = u[i];
118 count = nxl;
119 }
120 for (int i = 0; i < count; i++)
121 print_cell(buf[i]);
122 }
123 print_cell(rbound); // right boundary
124 }
125}
126
127/* exchange_ghost_cells: updates ghost cells using MPI communication */
128void exchange_ghost_cells() {
129 MPI_Sendrecv(&u[1], 1, MPI_DOUBLE, left, 0,
130 &u[nxl+1], 1, MPI_DOUBLE, right, 0,
131 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
132 MPI_Sendrecv(&u[nxl], 1, MPI_DOUBLE, right, 0,
133 &u[0], 1, MPI_DOUBLE, left, 0,
134 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
135}
136
137/* Updates u_new using u, then swaps u and u_new.
138 * Reads the ghost cells in u. Purely local operation. */
139void update() {
140 for (int i = 1; i <= nxl; i++)
141 u_new[i] = u[i] + k*(u[i+1] + u[i-1] - 2*u[i]);
142 double * tmp = u_new; u_new=u; u=tmp;
143}
144
145/* Executes the simulation. */
146int main(int argc, char *argv[]) {
147 MPI_Init(&argc, &argv);
148 initialize();
149 write_frame();
150 printf("nx = %d", nx);
151 for (time=1; time < nsteps; time++) {
152 exchange_ghost_cells();
153 update();
154 if (time%wstep==0)
155 write_frame();
156 }
157 MPI_Finalize();
158 free(u);
159 free(u_new);
160 if (rank == 0)
161 free(buf);
162 return 0;
163}
Note: See TracBrowser for help on using the repository browser.