-
Notifications
You must be signed in to change notification settings - Fork 0
/
nse3dMpiOpenMp.cpp
330 lines (260 loc) · 11 KB
/
nse3dMpiOpenMp.cpp
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
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
#include <fftw3-mpi.h>
#include <complex>
#include <vector>
#include <math.h>
#include <omp.h>
using namespace std;
//mpic++ -std=c++11 -O3 nse3dMpiOpenMp.cpp -o nse3dMpiOpenMp -lfftw3_mpi -lfftw3 -fopenmp -lm
//mpirun -np 4 ./nse3dMpiOpenMp
#define vd vector<double>
#define vcd vector<complex<double>>
// pseudospectral Fourier-Galerkin method
// uing 1d slab decompositon
int main() {
int rank, N, Nf;
double L, dx, nu, dt, T;
double t1, start_time;
MPI::Init();
fftw_mpi_init();
rank = MPI::COMM_WORLD.Get_rank();
double pi = atan(1.0)*4.0;;
long alloc_local, local_n0, local_0_start, local_n1, local_1_start, i, j, k;
vd s_in(1), s_out(1);
omp_set_num_threads(omp_get_max_threads());
N = pow(2,6); // --> change N here
nu = 0.000625;
T = 0.1;
dt = 0.001;
L = 2*pi;
dx = L / N;
if (rank==0)
std::cout << "N = " << N << std::endl;
vd a(4), b(3);
//coefficients for RK4 method
a[0] = 1./6.; a[1] = 1./3.; a[2] = 1./3.; a[3] = 1./6.;
b[0] = 0.5; b[1] = 0.5; b[2] = 1.0;
int tot = N*N*N;
Nf = N/2+1;
alloc_local = fftw_mpi_local_size_3d_transposed(N, N, Nf, MPI::COMM_WORLD,
&local_n0, &local_0_start,
&local_n1, &local_1_start);
vd U(2*alloc_local), V(2*alloc_local), W(2*alloc_local), U_tmp(2*alloc_local), V_tmp(2*alloc_local),
W_tmp(2*alloc_local), CU(2*alloc_local), CV(2*alloc_local), CW(2*alloc_local),
Fx(2*alloc_local),Fy(2*alloc_local),Fz(2*alloc_local);
vector<int> dealias(2*alloc_local); // for dealiasing
vd kk(2*alloc_local);
vcd U_hat(alloc_local), V_hat(alloc_local), W_hat(alloc_local), P_hat(alloc_local), U_hat0(alloc_local),
V_hat0(alloc_local), W_hat0(alloc_local), U_hat1(alloc_local), V_hat1(alloc_local), W_hat1(alloc_local),
dU(alloc_local), dV(alloc_local), dW(alloc_local), curlX(alloc_local), curlY(alloc_local), curlZ(alloc_local);
// Starting time
MPI::COMM_WORLD.Barrier();
start_time = MPI::Wtime();
vd kx(N), kz(Nf);
//fourier grid
for (i=0; i<N/2; i++)
{
kx[i] = i;
kz[i] = i;
}
kz[N/2] = N/2;
for (i=-N/2; i<0; i++)
kx[i+N] = i;
//fftw_plan plan_backward;
fftw_plan rfftn, irfftn;
rfftn = fftw_mpi_plan_dft_r2c_3d(N, N, N, U.data(), reinterpret_cast<fftw_complex*>(U_hat.data()),
MPI::COMM_WORLD, FFTW_MPI_TRANSPOSED_OUT);
irfftn = fftw_mpi_plan_dft_c2r_3d(N, N, N, reinterpret_cast<fftw_complex*>(U_hat.data()), U.data(),
MPI::COMM_WORLD, FFTW_MPI_TRANSPOSED_IN);
//Constant Force values
#pragma omp for
for (i=0; i<local_n0; i++)
for (j=0; j<N; j++)
for (k=0; k<N; k++) {
const int z = (i*N+j)*2*Nf+k;
Fx[z] = 0.01;
Fy[z] = 0.001;
Fz[z] = 0.0001;
}
// boundary conditions
#pragma omp for
for (i=0; i<local_n0; i++)
for (j=0; j<N; j++)
for (k=0; k<N; k++) {
const int z = (i*N+j)*2*Nf+k;
U[z] = sin(dx*(i+local_0_start))*cos(dx*j)*cos(dx*k); //---> initial boundary conditions
V[z] = -cos(dx*(i+local_0_start))*sin(dx*j)*cos(dx*k);
W[z] = 0.0;
}
//taking velocity field from physical to fourier
fftw_mpi_execute_dft_r2c( rfftn, U.data(), reinterpret_cast<fftw_complex*>(U_hat.data()));
fftw_mpi_execute_dft_r2c( rfftn, V.data(), reinterpret_cast<fftw_complex*>(V_hat.data()));
fftw_mpi_execute_dft_r2c( rfftn, W.data(), reinterpret_cast<fftw_complex*>(W_hat.data()));
double kmax = (2.0/3)*Nf;
//dealiasing using 2/3 rule
#pragma omp for
for (i=0; i<local_n1; i++)
for (j=0; j<N; j++)
for (k=0; k<Nf; k++)
{
const int z = (i*N+j)*Nf+k;
dealias[z] = (abs(kx[i+local_1_start])<kmax)*(abs(kx[j])<kmax)*(abs(kx[k])<kmax);
}
#pragma omp for
for (i=0; i<local_n1; i++)
for (j=0; j<N; j++)
for (k=0; k<Nf; k++)
{
const int z = (i*N+j)*Nf+k;
int m = kx[i+local_1_start]*kx[i+local_1_start] + kx[j]*kx[j] + kx[k]*kx[k];
kk[z] = m > 0 ? m : 1;
}
complex<double> one(0,1);
double t=0.0;
int tstep = 0;
//time stepping
while (t < T-1e-8) {
t += dt;
tstep++;
//copying fourier space velocity to mpi work arrays
std::copy(U_hat.begin(), U_hat.end(), U_hat0.begin());
std::copy(V_hat.begin(), V_hat.end(), V_hat0.begin());
std::copy(W_hat.begin(), W_hat.end(), W_hat0.begin());
std::copy(U_hat.begin(), U_hat.end(), U_hat1.begin());
std::copy(U_hat.begin(), U_hat.end(), U_hat1.begin());
std::copy(U_hat.begin(), U_hat.end(), U_hat1.begin());
//RK-4 Method
for (int rk=0; rk<4; rk++) {
if (rk > 0) {
fftw_mpi_execute_dft_c2r(irfftn, reinterpret_cast<fftw_complex*>(U_hat.data()), U.data());
fftw_mpi_execute_dft_c2r(irfftn, reinterpret_cast<fftw_complex*>(V_hat.data()), V.data());
fftw_mpi_execute_dft_c2r(irfftn, reinterpret_cast<fftw_complex*>(W_hat.data()), W.data());
for (k=0; k<U.size(); k++) {
U[k] /= tot;
V[k] /= tot;
W[k] /= tot;
}
}
// Compute curl
#pragma omp for
for (i=0; i<local_n1; i++)
for (j=0; j<N; j++)
for (k=0; k<Nf; k++) {
const int z = (i*N+j)*Nf+k;
curlZ[z] = one*(kx[i+local_1_start]*V_hat[z]-kx[j]*U_hat[z]);
curlY[z] = one*(kz[k]*U_hat[z]-kx[i+local_1_start]*W_hat[z]);
curlX[z] = one*(kx[j]*W_hat[z]-kz[k]*V_hat[z]);
}
// taking curl into fourier space
fftw_mpi_execute_dft_c2r(irfftn, reinterpret_cast<fftw_complex*>(curlX.data()), CU.data());
fftw_mpi_execute_dft_c2r(irfftn, reinterpret_cast<fftw_complex*>(curlY.data()), CV.data());
fftw_mpi_execute_dft_c2r(irfftn, reinterpret_cast<fftw_complex*>(curlZ.data()), CW.data());
for (k=0; k<CU.size(); k++)
{
CU[k] /= tot;
CV[k] /= tot;
CW[k] /= tot;
}
// Cross
#pragma omp for
for (i=0; i<local_n0; i++)
for (j=0; j<N; j++)
for (k=0; k<N; k++) {
const int z = (i*N+j)*2*Nf+k;
U_tmp[z] = V[z]*CW[z]-W[z]*CV[z]+Fx[z];
V_tmp[z] = W[z]*CU[z]-U[z]*CW[z]+Fy[z];
W_tmp[z] = U[z]*CV[z]-V[z]*CU[z]+Fz[z];
}
//taking cross into fourier space
fftw_mpi_execute_dft_r2c( rfftn, U_tmp.data(), reinterpret_cast<fftw_complex*>(dU.data()));
fftw_mpi_execute_dft_r2c( rfftn, V_tmp.data(), reinterpret_cast<fftw_complex*>(dV.data()));
fftw_mpi_execute_dft_r2c( rfftn, W_tmp.data(), reinterpret_cast<fftw_complex*>(dW.data()));
//dealiasing
#pragma omp for
for (i=0; i<local_n1; i++)
for (j=0; j<N; j++)
for (k=0; k<Nf; k++) {
const int z = (i*N+j)*Nf+k;
dU[z] *= (dealias[z]*dt);
dV[z] *= (dealias[z]*dt);
dW[z] *= (dealias[z]*dt);
}
//taking care of pressure and other linear terms
#pragma omp for
for (i=0; i<local_n1; i++)
for (j=0; j<N; j++)
for (k=0; k<Nf; k++) {
const int z = (i*N+j)*Nf+k;
P_hat[z] = (dU[z]*kx[i+local_1_start] + dV[z]*kx[j] + dW[z]*kz[k])/kk[z];
dU[z] -= (P_hat[z]*kx[i+local_1_start] + nu*dt*kk[z]*U_hat[z]);
dV[z] -= (P_hat[z]*kx[j] + nu*dt*kk[z]*V_hat[z]);
dW[z] -= (P_hat[z]*kz[k] + nu*dt*kk[z]*W_hat[z]);
}
if (rk < 3) {
#pragma omp for
for (i=0; i<local_n1; i++)
for (j=0; j<N; j++)
for (k=0; k<Nf; k++) {
const int z = (i*N+j)*Nf+k;
U_hat[z] = U_hat0[z] + b[rk]*dU[z];
V_hat[z] = V_hat0[z] + b[rk]*dV[z];
W_hat[z] = W_hat0[z] + b[rk]*dW[z];
}
}
#pragma omp for
for (i=0; i<local_n1; i++)
for (j=0; j<N; j++)
for (k=0; k<Nf; k++) {
const int z = (i*N+j)*Nf+k;
U_hat1[z] += a[rk]*dU[z];
V_hat1[z] += a[rk]*dV[z];
W_hat1[z] += a[rk]*dW[z];
}
}
//copying the temp fourier velocity into the hat arrays
#pragma omp for
for (i=0; i<local_n1; i++)
for (j=0; j<N; j++)
for (k=0; k<Nf; k++) {
const int z = (i*N+j)*Nf+k;
U_hat[z] = U_hat1[z];
V_hat[z] = V_hat1[z];
W_hat[z] = W_hat1[z];
}
//printing the total energy of the system at every even step
if (tstep % 2 == 0) {
s_in[0] = 0.0;
#pragma omp for
for (i=0; i<local_n0; i++)
for (j=0; j<N; j++)
for (k=0; k<N; k++) {
int z = (i*N+j)*2*Nf+k;
s_in[0] += (U[z]*U[z] + V[z]*V[z] + W[z]*W[z]);
}
s_in[0] *= (0.5*dx*dx*dx/L/L/L);
MPI::COMM_WORLD.Reduce(s_in.data(), s_out.data(), 1, MPI::DOUBLE, MPI::SUM, 0);
// if (rank==0)
// std::cout << " k = " << s_out[0] << std::endl;
}
}
s_in[0] = 0.0;
//get final velocities from U,V and W
//calculating total final energy
#pragma omp for
for (i=0; i<local_n0; i++)
for (j=0; j<N; j++)
for (k=0; k<N; k++) {
int z = (i*N+j)*2*Nf+k;
s_in[0] += (U[z]*U[z] + V[z]*V[z] + W[z]*W[z]);
}
s_in[0] *= (0.5*dx*dx*dx/L/L/L);
MPI::COMM_WORLD.Reduce(s_in.data(), s_out.data(), 1, MPI::DOUBLE, MPI::SUM, 0);
if (rank==0)
std::cout << "\nFinal Energy: \n k = " << s_out[0] << std::endl;
MPI::COMM_WORLD.Barrier();
t1 = MPI::Wtime();
if (rank == 0)
std::cout << "Total time taken = " << t1 - start_time << "s" << std::endl;
fftw_destroy_plan(rfftn);
fftw_destroy_plan(irfftn);
MPI::Finalize ( );
}