[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
rrk_ode_solver_base.h
1 #ifndef __RRK_ODE_SOLVER_BASE_H__
2 #define __RRK_ODE_SOLVER_BASE_H__
3 
4 #include "dg/dg_base.hpp"
5 #include "ode_solver/runge_kutta_ode_solver.h"
6 #include "ode_solver/relaxation_runge_kutta/runge_kutta_store_entropy.h"
7 
8 namespace PHiLiP {
9 namespace ODE {
10 
12 #if PHILIP_DIM==1
13 template <int dim, typename real, typename MeshType = dealii::Triangulation<dim>>
14 #else
15 template <int dim, typename real, typename MeshType = dealii::parallel::distributed::Triangulation<dim>>
16 #endif
17 class RRKODESolverBase: public RKNumEntropy<dim, real, MeshType>
18 {
19 public:
21  explicit RRKODESolverBase(
22  std::shared_ptr<RKTableauBase<dim,real,MeshType>> rk_tableau_input);
23 
25 
29 
30 protected:
31 
34  real update_relaxation_parameter (const real dt,
35  std::shared_ptr<DGBase<dim,real,MeshType>> dg,
36  const std::vector<dealii::LinearAlgebra::distributed::Vector<double>> &rk_stage,
37  const dealii::LinearAlgebra::distributed::Vector<double> & solution_update) override;
38 
39 
43  virtual real compute_relaxation_parameter(const real dt,
44  std::shared_ptr<DGBase<dim,real,MeshType>> dg,
45  const std::vector<dealii::LinearAlgebra::distributed::Vector<double>> &rk_stage,
46  const dealii::LinearAlgebra::distributed::Vector<double> &/*solution_update*/
47  ) = 0;
48 
49 };
50 
51 } // ODE namespace
52 } // PHiLiP namespace
53 
54 #endif
RRKODESolverBase(std::shared_ptr< RKTableauBase< dim, real, MeshType >> rk_tableau_input)
Constructor.
real relaxation_parameter
Relaxation Runge-Kutta parameter gamma^n.
Files for the baseline physics.
Definition: ADTypes.hpp:10
real update_relaxation_parameter(const real dt, std::shared_ptr< DGBase< dim, real, MeshType >> dg, const std::vector< dealii::LinearAlgebra::distributed::Vector< double >> &rk_stage, const dealii::LinearAlgebra::distributed::Vector< double > &solution_update) override
Base class for storing the RK method.
DGBase is independent of the number of state variables.
Definition: dg_base.hpp:82
Relaxation Runge-Kutta ODE solver base class.
virtual real compute_relaxation_parameter(const real dt, std::shared_ptr< DGBase< dim, real, MeshType >> dg, const std::vector< dealii::LinearAlgebra::distributed::Vector< double >> &rk_stage, const dealii::LinearAlgebra::distributed::Vector< double > &)=0