1 #include "flow_solver_case_base.h" 6 template<
int dim,
int nstate>
8 : initial_condition_function(
InitialConditionFactory<dim, nstate, double>::create_InitialConditionFunction(parameters_input))
9 , all_param(*parameters_input)
10 , mpi_communicator(MPI_COMM_WORLD)
11 , mpi_rank(dealii::Utilities::MPI::this_mpi_process(MPI_COMM_WORLD))
12 , n_mpi(dealii::Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD))
13 , pcout(
std::cout, mpi_rank==0)
16 template<
int dim,
int nstate>
25 std::string pde_string;
26 if (pde_type == PDE_enum::advection) {pde_string =
"advection";}
27 if (pde_type == PDE_enum::advection_vector) {pde_string =
"advection_vector";}
28 if (pde_type == PDE_enum::diffusion) {pde_string =
"diffusion";}
29 if (pde_type == PDE_enum::convection_diffusion) {pde_string =
"convection_diffusion";}
30 if (pde_type == PDE_enum::burgers_inviscid) {pde_string =
"burgers_inviscid";}
31 if (pde_type == PDE_enum::burgers_viscous) {pde_string =
"burgers_viscous";}
32 if (pde_type == PDE_enum::burgers_rewienski) {pde_string =
"burgers_rewienski";}
33 if (pde_type == PDE_enum::mhd) {pde_string =
"mhd";}
34 if (pde_type == PDE_enum::euler) {pde_string =
"euler";}
35 if (pde_type == PDE_enum::navier_stokes) {pde_string =
"navier_stokes";}
36 if (pde_type == PDE_enum::physics_model) {
37 pde_string =
"physics_model";
40 std::string model_string =
"WARNING: invalid model";
41 if(model == Model_enum::large_eddy_simulation) {
43 model_string =
"large_eddy_simulation";
46 std::string sgs_model_string =
"WARNING: invalid SGS model";
48 if (sgs_model==SGSModel_enum::smagorinsky) sgs_model_string =
"smagorinsky";
49 else if(sgs_model==SGSModel_enum::wall_adaptive_local_eddy_viscosity) sgs_model_string =
"wall_adaptive_local_eddy_viscosity";
50 else if(sgs_model==SGSModel_enum::vreman) sgs_model_string =
"vreman";
51 pde_string += std::string(
" (Model: ") + model_string + std::string(
", SGS Model: ") + sgs_model_string + std::string(
")");
53 else if(model == Model_enum::reynolds_averaged_navier_stokes) {
55 model_string =
"reynolds_averaged_navier_stokes";
58 std::string rans_model_string =
"WARNING: invalid RANS model";
60 if (rans_model==RANSModel_enum::SA_negative) rans_model_string =
"SA_negative";
61 pde_string += std::string(
" (Model: ") + model_string + std::string(
", RANS Model: ") + rans_model_string + std::string(
")");
63 if(pde_string ==
"physics_model") pde_string += std::string(
" (Model: ") + model_string + std::string(
")");
69 template<
int dim,
int nstate>
76 std::string flow_case_string;
77 if (flow_case_type == FlowCaseEnum::taylor_green_vortex) {flow_case_string =
"taylor_green_vortex";}
78 if (flow_case_type == FlowCaseEnum::burgers_viscous_snapshot) {flow_case_string =
"burgers_viscous_snapshot";}
79 if (flow_case_type == FlowCaseEnum::burgers_rewienski_snapshot) {flow_case_string =
"burgers_rewienski_snapshot";}
80 if (flow_case_type == FlowCaseEnum::naca0012) {flow_case_string =
"naca0012";}
81 if (flow_case_type == FlowCaseEnum::periodic_1D_unsteady) {flow_case_string =
"periodic_1D_unsteady";}
82 if (flow_case_type == FlowCaseEnum::gaussian_bump) {flow_case_string =
"gaussian_bump";}
83 if (flow_case_type == FlowCaseEnum::advection_limiter) {flow_case_string =
"advection_limiter";}
85 return flow_case_string;
88 template <
int dim,
int nstate>
92 pcout <<
"- PDE Type: " << pde_string <<
" " <<
"(dim=" << dim <<
", nstate=" << nstate <<
")" << std::endl;
97 const unsigned int number_of_degrees_of_freedom_per_state = dg->dof_handler.n_dofs()/nstate;
98 const double number_of_degrees_of_freedom_per_dim = pow(number_of_degrees_of_freedom_per_state,(1.0/dim));
99 pcout <<
"- Degrees of freedom (per state): " << number_of_degrees_of_freedom_per_state <<
" " <<
"(" << number_of_degrees_of_freedom_per_dim <<
" per state per dim)" << std::endl;
100 pcout <<
"- Number of active cells: " << dg->triangulation->n_global_active_cells() << std::endl;
104 if (use_weak_form ==
false){
105 this->
pcout <<
"- Using strong DG" << std::endl;
108 std::string c_parameter_string;
111 if (fr_type == FREnum::cDG) c_parameter_string =
"cDG";
112 else if (fr_type == FREnum::cSD) c_parameter_string =
"cSD";
113 else if (fr_type == FREnum::cHU) c_parameter_string =
"cHU";
114 else if (fr_type == FREnum::cNegative) c_parameter_string =
"cNegative";
115 else if (fr_type == FREnum::cNegative2) c_parameter_string =
"cNegative2";
116 else if (fr_type == FREnum::cPlus) c_parameter_string =
"cPlus";
117 else if (fr_type == FREnum::c10Thousand) c_parameter_string =
"c10Thousand";
118 else if (fr_type == FREnum::cHULumped) c_parameter_string =
"cHULumped";
120 if (c_parameter_string ==
"cDG" ) {
123 this->
pcout <<
"- - Using flux reconstruction c parameter: " << c_parameter_string << std::endl;
128 this->
pcout <<
"- - Using split form " << std::endl;
132 this->
pcout <<
"- Using weak DG" << std::endl;
137 pcout <<
"- Flow case: " << flow_case_string <<
" " << std::flush;
139 pcout <<
"(Steady state)" << std::endl;
141 pcout <<
"(Unsteady)" << std::endl;
148 template <
int dim,
int nstate>
154 template <
int dim,
int nstate>
166 template <
int dim,
int nstate>
169 pcout <<
"ERROR: Base definition for get_adaptive_time_step() has not yet been implemented. " <<std::flush;
174 template <
int dim,
int nstate>
177 pcout <<
"ERROR: Base definition for get_adaptive_time_step_initial() has not yet been implemented. " <<std::flush;
182 template <
int dim,
int nstate>
188 template <
int dim,
int nstate>
192 const std::shared_ptr <dealii::TableHandler> unsteady_data_table)
197 template <
int dim,
int nstate>
202 const std::shared_ptr <dealii::TableHandler> )
207 template <
int dim,
int nstate>
210 const std::string value_string,
211 const std::shared_ptr <dealii::TableHandler> data_table)
const 213 data_table->add_value(value_string, value);
214 data_table->set_precision(value_string, 16);
215 data_table->set_scientific(value_string,
true);
218 template <
int dim,
int nstate>
220 const double time_step_input)
225 template <
int dim,
int nstate>
const Parameters::AllParameters all_param
All parameters.
virtual double get_adaptive_time_step(std::shared_ptr< DGBase< dim, double >> dg) const
Virtual function to compute the adaptive time step.
FlowCaseType
Selects the flow case to be simulated.
double final_time
Final solution time.
virtual void set_higher_order_grid(std::shared_ptr< DGBase< dim, double >> dg) const
Set higher order grid.
PartialDifferentialEquation pde_type
Store the PDE type to be solved.
FlowCaseType flow_case_type
Selected FlowCaseType from the input file.
bool steady_state
Flag for solving steady state solution.
FlowSolverParam flow_solver_param
Contains the parameters for simulation cases (flow solver test)
virtual void steady_state_postprocessing(std::shared_ptr< DGBase< dim, double >> dg) const
Virtual function for postprocessing when solving for steady state.
double constant_time_step
Constant time step.
double get_time_step() const
Getter for time step.
PartialDifferentialEquation
Possible Partial Differential Equations to solve.
Files for the baseline physics.
bool use_weak_form
Flag to use weak or strong form of DG.
std::string get_flow_case_string() const
Returns the flow case type string from the all_param class member.
Flux_Reconstruction
Type of correction in Flux Reconstruction.
ModelType
Types of models available.
unsigned int poly_degree
Polynomial order (P) of the basis functions for DG.
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.
SubGridScaleModel SGS_model_type
Store the SubGridScale (SGS) model type.
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.
Flux_Reconstruction flux_reconstruction_type
Store flux reconstruction type.
ODESolverParam ode_solver_param
Contains parameters for ODE solver.
double initial_time_step
Time step used in ODE solver.
virtual double get_adaptive_time_step_initial(std::shared_ptr< DGBase< dim, double >> dg)
Virtual function to compute the initial adaptive time step.
Initial condition function factory.
virtual void display_additional_flow_case_specific_parameters() const =0
Display additional more specific flow case parameters.
bool use_split_form
Flag to use split form.
ReynoldsAveragedNavierStokesModel
Types of Reynolds-averaged Navier-Stokes (RANS) models that can be used.
virtual double get_constant_time_step(std::shared_ptr< DGBase< dim, double >> dg) const
Virtual function to compute the constant time step.
ReynoldsAveragedNavierStokesModel RANS_model_type
Store the Reynolds-averaged Navier-Stokes (RANS) model type.
std::string get_pde_string() const
Returns the pde type string from the all_param class member.
FlowSolverCaseBase(const Parameters::AllParameters *const parameters_input)
Constructor.
SubGridScaleModel
Types of sub-grid scale (SGS) models that can be used.
ModelType model_type
Store the model type.
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.
PhysicsModelParam physics_model_param
Contains parameters for Physics Model.
unsigned int max_poly_degree_for_adaptation
Maximum polynomial order of the DG basis functions for adaptation.
double time_step
Current time step.