[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
low_storage_rk_tableau_base.cpp
1 #include "low_storage_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  : rk_method_string(rk_method_string_input)
11  , pcout(std::cout, dealii::Utilities::MPI::this_mpi_process(MPI_COMM_WORLD)==0)
12 {
13  this->butcher_tableau_gamma.reinit(n_rk_stages+1,3);
14  this->butcher_tableau_beta.reinit(n_rk_stages+1);
15  this->butcher_tableau_delta.reinit(num_delta);
16  this->butcher_tableau_b_hat.reinit(n_rk_stages+1);
17 }
18 
19 template <int dim, typename real, typename MeshType>
21 {
22  set_gamma();
23  set_beta();
24  set_delta();
25  set_b_hat();
26  pcout << "Assigned RK method: " << rk_method_string << std::endl;
27 }
28 
29 template <int dim, typename real, typename MeshType>
30 double LowStorageRKTableauBase<dim,real, MeshType> :: get_gamma (const int i, const int j) const
31 {
32  return butcher_tableau_gamma[i][j];
33 }
34 
35 template <int dim, typename real, typename MeshType>
37 {
38  return butcher_tableau_beta[i];
39 }
40 
41 template <int dim, typename real, typename MeshType>
43 {
44  return butcher_tableau_delta[i];
45 }
46 
47 template <int dim, typename real, typename MeshType>
49 {
50  return butcher_tableau_b_hat[i];
51 }
52 
55 #if PHILIP_DIM != 1
57 #endif
58 
59 } // ODE namespace
60 } // PHiLiP namespace
LowStorageRKTableauBase(const int n_rk_stages, const int num_delta, const std::string rk_method_string_input)
Default constructor that will set the constants.
double get_b_hat(const int i) const
Returns Butcher tableau "b hat" coefficient at position [i].
double get_beta(const int i) const
Returns Butcher tableau "beta" coefficient at position [i].
dealii::Table< 1, double > butcher_tableau_b_hat
Butcher tableau "b hat".
double get_gamma(const int i, const int j) const
Returns Butcher tableau "gamma" coefficient at position [i][j].
dealii::Table< 2, double > butcher_tableau_gamma
Butcher tableau "gamma".
Files for the baseline physics.
Definition: ADTypes.hpp:10
Base class for storing the RK method.
void set_tableau()
Calls setters for butcher tableau.
dealii::ConditionalOStream pcout
Parallel std::cout that only outputs on mpi_rank==0.
virtual void set_gamma()=0
Setter for gamma.
virtual void set_b_hat()=0
Setter for b hat.
dealii::Table< 1, double > butcher_tableau_beta
Butcher tableau "beta".
virtual void set_beta()=0
Setter for beta.
const std::string rk_method_string
String identifying the RK method.
dealii::Table< 1, double > butcher_tableau_delta
Butcher tableau "delta".
virtual void set_delta()=0
Setter for delta.
double get_delta(const int i) const
Returns Butcher tableau "delta" coefficient at position [i].