MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
ex16.cpp
Go to the documentation of this file.
1// MFEM Example 16
2//
3// Compile with: make ex16
4//
5// Sample runs: ex16
6// ex16 -m ../data/inline-tri.mesh
7// ex16 -m ../data/disc-nurbs.mesh -tf 2
8// ex16 -s 21 -a 0.0 -k 1.0
9// ex16 -s 22 -a 1.0 -k 0.0
10// ex16 -s 23 -a 0.5 -k 0.5 -o 4
11// ex16 -s 4 -dt 1.0e-4 -tf 4.0e-2 -vs 40
12// ex16 -m ../data/fichera-q2.mesh
13// ex16 -m ../data/fichera-mixed.mesh
14// ex16 -m ../data/escher.mesh
15// ex16 -m ../data/beam-tet.mesh -tf 10 -dt 0.1
16// ex16 -m ../data/amr-quad.mesh -o 4 -r 0
17// ex16 -m ../data/amr-hex.mesh -o 2 -r 0
18//
19// Description: This example solves a time dependent nonlinear heat equation
20// problem of the form du/dt = C(u), with a non-linear diffusion
21// operator C(u) = \nabla \cdot (\kappa + \alpha u) \nabla u.
22//
23// The example demonstrates the use of nonlinear operators (the
24// class ConductionOperator defining C(u)), as well as their
25// implicit time integration. Note that implementing the method
26// ConductionOperator::ImplicitSolve is the only requirement for
27// high-order implicit (SDIRK) time integration. In this example,
28// the diffusion operator is linearized by evaluating with the
29// lagged solution from the previous timestep, so there is only
30// a linear solve.
31//
32// We recommend viewing examples 2, 9 and 10 before viewing this
33// example.
34
35#include "mfem.hpp"
36#include <fstream>
37#include <iostream>
38
39using namespace std;
40using namespace mfem;
41
42/** After spatial discretization, the conduction model can be written as:
43 *
44 * du/dt = M^{-1}(-Ku)
45 *
46 * where u is the vector representing the temperature, M is the mass matrix,
47 * and K is the diffusion operator with diffusivity depending on u:
48 * (\kappa + \alpha u).
49 *
50 * Class ConductionOperator represents the right-hand side of the above ODE.
51 */
52class ConductionOperator : public TimeDependentOperator
53{
54protected:
55 FiniteElementSpace &fespace;
56 Array<int> ess_tdof_list; // this list remains empty for pure Neumann b.c.
57
58 BilinearForm *M;
59 BilinearForm *K;
60
61 SparseMatrix Mmat, Kmat;
62 SparseMatrix *T; // T = M + dt K
63 real_t current_dt;
64
65 CGSolver M_solver; // Krylov solver for inverting the mass matrix M
66 DSmoother M_prec; // Preconditioner for the mass matrix M
67
68 CGSolver T_solver; // Implicit solver for T = M + dt K
69 DSmoother T_prec; // Preconditioner for the implicit solver
70
71 real_t alpha, kappa;
72
73 mutable Vector z; // auxiliary vector
74
75public:
76 ConductionOperator(FiniteElementSpace &f, real_t alpha, real_t kappa,
77 const Vector &u);
78
79 void Mult(const Vector &u, Vector &du_dt) const override;
80 /** Solve the Backward-Euler equation: k = f(u + dt*k, t), for the unknown k.
81 This is the only requirement for high-order SDIRK implicit integration.*/
82 void ImplicitSolve(const real_t dt, const Vector &u, Vector &k) override;
83
84 /// Update the diffusion BilinearForm K using the given true-dof vector `u`.
85 void SetParameters(const Vector &u);
86
87 ~ConductionOperator() override;
88};
89
91
92int main(int argc, char *argv[])
93{
94 // 1. Parse command-line options.
95 const char *mesh_file = "../data/star.mesh";
96 int ref_levels = 2;
97 int order = 2;
98
99 int ode_solver_type = 23; // SDIRK33Solver
100 real_t t_final = 0.5;
101 real_t dt = 1.0e-2;
102 real_t alpha = 1.0e-2;
103 real_t kappa = 0.5;
104
105 bool visualization = true;
106 bool visit = false;
107 int vis_steps = 5;
108
109 int precision = 8;
110 cout.precision(precision);
111
112 OptionsParser args(argc, argv);
113 args.AddOption(&mesh_file, "-m", "--mesh",
114 "Mesh file to use.");
115 args.AddOption(&ref_levels, "-r", "--refine",
116 "Number of times to refine the mesh uniformly.");
117 args.AddOption(&order, "-o", "--order",
118 "Order (degree) of the finite elements.");
119 args.AddOption(&ode_solver_type, "-s", "--ode-solver",
120 ODESolver::Types.c_str());
121 args.AddOption(&t_final, "-tf", "--t-final",
122 "Final time; start time is 0.");
123 args.AddOption(&dt, "-dt", "--time-step",
124 "Time step.");
125 args.AddOption(&alpha, "-a", "--alpha",
126 "Alpha coefficient.");
127 args.AddOption(&kappa, "-k", "--kappa",
128 "Kappa coefficient offset.");
129 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
130 "--no-visualization",
131 "Enable or disable GLVis visualization.");
132 args.AddOption(&visit, "-visit", "--visit-datafiles", "-no-visit",
133 "--no-visit-datafiles",
134 "Save data files for VisIt (visit.llnl.gov) visualization.");
135 args.AddOption(&vis_steps, "-vs", "--visualization-steps",
136 "Visualize every n-th timestep.");
137 args.Parse();
138 if (!args.Good())
139 {
140 args.PrintUsage(cout);
141 return 1;
142 }
143 args.PrintOptions(cout);
144
145 // 2. Read the mesh from the given mesh file. We can handle triangular,
146 // quadrilateral, tetrahedral and hexahedral meshes with the same code.
147 Mesh *mesh = new Mesh(mesh_file, 1, 1);
148 int dim = mesh->Dimension();
149
150 // 3. Define the ODE solver used for time integration. Several implicit
151 // singly diagonal implicit Runge-Kutta (SDIRK) methods, as well as
152 // explicit Runge-Kutta methods are available.
153 unique_ptr<ODESolver> ode_solver = ODESolver::Select(ode_solver_type);
154
155 // 4. Refine the mesh to increase the resolution. In this example we do
156 // 'ref_levels' of uniform refinement, where 'ref_levels' is a
157 // command-line parameter.
158 for (int lev = 0; lev < ref_levels; lev++)
159 {
160 mesh->UniformRefinement();
161 }
162
163 // 5. Define the vector finite element space representing the current and the
164 // initial temperature, u_ref.
165 H1_FECollection fe_coll(order, dim);
166 FiniteElementSpace fespace(mesh, &fe_coll);
167
168 int fe_size = fespace.GetTrueVSize();
169 cout << "Number of temperature unknowns: " << fe_size << endl;
170
171 GridFunction u_gf(&fespace);
172
173 // 6. Set the initial conditions for u. All boundaries are considered
174 // natural.
176 u_gf.ProjectCoefficient(u_0);
177 Vector u;
178 u_gf.GetTrueDofs(u);
179
180 // 7. Initialize the conduction operator and the visualization.
181 ConductionOperator oper(fespace, alpha, kappa, u);
182
183 u_gf.SetFromTrueDofs(u);
184 {
185 ofstream omesh("ex16.mesh");
186 omesh.precision(precision);
187 mesh->Print(omesh);
188 ofstream osol("ex16-init.gf");
189 osol.precision(precision);
190 u_gf.Save(osol);
191 }
192
193 VisItDataCollection visit_dc("Example16", mesh);
194 visit_dc.RegisterField("temperature", &u_gf);
195 if (visit)
196 {
197 visit_dc.SetCycle(0);
198 visit_dc.SetTime(0.0);
199 visit_dc.Save();
200 }
201
202 socketstream sout;
203 if (visualization)
204 {
205 char vishost[] = "localhost";
206 int visport = 19916;
207 sout.open(vishost, visport);
208 if (!sout)
209 {
210 cout << "Unable to connect to GLVis server at "
211 << vishost << ':' << visport << endl;
212 visualization = false;
213 cout << "GLVis visualization disabled.\n";
214 }
215 else
216 {
217 sout.precision(precision);
218 sout << "solution\n" << *mesh << u_gf;
219 sout << "pause\n";
220 sout << flush;
221 cout << "GLVis visualization paused."
222 << " Press space (in the GLVis window) to resume it.\n";
223 }
224 }
225
226 // 8. Perform time-integration (looping over the time iterations, ti, with a
227 // time-step dt).
228 ode_solver->Init(oper);
229 real_t t = 0.0;
230
231 bool last_step = false;
232 for (int ti = 1; !last_step; ti++)
233 {
234 if (t + dt >= t_final - dt/2)
235 {
236 last_step = true;
237 }
238
239 ode_solver->Step(u, t, dt);
240
241 if (last_step || (ti % vis_steps) == 0)
242 {
243 cout << "step " << ti << ", t = " << t << endl;
244
245 u_gf.SetFromTrueDofs(u);
246 if (visualization)
247 {
248 sout << "solution\n" << *mesh << u_gf << flush;
249 }
250
251 if (visit)
252 {
253 visit_dc.SetCycle(ti);
254 visit_dc.SetTime(t);
255 visit_dc.Save();
256 }
257 }
258 oper.SetParameters(u);
259 }
260
261 // 9. Save the final solution. This output can be viewed later using GLVis:
262 // "glvis -m ex16.mesh -g ex16-final.gf".
263 {
264 ofstream osol("ex16-final.gf");
265 osol.precision(precision);
266 u_gf.Save(osol);
267 }
268
269 // 10. Free the used memory.
270 delete mesh;
271
272 return 0;
273}
274
275ConductionOperator::ConductionOperator(FiniteElementSpace &f, real_t al,
276 real_t kap, const Vector &u)
277 : TimeDependentOperator(f.GetTrueVSize(), (real_t) 0.0), fespace(f),
278 M(NULL), K(NULL), T(NULL), current_dt(0.0), z(height)
279{
280 const real_t rel_tol = 1e-8;
281
282 M = new BilinearForm(&fespace);
283 M->AddDomainIntegrator(new MassIntegrator());
284 M->Assemble();
285 M->FormSystemMatrix(ess_tdof_list, Mmat);
286
287 M_solver.iterative_mode = false;
288 M_solver.SetRelTol(rel_tol);
289 M_solver.SetAbsTol(0.0);
290 M_solver.SetMaxIter(30);
291 M_solver.SetPrintLevel(0);
292 M_solver.SetPreconditioner(M_prec);
293 M_solver.SetOperator(Mmat);
294
295 alpha = al;
296 kappa = kap;
297
298 T_solver.iterative_mode = false;
299 T_solver.SetRelTol(rel_tol);
300 T_solver.SetAbsTol(0.0);
301 T_solver.SetMaxIter(100);
302 T_solver.SetPrintLevel(0);
303 T_solver.SetPreconditioner(T_prec);
304
305 SetParameters(u);
306}
307
308void ConductionOperator::Mult(const Vector &u, Vector &du_dt) const
309{
310 // Compute:
311 // du_dt = M^{-1}*-Ku
312 // for du_dt, where K is linearized by using u from the previous timestep
313 Kmat.Mult(u, z);
314 z.Neg(); // z = -z
315 M_solver.Mult(z, du_dt);
316}
317
318void ConductionOperator::ImplicitSolve(const real_t dt,
319 const Vector &u, Vector &du_dt)
320{
321 // Solve the equation:
322 // du_dt = M^{-1}*[-K(u + dt*du_dt)]
323 // for du_dt, where K is linearized by using u from the previous timestep
324 if (!T)
325 {
326 T = Add(1.0, Mmat, dt, Kmat);
327 current_dt = dt;
328 T_solver.SetOperator(*T);
329 }
330 MFEM_VERIFY(dt == current_dt, ""); // SDIRK methods use the same dt
331 Kmat.Mult(u, z);
332 z.Neg();
333 T_solver.Mult(z, du_dt);
334}
335
336void ConductionOperator::SetParameters(const Vector &u)
337{
338 GridFunction u_alpha_gf(&fespace);
339 u_alpha_gf.SetFromTrueDofs(u);
340 for (int i = 0; i < u_alpha_gf.Size(); i++)
341 {
342 u_alpha_gf(i) = kappa + alpha*u_alpha_gf(i);
343 }
344
345 delete K;
346 K = new BilinearForm(&fespace);
347
348 GridFunctionCoefficient u_coeff(&u_alpha_gf);
349
351 K->Assemble();
352 K->FormSystemMatrix(ess_tdof_list, Kmat);
353 delete T;
354 T = NULL; // re-compute T on the next ImplicitSolve
355}
356
357ConductionOperator::~ConductionOperator()
358{
359 delete T;
360 delete M;
361 delete K;
362}
363
365{
366 if (x.Norml2() < 0.5)
367 {
368 return 2.0;
369 }
370 else
371 {
372 return 1.0;
373 }
374}
A "square matrix" operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
Conjugate gradient method.
Definition solvers.hpp:538
void Mult(const Vector &b, Vector &x) const override
Iterative solution of the linear system using the Conjugate Gradient method.
Definition solvers.cpp:751
void SetOperator(const Operator &op) override
Set/update the solver for the given operator.
Definition solvers.hpp:551
Data type for scaled Jacobi-type smoother of sparse matrix.
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
void SetTime(real_t t)
Set physical time (for time-dependent simulations)
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:244
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition fespace.hpp:851
A general function coefficient.
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:31
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
Definition gridfunc.cpp:378
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
void GetTrueDofs(Vector &tv) const
Extract the true-dofs from the GridFunction.
Definition gridfunc.cpp:363
Arbitrary order H1-conforming (continuous) finite elements.
Definition fe_coll.hpp:275
Mesh data type.
Definition mesh.hpp:64
virtual void Print(std::ostream &os=mfem::out, const std::string &comments="") const
Definition mesh.hpp:2433
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1216
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:11295
static MFEM_EXPORT std::string Types
Definition ode.hpp:187
static MFEM_EXPORT std::unique_ptr< ODESolver > Select(const int ode_solver_type)
Definition ode.cpp:34
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void PrintUsage(std::ostream &out) const
Print the usage message.
void PrintOptions(std::ostream &out) const
Print the options.
void AddOption(bool *var, const char *enable_short_name, const char *enable_long_name, const char *disable_short_name, const char *disable_long_name, const char *description, bool required=false)
Add a boolean option and set 'var' to receive the value. Enable/disable tags are used to set the bool...
Definition optparser.hpp:82
bool Good() const
Return true if the command line options were parsed successfully.
Data type sparse matrix.
Definition sparsemat.hpp:51
void Mult(const Vector &x, Vector &y) const override
Matrix vector multiplication.
Base abstract class for first order time dependent operators.
Definition operator.hpp:332
Vector data type.
Definition vector.hpp:82
void Neg()
(*this) = -(*this)
Definition vector.cpp:375
real_t Norml2() const
Returns the l2 norm of the vector.
Definition vector.cpp:931
Data collection with VisIt I/O routines.
void Save() override
Save the collection and a VisIt root file.
void RegisterField(const std::string &field_name, GridFunction *gf) override
Add a grid function to the collection and update the root file.
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
const real_t alpha
Definition ex15.cpp:369
real_t InitialTemperature(const Vector &x)
Definition ex16.cpp:364
real_t kappa
Definition ex24.cpp:54
int dim
Definition ex24.cpp:53
int main()
real_t u(const Vector &xvec)
Definition lor_mms.hpp:22
float real_t
Definition config.hpp:43
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30
void Add(const DenseMatrix &A, const DenseMatrix &B, real_t alpha, DenseMatrix &C)
C = A + alpha*B.
const char vishost[]