ryujin 2.1.1 revision 0348cbb53a3e4b1da2a4c037e81f88f2d21ce219
|
#include <source/navier_stokes/parabolic_solver.h>
Public Member Functions | |
Constructor and setup | |
ParabolicSolver (const MPIEnsemble &mpi_ensemble, std::map< std::string, dealii::Timer > &computing_timer, const HyperbolicSystem &hyperbolic_system, const ParabolicSystem ¶bolic_system, const OfflineData< dim, Number > &offline_data, const InitialValues< Description, dim, Number > &initial_values, const std::string &subsection="ParabolicSolver") | |
void | prepare () |
Functions for performing implicit time steps | |
void | backward_euler_step (const StateVector &old_state_vector, const Number old_t, StateVector &new_state_vector, Number tau, const IDViolationStrategy id_violation_strategy, const bool reinitialize_gmg) const |
void | print_solver_statistics (std::ostream &output) const |
Typedefs and constexpr constants | |
using | HyperbolicSystem = typename Description::HyperbolicSystem |
using | View = typename Description::template HyperbolicSystemView< dim, Number > |
using | ParabolicSystem = typename Description::ParabolicSystem |
using | ScalarNumber = typename View::ScalarNumber |
using | state_type = typename View::state_type |
using | StateVector = typename View::StateVector |
using | ScalarVector = Vectors::ScalarVector< Number > |
using | BlockVector = Vectors::BlockVector< Number > |
static constexpr auto | problem_dimension = View::problem_dimension |
Implicit backward-Euler time stepping for the parabolic limiting equation [11], Eq. 3.3:
\begin{align} \newcommand{\bbm}{{\boldsymbol m}} \newcommand{\bef}{{\boldsymbol f}} \newcommand{\bk}{{\boldsymbol k}} \newcommand{\bu}{{\boldsymbol u}} \newcommand{\bv}{{\boldsymbol v}} \newcommand{\bn}{{\boldsymbol n}} \newcommand{\pols}{{\mathbb s}} \newcommand{\Hflux}{\bk} &\partial_t \rho = 0, \\ &\partial_t \bbm - \nabla\cdot(\pols(\bv)) = \bef, \\ &\partial_t E + \nabla\cdot(\Hflux(\bu)- \pols(\bv) \bv) = \bef\cdot\bv, \\ &\bv_{|\partial D}=\boldsymbol 0, \qquad \Hflux(\bu)\cdot\bn_{|\partial D}=0 . \end{align}
Internally, the module first performs an implicit backward Euler step updating the velocity (see [11], Eq. 5.5):
\begin{align} \begin{cases} \newcommand\bsfV{{\textbf V}} \newcommand{\polB}{{\mathbb B}} \newcommand{\calI}{{\mathcal I}} \newcommand\bsfF{{\textbf F}} \newcommand\bsfM{{\textbf M}} \newcommand{\upint}{^\circ} \newcommand{\upbnd}{^\partial} \newcommand{\dt}{{\tau}} \newcommand{\calV}{{\mathcal V}} \varrho^{n}_i m_i \bsfV^{n+1} + \dt\sum_{j\in\calI(i)} \polB_{ij} \bsfV^{n+1} = m_i \bsfM_i^{n} + \dt m_i \bsfF_i^{n+1}, & \forall i\in \calV\upint \\[0.3em] \bsfV_i^{n+1} = \boldsymbol 0, & \forall i\in \calV\upbnd, \end{cases} \end{align}
We then postprocess and compute an internal energy update with an additional backward Euler step, (cf. [11], Eq. 5.13)
\begin{align} \newcommand\bsfV{{\textbf V}} \newcommand\sfe{{\mathsf e}} \newcommand{\upHnph}{^{\text{H},n+1}} \newcommand{\calI}{{\mathcal I}} \newcommand\sfK{{\mathsf K}} \newcommand{\calV}{{\mathcal V}} m_i \varrho_i^{n}(\sfe_i{\upHnph} - \sfe_i^{n})+\dt \sum_{j\in\calI(i)} \beta_{ij}\sfe_i{\upHnph} = \tfrac12 m_i\|\bsfV^{n+1}-\bsfV^{n}\|^2 + \dt m_i\sfK_i^{n+1}, \qquad \forall i\in \calV. \end{align}
The result is then transformed back into conserved quantities and written to the output vector.
Definition at line 114 of file parabolic_solver.h.
using ryujin::NavierStokes::ParabolicSolver< Description, dim, Number >::HyperbolicSystem = typename Description::HyperbolicSystem |
Definition at line 122 of file parabolic_solver.h.
using ryujin::NavierStokes::ParabolicSolver< Description, dim, Number >::View = typename Description::template HyperbolicSystemView<dim, Number> |
Definition at line 124 of file parabolic_solver.h.
using ryujin::NavierStokes::ParabolicSolver< Description, dim, Number >::ParabolicSystem = typename Description::ParabolicSystem |
Definition at line 127 of file parabolic_solver.h.
using ryujin::NavierStokes::ParabolicSolver< Description, dim, Number >::ScalarNumber = typename View::ScalarNumber |
Definition at line 129 of file parabolic_solver.h.
using ryujin::NavierStokes::ParabolicSolver< Description, dim, Number >::state_type = typename View::state_type |
Definition at line 133 of file parabolic_solver.h.
using ryujin::NavierStokes::ParabolicSolver< Description, dim, Number >::StateVector = typename View::StateVector |
Definition at line 135 of file parabolic_solver.h.
using ryujin::NavierStokes::ParabolicSolver< Description, dim, Number >::ScalarVector = Vectors::ScalarVector<Number> |
Definition at line 137 of file parabolic_solver.h.
using ryujin::NavierStokes::ParabolicSolver< Description, dim, Number >::BlockVector = Vectors::BlockVector<Number> |
Definition at line 139 of file parabolic_solver.h.
ryujin::NavierStokes::ParabolicSolver< Description, dim, Number >::ParabolicSolver | ( | const MPIEnsemble & | mpi_ensemble, |
std::map< std::string, dealii::Timer > & | computing_timer, | ||
const HyperbolicSystem & | hyperbolic_system, | ||
const ParabolicSystem & | parabolic_system, | ||
const OfflineData< dim, Number > & | offline_data, | ||
const InitialValues< Description, dim, Number > & | initial_values, | ||
const std::string & | subsection = "ParabolicSolver< Description, dim, Number >" |
||
) |
Constructor.
Definition at line 34 of file parabolic_solver.template.h.
void ryujin::NavierStokes::ParabolicSolver< Description, dim, Number >::prepare |
Prepare time stepping. A call to prepare() allocates temporary storage and is necessary before any of the following time-stepping functions can be called.
Definition at line 125 of file parabolic_solver.template.h.
References ryujin::cg_q1, ryujin::dirichlet, and ryujin::no_slip.
void ryujin::NavierStokes::ParabolicSolver< Description, dim, Number >::backward_euler_step | ( | const StateVector & | old_state_vector, |
const Number | old_t, | ||
StateVector & | new_state_vector, | ||
Number | tau, | ||
const IDViolationStrategy | id_violation_strategy, | ||
const bool | reinitialize_gmg | ||
) | const |
Given a reference to a previous state vector old_state_vector
at time old_t
and a time-step size tau
perform an implicit backward Euler step (and store the result in new_state_vector
).
Definition at line 212 of file parabolic_solver.template.h.
References CALLGRIND_START_INSTRUMENTATION, CALLGRIND_STOP_INSTRUMENTATION, ryujin::dirichlet, ryujin::NavierStokes::EnergyMatrix< dim, Number, Number2 >::initialize(), ryujin::NavierStokes::VelocityMatrix< dim, Number, Number2 >::initialize(), LIKWID_MARKER_START, LIKWID_MARKER_STOP, ryujin::negative_part(), ryujin::no_slip, ryujin::raise_exception, ryujin::NavierStokes::DiagonalMatrix< dim, Number >::reinit(), RYUJIN_OMP_FOR, RYUJIN_PARALLEL_REGION_BEGIN, RYUJIN_PARALLEL_REGION_END, ryujin::slip, and ryujin::warn.
void ryujin::NavierStokes::ParabolicSolver< Description, dim, Number >::print_solver_statistics | ( | std::ostream & | output | ) | const |
Print a status line with solver statistics. This function is used for constructing the status message displayed periodically in the TimeLoop.
Definition at line 871 of file parabolic_solver.template.h.
|
staticconstexpr |
Definition at line 131 of file parabolic_solver.h.