1 #ifndef __FLOW_SOLVER_CASE_BASE__ 2 #define __FLOW_SOLVER_CASE_BASE__ 5 #include <deal.II/base/table_handler.h> 6 #include <deal.II/distributed/shared_tria.h> 7 #include <deal.II/distributed/tria.h> 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" 15 namespace FlowSolver {
18 using Triangulation = dealii::Triangulation<PHILIP_DIM>;
20 using Triangulation = dealii::parallel::distributed::Triangulation<PHILIP_DIM>;
23 template <
int dim,
int nstate>
39 virtual std::shared_ptr<Triangulation>
generate_grid()
const = 0;
49 const std::shared_ptr<dealii::TableHandler> unsteady_data_table);
53 const unsigned int current_iteration,
54 const double current_time,
56 const std::shared_ptr<dealii::TableHandler> unsteady_data_table);
82 dealii::ConditionalOStream
pcout;
87 const std::string value_string,
88 const std::shared_ptr <dealii::TableHandler> data_table)
const;
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.
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.
const int mpi_rank
MPI rank.
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.
dealii::ConditionalOStream pcout
ConditionalOStream.
void set_time_step(const double time_step_input)
Setter for time step.
double time_step
Current time step.