spatial
Definitions of the spatial discretizations for the scalar transport equation.
- class DG(mesh, root=None, **default)
Class containing the tools that define a standard discontinuous Galerkin (DG) finite element method for the scalar transport equation.
This may be compactly written as: Find \(u_h, v_h \in U_h\), such that
\[\bm{M} \partial_t u_h + B_c(u_h, v_h) + B_d(u_h, v_h) = l_c(v_h) + l_d(v_h),\]where the first term, involving the mass matrix \(\bm{M} = (u_h, v_h)\) is treated separately, depending on the chosen temporal scheme, see
time
. The discrete space is defined as\[U_h := L^2\left(\mathbb{P}^k(\Omega, \mathbb{R}^{d}) \right).\]The above bilinear and linear forms correspond to
\(B_d\): bilinear form for the diffusion/elliptic terms, see
add_diffusion_form()
.\(B_c\): bilinear form for the convective terms, see
add_convection_form()
.\(l_d\): linear form for the diffusion/elliptic terms, see
add_boundary_conditions()
.\(l_c\): linear form for the convective terms, see
add_boundary_conditions()
.
- Note:
For readability, the \((\cdot)_h\) is dropped throughout the documentation.
- property scheme: ImplicitEuler | BDF2 | SDIRK22 | SDIRK33 | IMEXRK_ARS443 | ExplicitEuler | SSPRK3 | CRK4
Time scheme for the DG method depending on the choosen time routine.
- Getter:
Returns the time scheme
- Setter:
Sets the time scheme
- add_finite_element_spaces(fes: dict[str, FESpace]) None
Defines the spaces for the test and trial functions based on a DG discretization.
- add_convection_form(blf: dict[str, dict[str, SumOfIntegrals]], lf: dict[str, dict[str, SumOfIntegrals]])
Implementation of the spatial discretization, in the internal domain, of the convective terms using a standard Riemann solver.
\[B_c\big(u,v\big) = -\int\limits_{D} \bm{f}(u) \cdot \nabla v\, d\bm{x}\, + \hspace{-4mm} \int\limits_{\partial D \backslash \partial \Omega} f^*(u,\hat{u})\, v\, d\bm{s},\]where the numerical flux \(f^*\) is defined based on a Riemann solver, see
get_convective_numerical_flux_dg()
.
- add_diffusion_form(blf: dict[str, dict[str, SumOfIntegrals]], lf: dict[str, dict[str, SumOfIntegrals]])
DG discretization (in the internal domain) of the elliptic terms using a symmetric interior penalty method. The diffusion bilinear form over an internal element is:
\[\begin{split}\begin{aligned} B_d\big(u,v\big) &= \int\limits_D \kappa \nabla u \cdot \nabla v\, d\bm{x}\, -\hspace{-4mm} \int\limits_{\partial D \setminus \partial \Omega} \frac{\kappa}{2} (\nabla u_i + \nabla u_j) \cdot \bm{n}\, v\, d\bm{s}\\ &\qquad\qquad -\hspace{-4mm} \int\limits_{\partial D \setminus \partial \Omega} \frac{\kappa}{2} (u_i - u_j) (\nabla v_i + \nabla v_j) \cdot \bm{n}\, d\bm{s}\, +\hspace{-4mm} \int\limits_{\partial D \setminus \partial \Omega} \kappa \tau (u_i - u_j)(v_i - v_j)\, d\bm{s}, \end{aligned}\end{split}\]where the ith and jth subscripts correspond to the indices of the local solution and its neighobring solution, respectively. Additionally, \(\tau\) is the interior penalty coefficient, see
get_penalty_coefficient()
.
- add_farfield_formulation(blf: dict[str, dict[str, SumOfIntegrals]], lf: dict[str, dict[str, SumOfIntegrals]], bc: FarField, bnd: str)
Imposes a farfield boundary condition using both convective and diffusive (if turned on) fluxes in a DG formulation. The bilinear and linear forms associated with a farfield BC are:
\[B^{\partial} = B_c^{\partial} + B_d^{\partial}, \qquad l^{\partial} = l_c^{\partial} + l_d^{\partial}.\]The convection forms are defined as
\[\begin{split}B_c^{\partial}\big(u, v\big) &= \hspace{-4mm}\int\limits_{\partial D \cap \partial \Omega} \frac{1}{2} (b_n + |b_n|) u\, v\, d\bm{s},\\ l_c^{\partial}(v) &= -\hspace{-4mm}\int\limits_{\partial D \cap \partial \Omega} \frac{1}{2} (b_n - |b_n|) \hat{u}\, v\, d\bm{s}.\end{split}\]The diffusion forms are defined as
\[\begin{split}B_d^{\partial}\big(u, v\big) &= \hspace{-4mm}\int\limits_{\partial D \cap \partial \Omega} \tau \kappa u\, v\, d\bm{s},\\ l_d^{\partial}(v) &= \hspace{-4mm}\int\limits_{\partial D \cap \partial \Omega} \tau \kappa \bar{u}\, v\, d\bm{s},\end{split}\]with \(\bar{u}\) denoting the farfield boundary value we impose.
- class HDG(mesh, root=None, **default)
Class containing the tools that define a hybridized discontinuous Galerkin (HDG) finite element method for the scalar transport equation.
This may be compactly written as: Find \(u_h, v_h \in U_h\) and \(\hat{u}_{h}, \hat{v}_{h} \in \hat{U}_h\), such that
\[\bm{M} \partial_t u_h + B_c\Big(\{u_h, \hat{u}_{h}\}, \{v_h, \hat{v}_{h}\}\Big) + B_d\Big(\{u_h, \hat{u}_{h}\}, \{v_h, \hat{v}_{h}\}\Big) = l_c\Big(\{v_h, \hat{v}_{h}\}\Big) + l_d\Big(\{v_h, \hat{v}_{h}\}\Big),\]where the first term, involving the mass matrix \(\bm{M} = (u_h, v_h)\) is treated separately, depending on the chosen temporal scheme, see
time
. The discrete spaces are defined as\[\begin{split}U_h &:= L^2\left(\mathbb{P}^k(\Omega, \mathbb{R}^{d}) \right),\\ \hat{U}_h &:= L^2\left(\mathbb{P}^k(\facets, \mathbb{R}^{d}) \right). \end{split}\]The above bilinear and linear forms correspond to
\(B_d\): bilinear form for the diffusion/elliptic terms, see
add_diffusion_form()
.\(B_c\): bilinear form for the convective terms, see
add_convection_form()
.\(l_d\): linear form for the diffusion/elliptic terms, see
add_boundary_conditions()
.\(l_c\): linear form for the convective terms, see
add_boundary_conditions()
.
- Note:
For readability, the \((\cdot)_h\) is dropped throughout the documentation.
- property scheme: ImplicitEuler | BDF2 | SDIRK22 | SDIRK33
Time scheme for the HDG method depending on the choosen time routine.
- Getter:
Returns the time scheme
- Setter:
Sets the time scheme
- add_finite_element_spaces(fes: dict[str, FESpace]) None
Defines the spaces for the test and trial functions based on an HDG discretization. This also sets periodic boundaries on the facet space, in case periodic boundary conditions are prescribed.
- add_convection_form(blf: dict[str, dict[str, SumOfIntegrals]], lf: dict[str, dict[str, SumOfIntegrals]])
Discretization (in the internal domain) of the convection terms using a standard Riemann solver. The convective bilinear form over an internal element is:
\[\begin{split}B_c\Big(\{u,\hat{u}\},v\Big) &= -\int\limits_{D} \bm{f}(u) \cdot \nabla v\, d\bm{x}\, + \hspace{-4mm} \int\limits_{\partial D \backslash \partial \Omega} f^*(u,\hat{u})\, v\, d\bm{s},\\ B_c\Big(\{u,\hat{u}\},\hat{v}\Big) &= -\hspace{-4mm} \int\limits_{\partial D \backslash \partial \Omega} f^*(u, \hat{u})\, \hat{v}\, d\bm{s},\end{split}\]where the physical (inviscid) flux is \(\bm{f}(u) = \bm{b} u\) and the numerical flux \(f^* = f^*(u,\hat{u})\), over an element boundary, is defined based on the local solution \(u\) and its neighboring facet variable \(\hat{u}\). See
get_convective_numerical_flux_hdg()
for details.
- add_diffusion_form(blf: dict[str, dict[str, SumOfIntegrals]], lf: dict[str, dict[str, SumOfIntegrals]])
HDG discretization (in the internal domain) of the elliptic terms using a symmetric interior penalty method. The diffusion bilinear form over an internal element is:
\[\begin{split}\begin{aligned} B_d\Big(\{u,\hat{u}\},v\Big) &= \int\limits_D \kappa \nabla u \cdot \nabla v \, d\bm{x}\, -\hspace{-4mm} \int\limits_{\partial D \setminus \partial \Omega} \kappa (\nabla u \cdot \bm{n}) v\, d\bm{s}\\ &\qquad\qquad -\hspace{-4mm} \int\limits_{\partial D \setminus \partial \Omega} \kappa (u - \hat{u}) (\bm{n} \cdot \nabla v)\, d\bm{s}\, +\hspace{-4mm} \int\limits_{\partial D \setminus \partial \Omega} \kappa \tau (u - \hat{u}) v \, d\bm{s}, \\[1ex] B_d\Big(\{u,\hat{u}\},\hat{v}\Big) &= \hspace{-4mm} \int\limits_{\partial D \setminus \partial \Omega} \kappa (\nabla u \cdot \bm{n}) \hat{v}\, d\bm{s}\, -\hspace{-4mm} \int\limits_{\partial D \setminus \partial \Omega} \kappa \tau (u - \hat{u}) \hat{v}\, d\bm{s}, \end{aligned}\end{split}\]where \(\tau\) is the interior penalty coefficient, see
get_penalty_coefficient()
.
- set_initial_conditions(U: CoefficientFunction = None)
Initializes the solution by projecting the initial condition on both the volume elements and their facets.
- add_farfield_formulation(blf: dict[str, dict[str, SumOfIntegrals]], lf: dict[str, dict[str, SumOfIntegrals]], bc: FarField, bnd: str)
Imposes a farfield boundary condition using both convective and diffusive (if turned on) fluxes in an HDG formulation. The bilinear and linear forms associated with a farfield BC are:
\[B^{\partial} = B_c^{\partial} + B_d^{\partial}, \qquad l^{\partial} = l_c^{\partial} + l_d^{\partial}.\]The convection forms are defined as
\[\begin{split}B_c^{\partial}\Big(\{u,\hat{u}\}, \hat{v}\Big) &= -\hspace{-4mm}\int\limits_{\partial D \cap \partial \Omega} \Big( (b_n^{+} + b_n^{-}) \hat{u} + b_n^{+} u\Big)\, \hat{v}\, d\bm{s},\\ l_c^{\partial}(\hat{v}) &= -\hspace{-4mm}\int\limits_{\partial D \cap \partial \Omega} b_n^{-} \bar{u}\, \hat{v}\, d\bm{s},\end{split}\]with \(b_n^{+} = \max(b_n, 0)\), \(b_n^{-} = \min(b_n, 0)\) and \(b_n = \bm{b} \cdot \bm{n}\).
The diffusion forms are defined as
\[\begin{split}B_d^{\partial}\Big(\{u,\hat{u}\}, \hat{v}\Big) &= -\hspace{-4mm}\int\limits_{\partial D \cap \partial \Omega} \tau \kappa \hat{u}\, \hat{v}\, d\bm{s},\\ l_d^{\partial}(\hat{v}) &= -\hspace{-4mm}\int\limits_{\partial D \cap \partial \Omega} \tau \kappa \bar{u}\, \hat{v}\, d\bm{s},\end{split}\]with \(\bar{u}\) denoting the farfield boundary value we impose.
- class ScalarTransportFiniteElementMethod(mesh, root=None, **default)
- property interior_penalty_coefficient: float
Sets the interior penalty constant.
- Getter:
Returns the interior penalty constant
- Setter:
Sets the interior penalty constant, defaults to 1.0
- get_domain_boundary_mask() GridFunction
Returns a Gridfunction that is 0 on the domain boundaries and 1 on the domain interior.
- get_penalty_coefficient() CoefficientFunction
Returns the (dimensionally-consistent) interior penalty coefficient. Its definition is
\(\tau = \alpha \frac{(N+1)^2}{h}\),
where N denotes the polynomial order and \(\alpha\) is (user-specified) penalty constant from
interior_penalty_coefficient()
- add_symbolic_spatial_forms(blf: dict[str, dict[str, SumOfIntegrals]], lf: dict[str, dict[str, SumOfIntegrals]])
Defines the bilinear forms and linear functionals associated with an spatial formulation.
See
add_convection_form()
andadd_convection_form()
for the definition of the convection terms for the HDG and DG formulation, respectively.See
add_diffusion_form()
andadd_diffusion_form()
for the implementation of the diffusion terms for the HDG and DG formulation, respectively. Note,is_inviscid()
must be False.See
add_boundary_conditions()
for physical boundary condition imposition.See
add_domain_conditions()
for imposing domain conditions (e.g. initial conditions).
- add_boundary_conditions(blf: dict[str, dict[str, SumOfIntegrals]], lf: dict[str, dict[str, SumOfIntegrals]])
Adds specific boundary conditions to the blinear/linear forms.
The available options are:
FarField
andPeriodic
- add_domain_conditions(blf: dict[str, dict[str, SumOfIntegrals]], lf: dict[str, dict[str, SumOfIntegrals]])
Adds domain conditions.
- set_boundary_conditions() None
Boundary conditions for the scalar transport equation are set weakly. Therefore we do nothing here.