MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
ex18.hpp
Go to the documentation of this file.
1// MFEM Example 18 - Serial/Parallel Shared Code
2// (Implementation of Time-dependent DG Operator)
3//
4// This code provide example problems for the Euler equations and implements
5// the time-dependent DG operator given by the equation:
6//
7// (u_t, v)_T - (F(u), ∇ v)_T + <F̂(u, n), [[v]]>_F = 0.
8//
9// This operator is designed for explicit time stepping methods. Specifically,
10// the function DGHyperbolicConservationLaws::Mult implements the following
11// transformation:
12//
13// u ↦ M⁻¹(-DF(u) + NF(u))
14//
15// where M is the mass matrix, DF is the weak divergence of flux, and NF is the
16// interface flux. The inverse of the mass matrix is computed element-wise by
17// leveraging the block-diagonal structure of the DG mass matrix. Additionally,
18// the flux-related terms are computed using the HyperbolicFormIntegrator.
19//
20// The maximum characteristic speed is determined for each time step. For more
21// details, refer to the documentation of DGHyperbolicConservationLaws::Mult.
22//
23
24#include <functional>
25#include "mfem.hpp"
26
27namespace mfem
28{
29
30/// @brief Time dependent DG operator for hyperbolic conservation laws
32{
33private:
34 const int num_equations; // the number of equations
35 const int dim;
36 FiniteElementSpace &vfes; // vector finite element space
37 // Element integration form. Should contain ComputeFlux
38 std::unique_ptr<HyperbolicFormIntegrator> formIntegrator;
39 // Base Nonlinear Form
40 std::unique_ptr<NonlinearForm> nonlinearForm;
41 // element-wise inverse mass matrix
42 std::vector<DenseMatrix> invmass; // local scalar inverse mass
43 std::vector<DenseMatrix> weakdiv; // local weak divergence (trial space ByDim)
44 // global maximum characteristic speed. Updated by form integrators
45 mutable real_t max_char_speed;
46 // auxiliary variable used in Mult
47 mutable Vector z;
48
49 // Compute element-wise inverse mass matrix
50 void ComputeInvMass();
51 // Compute element-wise weak-divergence matrix
52 void ComputeWeakDivergence();
53
54public:
55 /**
56 * @brief Construct a new DGHyperbolicConservationLaws object
57 *
58 * @param vfes_ vector finite element space. Only tested for DG [Pₚ]ⁿ
59 * @param formIntegrator_ integrator (F(u,x), grad v)
60 * @param preassembleWeakDivergence preassemble weak divergence for faster
61 * assembly
62 */
64 FiniteElementSpace &vfes_,
65 std::unique_ptr<HyperbolicFormIntegrator> formIntegrator_,
66 bool preassembleWeakDivergence=true);
67 /**
68 * @brief Apply nonlinear form to obtain M⁻¹(DIVF + JUMP HAT(F))
69 *
70 * @param x current solution vector
71 * @param y resulting dual vector to be used in an EXPLICIT solver
72 */
73 void Mult(const Vector &x, Vector &y) const override;
74 // get global maximum characteristic speed to be used in CFL condition
75 // where max_char_speed is updated during Mult.
76 real_t GetMaxCharSpeed() { return max_char_speed; }
77 void Update();
78
79};
80
81//////////////////////////////////////////////////////////////////
82/// HYPERBOLIC CONSERVATION LAWS IMPLEMENTATION ///
83//////////////////////////////////////////////////////////////////
84
85// Implementation of class DGHyperbolicConservationLaws
87 FiniteElementSpace &vfes_,
88 std::unique_ptr<HyperbolicFormIntegrator> formIntegrator_,
89 bool preassembleWeakDivergence)
90 : TimeDependentOperator(vfes_.GetTrueVSize()),
91 num_equations(formIntegrator_->num_equations),
92 dim(vfes_.GetMesh()->SpaceDimension()),
93 vfes(vfes_),
94 formIntegrator(std::move(formIntegrator_)),
95 z(vfes_.GetTrueVSize())
96{
97 // Standard local assembly and inversion for energy mass matrices.
98 ComputeInvMass();
99#ifndef MFEM_USE_MPI
100 nonlinearForm.reset(new NonlinearForm(&vfes));
101#else
102 ParFiniteElementSpace *pvfes = dynamic_cast<ParFiniteElementSpace *>(&vfes);
103 if (pvfes)
104 {
105 nonlinearForm.reset(new ParNonlinearForm(pvfes));
106 }
107 else
108 {
109 nonlinearForm.reset(new NonlinearForm(&vfes));
110 }
111#endif
112 if (preassembleWeakDivergence)
113 {
114 ComputeWeakDivergence();
115 }
116 else
117 {
118 nonlinearForm->AddDomainIntegrator(formIntegrator.get());
119 }
120 nonlinearForm->AddInteriorFaceIntegrator(formIntegrator.get());
121 nonlinearForm->UseExternalIntegrators();
122
123}
124
125void DGHyperbolicConservationLaws::ComputeInvMass()
126{
127 InverseIntegrator inv_mass(new MassIntegrator());
128
129 invmass.resize(vfes.GetNE());
130 for (int i=0; i<vfes.GetNE(); i++)
131 {
132 int dof = vfes.GetFE(i)->GetDof();
133 invmass[i].SetSize(dof);
134 inv_mass.AssembleElementMatrix(*vfes.GetFE(i),
136 invmass[i]);
137 }
138}
139
140void DGHyperbolicConservationLaws::ComputeWeakDivergence()
141{
142 TransposeIntegrator weak_div(new GradientIntegrator());
143 DenseMatrix weakdiv_bynodes;
144
145 weakdiv.resize(vfes.GetNE());
146 for (int i=0; i<vfes.GetNE(); i++)
147 {
148 int dof = vfes.GetFE(i)->GetDof();
149 weakdiv_bynodes.SetSize(dof, dof*dim);
150 weak_div.AssembleElementMatrix2(*vfes.GetFE(i), *vfes.GetFE(i),
152 weakdiv_bynodes);
153 weakdiv[i].SetSize(dof, dof*dim);
154 // Reorder so that trial space is ByDim.
155 // This makes applying weak divergence to flux value simpler.
156 for (int j=0; j<dof; j++)
157 {
158 for (int d=0; d<dim; d++)
159 {
160 weakdiv[i].SetCol(j*dim + d, weakdiv_bynodes.GetColumn(d*dof + j));
161 }
162 }
163
164 }
165}
166
167
169{
170 // 0. Reset wavespeed computation before operator application.
171 formIntegrator->ResetMaxCharSpeed();
172 // 1. Apply Nonlinear form to obtain an auxiliary result
173 // z = - <F̂(u_h,n), [[v]]>_e
174 // If weak-divergence is not preassembled, we also have weak-divergence
175 // z = - <F̂(u_h,n), [[v]]>_e + (F(u_h), ∇v)
176 nonlinearForm->Mult(x, z);
177 if (!weakdiv.empty()) // if weak divergence is pre-assembled
178 {
179 // Apply weak divergence to F(u_h), and inverse mass to z_loc + weakdiv_loc
180 Vector current_state; // view of current state at a node
181 DenseMatrix current_flux; // flux of current state
182 DenseMatrix flux; // element flux value. Whose column is ordered by dim.
183 DenseMatrix current_xmat; // view of current states in an element, dof x num_eq
184 DenseMatrix current_zmat; // view of element auxiliary result, dof x num_eq
185 DenseMatrix current_ymat; // view of element result, dof x num_eq
186 const FluxFunction &fluxFunction = formIntegrator->GetFluxFunction();
187 Array<int> vdofs;
188 Vector xval, zval;
189 for (int i=0; i<vfes.GetNE(); i++)
190 {
192 int dof = vfes.GetFE(i)->GetDof();
193 vfes.GetElementVDofs(i, vdofs);
194 x.GetSubVector(vdofs, xval);
195 current_xmat.UseExternalData(xval.GetData(), dof, num_equations);
196 flux.SetSize(num_equations, dim*dof);
197 for (int j=0; j<dof; j++) // compute flux for all nodes in the element
198 {
199 current_xmat.GetRow(j, current_state);
200 current_flux.UseExternalData(flux.GetData() + num_equations*dim*j,
201 num_equations, dof);
202 fluxFunction.ComputeFlux(current_state, *Tr, current_flux);
203 }
204 // Compute weak-divergence and add it to auxiliary result, z
205 // Recalling that weakdiv is reordered by dim, we can apply
206 // weak-divergence to the transpose of flux.
207 z.GetSubVector(vdofs, zval);
208 current_zmat.UseExternalData(zval.GetData(), dof, num_equations);
209 mfem::AddMult_a_ABt(1.0, weakdiv[i], flux, current_zmat);
210 // Apply inverse mass to auxiliary result to obtain the final result
211 current_ymat.SetSize(dof, num_equations);
212 mfem::Mult(invmass[i], current_zmat, current_ymat);
213 y.SetSubVector(vdofs, current_ymat.GetData());
214 }
215 }
216 else
217 {
218 // Apply block inverse mass
219 Vector zval; // z_loc, dof*num_eq
220
221 DenseMatrix current_zmat; // view of element auxiliary result, dof x num_eq
222 DenseMatrix current_ymat; // view of element result, dof x num_eq
223 Array<int> vdofs;
224 for (int i=0; i<vfes.GetNE(); i++)
225 {
226 int dof = vfes.GetFE(i)->GetDof();
227 vfes.GetElementVDofs(i, vdofs);
228 z.GetSubVector(vdofs, zval);
229 current_zmat.UseExternalData(zval.GetData(), dof, num_equations);
230 current_ymat.SetSize(dof, num_equations);
231 mfem::Mult(invmass[i], current_zmat, current_ymat);
232 y.SetSubVector(vdofs, current_ymat.GetData());
233 }
234 }
235 max_char_speed = formIntegrator->GetMaxCharSpeed();
236}
237
239{
240 nonlinearForm->Update();
241 height = nonlinearForm->Height();
242 width = height;
243 z.SetSize(height);
244
245 ComputeInvMass();
246 if (!weakdiv.empty()) {ComputeWeakDivergence();}
247}
248
249std::function<void(const Vector&, Vector&)> GetMovingVortexInit(
250 const real_t radius, const real_t Minf, const real_t beta,
251 const real_t gas_constant, const real_t specific_heat_ratio)
252{
253 return [specific_heat_ratio,
254 gas_constant, Minf, radius, beta](const Vector &x, Vector &y)
255 {
256 MFEM_ASSERT(x.Size() == 2, "");
257
258 const real_t xc = 0.0, yc = 0.0;
259
260 // Nice units
261 const real_t vel_inf = 1.;
262 const real_t den_inf = 1.;
263
264 // Derive remainder of background state from this and Minf
265 const real_t pres_inf = (den_inf / specific_heat_ratio) *
266 (vel_inf / Minf) * (vel_inf / Minf);
267 const real_t temp_inf = pres_inf / (den_inf * gas_constant);
268
269 real_t r2rad = 0.0;
270 r2rad += (x(0) - xc) * (x(0) - xc);
271 r2rad += (x(1) - yc) * (x(1) - yc);
272 r2rad /= (radius * radius);
273
274 const real_t shrinv1 = 1.0 / (specific_heat_ratio - 1.);
275
276 const real_t velX =
277 vel_inf * (1 - beta * (x(1) - yc) / radius * std::exp(-0.5 * r2rad));
278 const real_t velY =
279 vel_inf * beta * (x(0) - xc) / radius * std::exp(-0.5 * r2rad);
280 const real_t vel2 = velX * velX + velY * velY;
281
282 const real_t specific_heat =
283 gas_constant * specific_heat_ratio * shrinv1;
284 const real_t temp = temp_inf - 0.5 * (vel_inf * beta) *
285 (vel_inf * beta) / specific_heat *
286 std::exp(-r2rad);
287
288 const real_t den = den_inf * std::pow(temp / temp_inf, shrinv1);
289 const real_t pres = den * gas_constant * temp;
290 const real_t energy = shrinv1 * pres / den + 0.5 * vel2;
291
292 y(0) = den;
293 y(1) = den * velX;
294 y(2) = den * velY;
295 y(3) = den * energy;
296 };
297}
298
300{
301 switch (problem)
302 {
303 case 1:
304 case 2:
305 case 3:
306 return Mesh("../data/periodic-square.mesh");
307 break;
308 case 4:
309 return Mesh("../data/periodic-segment.mesh");
310 break;
311 default:
312 MFEM_ABORT("Problem Undefined");
313 }
314}
315
316// Initial condition
318 const real_t specific_heat_ratio,
319 const real_t gas_constant)
320{
321 switch (problem)
322 {
323 case 1: // fast moving vortex
325 4, GetMovingVortexInit(0.2, 0.5, 1. / 5., gas_constant,
326 specific_heat_ratio));
327 case 2: // slow moving vortex
329 4, GetMovingVortexInit(0.2, 0.05, 1. / 50., gas_constant,
330 specific_heat_ratio));
331 case 3: // moving sine wave
332 return VectorFunctionCoefficient(4, [](const Vector &x, Vector &y)
333 {
334 MFEM_ASSERT(x.Size() == 2, "");
335 const real_t density = 1.0 + 0.2 * std::sin(M_PI*(x(0) + x(1)));
336 const real_t velocity_x = 0.7;
337 const real_t velocity_y = 0.3;
338 const real_t pressure = 1.0;
339 const real_t energy =
340 pressure / (1.4 - 1.0) +
341 density * 0.5 * (velocity_x * velocity_x + velocity_y * velocity_y);
342
343 y(0) = density;
344 y(1) = density * velocity_x;
345 y(2) = density * velocity_y;
346 y(3) = energy;
347 });
348 case 4:
349 return VectorFunctionCoefficient(3, [](const Vector &x, Vector &y)
350 {
351 MFEM_ASSERT(x.Size() == 1, "");
352 const real_t density = 1.0 + 0.2 * std::sin(M_PI * 2 * x(0));
353 const real_t velocity_x = 1.0;
354 const real_t pressure = 1.0;
355 const real_t energy =
356 pressure / (1.4 - 1.0) + density * 0.5 * (velocity_x * velocity_x);
357
358 y(0) = density;
359 y(1) = density * velocity_x;
360 y(2) = energy;
361 });
362 default:
363 MFEM_ABORT("Problem Undefined");
364 }
365}
366
367} // namespace mfem
Time dependent DG operator for hyperbolic conservation laws.
Definition ex18.hpp:32
void Mult(const Vector &x, Vector &y) const override
Apply nonlinear form to obtain M⁻¹(DIVF + JUMP HAT(F))
Definition ex18.hpp:168
DGHyperbolicConservationLaws(FiniteElementSpace &vfes_, std::unique_ptr< HyperbolicFormIntegrator > formIntegrator_, bool preassembleWeakDivergence=true)
Construct a new DGHyperbolicConservationLaws object.
Definition ex18.hpp:86
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
real_t * GetData() const
Returns the matrix data array.
Definition densemat.hpp:115
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition densemat.hpp:105
void UseExternalData(real_t *d, int h, int w)
Change the data array and the size of the DenseMatrix.
Definition densemat.hpp:80
void GetRow(int r, Vector &row) const
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Definition fespace.hpp:773
DofTransformation * GetElementVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom for the i'th element. The returned indices are offsets into an ...
Definition fespace.cpp:287
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
Definition fespace.cpp:3168
int GetNE() const
Returns number of elements in the mesh.
Definition fespace.hpp:740
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition fe_base.hpp:329
Abstract class for hyperbolic flux for a system of hyperbolic conservation laws.
virtual real_t ComputeFlux(const Vector &state, ElementTransformation &Tr, DenseMatrix &flux) const =0
Compute flux F(u, x) for given state u and physical point x.
Integrator that inverts the matrix assembled by another integrator.
Mesh data type.
Definition mesh.hpp:56
int width
Dimension of the input / number of columns in the matrix.
Definition operator.hpp:28
int height
Dimension of the output / number of rows in the matrix.
Definition operator.hpp:27
Abstract parallel finite element space.
Definition pfespace.hpp:29
Parallel non-linear operator on the true dofs.
Base abstract class for first order time dependent operators.
Definition operator.hpp:317
A general vector function coefficient.
Vector data type.
Definition vector.hpp:80
void SetSubVector(const Array< int > &dofs, const real_t value)
Set the entries listed in dofs to the given value.
Definition vector.cpp:604
int Size() const
Returns the size of the vector.
Definition vector.hpp:218
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:538
real_t * GetData() const
Return a pointer to the beginning of the Vector data.
Definition vector.hpp:227
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
Definition vector.cpp:578
Vector beta
const real_t radius
Definition distance.cpp:107
int problem
Definition ex15.cpp:62
int dim
Definition ex24.cpp:53
Mesh * GetMesh(int type)
Definition ex29.cpp:218
void AddMult_a_ABt(real_t a, const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
ABt += a * A * B^t.
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition table.cpp:476
Mesh EulerMesh(const int problem)
Definition ex18.hpp:299
std::function< void(const Vector &, Vector &)> GetMovingVortexInit(const real_t radius, const real_t Minf, const real_t beta, const real_t gas_constant, const real_t specific_heat_ratio)
Definition ex18.hpp:249
VectorFunctionCoefficient EulerInitialCondition(const int problem, const real_t specific_heat_ratio, const real_t gas_constant)
Definition ex18.hpp:317
float real_t
Definition config.hpp:43