-
Notifications
You must be signed in to change notification settings - Fork 15
/
Copy pathnagumo.cpp
294 lines (247 loc) · 10.4 KB
/
nagumo.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
// Copyright 2018-2024 the samurai's authors
// SPDX-License-Identifier: BSD-3-Clause
#include <samurai/hdf5.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_);
}
int main(int argc, char* argv[])
{
auto& app = samurai::initialize("Finite volume example for the Nagumo equation", argc, argv);
static constexpr std::size_t dim = 1;
static constexpr std::size_t field_size = 1;
using Config = samurai::MRConfig<dim>;
using Box = samurai::Box<double, dim>;
using point_t = typename Box::point_t;
std::cout << "------------------------- Nagumo -------------------------" << std::endl;
/**
* Nagumo, or Fisher-KPP equation:
*
* du/dt -D*Lap(u) = k u^2 (1-u)
*/
//--------------------//
// Program parameters //
//--------------------//
// Simulation parameters
double left_box = -10;
double right_box = 10;
double D = 1;
double k = 10;
bool explicit_diffusion = false;
bool explicit_reaction = false;
// Time integration
double Tf = 1.;
double dt = 0;
double cfl = 0.95;
// Multiresolution parameters
std::size_t min_level = 4;
std::size_t max_level = 8;
double mr_epsilon = 1e-5; // Threshold used by multiresolution
double mr_regularity = 1.; // Regularity guess for multiresolution
// Output parameters
fs::path path = fs::current_path();
std::string filename = "nagumo";
bool save_final_state_only = false;
app.add_option("--left", left_box, "The left border of the box")->capture_default_str()->group("Simulation parameters");
app.add_option("--right", right_box, "The right border of the box")->capture_default_str()->group("Simulation parameters");
app.add_option("--D", D, "Diffusion coefficient")->capture_default_str()->group("Simulation parameters");
app.add_option("--k", k, "Parameter of the reaction operator")->capture_default_str()->group("Simulation parameters");
app.add_option("--Tf", Tf, "Final time")->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_flag("--explicit-reaction", explicit_reaction, "Explicit the reaction term")->capture_default_str()->group("Simulation parameters");
app.add_flag("--explicit-diffusion", explicit_diffusion, "Explicit the diffusion term")->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 //
//------------------//
samurai::times::timers.start("petsc init");
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
samurai::times::timers.stop("petsc init");
//--------------------//
// 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{box, min_level, max_level};
auto u = samurai::make_field<field_size>("u", mesh);
double z0 = left_box / 5; // wave initial position
double c = sqrt(k * D / 2); // wave velocity
auto beta = [&](double z)
{
double e = exp(-sqrt(k / (2 * D)) * (z - z0));
return e / (1 + e);
};
auto exact_solution = [&](double x, double t)
{
return beta(x - c * t);
};
// Initial solution
samurai::for_each_cell(mesh,
[&](auto& cell)
{
u[cell] = exact_solution(cell.center(0), 0);
});
auto unp1 = samurai::make_field<field_size>("unp1", mesh);
samurai::make_bc<samurai::Neumann<1>>(u);
samurai::make_bc<samurai::Neumann<1>>(unp1);
auto diff = samurai::make_diffusion_order2<decltype(u)>(D);
auto id = samurai::make_identity<decltype(u)>();
// Reaction operator
using cfg = samurai::LocalCellSchemeConfig<samurai::SchemeType::NonLinear, field_size, decltype(u)>;
auto react = samurai::make_cell_based_scheme<cfg>();
react.set_name("Reaction");
react.set_scheme_function(
[&](const auto& cell, const auto& field) -> samurai::SchemeValue<cfg>
{
auto v = field[cell];
return k * v * v * (1 - v);
});
react.set_jacobian_function(
[&](const auto& cell, const auto& field) -> samurai::JacobianMatrix<cfg>
{
auto v = field[cell];
return k * (2 * v * (1 - v) - v * v);
});
//--------------------//
// Time iteration //
//--------------------//
if (dt == 0)
{
dt = Tf / 100;
if (explicit_diffusion)
{
double dx = mesh.cell_length(max_level);
dt = cfl * (dx * dx) / (pow(2, dim) * D);
}
}
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++));
}
auto rhs = samurai::make_field<field_size>("rhs", mesh);
double t = 0;
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();
rhs.resize();
if (explicit_diffusion && explicit_reaction)
{
unp1 = u + dt * react(u) - dt * diff(u);
}
else if (!explicit_diffusion && explicit_reaction)
{
// u_np1 + dt*diff(u_np1) = u + dt*react(u)
auto implicit_operator = id + dt * diff;
rhs = u + dt * react(u);
// Solve the linear equation [Id + dt*Diff](unp1) = rhs
samurai::petsc::solve(implicit_operator, unp1, rhs);
}
else if (explicit_diffusion && !explicit_reaction)
{
// u_np1 - dt*react(u_np1) = u - dt*diff(u)
auto implicit_operator = id - dt * react;
rhs = u - dt * diff(u);
// Set initial guess for the Newton algorithm
unp1 = u;
// Solve the non-linear equation [Id - dt*React](unp1) = u - dt*Diff(u)
// Here, small independent local Newton methods are used.
samurai::petsc::solve(implicit_operator, unp1, rhs);
}
else
{
// u_np1 + dt*diff(u_np1) - dt*react(u_np1) = u
auto implicit_operator = id + dt * diff - dt * react;
// Set initial guess for the Newton algorithm
unp1 = u;
// Solve the non-linear equation [Id + dt*Diff - dt*React](unp1) = u
samurai::petsc::solve(implicit_operator, unp1, u);
}
// u <-- unp1
std::swap(u.array(), unp1.array());
// Compute error
double error = samurai::L2_error(u,
[&](const auto& coord)
{
return exact_solution(coord(0), t);
});
std::cout.precision(2);
std::cout << ", L2-error: " << std::scientific << error;
// 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 && dim == 1)
{
std::cout << std::endl;
std::cout << "Run the following command to view the results:" << std::endl;
std::cout << "python <<path to samurai>>/python/read_mesh.py " << filename << "_ite_ --field u level --start 1 --end " << nsave
<< std::endl;
}
if (save_final_state_only)
{
save(path, filename, u);
}
samurai::times::timers.start("petsc finalize");
PetscFinalize();
samurai::times::timers.stop("petsc finalize");
samurai::finalize();
return 0;
}