forked from mfem/mfem
-
Notifications
You must be signed in to change notification settings - Fork 0
/
ex16.cpp
374 lines (321 loc) · 11 KB
/
ex16.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
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
// MFEM Example 16
//
// Compile with: make ex16
//
// Sample runs: ex16
// ex16 -m ../data/inline-tri.mesh
// ex16 -m ../data/disc-nurbs.mesh -tf 2
// ex16 -s 21 -a 0.0 -k 1.0
// ex16 -s 22 -a 1.0 -k 0.0
// ex16 -s 23 -a 0.5 -k 0.5 -o 4
// ex16 -s 4 -dt 1.0e-4 -tf 4.0e-2 -vs 40
// ex16 -m ../data/fichera-q2.mesh
// ex16 -m ../data/fichera-mixed.mesh
// ex16 -m ../data/escher.mesh
// ex16 -m ../data/beam-tet.mesh -tf 10 -dt 0.1
// ex16 -m ../data/amr-quad.mesh -o 4 -r 0
// ex16 -m ../data/amr-hex.mesh -o 2 -r 0
//
// Description: This example solves a time dependent nonlinear heat equation
// problem of the form du/dt = C(u), with a non-linear diffusion
// operator C(u) = \nabla \cdot (\kappa + \alpha u) \nabla u.
//
// The example demonstrates the use of nonlinear operators (the
// class ConductionOperator defining C(u)), as well as their
// implicit time integration. Note that implementing the method
// ConductionOperator::ImplicitSolve is the only requirement for
// high-order implicit (SDIRK) time integration. In this example,
// the diffusion operator is linearized by evaluating with the
// lagged solution from the previous timestep, so there is only
// a linear solve.
//
// We recommend viewing examples 2, 9 and 10 before viewing this
// example.
#include "mfem.hpp"
#include <fstream>
#include <iostream>
using namespace std;
using namespace mfem;
/** After spatial discretization, the conduction model can be written as:
*
* du/dt = M^{-1}(-Ku)
*
* where u is the vector representing the temperature, M is the mass matrix,
* and K is the diffusion operator with diffusivity depending on u:
* (\kappa + \alpha u).
*
* Class ConductionOperator represents the right-hand side of the above ODE.
*/
class ConductionOperator : public TimeDependentOperator
{
protected:
FiniteElementSpace &fespace;
Array<int> ess_tdof_list; // this list remains empty for pure Neumann b.c.
BilinearForm *M;
BilinearForm *K;
SparseMatrix Mmat, Kmat;
SparseMatrix *T; // T = M + dt K
real_t current_dt;
CGSolver M_solver; // Krylov solver for inverting the mass matrix M
DSmoother M_prec; // Preconditioner for the mass matrix M
CGSolver T_solver; // Implicit solver for T = M + dt K
DSmoother T_prec; // Preconditioner for the implicit solver
real_t alpha, kappa;
mutable Vector z; // auxiliary vector
public:
ConductionOperator(FiniteElementSpace &f, real_t alpha, real_t kappa,
const Vector &u);
void Mult(const Vector &u, Vector &du_dt) const override;
/** Solve the Backward-Euler equation: k = f(u + dt*k, t), for the unknown k.
This is the only requirement for high-order SDIRK implicit integration.*/
void ImplicitSolve(const real_t dt, const Vector &u, Vector &k) override;
/// Update the diffusion BilinearForm K using the given true-dof vector `u`.
void SetParameters(const Vector &u);
~ConductionOperator() override;
};
real_t InitialTemperature(const Vector &x);
int main(int argc, char *argv[])
{
// 1. Parse command-line options.
const char *mesh_file = "../data/star.mesh";
int ref_levels = 2;
int order = 2;
int ode_solver_type = 23; // SDIRK33Solver
real_t t_final = 0.5;
real_t dt = 1.0e-2;
real_t alpha = 1.0e-2;
real_t kappa = 0.5;
bool visualization = true;
bool visit = false;
int vis_steps = 5;
int precision = 8;
cout.precision(precision);
OptionsParser args(argc, argv);
args.AddOption(&mesh_file, "-m", "--mesh",
"Mesh file to use.");
args.AddOption(&ref_levels, "-r", "--refine",
"Number of times to refine the mesh uniformly.");
args.AddOption(&order, "-o", "--order",
"Order (degree) of the finite elements.");
args.AddOption(&ode_solver_type, "-s", "--ode-solver",
ODESolver::Types.c_str());
args.AddOption(&t_final, "-tf", "--t-final",
"Final time; start time is 0.");
args.AddOption(&dt, "-dt", "--time-step",
"Time step.");
args.AddOption(&alpha, "-a", "--alpha",
"Alpha coefficient.");
args.AddOption(&kappa, "-k", "--kappa",
"Kappa coefficient offset.");
args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
"--no-visualization",
"Enable or disable GLVis visualization.");
args.AddOption(&visit, "-visit", "--visit-datafiles", "-no-visit",
"--no-visit-datafiles",
"Save data files for VisIt (visit.llnl.gov) visualization.");
args.AddOption(&vis_steps, "-vs", "--visualization-steps",
"Visualize every n-th timestep.");
args.Parse();
if (!args.Good())
{
args.PrintUsage(cout);
return 1;
}
args.PrintOptions(cout);
// 2. Read the mesh from the given mesh file. We can handle triangular,
// quadrilateral, tetrahedral and hexahedral meshes with the same code.
Mesh *mesh = new Mesh(mesh_file, 1, 1);
int dim = mesh->Dimension();
// 3. Define the ODE solver used for time integration. Several implicit
// singly diagonal implicit Runge-Kutta (SDIRK) methods, as well as
// explicit Runge-Kutta methods are available.
unique_ptr<ODESolver> ode_solver = ODESolver::Select(ode_solver_type);
// 4. Refine the mesh to increase the resolution. In this example we do
// 'ref_levels' of uniform refinement, where 'ref_levels' is a
// command-line parameter.
for (int lev = 0; lev < ref_levels; lev++)
{
mesh->UniformRefinement();
}
// 5. Define the vector finite element space representing the current and the
// initial temperature, u_ref.
H1_FECollection fe_coll(order, dim);
FiniteElementSpace fespace(mesh, &fe_coll);
int fe_size = fespace.GetTrueVSize();
cout << "Number of temperature unknowns: " << fe_size << endl;
GridFunction u_gf(&fespace);
// 6. Set the initial conditions for u. All boundaries are considered
// natural.
FunctionCoefficient u_0(InitialTemperature);
u_gf.ProjectCoefficient(u_0);
Vector u;
u_gf.GetTrueDofs(u);
// 7. Initialize the conduction operator and the visualization.
ConductionOperator oper(fespace, alpha, kappa, u);
u_gf.SetFromTrueDofs(u);
{
ofstream omesh("ex16.mesh");
omesh.precision(precision);
mesh->Print(omesh);
ofstream osol("ex16-init.gf");
osol.precision(precision);
u_gf.Save(osol);
}
VisItDataCollection visit_dc("Example16", mesh);
visit_dc.RegisterField("temperature", &u_gf);
if (visit)
{
visit_dc.SetCycle(0);
visit_dc.SetTime(0.0);
visit_dc.Save();
}
socketstream sout;
if (visualization)
{
char vishost[] = "localhost";
int visport = 19916;
sout.open(vishost, visport);
if (!sout)
{
cout << "Unable to connect to GLVis server at "
<< vishost << ':' << visport << endl;
visualization = false;
cout << "GLVis visualization disabled.\n";
}
else
{
sout.precision(precision);
sout << "solution\n" << *mesh << u_gf;
sout << "pause\n";
sout << flush;
cout << "GLVis visualization paused."
<< " Press space (in the GLVis window) to resume it.\n";
}
}
// 8. Perform time-integration (looping over the time iterations, ti, with a
// time-step dt).
ode_solver->Init(oper);
real_t t = 0.0;
bool last_step = false;
for (int ti = 1; !last_step; ti++)
{
if (t + dt >= t_final - dt/2)
{
last_step = true;
}
ode_solver->Step(u, t, dt);
if (last_step || (ti % vis_steps) == 0)
{
cout << "step " << ti << ", t = " << t << endl;
u_gf.SetFromTrueDofs(u);
if (visualization)
{
sout << "solution\n" << *mesh << u_gf << flush;
}
if (visit)
{
visit_dc.SetCycle(ti);
visit_dc.SetTime(t);
visit_dc.Save();
}
}
oper.SetParameters(u);
}
// 9. Save the final solution. This output can be viewed later using GLVis:
// "glvis -m ex16.mesh -g ex16-final.gf".
{
ofstream osol("ex16-final.gf");
osol.precision(precision);
u_gf.Save(osol);
}
// 10. Free the used memory.
delete mesh;
return 0;
}
ConductionOperator::ConductionOperator(FiniteElementSpace &f, real_t al,
real_t kap, const Vector &u)
: TimeDependentOperator(f.GetTrueVSize(), (real_t) 0.0), fespace(f),
M(NULL), K(NULL), T(NULL), current_dt(0.0), z(height)
{
const real_t rel_tol = 1e-8;
M = new BilinearForm(&fespace);
M->AddDomainIntegrator(new MassIntegrator());
M->Assemble();
M->FormSystemMatrix(ess_tdof_list, Mmat);
M_solver.iterative_mode = false;
M_solver.SetRelTol(rel_tol);
M_solver.SetAbsTol(0.0);
M_solver.SetMaxIter(30);
M_solver.SetPrintLevel(0);
M_solver.SetPreconditioner(M_prec);
M_solver.SetOperator(Mmat);
alpha = al;
kappa = kap;
T_solver.iterative_mode = false;
T_solver.SetRelTol(rel_tol);
T_solver.SetAbsTol(0.0);
T_solver.SetMaxIter(100);
T_solver.SetPrintLevel(0);
T_solver.SetPreconditioner(T_prec);
SetParameters(u);
}
void ConductionOperator::Mult(const Vector &u, Vector &du_dt) const
{
// Compute:
// du_dt = M^{-1}*-Ku
// for du_dt, where K is linearized by using u from the previous timestep
Kmat.Mult(u, z);
z.Neg(); // z = -z
M_solver.Mult(z, du_dt);
}
void ConductionOperator::ImplicitSolve(const real_t dt,
const Vector &u, Vector &du_dt)
{
// Solve the equation:
// du_dt = M^{-1}*[-K(u + dt*du_dt)]
// for du_dt, where K is linearized by using u from the previous timestep
if (!T)
{
T = Add(1.0, Mmat, dt, Kmat);
current_dt = dt;
T_solver.SetOperator(*T);
}
MFEM_VERIFY(dt == current_dt, ""); // SDIRK methods use the same dt
Kmat.Mult(u, z);
z.Neg();
T_solver.Mult(z, du_dt);
}
void ConductionOperator::SetParameters(const Vector &u)
{
GridFunction u_alpha_gf(&fespace);
u_alpha_gf.SetFromTrueDofs(u);
for (int i = 0; i < u_alpha_gf.Size(); i++)
{
u_alpha_gf(i) = kappa + alpha*u_alpha_gf(i);
}
delete K;
K = new BilinearForm(&fespace);
GridFunctionCoefficient u_coeff(&u_alpha_gf);
K->AddDomainIntegrator(new DiffusionIntegrator(u_coeff));
K->Assemble();
K->FormSystemMatrix(ess_tdof_list, Kmat);
delete T;
T = NULL; // re-compute T on the next ImplicitSolve
}
ConductionOperator::~ConductionOperator()
{
delete T;
delete M;
delete K;
}
real_t InitialTemperature(const Vector &x)
{
if (x.Norml2() < 0.5)
{
return 2.0;
}
else
{
return 1.0;
}
}