[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.h
1 #ifndef __RK_TABLEAU_BASE__
2 #define __RK_TABLEAU_BASE__
3 
4 #include <deal.II/base/conditional_ostream.h>
5 
6 #include <deal.II/grid/tria.h>
7 #include <deal.II/distributed/shared_tria.h>
8 #include <deal.II/distributed/tria.h>
9 
10 namespace PHiLiP {
11 namespace ODE {
12 
14 #if PHILIP_DIM==1
15 template <int dim, typename real, typename MeshType = dealii::Triangulation<dim>>
16 #else
17 template <int dim, typename real, typename MeshType = dealii::parallel::distributed::Triangulation<dim>>
18 #endif
20 {
21 public:
23  RKTableauBase(const int n_rk_stages, const std::string rk_method_string_input);
24 
26  virtual ~RKTableauBase() = default;
27 
29  double get_a(const int i, const int j) const;
30 
32  double get_b(const int i) const;
33 
35  double get_c(const int i) const;
36 
38  void set_tableau();
39 
41  const int n_rk_stages;
42 
43 protected:
44 
45  dealii::ConditionalOStream pcout;
46 
48  const std::string rk_method_string;
49 
51  dealii::Table<2,double> butcher_tableau_a;
52 
54  dealii::Table<1,double> butcher_tableau_b;
55 
57  dealii::Table<1,double> butcher_tableau_c;
58 
60  virtual void set_a() = 0;
61 
63  virtual void set_b() = 0;
64 
66  virtual void set_c() = 0;
67 
68 
69 };
70 
71 } // ODE namespace
72 } // PHiLiP namespace
73 
74 #endif
virtual void set_b()=0
Setter for butcher_tableau_b.
virtual ~RKTableauBase()=default
Destructor.
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].