[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
rk_tableau_base.cpp
1 #include "rk_tableau_base.h"
2 
3 
4 namespace PHiLiP {
5 namespace ODE {
6 
7 template <int dim, typename real, typename MeshType>
9  const std::string rk_method_string_input)
10  : n_rk_stages(n_rk_stages_input)
11  , pcout(std::cout, dealii::Utilities::MPI::this_mpi_process(MPI_COMM_WORLD)==0)
12  , rk_method_string(rk_method_string_input)
13 {
15  this->butcher_tableau_b.reinit(n_rk_stages);
16  this->butcher_tableau_c.reinit(n_rk_stages);
17 }
18 
19 template <int dim, typename real, typename MeshType>
21 {
22  set_a();
23  set_b();
24  set_c();
25  pcout << "Assigned RK method: " << rk_method_string << std::endl;
26 }
27 
28 template <int dim, typename real, typename MeshType>
29 double RKTableauBase<dim,real, MeshType> :: get_a (const int i, const int j) const
30 {
31  return butcher_tableau_a[i][j];
32 }
33 
34 template <int dim, typename real, typename MeshType>
36 {
37  return butcher_tableau_b[i];
38 }
39 
40 template <int dim, typename real, typename MeshType>
42 {
43  return butcher_tableau_c[i];
44 }
45 
48 #if PHILIP_DIM != 1
50 #endif
51 
52 } // ODE namespace
53 } // PHiLiP namespace
virtual void set_b()=0
Setter for butcher_tableau_b.
const int n_rk_stages
Store number of stages.
const std::string rk_method_string
String identifying the RK method.
double get_a(const int i, const int j) const
Returns Butcher tableau "a" coefficient at position [i][j].
dealii::ConditionalOStream pcout
Parallel std::cout that only outputs on mpi_rank==0.
dealii::Table< 1, double > butcher_tableau_c
Butcher tableau "c".
Files for the baseline physics.
Definition: ADTypes.hpp:10
dealii::Table< 1, double > butcher_tableau_b
Butcher tableau "b".
RKTableauBase(const int n_rk_stages, const std::string rk_method_string_input)
Default constructor that will set the constants.
virtual void set_a()=0
Setter for butcher_tableau_a.
dealii::Table< 2, double > butcher_tableau_a
Butcher tableau "a".
void set_tableau()
Calls setters for butcher tableau.
virtual void set_c()=0
Setter for butcher_tableau_c.
Base class for storing the RK method.
double get_c(const int i) const
Returns Butcher tableau "c" coefficient at position [i].
double get_b(const int i) const
Returns Butcher tableau "b" coefficient at position [i].