[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
flow_solver_case_base.h
1 #ifndef __FLOW_SOLVER_CASE_BASE__
2 #define __FLOW_SOLVER_CASE_BASE__
3 
4 // for FlowSolver class:
5 #include <deal.II/base/table_handler.h>
6 #include <deal.II/distributed/shared_tria.h>
7 #include <deal.II/distributed/tria.h>
8 
9 #include "dg/dg_base.hpp"
10 #include "parameters/all_parameters.h"
11 #include "physics/initial_conditions/initial_condition_function.h"
12 #include "ode_solver/ode_solver_base.h"
13 
14 namespace PHiLiP {
15 namespace FlowSolver {
16 
17 #if PHILIP_DIM==1
18 using Triangulation = dealii::Triangulation<PHILIP_DIM>;
19 #else
20 using Triangulation = dealii::parallel::distributed::Triangulation<PHILIP_DIM>;
21 #endif
22 
23 template <int dim, int nstate>
25 {
26 public:
28  explicit FlowSolverCaseBase(const Parameters::AllParameters *const parameters_input);
29 
30  std::shared_ptr<InitialConditionFunction<dim,nstate,double>> initial_condition_function;
31 
33  virtual ~FlowSolverCaseBase() = default;
34 
36  void display_flow_solver_setup(std::shared_ptr<DGBase<dim,double>> dg) const;
37 
39  virtual std::shared_ptr<Triangulation> generate_grid() const = 0;
40 
42  virtual void set_higher_order_grid(std::shared_ptr <DGBase<dim, double>> dg) const;
43 
47  const std::shared_ptr<ODE::ODESolverBase<dim, double>> ode_solver,
48  const std::shared_ptr <DGBase<dim, double>> dg,
49  const std::shared_ptr<dealii::TableHandler> unsteady_data_table);
50 
53  const unsigned int current_iteration,
54  const double current_time,
55  const std::shared_ptr <DGBase<dim, double>> dg,
56  const std::shared_ptr<dealii::TableHandler> unsteady_data_table);
57 
59  virtual double get_constant_time_step(std::shared_ptr <DGBase<dim, double>> dg) const;
60 
62  virtual double get_adaptive_time_step(std::shared_ptr <DGBase<dim, double>> dg) const;
63 
65  virtual double get_adaptive_time_step_initial(std::shared_ptr <DGBase<dim, double>> dg);
66 
68  virtual void steady_state_postprocessing(std::shared_ptr <DGBase<dim, double>> dg) const;
69 
71  void set_time_step(const double time_step_input);
72 
73 protected:
75  const MPI_Comm mpi_communicator;
76  const int mpi_rank;
77  const int n_mpi;
78 
80 
82  dealii::ConditionalOStream pcout;
83 
86  const double value,
87  const std::string value_string,
88  const std::shared_ptr <dealii::TableHandler> data_table) const;
89 
91  virtual void display_additional_flow_case_specific_parameters() const = 0;
92 
94  double get_time_step() const;
95 
96 private:
98  std::string get_pde_string() const;
99 
101  std::string get_flow_case_string() const;
102 
104  double time_step;
105 };
106 
107 } // FlowSolver namespace
108 } // PHiLiP namespace
109 
110 #endif
const Parameters::AllParameters all_param
All parameters.
const MPI_Comm mpi_communicator
MPI communicator.
virtual double get_adaptive_time_step(std::shared_ptr< DGBase< dim, double >> dg) const
Virtual function to compute the adaptive time step.
virtual void set_higher_order_grid(std::shared_ptr< DGBase< dim, double >> dg) const
Set higher order grid.
std::shared_ptr< InitialConditionFunction< dim, nstate, double > > initial_condition_function
Initial condition function.
virtual void steady_state_postprocessing(std::shared_ptr< DGBase< dim, double >> dg) const
Virtual function for postprocessing when solving for steady state.
double get_time_step() const
Getter for time step.
Files for the baseline physics.
Definition: ADTypes.hpp:10
Base class ODE solver.
std::string get_flow_case_string() const
Returns the flow case type string from the all_param class member.
void display_flow_solver_setup(std::shared_ptr< DGBase< dim, double >> dg) const
Displays the flow setup parameters.
Main parameter class that contains the various other sub-parameter classes.
const int n_mpi
Number of MPI processes.
void add_value_to_data_table(const double value, const std::string value_string, const std::shared_ptr< dealii::TableHandler > data_table) const
Add a value to a given data table with scientific format.
virtual double get_adaptive_time_step_initial(std::shared_ptr< DGBase< dim, double >> dg)
Virtual function to compute the initial adaptive time step.
virtual void display_additional_flow_case_specific_parameters() const =0
Display additional more specific flow case parameters.
virtual double get_constant_time_step(std::shared_ptr< DGBase< dim, double >> dg) const
Virtual function to compute the constant time step.
std::string get_pde_string() const
Returns the pde type string from the all_param class member.
virtual ~FlowSolverCaseBase()=default
Destructor.
virtual std::shared_ptr< Triangulation > generate_grid() const =0
Pure Virtual function to generate the grid.
FlowSolverCaseBase(const Parameters::AllParameters *const parameters_input)
Constructor.
virtual void compute_unsteady_data_and_write_to_table(const std::shared_ptr< ODE::ODESolverBase< dim, double >> ode_solver, const std::shared_ptr< DGBase< dim, double >> dg, const std::shared_ptr< dealii::TableHandler > unsteady_data_table)
DGBase is independent of the number of state variables.
Definition: dg_base.hpp:82
dealii::ConditionalOStream pcout
ConditionalOStream.
void set_time_step(const double time_step_input)
Setter for time step.