-
Notifications
You must be signed in to change notification settings - Fork 16
/
Copy pathadvection_2d.cpp
285 lines (239 loc) · 10.6 KB
/
advection_2d.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
// Copyright 2018-2025 the samurai's authors
// SPDX-License-Identifier: BSD-3-Clause
#include <array>
#include <xtensor/xfixed.hpp>
#include <samurai/algorithm.hpp>
#include <samurai/bc.hpp>
#include <samurai/field.hpp>
#include <samurai/io/hdf5.hpp>
#include <samurai/io/restart.hpp>
#include <samurai/mr/adapt.hpp>
#include <samurai/mr/mesh.hpp>
#include <samurai/samurai.hpp>
#include <samurai/stencil_field.hpp>
#include <samurai/subset/node.hpp>
#include <filesystem>
namespace fs = std::filesystem;
template <class Field>
void init(Field& u)
{
auto& mesh = u.mesh();
u.resize();
samurai::for_each_cell(
mesh,
[&](auto& cell)
{
auto center = cell.center();
const double radius = .2;
const double x_center = 0.3;
const double y_center = 0.3;
if (((center[0] - x_center) * (center[0] - x_center) + (center[1] - y_center) * (center[1] - y_center)) <= radius * radius)
{
u[cell] = 1;
}
else
{
u[cell] = 0;
}
});
}
template <class Field>
void flux_correction(double dt, const std::array<double, 2>& a, const Field& u, Field& unp1)
{
using mesh_t = typename Field::mesh_t;
using mesh_id_t = typename mesh_t::mesh_id_t;
using interval_t = typename mesh_t::interval_t;
constexpr std::size_t dim = Field::dim;
auto mesh = u.mesh();
for (std::size_t level = mesh.min_level(); level < mesh.max_level(); ++level)
{
xt::xtensor_fixed<int, xt::xshape<dim>> stencil;
stencil = {
{-1, 0}
};
auto subset_right = samurai::intersection(samurai::translate(mesh[mesh_id_t::cells][level + 1], stencil),
mesh[mesh_id_t::cells][level])
.on(level);
subset_right(
[&](const auto& i, const auto& index)
{
auto j = index[0];
const double dx = mesh.cell_length(level);
unp1(level, i, j) = unp1(level, i, j)
+ dt / dx
* (samurai::upwind_op<dim, interval_t>(level, i, j).right_flux(a, u)
- .5 * samurai::upwind_op<dim, interval_t>(level + 1, 2 * i + 1, 2 * j).right_flux(a, u)
- .5 * samurai::upwind_op<dim, interval_t>(level + 1, 2 * i + 1, 2 * j + 1).right_flux(a, u));
});
stencil = {
{1, 0}
};
auto subset_left = samurai::intersection(samurai::translate(mesh[mesh_id_t::cells][level + 1], stencil),
mesh[mesh_id_t::cells][level])
.on(level);
subset_left(
[&](const auto& i, const auto& index)
{
auto j = index[0];
const double dx = mesh.cell_length(level);
unp1(level, i, j) = unp1(level, i, j)
- dt / dx
* (samurai::upwind_op<dim, interval_t>(level, i, j).left_flux(a, u)
- .5 * samurai::upwind_op<dim, interval_t>(level + 1, 2 * i, 2 * j).left_flux(a, u)
- .5 * samurai::upwind_op<dim, interval_t>(level + 1, 2 * i, 2 * j + 1).left_flux(a, u));
});
stencil = {
{0, -1}
};
auto subset_up = samurai::intersection(samurai::translate(mesh[mesh_id_t::cells][level + 1], stencil), mesh[mesh_id_t::cells][level])
.on(level);
subset_up(
[&](const auto& i, const auto& index)
{
auto j = index[0];
const double dx = mesh.cell_length(level);
unp1(level, i, j) = unp1(level, i, j)
+ dt / dx
* (samurai::upwind_op<dim, interval_t>(level, i, j).up_flux(a, u)
- .5 * samurai::upwind_op<dim, interval_t>(level + 1, 2 * i, 2 * j + 1).up_flux(a, u)
- .5 * samurai::upwind_op<dim, interval_t>(level + 1, 2 * i + 1, 2 * j + 1).up_flux(a, u));
});
stencil = {
{0, 1}
};
auto subset_down = samurai::intersection(samurai::translate(mesh[mesh_id_t::cells][level + 1], stencil),
mesh[mesh_id_t::cells][level])
.on(level);
subset_down(
[&](const auto& i, const auto& index)
{
auto j = index[0];
const double dx = mesh.cell_length(level);
unp1(level, i, j) = unp1(level, i, j)
- dt / dx
* (samurai::upwind_op<dim, interval_t>(level, i, j).down_flux(a, u)
- .5 * samurai::upwind_op<dim, interval_t>(level + 1, 2 * i, 2 * j).down_flux(a, u)
- .5 * samurai::upwind_op<dim, interval_t>(level + 1, 2 * i + 1, 2 * j).down_flux(a, u));
});
}
}
template <class Field>
void save(const fs::path& path, const std::string& filename, const Field& u, const std::string& suffix = "")
{
auto mesh = u.mesh();
auto level_ = samurai::make_field<std::size_t, 1>("level", mesh);
if (!fs::exists(path))
{
fs::create_directory(path);
}
samurai::for_each_cell(mesh,
[&](const auto& cell)
{
level_[cell] = cell.level;
});
#ifdef SAMURAI_WITH_MPI
mpi::communicator world;
samurai::save(path, fmt::format("{}_size_{}{}", filename, world.size(), suffix), mesh, u, level_);
#else
samurai::save(path, fmt::format("{}{}", filename, suffix), mesh, u, level_);
samurai::dump(path, fmt::format("{}_restart{}", filename, suffix), mesh, u);
#endif
}
int main(int argc, char* argv[])
{
auto& app = samurai::initialize("Finite volume example for the advection equation in 2d using multiresolution", argc, argv);
constexpr std::size_t dim = 2;
using Config = samurai::MRConfig<dim>;
// Simulation parameters
xt::xtensor_fixed<double, xt::xshape<dim>> min_corner = {0., 0.};
xt::xtensor_fixed<double, xt::xshape<dim>> max_corner = {1., 1.};
std::array<double, dim> a{
{1, 1}
};
double Tf = .1;
double cfl = 0.5;
double t = 0.;
std::string restart_file;
// Multiresolution parameters
std::size_t min_level = 4;
std::size_t max_level = 10;
double mr_epsilon = 2.e-4; // Threshold used by multiresolution
double mr_regularity = 1.; // Regularity guess for multiresolution
bool correction = false;
// Output parameters
fs::path path = fs::current_path();
std::string filename = "FV_advection_2d";
std::size_t nfiles = 1;
app.add_option("--min-corner", min_corner, "The min corner of the box")->capture_default_str()->group("Simulation parameters");
app.add_option("--max-corner", max_corner, "The max corner of the box")->capture_default_str()->group("Simulation parameters");
app.add_option("--velocity", a, "The velocity of the advection equation")->capture_default_str()->group("Simulation parameters");
app.add_option("--cfl", cfl, "The CFL")->capture_default_str()->group("Simulation parameters");
app.add_option("--Ti", t, "Initial time")->capture_default_str()->group("Simulation parameters");
app.add_option("--Tf", Tf, "Final time")->capture_default_str()->group("Simulation parameters");
app.add_option("--restart-file", restart_file, "Restart file")->capture_default_str()->group("Simulation parameters");
app.add_option("--min-level", min_level, "Minimum level of the multiresolution")->capture_default_str()->group("Multiresolution");
app.add_option("--max-level", max_level, "Maximum level of the multiresolution")->capture_default_str()->group("Multiresolution");
app.add_option("--mr-eps", mr_epsilon, "The epsilon used by the multiresolution to adapt the mesh")
->capture_default_str()
->group("Multiresolution");
app.add_option("--mr-reg",
mr_regularity,
"The regularity criteria used by the multiresolution to "
"adapt the mesh")
->capture_default_str()
->group("Multiresolution");
app.add_option("--with-correction", correction, "Apply flux correction at the interface of two refinement levels")
->capture_default_str()
->group("Multiresolution");
app.add_option("--path", path, "Output path")->capture_default_str()->group("Output");
app.add_option("--filename", filename, "File name prefix")->capture_default_str()->group("Output");
app.add_option("--nfiles", nfiles, "Number of output files")->capture_default_str()->group("Output");
SAMURAI_PARSE(argc, argv);
const samurai::Box<double, dim> box(min_corner, max_corner);
samurai::MRMesh<Config> mesh;
auto u = samurai::make_field<double, 1>("u", mesh);
if (restart_file.empty())
{
mesh = {box, min_level, max_level};
init(u);
}
else
{
samurai::load(restart_file, mesh, u);
}
samurai::make_bc<samurai::Dirichlet<1>>(u, 0.);
double dt = cfl * mesh.cell_length(max_level);
const double dt_save = Tf / static_cast<double>(nfiles);
auto unp1 = samurai::make_field<double, 1>("unp1", mesh);
auto MRadaptation = samurai::make_MRAdapt(u);
MRadaptation(mr_epsilon, mr_regularity);
save(path, filename, u, "_init");
std::size_t nsave = 1;
std::size_t nt = 0;
while (t != Tf)
{
MRadaptation(mr_epsilon, mr_regularity);
t += dt;
if (t > Tf)
{
dt += Tf - t;
t = Tf;
}
std::cout << fmt::format("iteration {}: t = {}, dt = {}", nt++, t, dt) << std::endl;
samurai::update_ghost_mr(u);
unp1.resize();
unp1 = u - dt * samurai::upwind(a, u);
if (correction)
{
flux_correction(dt, a, u, unp1);
}
std::swap(u.array(), unp1.array());
if (t >= static_cast<double>(nsave + 1) * dt_save || t == Tf)
{
const std::string suffix = (nfiles != 1) ? fmt::format("_ite_{}", nsave++) : "";
save(path, filename, u, suffix);
}
}
samurai::finalize();
return 0;
}