[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 "stability_fr_parameter_range.h"
48 #include "bound_preserving_limiter_tests.h"
49 #include "naca0012_unsteady_check_quick.h"
50 #include "build_NNLS_problem.h"
51 #include "hyper_reduction_comparison.h"
52 #include "hyper_adaptive_sampling_run.h"
53 #include "hyper_reduction_post_sampling.h"
54 #include "ROM_error_post_sampling.h"
55 #include "HROM_error_post_sampling.h"
56 #include "hyper_adaptive_sampling_new_error.h"
57 #include "halton_sampling_run.h"
58 
59 namespace PHiLiP {
60 namespace Tests {
61 
62 using AllParam = Parameters::AllParameters;
63 
64 TestsBase::TestsBase(Parameters::AllParameters const *const parameters_input)
65  : all_parameters(parameters_input)
66  , mpi_communicator(MPI_COMM_WORLD)
67  , mpi_rank(dealii::Utilities::MPI::this_mpi_process(MPI_COMM_WORLD))
68  , n_mpi(dealii::Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD))
69  , pcout(std::cout, mpi_rank==0)
70 {}
71 
72 std::vector<int> TestsBase::get_number_1d_cells(const int n_grids) const
73 {
74  std::vector<int> n_1d_cells(n_grids);
76  n_1d_cells[0] = param.initial_grid_size;
77  for (int igrid=1;igrid<n_grids;++igrid) {
78  n_1d_cells[igrid] = static_cast<int>(n_1d_cells[igrid-1]*param.grid_progression) + param.grid_progression_add;
79  }
80  return n_1d_cells;
81 }
82 
83 std::string TestsBase::get_pde_string(const Parameters::AllParameters *const param) const
84 {
89 
90  const PDE_enum pde_type = param->pde_type;
91  std::string pde_string;
92  if (pde_type == PDE_enum::advection) {pde_string = "advection";}
93  if (pde_type == PDE_enum::advection_vector) {pde_string = "advection_vector";}
94  if (pde_type == PDE_enum::diffusion) {pde_string = "diffusion";}
95  if (pde_type == PDE_enum::convection_diffusion) {pde_string = "convection_diffusion";}
96  if (pde_type == PDE_enum::burgers_inviscid) {pde_string = "burgers_inviscid";}
97  if (pde_type == PDE_enum::burgers_viscous) {pde_string = "burgers_viscous";}
98  if (pde_type == PDE_enum::burgers_rewienski) {pde_string = "burgers_rewienski";}
99  if (pde_type == PDE_enum::euler) {pde_string = "euler";}
100  if (pde_type == PDE_enum::navier_stokes) {pde_string = "navier_stokes";}
101  if (pde_type == PDE_enum::physics_model) {
102  pde_string = "physics_model";
103  // add the model name + sub model name (if applicable)
104  const Model_enum model = param->model_type;
105  std::string model_string = "WARNING: invalid model";
106  if(model == Model_enum::large_eddy_simulation) {
107  // assign model string
108  model_string = "large_eddy_simulation";
109  // sub-grid scale (SGS)
110  const SGSModel_enum sgs_model = param->physics_model_param.SGS_model_type;
111  std::string sgs_model_string = "WARNING: invalid SGS model";
112  // assign SGS model string
113  if (sgs_model==SGSModel_enum::smagorinsky) sgs_model_string = "smagorinsky";
114  else if(sgs_model==SGSModel_enum::wall_adaptive_local_eddy_viscosity) sgs_model_string = "wall_adaptive_local_eddy_viscosity";
115  else if(sgs_model==SGSModel_enum::vreman) sgs_model_string = "vreman";
116  pde_string += std::string(" (Model: ") + model_string + std::string(", SGS Model: ") + sgs_model_string + std::string(")");
117  }
118  else if(model == Model_enum::reynolds_averaged_navier_stokes) {
119  // assign model string
120  model_string = "reynolds_averaged_navier_stokes";
121  // reynolds-averaged navier-stokes (RANS)
122  const RANSModel_enum rans_model = param->physics_model_param.RANS_model_type;
123  std::string rans_model_string = "WARNING: invalid RANS model";
124  // assign RANS model string
125  if (rans_model==RANSModel_enum::SA_negative) rans_model_string = "SA_negative";
126  pde_string += std::string(" (Model: ") + model_string + std::string(", RANS Model: ") + rans_model_string + std::string(")");
127  }
128  if(pde_string == "physics_model") pde_string += std::string(" (Model: ") + model_string + std::string(")");
129  }
130  return pde_string;
131 }
132 
134 {
136  const CNF_enum CNF_type = param->conv_num_flux_type;
137  std::string conv_num_flux_string;
138  if (CNF_type == CNF_enum::lax_friedrichs) {conv_num_flux_string = "lax_friedrichs";}
139  if (CNF_type == CNF_enum::roe) {conv_num_flux_string = "roe";}
140  if (CNF_type == CNF_enum::l2roe) {conv_num_flux_string = "l2roe";}
141  if (CNF_type == CNF_enum::central_flux) {conv_num_flux_string = "central_flux";}
142  if (CNF_type == CNF_enum::two_point_flux) {conv_num_flux_string = "two_point_flux";}
143  if (CNF_type == CNF_enum::two_point_flux_with_lax_friedrichs_dissipation) {
144  conv_num_flux_string = "two_point_flux_with_lax_friedrichs_dissipation";
145  } if (CNF_type == CNF_enum::two_point_flux_with_roe_dissipation) {
146  conv_num_flux_string = "two_point_flux_with_roe_dissipation";
147  } if (CNF_type == CNF_enum::two_point_flux_with_l2roe_dissipation) {
148  conv_num_flux_string = "two_point_flux_with_l2roe_dissipation";
149  }
150 
151  return conv_num_flux_string;
152 }
153 
155 {
157  const DNF_enum DNF_type = param->diss_num_flux_type;
158  std::string diss_num_flux_string;
159  if (DNF_type == DNF_enum::symm_internal_penalty) {diss_num_flux_string = "symm_internal_penalty";}
160  if (DNF_type == DNF_enum::bassi_rebay_2) {diss_num_flux_string = "bassi_rebay_2";}
161  return diss_num_flux_string;
162 }
163 
165 {
167  ManParam manu_grid_conv_param = param->manufactured_convergence_study_param;
169  const ManufacturedSolutionEnum MS_type = manu_grid_conv_param.manufactured_solution_param.manufactured_solution_type;
170  std::string manufactured_solution_string;
171  if (MS_type == ManufacturedSolutionEnum::sine_solution) {manufactured_solution_string = "sine_solution";}
172  if (MS_type == ManufacturedSolutionEnum::cosine_solution) {manufactured_solution_string = "cosine_solution";}
173  if (MS_type == ManufacturedSolutionEnum::additive_solution) {manufactured_solution_string = "additive_solution";}
174  if (MS_type == ManufacturedSolutionEnum::exp_solution) {manufactured_solution_string = "exp_solution";}
175  if (MS_type == ManufacturedSolutionEnum::poly_solution) {manufactured_solution_string = "poly_solution";}
176  if (MS_type == ManufacturedSolutionEnum::even_poly_solution) {manufactured_solution_string = "even_poly_solution";}
177  if (MS_type == ManufacturedSolutionEnum::atan_solution) {manufactured_solution_string = "atan_solution";}
178  if (MS_type == ManufacturedSolutionEnum::boundary_layer_solution) {manufactured_solution_string = "boundary_layer_solution";}
179  if (MS_type == ManufacturedSolutionEnum::s_shock_solution) {manufactured_solution_string = "s_shock_solution";}
180  if (MS_type == ManufacturedSolutionEnum::quadratic_solution) {manufactured_solution_string = "quadratic_solution";}
181  if (MS_type == ManufacturedSolutionEnum::example_solution) {manufactured_solution_string = "example_solution";}
182  if (MS_type == ManufacturedSolutionEnum::navah_solution_1) {manufactured_solution_string = "navah_solution_1";}
183  if (MS_type == ManufacturedSolutionEnum::navah_solution_2) {manufactured_solution_string = "navah_solution_2";}
184  if (MS_type == ManufacturedSolutionEnum::navah_solution_3) {manufactured_solution_string = "navah_solution_3";}
185  if (MS_type == ManufacturedSolutionEnum::navah_solution_4) {manufactured_solution_string = "navah_solution_4";}
186  if (MS_type == ManufacturedSolutionEnum::navah_solution_5) {manufactured_solution_string = "navah_solution_5";}
187  return manufactured_solution_string;
188 }
189 
190 //template<int dim, int nstate>
191 // void TestsBase::globally_refine_and_interpolate(DGBase<dim, double> &dg) const
192 //{
193 // dealii::LinearAlgebra::distributed::Vector<double> old_solution(dg->solution);
194 // dealii::parallel::distributed::SolutionTransfer<dim, dealii::LinearAlgebra::distributed::Vector<double>, dealii::hp::DoFHandler<dim>> solution_transfer(dg->dof_handler);
195 // solution_transfer.prepare_for_coarsening_and_refinement(old_solution);
196 // grid.refine_global (1);
197 // dg->allocate_system ();
198 // solution_transfer.interpolate(old_solution, dg->solution);
199 // solution_transfer.clear();
200 //}
201 
202 template<int dim, int nstate, typename MeshType>
203 std::unique_ptr< TestsBase > TestsFactory<dim,nstate,MeshType>
204 ::select_mesh(const AllParam *const parameters_input,
205  dealii::ParameterHandler &parameter_handler_input) {
206  using Mesh_enum = AllParam::MeshType;
207  Mesh_enum mesh_type = parameters_input->mesh_type;
208 
209  if(mesh_type == Mesh_enum::default_triangulation) {
210  #if PHILIP_DIM == 1
211  return TestsFactory<dim,nstate,dealii::Triangulation<dim>>::select_test(parameters_input,parameter_handler_input);
212  #else
213  return TestsFactory<dim,nstate,dealii::parallel::distributed::Triangulation<dim>>::select_test(parameters_input,parameter_handler_input);
214  #endif
215  } else if(mesh_type == Mesh_enum::triangulation) {
216  return TestsFactory<dim,nstate,dealii::Triangulation<dim>>::select_test(parameters_input,parameter_handler_input);
217  } else if(mesh_type == Mesh_enum::parallel_shared_triangulation) {
218  return TestsFactory<dim,nstate,dealii::parallel::shared::Triangulation<dim>>::select_test(parameters_input,parameter_handler_input);
219  } else if(mesh_type == Mesh_enum::parallel_distributed_triangulation) {
220  #if PHILIP_DIM == 1
221  std::cout << "dealii::parallel::distributed::Triangulation is unavailible in 1D." << std::endl;
222  #else
223  return TestsFactory<dim,nstate,dealii::parallel::distributed::Triangulation<dim>>::select_test(parameters_input,parameter_handler_input);
224  #endif
225  } else {
226  std::cout << "Invalid mesh type." << std::endl;
227  }
228 
229  return nullptr;
230 }
231 
232 template<int dim, int nstate, typename MeshType>
233 std::unique_ptr< TestsBase > TestsFactory<dim,nstate,MeshType>
234 ::select_test(const AllParam *const parameters_input,
235  dealii::ParameterHandler &parameter_handler_input) {
236  using Test_enum = AllParam::TestType;
237  const Test_enum test_type = parameters_input->test_type;
238 
239  // prevent warnings for when a create_FlowSolver is not being called (explicit and implicit cases)
240  if((test_type != Test_enum::finite_difference_sensitivity) &&
241  (test_type != Test_enum::taylor_green_vortex_energy_check) &&
242  (test_type != Test_enum::taylor_green_vortex_restart_check)) {
243  (void) parameter_handler_input;
244  } else if (!((dim==3 && nstate==dim+2) || (dim==1 && nstate==1))) {
245  (void) parameter_handler_input;
246  }
247 
248  if(test_type == Test_enum::run_control) { // TO DO: rename to grid_study
249  return std::make_unique<GridStudy<dim,nstate>>(parameters_input);
250  } else if(test_type == Test_enum::grid_refinement_study) {
251  return std::make_unique<GridRefinementStudy<dim,nstate,MeshType>>(parameters_input);
252  } else if(test_type == Test_enum::stability_fr_parameter_range) {
253  return std::make_unique<StabilityFRParametersRange<dim,nstate>>(parameters_input);
254  } else if(test_type == Test_enum::burgers_energy_stability) {
255  if constexpr (dim==1 && nstate==1) return std::make_unique<BurgersEnergyStability<dim,nstate>>(parameters_input);
256  } else if(test_type == Test_enum::diffusion_exact_adjoint) {
257  if constexpr (dim>=1 && nstate==1) return std::make_unique<DiffusionExactAdjoint<dim,nstate>>(parameters_input);
258  } else if (test_type == Test_enum::advection_periodicity){
259  if constexpr (nstate == 1) return std::make_unique<AdvectionPeriodic<dim,nstate>> (parameters_input);
260  } else if (test_type == Test_enum::convection_diffusion_periodicity){
261  if constexpr (nstate == 1) return std::make_unique<ConvectionDiffusionPeriodic<dim,nstate>> (parameters_input);
262  } else if(test_type == Test_enum::euler_gaussian_bump) {
263  if constexpr (dim==2 && nstate==dim+2) return std::make_unique<EulerGaussianBump<dim,nstate>>(parameters_input,parameter_handler_input);
264  } else if(test_type == Test_enum::euler_gaussian_bump_enthalpy) {
265  if constexpr (dim==2 && nstate==dim+2) return std::make_unique<EulerGaussianBumpEnthalpyCheck<dim,nstate>>(parameters_input, parameter_handler_input);
266  //} else if(test_type == Test_enum::euler_gaussian_bump_adjoint){
267  // if constexpr (dim==2 && nstate==dim+2) return std::make_unique<EulerGaussianBumpAdjoint<dim,nstate>>(parameters_input);
268  } else if(test_type == Test_enum::euler_cylinder) {
269  if constexpr (dim==2 && nstate==dim+2) return std::make_unique<EulerCylinder<dim,nstate>>(parameters_input);
270  } else if(test_type == Test_enum::euler_cylinder_adjoint) {
271  if constexpr (dim==2 && nstate==dim+2) return std::make_unique<EulerCylinderAdjoint<dim,nstate>>(parameters_input);
272  } else if(test_type == Test_enum::euler_vortex) {
273  if constexpr (dim==2 && nstate==dim+2) return std::make_unique<EulerVortex<dim,nstate>>(parameters_input);
274  } else if(test_type == Test_enum::euler_entropy_waves) {
275  if constexpr (dim>=2 && nstate==PHILIP_DIM+2) return std::make_unique<EulerEntropyWaves<dim,nstate>>(parameters_input);
276  } else if(test_type == Test_enum::euler_split_taylor_green) {
277  if constexpr (dim==3 && nstate == dim+2) return std::make_unique<EulerTaylorGreen<dim,nstate>>(parameters_input);
278  } else if(test_type == Test_enum::taylor_green_scaling) {
279  if constexpr (dim==3 && nstate == dim+2) return std::make_unique<EulerTaylorGreenScaling<dim,nstate>>(parameters_input);
280  } else if(test_type == Test_enum::optimization_inverse_manufactured) {
281  return std::make_unique<OptimizationInverseManufactured<dim,nstate>>(parameters_input);
282  } else if(test_type == Test_enum::euler_bump_optimization) {
283  if constexpr (dim==2 && nstate==dim+2) return std::make_unique<EulerBumpOptimization<dim,nstate>>(parameters_input);
284  } else if(test_type == Test_enum::euler_naca_optimization) {
285  if constexpr (dim==2 && nstate==dim+2) return std::make_unique<EulerNACAOptimization<dim,nstate>>(parameters_input);
286  } else if(test_type == Test_enum::shock_1d) {
287  if constexpr (dim==1 && nstate==1) return std::make_unique<Shock1D<dim,nstate>>(parameters_input);
288  } else if(test_type == Test_enum::reduced_order) {
289  if constexpr ((dim==2 && nstate==dim+2) || (dim==1 && nstate==1)) return std::make_unique<ReducedOrder<dim,nstate>>(parameters_input, parameter_handler_input);
290  } else if(test_type == Test_enum::unsteady_reduced_order) {
291  if constexpr (dim==2 && nstate==dim+2) return std::make_unique<UnsteadyReducedOrder<dim,nstate>>(parameters_input, parameter_handler_input);
292  } else if(test_type == Test_enum::POD_adaptive_sampling_run) {
293  if constexpr ((dim==2 && nstate==dim+2) || (dim==1 && nstate==1)) return std::make_unique<AdaptiveSamplingRun<dim,nstate>>(parameters_input,parameter_handler_input);
294  } else if(test_type == Test_enum::adaptive_sampling_testing) {
295  if constexpr ((dim==2 && nstate==dim+2) || (dim==1 && nstate==1)) return std::make_unique<AdaptiveSamplingTesting<dim,nstate>>(parameters_input,parameter_handler_input);
296  } else if(test_type == Test_enum::euler_naca0012) {
297  if constexpr (dim==2 && nstate==dim+2) return std::make_unique<EulerNACA0012<dim,nstate>>(parameters_input,parameter_handler_input);
298  } else if(test_type == Test_enum::dual_weighted_residual_mesh_adaptation) {
299  if constexpr (dim==2 && nstate==1) return std::make_unique<DualWeightedResidualMeshAdaptation<dim, nstate>>(parameters_input,parameter_handler_input);
300  } else if(test_type == Test_enum::anisotropic_mesh_adaptation) {
301  if constexpr( (dim==2 && nstate==1) || (dim==2 && nstate==dim+2)) return std::make_unique<AnisotropicMeshAdaptationCases<dim, nstate>>(parameters_input,parameter_handler_input);
302  } else if(test_type == Test_enum::taylor_green_vortex_energy_check) {
303  if constexpr (dim==3 && nstate==dim+2) return std::make_unique<TaylorGreenVortexEnergyCheck<dim,nstate>>(parameters_input,parameter_handler_input);
304  } else if(test_type == Test_enum::taylor_green_vortex_restart_check) {
305  if constexpr (dim==3 && nstate==dim+2) return std::make_unique<TaylorGreenVortexRestartCheck<dim,nstate>>(parameters_input,parameter_handler_input);
306  } else if(test_type == Test_enum::homogeneous_isotropic_turbulence_initialization_check){
307  if constexpr (dim==3 && nstate==dim+2) return std::make_unique<HomogeneousIsotropicTurbulenceInitializationCheck<dim,nstate>>(parameters_input,parameter_handler_input);
308  } else if(test_type == Test_enum::time_refinement_study) {
309  if constexpr (dim==1 && nstate==1) return std::make_unique<TimeRefinementStudy<dim, nstate>>(parameters_input, parameter_handler_input);
310  } else if(test_type == Test_enum::h_refinement_study_isentropic_vortex) {
311  if constexpr (dim+2==nstate && dim!=1) return std::make_unique<HRefinementStudyIsentropicVortex<dim, nstate>>(parameters_input, parameter_handler_input);
312  } else if(test_type == Test_enum::time_refinement_study_reference) {
313  if constexpr (dim==1 && nstate==1) return std::make_unique<TimeRefinementStudyReference<dim, nstate>>(parameters_input, parameter_handler_input);
314  } else if(test_type == Test_enum::rrk_numerical_entropy_conservation_check) {
315  if constexpr ((dim==1 && nstate==1) || (dim==3 && nstate==dim+2)) return std::make_unique<RRKNumericalEntropyConservationCheck<dim, nstate>>(parameters_input, parameter_handler_input);
316  } else if(test_type == Test_enum::euler_entropy_conserving_split_forms_check) {
317  if constexpr (dim==3 && nstate==dim+2) return std::make_unique<EulerSplitEntropyCheck<dim, nstate>>(parameters_input, parameter_handler_input);
318  } else if(test_type == Test_enum::khi_robustness) {
319  if constexpr (dim==2 && nstate==dim+2) return std::make_unique<KHIRobustness<dim, nstate>>(parameters_input, parameter_handler_input);
320  } else if(test_type == Test_enum::build_NNLS_problem) {
321  if constexpr (dim==1 && nstate==1) return std::make_unique<BuildNNLSProblem<dim,nstate>>(parameters_input, parameter_handler_input);
322  } else if(test_type == Test_enum::hyper_reduction_comparison) {
323  if constexpr (dim==1 && nstate==1) return std::make_unique<HyperReductionComparison<dim,nstate>>(parameters_input, parameter_handler_input);
324  } else if(test_type == Test_enum::hyper_adaptive_sampling_run) {
325  if constexpr ((dim==2 && nstate==dim+2) || (dim==1 && nstate==1)) return std::make_unique<HyperAdaptiveSamplingRun<dim,nstate>>(parameters_input, parameter_handler_input);
326  } else if(test_type == Test_enum::hyper_reduction_post_sampling) {
327  if constexpr ((dim==2 && nstate==dim+2) || (dim==1 && nstate==1)) return std::make_unique<HyperReductionPostSampling<dim,nstate>>(parameters_input, parameter_handler_input);
328  } else if(test_type == Test_enum::ROM_error_post_sampling) {
329  if constexpr ((dim==2 && nstate==dim+2) || (dim==1 && nstate==1)) return std::make_unique<ROMErrorPostSampling<dim,nstate>>(parameters_input, parameter_handler_input);
330  } else if(test_type == Test_enum::HROM_error_post_sampling) {
331  if constexpr ((dim==2 && nstate==dim+2) || (dim==1 && nstate==1)) return std::make_unique<HROMErrorPostSampling<dim,nstate>>(parameters_input, parameter_handler_input);
332  } else if(test_type == Test_enum::hyper_adaptive_sampling_new_error) {
333  if constexpr ((dim==2 && nstate==dim+2) || (dim==1 && nstate==1)) return std::make_unique<HyperAdaptiveSamplingNewError<dim,nstate>>(parameters_input, parameter_handler_input);
334  } else if(test_type == Test_enum::halton_sampling_run) {
335  if constexpr ((dim==2 && nstate==dim+2) || (dim==1 && nstate==1)) return std::make_unique<HaltonSamplingRun<dim,nstate>>(parameters_input, parameter_handler_input);
336  } else if (test_type == Test_enum::advection_limiter) {
337  if constexpr (nstate == 1 && dim < 3) return std::make_unique<BoundPreservingLimiterTests<dim, nstate>>(parameters_input, parameter_handler_input);
338  } else if (test_type == Test_enum::burgers_limiter) {
339  if constexpr (nstate == dim && dim < 3) return std::make_unique<BoundPreservingLimiterTests<dim, nstate>>(parameters_input, parameter_handler_input);
340  } else if(test_type == Test_enum::low_density) {
341  if constexpr (dim<3 && nstate==dim+2) return std::make_unique<BoundPreservingLimiterTests<dim, nstate>>(parameters_input, parameter_handler_input);
342  } else if(test_type == Test_enum::naca0012_unsteady_check_quick){
343  if constexpr (dim==2 && nstate==dim+2) return std::make_unique<NACA0012UnsteadyCheckQuick<dim, nstate>>(parameters_input, parameter_handler_input);
344  } else {
345  std::cout << "Invalid test. You probably forgot to add it to the list of tests in tests.cpp" << std::endl;
346  std::abort();
347  }
348 
349  return nullptr;
350 }
351 
352 template<int dim, int nstate, typename MeshType>
353 std::unique_ptr< TestsBase > TestsFactory<dim,nstate,MeshType>
354 ::create_test(AllParam const *const parameters_input,
355  dealii::ParameterHandler &parameter_handler_input)
356 {
357  // Recursive templating required because template parameters must be compile time constants
358  // As a results, this recursive template initializes all possible dimensions with all possible nstate
359  // without having 15 different if-else statements
360  if(dim == parameters_input->dimension)
361  {
362  // This template parameters dim and nstate match the runtime parameters
363  // then create the selected test with template parameters dim and nstate
364  // Otherwise, keep decreasing nstate and dim until it matches
365  if(nstate == parameters_input->nstate)
366  return TestsFactory<dim,nstate>::select_mesh(parameters_input,parameter_handler_input);
367  else if constexpr (nstate > 1)
368  return TestsFactory<dim,nstate-1>::create_test(parameters_input,parameter_handler_input);
369  else
370  return nullptr;
371  }
372  else if constexpr (dim > 1)
373  {
374  //return TestsFactory<dim-1,nstate>::create_test(parameters_input);
375  return nullptr;
376  }
377  else
378  {
379  return nullptr;
380  }
381 }
382 
383 // Will recursively create all the possible test sizes
384 //template class TestsFactory <PHILIP_DIM,1>;
385 //template class TestsFactory <PHILIP_DIM,2>;
386 //template class TestsFactory <PHILIP_DIM,3>;
387 //template class TestsFactory <PHILIP_DIM,4>;
388 //template class TestsFactory <PHILIP_DIM,5>;
389 
392 #if PHILIP_DIM!=1
394 #endif
395 
396 } // Tests namespace
397 } // 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:354
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:133
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:83
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:154
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:164
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:204
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:234
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:72
PhysicsModelParam physics_model_param
Contains parameters for Physics Model.
MeshType mesh_type
Store selected MeshType from the input file.