-
Notifications
You must be signed in to change notification settings - Fork 16
/
Copy pathheat_heterogeneous.cpp
237 lines (199 loc) · 8.19 KB
/
heat_heterogeneous.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
// Copyright 2018-2025 the samurai's authors
// SPDX-License-Identifier: BSD-3-Clause
#include <samurai/io/hdf5.hpp>
#include <samurai/io/restart.hpp>
#include <samurai/mr/adapt.hpp>
#include <samurai/mr/mesh.hpp>
#include <samurai/petsc.hpp>
#include <samurai/samurai.hpp>
#include <filesystem>
namespace fs = std::filesystem;
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;
});
samurai::save(path, fmt::format("{}{}", filename, suffix), mesh, u, level_);
samurai::dump(path, fmt::format("{}_restart{}", filename, suffix), mesh, u);
}
int main(int argc, char* argv[])
{
auto& app = samurai::initialize("Finite volume example for the heat equation", argc, argv);
static constexpr std::size_t dim = 2;
using Config = samurai::MRConfig<dim>;
using Box = samurai::Box<double, dim>;
using point_t = typename Box::point_t;
std::cout << "------------------------- Heat -------------------------" << std::endl;
//--------------------//
// Program parameters //
//--------------------//
// Simulation parameters
double left_box = -4;
double right_box = 4;
// Time integration
double Tf = 1.;
double dt = Tf / 100;
bool explicit_scheme = false;
double cfl = 0.95;
double t = 0.;
std::string restart_file;
// Multiresolution parameters
std::size_t min_level = 3;
std::size_t max_level = 6;
double mr_epsilon = 1e-4; // Threshold used by multiresolution
double mr_regularity = 1.; // Regularity guess for multiresolution
// Output parameters
fs::path path = fs::current_path();
std::string filename = "heat_heterog_" + std::to_string(dim) + "D";
bool save_final_state_only = false;
app.add_flag("--explicit", explicit_scheme, "Explicit scheme instead of implicit")->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("--dt", dt, "Time step")->capture_default_str()->group("Simulation parameters");
app.add_option("--cfl", cfl, "The CFL")->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("--path", path, "Output path")->capture_default_str()->group("Output");
app.add_option("--filename", filename, "File name prefix")->capture_default_str()->group("Output");
app.add_flag("--save-final-state-only", save_final_state_only, "Save final state only")->group("Output");
app.allow_extras();
SAMURAI_PARSE(argc, argv);
//------------------//
// Petsc initialize //
//------------------//
PetscInitialize(&argc, &argv, 0, nullptr);
PetscMPIInt size;
PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
PetscOptionsSetValue(NULL, "-options_left", "off"); // If on, Petsc will issue warnings saying that
// the options managed by CLI are unused
//--------------------//
// Problem definition //
//--------------------//
point_t box_corner1, box_corner2;
box_corner1.fill(left_box);
box_corner2.fill(right_box);
Box box(box_corner1, box_corner2);
samurai::MRMesh<Config> mesh;
auto u = samurai::make_field<1>("u", mesh);
auto unp1 = samurai::make_field<1>("unp1", mesh);
if (restart_file.empty())
{
mesh = {box, min_level, max_level};
u.resize();
// Initial solution: crenel
samurai::for_each_cell(mesh,
[&](auto& cell)
{
bool is_in_crenel = true;
for (std::size_t d = 0; d < dim; ++d)
{
is_in_crenel = is_in_crenel && (abs(cell.center(d)) < right_box / 3);
}
u[cell] = is_in_crenel ? 1 : 0;
});
}
else
{
samurai::load(restart_file, mesh, u);
}
samurai::make_bc<samurai::Neumann<1>>(u, 0.);
samurai::make_bc<samurai::Neumann<1>>(unp1, 0.);
auto K = samurai::make_field<samurai::DiffCoeff<dim>, 1>("K", mesh);
auto set_K_values = [&]()
{
samurai::for_each_cell(mesh[decltype(mesh)::mesh_id_t::reference],
[&](auto& cell)
{
double x = cell.center(0);
double y = cell.center(1);
// Domain divided into 4 quadrants
if ((x < 0 && y > 0) || (x > 0 && y < 0))
{
K[cell] = {0.1, 4};
}
else
{
K[cell] = {4, 0.1};
}
});
};
auto diff = samurai::make_diffusion_order2<decltype(u)>(K);
auto id = samurai::make_identity<decltype(u)>();
//--------------------//
// Time iteration //
//--------------------//
if (explicit_scheme)
{
double diff_coeff = 4;
double dx = mesh.cell_length(max_level);
dt = cfl * (dx * dx) / (pow(2, dim) * diff_coeff);
}
auto MRadaptation = samurai::make_MRAdapt(u);
MRadaptation(mr_epsilon, mr_regularity);
std::size_t nsave = 0, nt = 0;
if (!save_final_state_only)
{
save(path, filename, u, fmt::format("_ite_{}", nsave++));
}
while (t != Tf)
{
// Move to next timestep
t += dt;
if (t > Tf)
{
dt += Tf - t;
t = Tf;
}
std::cout << fmt::format("iteration {}: t = {:.2f}, dt = {}", nt++, t, dt) << std::flush;
// Mesh adaptation
MRadaptation(mr_epsilon, mr_regularity);
samurai::update_ghost_mr(u);
unp1.resize();
K.resize();
set_K_values();
if (explicit_scheme)
{
unp1 = u - dt * diff(u);
}
else
{
auto back_euler = id + dt * diff;
samurai::petsc::solve(back_euler, unp1, u); // solves the linear equation [Id + dt*Diff](unp1) = u
}
// u <-- unp1
std::swap(u.array(), unp1.array());
// Save the result
if (!save_final_state_only)
{
save(path, filename, u, fmt::format("_ite_{}", nsave++));
}
std::cout << std::endl;
}
if (save_final_state_only)
{
save(path, filename, u);
}
PetscFinalize();
samurai::finalize();
return 0;
}