[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
algebraic_rrk_ode_solver.h
1 #ifndef __ENERGY_RRK_ODESOLVER_H__
2 #define __ENERGY_RRK_ODESOLVER_H__
3 
4 #include "dg/dg_base.hpp"
5 #include "rrk_ode_solver_base.h"
6 #include "ode_solver/ode_solver_base.h"
7 #include "ode_solver/runge_kutta_ode_solver.h"
8 
9 namespace PHiLiP {
10 namespace ODE {
11 
15 #if PHILIP_DIM==1
16 template <int dim, typename real, typename MeshType = dealii::Triangulation<dim>>
17 #else
18 template <int dim, typename real, typename MeshType = dealii::parallel::distributed::Triangulation<dim>>
19 #endif
20 class AlgebraicRRKODESolver: public RRKODESolverBase<dim, real, MeshType>
21 {
22 public:
24  explicit AlgebraicRRKODESolver(
25  std::shared_ptr<RKTableauBase<dim,real,MeshType>> rk_tableau_input);
26 
27 protected:
28 
31  real compute_relaxation_parameter(const real dt,
32  std::shared_ptr<DGBase<dim,real,MeshType>> dg,
33  const std::vector<dealii::LinearAlgebra::distributed::Vector<double>> &rk_stage,
34  const dealii::LinearAlgebra::distributed::Vector<double> &/*solution_update*/
35  ) override;
36 
39  const dealii::LinearAlgebra::distributed::Vector<double> &stage_i,
40  const dealii::LinearAlgebra::distributed::Vector<double> &stage_j,
41  std::shared_ptr<DGBase<dim,real,MeshType>> dg
42  ) const;
43 
44 };
45 
46 } // ODE namespace
47 } // PHiLiP namespace
48 
49 #endif
AlgebraicRRKODESolver(std::shared_ptr< RKTableauBase< dim, real, MeshType >> rk_tableau_input)
Default constructor that will set the constants.
Files for the baseline physics.
Definition: ADTypes.hpp:10
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 > &) override
real compute_inner_product(const dealii::LinearAlgebra::distributed::Vector< double > &stage_i, const dealii::LinearAlgebra::distributed::Vector< double > &stage_j, std::shared_ptr< DGBase< dim, real, MeshType >> dg) const
Compute inner product according to the nodes being used.
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.