[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
tests.cpp
1 #include <iostream>
2 
3 #include <deal.II/grid/grid_out.h>
4 
5 #include <deal.II/distributed/solution_transfer.h>
6 
7 #include <deal.II/grid/tria.h>
8 #include <deal.II/distributed/shared_tria.h>
9 #include <deal.II/distributed/tria.h>
10 
11 #include "tests.h"
12 #include "grid_study.h"
13 #include "grid_refinement_study.h"
14 #include "burgers_stability.h"
15 #include "diffusion_exact_adjoint.h"
16 #include "euler_gaussian_bump.h"
17 #include "euler_gaussian_bump_enthalpy_check.h"
18 #include "euler_gaussian_bump_adjoint.h"
19 #include "euler_cylinder.h"
20 #include "euler_cylinder_adjoint.h"
21 #include "euler_vortex.h"
22 #include "euler_entropy_waves.h"
23 #include "advection_explicit_periodic.h"
24 #include "euler_split_inviscid_taylor_green_vortex.h"
25 #include "TGV_scaling.h"
26 #include "optimization_inverse_manufactured/optimization_inverse_manufactured.h"
27 #include "euler_bump_optimization.h"
28 #include "euler_naca0012_optimization.hpp"
29 #include "shock_1d.h"
30 #include "euler_naca0012.hpp"
31 #include "reduced_order.h"
32 #include "unsteady_reduced_order.h"
33 #include "convection_diffusion_explicit_periodic.h"
34 #include "dual_weighted_residual_mesh_adaptation.h"
35 #include "anisotropic_mesh_adaptation_cases.h"
36 #include "pod_adaptive_sampling_run.h"
37 #include "pod_adaptive_sampling_testing.h"
38 #include "taylor_green_vortex_energy_check.h"
39 #include "taylor_green_vortex_restart_check.h"
40 #include "time_refinement_study.h"
41 #include "time_refinement_study_reference.h"
42 #include "h_refinement_study_isentropic_vortex.h"
43 #include "rrk_numerical_entropy_conservation_check.h"
44 #include "euler_entropy_conserving_split_forms_check.h"
45 #include "homogeneous_isotropic_turbulence_initialization_check.h"
46 #include "khi_robustness.h"
47 #include "bound_preserving_limiter_tests.h"
48 #include "naca0012_unsteady_check_quick.h"
49 #include "build_NNLS_problem.h"
50 #include "hyper_reduction_comparison.h"
51 #include "hyper_adaptive_sampling_run.h"
52 #include "hyper_reduction_post_sampling.h"
53 #include "ROM_error_post_sampling.h"
54 #include "HROM_error_post_sampling.h"
55 #include "hyper_adaptive_sampling_new_error.h"
56 #include "halton_sampling_run.h"
57 
58 namespace PHiLiP {
59 namespace Tests {
60 
61 using AllParam = Parameters::AllParameters;
62 
63 TestsBase::TestsBase(Parameters::AllParameters const *const parameters_input)
64  : all_parameters(parameters_input)
65  , mpi_communicator(MPI_COMM_WORLD)
66  , mpi_rank(dealii::Utilities::MPI::this_mpi_process(MPI_COMM_WORLD))
67  , n_mpi(dealii::Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD))
68  , pcout(std::cout, mpi_rank==0)
69 {}
70 
71 std::vector<int> TestsBase::get_number_1d_cells(const int n_grids) const
72 {
73  std::vector<int> n_1d_cells(n_grids);
75  n_1d_cells[0] = param.initial_grid_size;
76  for (int igrid=1;igrid<n_grids;++igrid) {
77  n_1d_cells[igrid] = static_cast<int>(n_1d_cells[igrid-1]*param.grid_progression) + param.grid_progression_add;
78  }
79  return n_1d_cells;
80 }
81 
82 std::string TestsBase::get_pde_string(const Parameters::AllParameters *const param) const
83 {
88 
89  const PDE_enum pde_type = param->pde_type;
90  std::string pde_string;
91  if (pde_type == PDE_enum::advection) {pde_string = "advection";}
92  if (pde_type == PDE_enum::advection_vector) {pde_string = "advection_vector";}
93  if (pde_type == PDE_enum::diffusion) {pde_string = "diffusion";}
94  if (pde_type == PDE_enum::convection_diffusion) {pde_string = "convection_diffusion";}
95  if (pde_type == PDE_enum::burgers_inviscid) {pde_string = "burgers_inviscid";}
96  if (pde_type == PDE_enum::burgers_viscous) {pde_string = "burgers_viscous";}
97  if (pde_type == PDE_enum::burgers_rewienski) {pde_string = "burgers_rewienski";}
98  if (pde_type == PDE_enum::euler) {pde_string = "euler";}
99  if (pde_type == PDE_enum::navier_stokes) {pde_string = "navier_stokes";}
100  if (pde_type == PDE_enum::physics_model) {
101  pde_string = "physics_model";
102  // add the model name + sub model name (if applicable)
103  const Model_enum model = param->model_type;
104  std::string model_string = "WARNING: invalid model";
105  if(model == Model_enum::large_eddy_simulation) {
106  // assign model string
107  model_string = "large_eddy_simulation";
108  // sub-grid scale (SGS)
109  const SGSModel_enum sgs_model = param->physics_model_param.SGS_model_type;
110  std::string sgs_model_string = "WARNING: invalid SGS model";
111  // assign SGS model string
112  if (sgs_model==SGSModel_enum::smagorinsky) sgs_model_string = "smagorinsky";
113  else if(sgs_model==SGSModel_enum::wall_adaptive_local_eddy_viscosity) sgs_model_string = "wall_adaptive_local_eddy_viscosity";
114  else if(sgs_model==SGSModel_enum::vreman) sgs_model_string = "vreman";
115  pde_string += std::string(" (Model: ") + model_string + std::string(", SGS Model: ") + sgs_model_string + std::string(")");
116  }
117  else if(model == Model_enum::reynolds_averaged_navier_stokes) {
118  // assign model string
119  model_string = "reynolds_averaged_navier_stokes";
120  // reynolds-averaged navier-stokes (RANS)
121  const RANSModel_enum rans_model = param->physics_model_param.RANS_model_type;
122  std::string rans_model_string = "WARNING: invalid RANS model";
123  // assign RANS model string
124  if (rans_model==RANSModel_enum::SA_negative) rans_model_string = "SA_negative";
125  pde_string += std::string(" (Model: ") + model_string + std::string(", RANS Model: ") + rans_model_string + std::string(")");
126  }
127  if(pde_string == "physics_model") pde_string += std::string(" (Model: ") + model_string + std::string(")");
128  }
129  return pde_string;
130 }
131 
133 {
135  const CNF_enum CNF_type = param->conv_num_flux_type;
136  std::string conv_num_flux_string;
137  if (CNF_type == CNF_enum::lax_friedrichs) {conv_num_flux_string = "lax_friedrichs";}
138  if (CNF_type == CNF_enum::roe) {conv_num_flux_string = "roe";}
139  if (CNF_type == CNF_enum::l2roe) {conv_num_flux_string = "l2roe";}
140  if (CNF_type == CNF_enum::central_flux) {conv_num_flux_string = "central_flux";}
141  if (CNF_type == CNF_enum::two_point_flux) {conv_num_flux_string = "two_point_flux";}
142  if (CNF_type == CNF_enum::two_point_flux_with_lax_friedrichs_dissipation) {
143  conv_num_flux_string = "two_point_flux_with_lax_friedrichs_dissipation";
144  } if (CNF_type == CNF_enum::two_point_flux_with_roe_dissipation) {
145  conv_num_flux_string = "two_point_flux_with_roe_dissipation";
146  } if (CNF_type == CNF_enum::two_point_flux_with_l2roe_dissipation) {
147  conv_num_flux_string = "two_point_flux_with_l2roe_dissipation";
148  }
149 
150  return conv_num_flux_string;
151 }
152 
154 {
156  const DNF_enum DNF_type = param->diss_num_flux_type;
157  std::string diss_num_flux_string;
158  if (DNF_type == DNF_enum::symm_internal_penalty) {diss_num_flux_string = "symm_internal_penalty";}
159  if (DNF_type == DNF_enum::bassi_rebay_2) {diss_num_flux_string = "bassi_rebay_2";}
160  return diss_num_flux_string;
161 }
162 
164 {
166  ManParam manu_grid_conv_param = param->manufactured_convergence_study_param;
168  const ManufacturedSolutionEnum MS_type = manu_grid_conv_param.manufactured_solution_param.manufactured_solution_type;
169  std::string manufactured_solution_string;
170  if (MS_type == ManufacturedSolutionEnum::sine_solution) {manufactured_solution_string = "sine_solution";}
171  if (MS_type == ManufacturedSolutionEnum::cosine_solution) {manufactured_solution_string = "cosine_solution";}
172  if (MS_type == ManufacturedSolutionEnum::additive_solution) {manufactured_solution_string = "additive_solution";}
173  if (MS_type == ManufacturedSolutionEnum::exp_solution) {manufactured_solution_string = "exp_solution";}
174  if (MS_type == ManufacturedSolutionEnum::poly_solution) {manufactured_solution_string = "poly_solution";}
175  if (MS_type == ManufacturedSolutionEnum::even_poly_solution) {manufactured_solution_string = "even_poly_solution";}
176  if (MS_type == ManufacturedSolutionEnum::atan_solution) {manufactured_solution_string = "atan_solution";}
177  if (MS_type == ManufacturedSolutionEnum::boundary_layer_solution) {manufactured_solution_string = "boundary_layer_solution";}
178  if (MS_type == ManufacturedSolutionEnum::s_shock_solution) {manufactured_solution_string = "s_shock_solution";}
179  if (MS_type == ManufacturedSolutionEnum::quadratic_solution) {manufactured_solution_string = "quadratic_solution";}
180  if (MS_type == ManufacturedSolutionEnum::example_solution) {manufactured_solution_string = "example_solution";}
181  if (MS_type == ManufacturedSolutionEnum::navah_solution_1) {manufactured_solution_string = "navah_solution_1";}
182  if (MS_type == ManufacturedSolutionEnum::navah_solution_2) {manufactured_solution_string = "navah_solution_2";}
183  if (MS_type == ManufacturedSolutionEnum::navah_solution_3) {manufactured_solution_string = "navah_solution_3";}
184  if (MS_type == ManufacturedSolutionEnum::navah_solution_4) {manufactured_solution_string = "navah_solution_4";}
185  if (MS_type == ManufacturedSolutionEnum::navah_solution_5) {manufactured_solution_string = "navah_solution_5";}
186  return manufactured_solution_string;
187 }
188 
189 //template<int dim, int nstate>
190 // void TestsBase::globally_refine_and_interpolate(DGBase<dim, double> &dg) const
191 //{
192 // dealii::LinearAlgebra::distributed::Vector<double> old_solution(dg->solution);
193 // dealii::parallel::distributed::SolutionTransfer<dim, dealii::LinearAlgebra::distributed::Vector<double>, dealii::hp::DoFHandler<dim>> solution_transfer(dg->dof_handler);
194 // solution_transfer.prepare_for_coarsening_and_refinement(old_solution);
195 // grid.refine_global (1);
196 // dg->allocate_system ();
197 // solution_transfer.interpolate(old_solution, dg->solution);
198 // solution_transfer.clear();
199 //}
200 
201 template<int dim, int nstate, typename MeshType>
202 std::unique_ptr< TestsBase > TestsFactory<dim,nstate,MeshType>
203 ::select_mesh(const AllParam *const parameters_input,
204  dealii::ParameterHandler &parameter_handler_input) {
205  using Mesh_enum = AllParam::MeshType;
206  Mesh_enum mesh_type = parameters_input->mesh_type;
207 
208  if(mesh_type == Mesh_enum::default_triangulation) {
209  #if PHILIP_DIM == 1
210  return TestsFactory<dim,nstate,dealii::Triangulation<dim>>::select_test(parameters_input,parameter_handler_input);
211  #else
212  return TestsFactory<dim,nstate,dealii::parallel::distributed::Triangulation<dim>>::select_test(parameters_input,parameter_handler_input);
213  #endif
214  } else if(mesh_type == Mesh_enum::triangulation) {
215  return TestsFactory<dim,nstate,dealii::Triangulation<dim>>::select_test(parameters_input,parameter_handler_input);
216  } else if(mesh_type == Mesh_enum::parallel_shared_triangulation) {
217  return TestsFactory<dim,nstate,dealii::parallel::shared::Triangulation<dim>>::select_test(parameters_input,parameter_handler_input);
218  } else if(mesh_type == Mesh_enum::parallel_distributed_triangulation) {
219  #if PHILIP_DIM == 1
220  std::cout << "dealii::parallel::distributed::Triangulation is unavailible in 1D." << std::endl;
221  #else
222  return TestsFactory<dim,nstate,dealii::parallel::distributed::Triangulation<dim>>::select_test(parameters_input,parameter_handler_input);
223  #endif
224  } else {
225  std::cout << "Invalid mesh type." << std::endl;
226  }
227 
228  return nullptr;
229 }
230 
231 template<int dim, int nstate, typename MeshType>
232 std::unique_ptr< TestsBase > TestsFactory<dim,nstate,MeshType>
233 ::select_test(const AllParam *const parameters_input,
234  dealii::ParameterHandler &parameter_handler_input) {
235  using Test_enum = AllParam::TestType;
236  const Test_enum test_type = parameters_input->test_type;
237 
238  // prevent warnings for when a create_FlowSolver is not being called (explicit and implicit cases)
239  if((test_type != Test_enum::finite_difference_sensitivity) &&
240  (test_type != Test_enum::taylor_green_vortex_energy_check) &&
241  (test_type != Test_enum::taylor_green_vortex_restart_check)) {
242  (void) parameter_handler_input;
243  } else if (!((dim==3 && nstate==dim+2) || (dim==1 && nstate==1))) {
244  (void) parameter_handler_input;
245  }
246 
247  if(test_type == Test_enum::run_control) { // TO DO: rename to grid_study
248  return std::make_unique<GridStudy<dim,nstate>>(parameters_input);
249  } else if(test_type == Test_enum::grid_refinement_study) {
250  return std::make_unique<GridRefinementStudy<dim,nstate,MeshType>>(parameters_input);
251  } else if(test_type == Test_enum::burgers_energy_stability) {
252  if constexpr (dim==1 && nstate==1) return std::make_unique<BurgersEnergyStability<dim,nstate>>(parameters_input);
253  } else if(test_type == Test_enum::diffusion_exact_adjoint) {
254  if constexpr (dim>=1 && nstate==1) return std::make_unique<DiffusionExactAdjoint<dim,nstate>>(parameters_input);
255  } else if (test_type == Test_enum::advection_periodicity){
256  if constexpr (nstate == 1) return std::make_unique<AdvectionPeriodic<dim,nstate>> (parameters_input);
257  } else if (test_type == Test_enum::convection_diffusion_periodicity){
258  if constexpr (nstate == 1) return std::make_unique<ConvectionDiffusionPeriodic<dim,nstate>> (parameters_input);
259  } else if(test_type == Test_enum::euler_gaussian_bump) {
260  if constexpr (dim==2 && nstate==dim+2) return std::make_unique<EulerGaussianBump<dim,nstate>>(parameters_input,parameter_handler_input);
261  } else if(test_type == Test_enum::euler_gaussian_bump_enthalpy) {
262  if constexpr (dim==2 && nstate==dim+2) return std::make_unique<EulerGaussianBumpEnthalpyCheck<dim,nstate>>(parameters_input, parameter_handler_input);
263  //} else if(test_type == Test_enum::euler_gaussian_bump_adjoint){
264  // if constexpr (dim==2 && nstate==dim+2) return std::make_unique<EulerGaussianBumpAdjoint<dim,nstate>>(parameters_input);
265  } else if(test_type == Test_enum::euler_cylinder) {
266  if constexpr (dim==2 && nstate==dim+2) return std::make_unique<EulerCylinder<dim,nstate>>(parameters_input);
267  } else if(test_type == Test_enum::euler_cylinder_adjoint) {
268  if constexpr (dim==2 && nstate==dim+2) return std::make_unique<EulerCylinderAdjoint<dim,nstate>>(parameters_input);
269  } else if(test_type == Test_enum::euler_vortex) {
270  if constexpr (dim==2 && nstate==dim+2) return std::make_unique<EulerVortex<dim,nstate>>(parameters_input);
271  } else if(test_type == Test_enum::euler_entropy_waves) {
272  if constexpr (dim>=2 && nstate==PHILIP_DIM+2) return std::make_unique<EulerEntropyWaves<dim,nstate>>(parameters_input);
273  } else if(test_type == Test_enum::euler_split_taylor_green) {
274  if constexpr (dim==3 && nstate == dim+2) return std::make_unique<EulerTaylorGreen<dim,nstate>>(parameters_input);
275  } else if(test_type == Test_enum::taylor_green_scaling) {
276  if constexpr (dim==3 && nstate == dim+2) return std::make_unique<EulerTaylorGreenScaling<dim,nstate>>(parameters_input);
277  } else if(test_type == Test_enum::optimization_inverse_manufactured) {
278  return std::make_unique<OptimizationInverseManufactured<dim,nstate>>(parameters_input);
279  } else if(test_type == Test_enum::euler_bump_optimization) {
280  if constexpr (dim==2 && nstate==dim+2) return std::make_unique<EulerBumpOptimization<dim,nstate>>(parameters_input);
281  } else if(test_type == Test_enum::euler_naca_optimization) {
282  if constexpr (dim==2 && nstate==dim+2) return std::make_unique<EulerNACAOptimization<dim,nstate>>(parameters_input);
283  } else if(test_type == Test_enum::shock_1d) {
284  if constexpr (dim==1 && nstate==1) return std::make_unique<Shock1D<dim,nstate>>(parameters_input);
285  } else if(test_type == Test_enum::reduced_order) {
286  if constexpr ((dim==2 && nstate==dim+2) || (dim==1 && nstate==1)) return std::make_unique<ReducedOrder<dim,nstate>>(parameters_input, parameter_handler_input);
287  } else if(test_type == Test_enum::unsteady_reduced_order) {
288  if constexpr (dim==2 && nstate==dim+2) return std::make_unique<UnsteadyReducedOrder<dim,nstate>>(parameters_input, parameter_handler_input);
289  } else if(test_type == Test_enum::POD_adaptive_sampling_run) {
290  if constexpr ((dim==2 && nstate==dim+2) || (dim==1 && nstate==1)) return std::make_unique<AdaptiveSamplingRun<dim,nstate>>(parameters_input,parameter_handler_input);
291  } else if(test_type == Test_enum::adaptive_sampling_testing) {
292  if constexpr ((dim==2 && nstate==dim+2) || (dim==1 && nstate==1)) return std::make_unique<AdaptiveSamplingTesting<dim,nstate>>(parameters_input,parameter_handler_input);
293  } else if(test_type == Test_enum::euler_naca0012) {
294  if constexpr (dim==2 && nstate==dim+2) return std::make_unique<EulerNACA0012<dim,nstate>>(parameters_input,parameter_handler_input);
295  } else if(test_type == Test_enum::dual_weighted_residual_mesh_adaptation) {
296  if constexpr (dim==2 && nstate==1) return std::make_unique<DualWeightedResidualMeshAdaptation<dim, nstate>>(parameters_input,parameter_handler_input);
297  } else if(test_type == Test_enum::anisotropic_mesh_adaptation) {
298  if constexpr( (dim==2 && nstate==1) || (dim==2 && nstate==dim+2)) return std::make_unique<AnisotropicMeshAdaptationCases<dim, nstate>>(parameters_input,parameter_handler_input);
299  } else if(test_type == Test_enum::taylor_green_vortex_energy_check) {
300  if constexpr (dim==3 && nstate==dim+2) return std::make_unique<TaylorGreenVortexEnergyCheck<dim,nstate>>(parameters_input,parameter_handler_input);
301  } else if(test_type == Test_enum::taylor_green_vortex_restart_check) {
302  if constexpr (dim==3 && nstate==dim+2) return std::make_unique<TaylorGreenVortexRestartCheck<dim,nstate>>(parameters_input,parameter_handler_input);
303  } else if(test_type == Test_enum::homogeneous_isotropic_turbulence_initialization_check){
304  if constexpr (dim==3 && nstate==dim+2) return std::make_unique<HomogeneousIsotropicTurbulenceInitializationCheck<dim,nstate>>(parameters_input,parameter_handler_input);
305  } else if(test_type == Test_enum::time_refinement_study) {
306  if constexpr (dim==1 && nstate==1) return std::make_unique<TimeRefinementStudy<dim, nstate>>(parameters_input, parameter_handler_input);
307  } else if(test_type == Test_enum::h_refinement_study_isentropic_vortex) {
308  if constexpr (dim+2==nstate && dim!=1) return std::make_unique<HRefinementStudyIsentropicVortex<dim, nstate>>(parameters_input, parameter_handler_input);
309  } else if(test_type == Test_enum::time_refinement_study_reference) {
310  if constexpr (dim==1 && nstate==1) return std::make_unique<TimeRefinementStudyReference<dim, nstate>>(parameters_input, parameter_handler_input);
311  } else if(test_type == Test_enum::rrk_numerical_entropy_conservation_check) {
312  if constexpr ((dim==1 && nstate==1) || (dim==3 && nstate==dim+2)) return std::make_unique<RRKNumericalEntropyConservationCheck<dim, nstate>>(parameters_input, parameter_handler_input);
313  } else if(test_type == Test_enum::euler_entropy_conserving_split_forms_check) {
314  if constexpr (dim==3 && nstate==dim+2) return std::make_unique<EulerSplitEntropyCheck<dim, nstate>>(parameters_input, parameter_handler_input);
315  } else if(test_type == Test_enum::khi_robustness) {
316  if constexpr (dim==2 && nstate==dim+2) return std::make_unique<KHIRobustness<dim, nstate>>(parameters_input, parameter_handler_input);
317  } else if(test_type == Test_enum::build_NNLS_problem) {
318  if constexpr (dim==1 && nstate==1) return std::make_unique<BuildNNLSProblem<dim,nstate>>(parameters_input, parameter_handler_input);
319  } else if(test_type == Test_enum::hyper_reduction_comparison) {
320  if constexpr (dim==1 && nstate==1) return std::make_unique<HyperReductionComparison<dim,nstate>>(parameters_input, parameter_handler_input);
321  } else if(test_type == Test_enum::hyper_adaptive_sampling_run) {
322  if constexpr ((dim==2 && nstate==dim+2) || (dim==1 && nstate==1)) return std::make_unique<HyperAdaptiveSamplingRun<dim,nstate>>(parameters_input, parameter_handler_input);
323  } else if(test_type == Test_enum::hyper_reduction_post_sampling) {
324  if constexpr ((dim==2 && nstate==dim+2) || (dim==1 && nstate==1)) return std::make_unique<HyperReductionPostSampling<dim,nstate>>(parameters_input, parameter_handler_input);
325  } else if(test_type == Test_enum::ROM_error_post_sampling) {
326  if constexpr ((dim==2 && nstate==dim+2) || (dim==1 && nstate==1)) return std::make_unique<ROMErrorPostSampling<dim,nstate>>(parameters_input, parameter_handler_input);
327  } else if(test_type == Test_enum::HROM_error_post_sampling) {
328  if constexpr ((dim==2 && nstate==dim+2) || (dim==1 && nstate==1)) return std::make_unique<HROMErrorPostSampling<dim,nstate>>(parameters_input, parameter_handler_input);
329  } else if(test_type == Test_enum::hyper_adaptive_sampling_new_error) {
330  if constexpr ((dim==2 && nstate==dim+2) || (dim==1 && nstate==1)) return std::make_unique<HyperAdaptiveSamplingNewError<dim,nstate>>(parameters_input, parameter_handler_input);
331  } else if(test_type == Test_enum::halton_sampling_run) {
332  if constexpr ((dim==2 && nstate==dim+2) || (dim==1 && nstate==1)) return std::make_unique<HaltonSamplingRun<dim,nstate>>(parameters_input, parameter_handler_input);
333  } else if (test_type == Test_enum::advection_limiter) {
334  if constexpr (nstate == 1 && dim < 3) return std::make_unique<BoundPreservingLimiterTests<dim, nstate>>(parameters_input, parameter_handler_input);
335  } else if (test_type == Test_enum::burgers_limiter) {
336  if constexpr (nstate == dim && dim < 3) return std::make_unique<BoundPreservingLimiterTests<dim, nstate>>(parameters_input, parameter_handler_input);
337  } else if(test_type == Test_enum::low_density) {
338  if constexpr (dim==2 && nstate==dim+2) return std::make_unique<BoundPreservingLimiterTests<dim, nstate>>(parameters_input, parameter_handler_input);
339  } else if(test_type == Test_enum::naca0012_unsteady_check_quick){
340  if constexpr (dim==2 && nstate==dim+2) return std::make_unique<NACA0012UnsteadyCheckQuick<dim, nstate>>(parameters_input, parameter_handler_input);
341  } else {
342  std::cout << "Invalid test. You probably forgot to add it to the list of tests in tests.cpp" << std::endl;
343  std::abort();
344  }
345 
346  return nullptr;
347 }
348 
349 template<int dim, int nstate, typename MeshType>
350 std::unique_ptr< TestsBase > TestsFactory<dim,nstate,MeshType>
351 ::create_test(AllParam const *const parameters_input,
352  dealii::ParameterHandler &parameter_handler_input)
353 {
354  // Recursive templating required because template parameters must be compile time constants
355  // As a results, this recursive template initializes all possible dimensions with all possible nstate
356  // without having 15 different if-else statements
357  if(dim == parameters_input->dimension)
358  {
359  // This template parameters dim and nstate match the runtime parameters
360  // then create the selected test with template parameters dim and nstate
361  // Otherwise, keep decreasing nstate and dim until it matches
362  if(nstate == parameters_input->nstate)
363  return TestsFactory<dim,nstate>::select_mesh(parameters_input,parameter_handler_input);
364  else if constexpr (nstate > 1)
365  return TestsFactory<dim,nstate-1>::create_test(parameters_input,parameter_handler_input);
366  else
367  return nullptr;
368  }
369  else if constexpr (dim > 1)
370  {
371  //return TestsFactory<dim-1,nstate>::create_test(parameters_input);
372  return nullptr;
373  }
374  else
375  {
376  return nullptr;
377  }
378 }
379 
380 // Will recursively create all the possible test sizes
381 //template class TestsFactory <PHILIP_DIM,1>;
382 //template class TestsFactory <PHILIP_DIM,2>;
383 //template class TestsFactory <PHILIP_DIM,3>;
384 //template class TestsFactory <PHILIP_DIM,4>;
385 //template class TestsFactory <PHILIP_DIM,5>;
386 
389 #if PHILIP_DIM!=1
391 #endif
392 
393 } // Tests namespace
394 } // PHiLiP namespace
ManufacturedSolutionType
Selects the manufactured solution to be used if use_manufactured_source_term=true.
PartialDifferentialEquation pde_type
Store the PDE type to be solved.
TestType
Possible integration tests to run.
Parameters related to the manufactured convergence study.
PartialDifferentialEquation
Possible Partial Differential Equations to solve.
Test factory, that will create the correct test with the right template parameters.
Definition: tests.h:71
Files for the baseline physics.
Definition: ADTypes.hpp:10
static std::unique_ptr< TestsBase > create_test(const Parameters::AllParameters *const parameters_input, dealii::ParameterHandler &parameter_handler_input)
Recursive factory that will create TestBase<int dim, int nstate>
Definition: tests.cpp:351
TestsBase()=delete
Constructor. Deleted the default constructor since it should not be used.
ModelType
Types of models available.
unsigned int dimension
Number of dimensions. Note that it has to match the executable PHiLiP_xD.
Main parameter class that contains the various other sub-parameter classes.
SubGridScaleModel SGS_model_type
Store the SubGridScale (SGS) model type.
ManufacturedConvergenceStudyParam manufactured_convergence_study_param
Contains parameters for manufactured convergence study.
const Parameters::AllParameters *const all_parameters
Pointer to all parameters.
Definition: tests.h:20
DissipativeNumericalFlux
Possible dissipative numerical flux types.
TestType test_type
Store selected TestType from the input file.
std::string get_conv_num_flux_string(const Parameters::AllParameters *const param) const
Returns a string describing which convective numerical flux is being used.
Definition: tests.cpp:132
ReynoldsAveragedNavierStokesModel
Types of Reynolds-averaged Navier-Stokes (RANS) models that can be used.
ConvectiveNumericalFlux
Possible convective numerical flux types.
MeshType
Mesh type to be used in defining the triangulation.
std::string get_pde_string(const Parameters::AllParameters *const param) const
Returns a string describing which PDE is being used.
Definition: tests.cpp:82
ReynoldsAveragedNavierStokesModel RANS_model_type
Store the Reynolds-averaged Navier-Stokes (RANS) model type.
ConvectiveNumericalFlux conv_num_flux_type
Store convective flux type.
std::string get_diss_num_flux_string(const Parameters::AllParameters *const param) const
Returns a string describing which dissipative numerical flux is being used.
Definition: tests.cpp:153
int nstate
Number of state variables. Will depend on PDE.
std::string get_manufactured_solution_string(const Parameters::AllParameters *const param) const
Returns a string describing which manufactured solution is being used.
Definition: tests.cpp:163
DissipativeNumericalFlux diss_num_flux_type
Store diffusive flux type.
static std::unique_ptr< TestsBase > select_mesh(const Parameters::AllParameters *const parameters_input, dealii::ParameterHandler &parameter_handler_input)
selects the mesh type to be used in the test
Definition: tests.cpp:203
SubGridScaleModel
Types of sub-grid scale (SGS) models that can be used.
double grid_progression
Multiplies the last grid size by this amount.
ModelType model_type
Store the model type.
static std::unique_ptr< TestsBase > select_test(const Parameters::AllParameters *const parameters_input, dealii::ParameterHandler &parameter_handler_input)
Selects the actual test such as grid convergence, numerical flux conversation, etc.
Definition: tests.cpp:233
std::vector< int > get_number_1d_cells(const int ngrids) const
Evaluates the number of cells to generate the grids for 1D grid based on input file.
Definition: tests.cpp:71
PhysicsModelParam physics_model_param
Contains parameters for Physics Model.
MeshType mesh_type
Store selected MeshType from the input file.