[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.cpp
1 #include "flow_solver_case_base.h"
2 
3 namespace PHiLiP {
4 namespace FlowSolver {
5 
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)
14  {}
15 
16 template<int dim, int nstate>
18 {
20  using Model_enum = Parameters::AllParameters::ModelType;
23 
24  const PDE_enum pde_type = this->all_param.pde_type;
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";
38  // add the model name + sub model name (if applicable)
39  const Model_enum model = this->all_param.model_type;
40  std::string model_string = "WARNING: invalid model";
41  if(model == Model_enum::large_eddy_simulation) {
42  // assign model string
43  model_string = "large_eddy_simulation";
44  // sub-grid scale (SGS)
45  const SGSModel_enum sgs_model = this->all_param.physics_model_param.SGS_model_type;
46  std::string sgs_model_string = "WARNING: invalid SGS model";
47  // assign SGS model string
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(")");
52  }
53  else if(model == Model_enum::reynolds_averaged_navier_stokes) {
54  // assign model string
55  model_string = "reynolds_averaged_navier_stokes";
56  // reynolds-averaged navier-stokes (RANS)
57  const RANSModel_enum rans_model = this->all_param.physics_model_param.RANS_model_type;
58  std::string rans_model_string = "WARNING: invalid RANS model";
59  // assign RANS model string
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(")");
62  }
63  if(pde_string == "physics_model") pde_string += std::string(" (Model: ") + model_string + std::string(")");
64  }
65 
66  return pde_string;
67 }
68 
69 template<int dim, int nstate>
71 {
72  // Get the flow case type
73  using FlowCaseEnum = Parameters::FlowSolverParam::FlowCaseType;
74  const FlowCaseEnum flow_case_type = this->all_param.flow_solver_param.flow_case_type;
75 
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";}
84 
85  return flow_case_string;
86 }
87 
88 template <int dim, int nstate>
90 {
91  const std::string pde_string = this->get_pde_string();
92  pcout << "- PDE Type: " << pde_string << " " << "(dim=" << dim << ", nstate=" << nstate << ")" << std::endl;
93 
94  pcout << "- Polynomial degree: " << this->all_param.flow_solver_param.poly_degree << std::endl;
95  pcout << "- Maximum polynomial degree for adaptation: " << this->all_param.flow_solver_param.max_poly_degree_for_adaptation << std::endl;
96 
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;
101 
102  const bool use_weak_form = this->all_param.use_weak_form;
103 
104  if (use_weak_form == false){
105  this->pcout << "- Using strong DG" << std::endl;
106 
107  // only print c param for strong DG as FR is implemented only for strong
108  std::string c_parameter_string;
110  FREnum fr_type = this->all_param.flux_reconstruction_type;
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";
119 
120  if (c_parameter_string == "cDG" ) {
121  // No additional output to indicate classical strong DG
122  } else {
123  this->pcout << "- - Using flux reconstruction c parameter: " << c_parameter_string << std::endl;
124  }
125 
126  const bool use_split_form = this->all_param.use_split_form;
127  if (use_split_form){
128  this->pcout << "- - Using split form " << std::endl;
129  }
130  }
131  else{
132  this->pcout << "- Using weak DG" << std::endl;
133 
134  }
135 
136  const std::string flow_case_string = this->get_flow_case_string();
137  pcout << "- Flow case: " << flow_case_string << " " << std::flush;
138  if(this->all_param.flow_solver_param.steady_state == true) {
139  pcout << "(Steady state)" << std::endl;
140  } else {
141  pcout << "(Unsteady)" << std::endl;
142  pcout << "- - Final time: " << this->all_param.flow_solver_param.final_time << std::endl;
143  }
144 
146 }
147 
148 template <int dim, int nstate>
150 {
151  // Do nothing
152 }
153 
154 template <int dim, int nstate>
156 {
158  // Using constant time step in FlowSolver parameters.
160  } else {
161  // Using initial time step in ODE parameters.
163  }
164 }
165 
166 template <int dim, int nstate>
168 {
169  pcout << "ERROR: Base definition for get_adaptive_time_step() has not yet been implemented. " <<std::flush;
170  std::abort();
171  return 0.0;
172 }
173 
174 template <int dim, int nstate>
176 {
177  pcout << "ERROR: Base definition for get_adaptive_time_step_initial() has not yet been implemented. " <<std::flush;
178  std::abort();
179  return 0.0;
180 }
181 
182 template <int dim, int nstate>
184 {
185  // do nothing by default
186 }
187 
188 template <int dim, int nstate>
190  const std::shared_ptr<ODE::ODESolverBase<dim, double>> ode_solver,
191  const std::shared_ptr <DGBase<dim, double>> dg,
192  const std::shared_ptr <dealii::TableHandler> unsteady_data_table)
193 {
194  this->compute_unsteady_data_and_write_to_table(ode_solver->current_iteration, ode_solver->current_time, dg, unsteady_data_table);
195 }
196 
197 template <int dim, int nstate>
199  const unsigned int /*current_iteration*/,
200  const double /*current_time*/,
201  const std::shared_ptr <DGBase<dim, double>> /*dg*/,
202  const std::shared_ptr <dealii::TableHandler> /*unsteady_data_table*/)
203 {
204  // do nothing by default
205 }
206 
207 template <int dim, int nstate>
209  const double value,
210  const std::string value_string,
211  const std::shared_ptr <dealii::TableHandler> data_table) const
212 {
213  data_table->add_value(value_string, value);
214  data_table->set_precision(value_string, 16);
215  data_table->set_scientific(value_string, true);
216 }
217 
218 template <int dim, int nstate>
220  const double time_step_input)
221 {
222  this->time_step = time_step_input;
223 }
224 
225 template <int dim, int nstate>
227 {
228  return this->time_step;
229 }
230 
231 template class FlowSolverCaseBase<PHILIP_DIM,1>;
232 template class FlowSolverCaseBase<PHILIP_DIM,2>;
233 template class FlowSolverCaseBase<PHILIP_DIM,3>;
234 template class FlowSolverCaseBase<PHILIP_DIM,4>;
235 template class FlowSolverCaseBase<PHILIP_DIM,5>;
236 template class FlowSolverCaseBase<PHILIP_DIM,6>;
237 
238 } // FlowSolver namespace
239 } // PHiLiP namespace
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.
Definition: ADTypes.hpp:10
Base class ODE solver.
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.
Definition: dg_base.hpp:82
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.