![]() |
ryujin 2.1.1 revision 9dcb748690310d6a540ebb8b066d1a0834fc7604
|
#include <source/parabolic_module.h>
Public Types | |
using | HyperbolicSystem = typename Description::HyperbolicSystem |
using | ParabolicSystem = typename Description::ParabolicSystem |
using | ParabolicSolver = typename Description::template ParabolicSolver< dim, Number > |
using | HyperbolicSystemView = typename HyperbolicSystem::template View< dim, Number > |
using | vector_type = typename HyperbolicSystemView::vector_type |
Public Member Functions | |
ParabolicModule (const MPI_Comm &mpi_communicator, std::map< std::string, dealii::Timer > &computing_timer, const OfflineData< dim, Number > &offline_data, const HyperbolicSystem &hyperbolic_system, const ParabolicSystem ¶bolic_system, const InitialValues< Description, dim, Number > &initial_values, const std::string &subsection="/ParabolicModule") | |
void | prepare () |
Functons for performing explicit time steps | |
template<int stages> | |
void | step (const vector_type &old_U, const Number old_t, std::array< std::reference_wrapper< const vector_type >, stages > stage_U, const std::array< Number, stages > stage_weights, vector_type &new_U, const Number tau) const |
void | crank_nicolson_step (const vector_type &old_U, const Number old_t, vector_type &new_U, const Number tau) const |
void | print_solver_statistics (std::ostream &output) const |
Accessors | |
IDViolationStrategy | id_violation_strategy_ |
auto & | n_restarts () const |
auto & | n_warnings () const |
Implicit backward Euler time-stepping and Crank-Nicolson time-stepping for the parabolic subsystem.
Definition at line 34 of file parabolic_module.h.
using ryujin::ParabolicModule< Description, dim, Number >::HyperbolicSystem = typename Description::HyperbolicSystem |
Definition at line 40 of file parabolic_module.h.
using ryujin::ParabolicModule< Description, dim, Number >::ParabolicSystem = typename Description::ParabolicSystem |
Definition at line 45 of file parabolic_module.h.
using ryujin::ParabolicModule< Description, dim, Number >::ParabolicSolver = typename Description::template ParabolicSolver<dim, Number> |
The ParabolicSolver
Definition at line 50 of file parabolic_module.h.
using ryujin::ParabolicModule< Description, dim, Number >::HyperbolicSystemView = typename HyperbolicSystem::template View<dim, Number> |
Definition at line 56 of file parabolic_module.h.
using ryujin::ParabolicModule< Description, dim, Number >::vector_type = typename HyperbolicSystemView::vector_type |
Definition at line 62 of file parabolic_module.h.
ryujin::ParabolicModule< Description, dim, Number >::ParabolicModule | ( | const MPI_Comm & | mpi_communicator, |
std::map< std::string, dealii::Timer > & | computing_timer, | ||
const OfflineData< dim, Number > & | offline_data, | ||
const HyperbolicSystem & | hyperbolic_system, | ||
const ParabolicSystem & | parabolic_system, | ||
const InitialValues< Description, dim, Number > & | initial_values, | ||
const std::string & | subsection = "/ParabolicModule< Description, dim, Number >" |
||
) |
void ryujin::ParabolicModule< 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 39 of file parabolic_module.template.h.
References ryujin::ScalarConservation::ParabolicSystem::is_identity.
void ryujin::ParabolicModule< Description, dim, Number >::step | ( | const vector_type & | old_U, |
const Number | old_t, | ||
std::array< std::reference_wrapper< const vector_type >, stages > | stage_U, | ||
const std::array< Number, stages > | stage_weights, | ||
vector_type & | new_U, | ||
const Number | tau | ||
) | const |
Given a reference to a previous state vector old_U
at time old_t
and a time-step size tau
perform an implicit backward euler step (and store the result in new_U
).
The function takes an optional array of states stage_U
together with a an array of weights stage_weights
to construct a modified high-order right-hand side / flux.
Definition at line 54 of file parabolic_module.template.h.
References ryujin::ScalarConservation::ParabolicSystem::is_identity.
void ryujin::ParabolicModule< Description, dim, Number >::crank_nicolson_step | ( | const vector_type & | old_U, |
const Number | old_t, | ||
vector_type & | new_U, | ||
const Number | tau | ||
) | const |
Given a reference to a previous state vector old_U
at time old_t
and a time-step size tau
perform an implicit Crank-Nicolson step (and store the result in new_U
).
Definition at line 78 of file parabolic_module.template.h.
void ryujin::ParabolicModule< 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 102 of file parabolic_module.template.h.
References ryujin::ScalarConservation::ParabolicSystem::is_identity.
|
inline |
The number of restarts issued by the step() function.
Definition at line 131 of file parabolic_module.h.
|
inline |
The number of ID violation warnings encounterd in the step() function.
Definition at line 137 of file parabolic_module.h.
|
mutable |
Definition at line 140 of file parabolic_module.h.