-
Notifications
You must be signed in to change notification settings - Fork 0
/
driver.c
160 lines (151 loc) · 5.93 KB
/
driver.c
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
/* Main solver routines for heat equation solver */
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <assert.h>
#include <mpi.h>
#include "heat.h"
/* Exchange the boundary values */
void exchange_init(field *temperature, parallel_data *parallel)
{
int ind, width;
width = temperature->ny + 2;
// Send to the up, receive from down
ind = idx(1, 0, width);
MPI_Isend(&temperature->data[ind], 1, parallel->rowtype,
parallel->nup, 11, parallel->comm, ¶llel->requests[0]);
ind = idx(temperature->nx + 1, 0, width);
MPI_Irecv(&temperature->data[ind], 1, parallel->rowtype,
parallel->ndown, 11, parallel->comm, ¶llel->requests[1]);
// Send to the down, receive from up
ind = idx(temperature->nx, 0, width);
MPI_Isend(&temperature->data[ind], 1, parallel->rowtype,
parallel->ndown, 12, parallel->comm, ¶llel->requests[2]);
ind = idx(0, 0, width);
MPI_Irecv(&temperature->data[ind], 1, parallel->rowtype,
parallel->nup, 12, parallel->comm, ¶llel->requests[3]);
// Send to the left, receive from right
ind = idx(0, 1, width);
MPI_Isend(&temperature->data[ind], 1, parallel->columntype,
parallel->nleft, 13, parallel->comm, ¶llel->requests[4]);
ind = idx(0, temperature->ny + 1, width);
MPI_Irecv(&temperature->data[ind], 1, parallel->columntype,
parallel->nright, 13, parallel->comm, ¶llel->requests[5]);
// Send to the right, receive from left
ind = idx(0, temperature->ny, width);
MPI_Isend(&temperature->data[ind], 1, parallel->columntype,
parallel->nright, 14, parallel->comm, ¶llel->requests[7]);
ind = 0;
MPI_Irecv(&temperature->data[ind], 1, parallel->columntype,
parallel->nleft, 14, parallel->comm, ¶llel->requests[6]);
}
/* complete the non-blocking communication */
void exchange_finalize(parallel_data *parallel)
{
MPI_Waitall(8, ¶llel->requests[0], MPI_STATUSES_IGNORE);
}
/* Update the temperature values using five-point stencil */
void evolve_interior(field *curr, field *prev, double a, double dt)
{
int i, j;
int ic, iu, id, il, ir; // indexes for center, up, down, left, right
int width;
width = curr->ny + 2;
double dx2, dy2;
/* Determine the temperature field at next time step
* As we have fixed boundary conditions, the outermost gridpoints
* are not updated. */
dx2 = prev->dx * prev->dx;
dy2 = prev->dy * prev->dy;
for (i = 2; i < curr->nx; i++) {
for (j = 2; j < curr->ny; j++) {
ic = idx(i, j, width);
iu = idx(i+1, j, width);
id = idx(i-1, j, width);
ir = idx(i, j+1, width);
il = idx(i, j-1, width);
curr->data[ic] = prev->data[ic] + a * dt *
((prev->data[iu] -
2.0 * prev->data[ic] +
prev->data[id]) / dx2 +
(prev->data[ir] -
2.0 * prev->data[ic] +
prev->data[il]) / dy2);
}
}
}
/* Update the temperature values using five-point stencil */
/* update only the border-dependent regions of the field */
void evolve_edges(field *curr, field *prev, double a, double dt)
{
int i, j;
int ic, iu, id, il, ir; // indexes for center, up, down, left, right
int width;
width = curr->ny + 2;
double dx2, dy2;
/* Determine the temperature field at next time step
* As we have fixed boundary conditions, the outermost gridpoints
* are not updated. */
dx2 = prev->dx * prev->dx;
dy2 = prev->dy * prev->dy;
i = 1;
for (j = 1; j < curr->ny + 1; j++) {
ic = idx(i, j, width);
iu = idx(i+1, j, width);
id = idx(i-1, j, width);
ir = idx(i, j+1, width);
il = idx(i, j-1, width);
curr->data[ic] = prev->data[ic] + a * dt *
((prev->data[iu] -
2.0 * prev->data[ic] +
prev->data[id]) / dx2 +
(prev->data[ir] -
2.0 * prev->data[ic] +
prev->data[il]) / dy2);
}
i = curr -> nx;
for (j = 1; j < curr->ny + 1; j++) {
ic = idx(i, j, width);
iu = idx(i+1, j, width);
id = idx(i-1, j, width);
ir = idx(i, j+1, width);
il = idx(i, j-1, width);
curr->data[ic] = prev->data[ic] + a * dt *
((prev->data[iu] -
2.0 * prev->data[ic] +
prev->data[id]) / dx2 +
(prev->data[ir] -
2.0 * prev->data[ic] +
prev->data[il]) / dy2);
}
j = 1;
for (i = 1; i < curr->nx + 1; i++) {
ic = idx(i, j, width);
iu = idx(i+1, j, width);
id = idx(i-1, j, width);
ir = idx(i, j+1, width);
il = idx(i, j-1, width);
curr->data[ic] = prev->data[ic] + a * dt *
((prev->data[iu] -
2.0 * prev->data[ic] +
prev->data[id]) / dx2 +
(prev->data[ir] -
2.0 * prev->data[ic] +
prev->data[il]) / dy2);
}
j = curr -> ny;
for (i = 1; i < curr->nx + 1; i++) {
ic = idx(i, j, width);
iu = idx(i+1, j, width);
id = idx(i-1, j, width);
ir = idx(i, j+1, width);
il = idx(i, j-1, width);
curr->data[ic] = prev->data[ic] + a * dt *
((prev->data[iu] -
2.0 * prev->data[ic] +
prev->data[id]) / dx2 +
(prev->data[ir] -
2.0 * prev->data[ic] +
prev->data[il]) / dy2);
}
}