[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_butcher_base.cpp
1 #include "rk_tableau_butcher_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  : RKTableauBase<dim,real,MeshType>(n_rk_stages_input,rk_method_string_input)
11 {
12  this->butcher_tableau_a.reinit(this->n_rk_stages,this->n_rk_stages);
13  this->butcher_tableau_b.reinit(this->n_rk_stages);
14  this->butcher_tableau_c.reinit(this->n_rk_stages);
15 }
16 
17 template <int dim, typename real, typename MeshType>
19 {
20  set_a();
21  set_b();
22  set_c();
23  this->pcout << "Assigned RK method: " << this->rk_method_string << std::endl;
24 }
25 
26 template <int dim, typename real, typename MeshType>
27 double RKTableauButcherBase<dim,real, MeshType> :: get_a (const int i, const int j) const
28 {
29  return butcher_tableau_a[i][j];
30 }
31 
32 template <int dim, typename real, typename MeshType>
34 {
35  return butcher_tableau_b[i];
36 }
37 
38 
39 template <int dim, typename real, typename MeshType>
41 {
42  return butcher_tableau_c[i];
43 }
44 
47 #if PHILIP_DIM != 1
49 #endif
50 
51 } // ODE namespace
52 } // PHiLiP namespace
void set_tableau() override
Calls setters for butcher tableau.
const int n_rk_stages
Store number of stages.
const std::string rk_method_string
String identifying the RK method.
dealii::ConditionalOStream pcout
Parallel std::cout that only outputs on mpi_rank==0.
Files for the baseline physics.
Definition: ADTypes.hpp:10
virtual void set_a()=0
Setter for butcher_tableau_a.
double get_b(const int i) const
Returns Butcher tableau "b" coefficient at position [i].
double get_c(const int i) const
Returns Butcher tableau "c" coefficient at position [i].
dealii::Table< 2, double > butcher_tableau_a
Butcher tableau "a".
double get_a(const int i, const int j) const
Returns Butcher tableau "a" coefficient at position [i][j].
RKTableauButcherBase(const int n_rk_stages, const std::string rk_method_string_input)
Default constructor that will set the constants.
virtual void set_b()=0
Setter for butcher_tableau_b.
Base class for storing the RK method.
Base class for storing the RK method.
virtual void set_c()=0
Setter for butcher_tableau_c.
dealii::Table< 1, double > butcher_tableau_c
Butcher tableau "c".
dealii::Table< 1, double > butcher_tableau_b
Butcher tableau "c".