[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
euler_naca0012.cpp
1 #include <stdlib.h> /* srand, rand */
2 #include <iostream>
3 #include <deal.II/grid/grid_refinement.h>
4 #include "physics/manufactured_solution.h"
5 #include "euler_naca0012.hpp"
6 #include "flow_solver/flow_solver_factory.h"
7 
8 
9 namespace PHiLiP {
10 namespace Tests {
11 
12 template <int dim, int nstate>
14  const dealii::ParameterHandler &parameter_handler_input)
15  :
16  TestsBase::TestsBase(parameters_input)
17  , parameter_handler(parameter_handler_input)
18 {}
19 
20 template<int dim, int nstate>
22 ::run_test () const
23 {
25 
26  const unsigned int p_start = param.manufactured_convergence_study_param.degree_start;
27  const unsigned int p_end = param.manufactured_convergence_study_param.degree_end;
28  const unsigned int n_grids_input = param.manufactured_convergence_study_param.number_of_grids;
29 
30  for (unsigned int poly_degree = p_start; poly_degree <= p_end; ++poly_degree) {
31  for (unsigned int igrid=0; igrid<n_grids_input; ++igrid) {
32  param.flow_solver_param.poly_degree = poly_degree;
35  std::unique_ptr<FlowSolver::FlowSolver<dim,nstate>> flow_solver = FlowSolver::FlowSolverFactory<dim,nstate>::select_flow_case(&param, parameter_handler);
36  flow_solver->run();
37  }
38  }
39  return 0;
40 }
41 
42 
43 #if PHILIP_DIM==2
45 #endif
46 
47 } // Tests namespace
48 } // PHiLiP namespace
49 
50 
unsigned int degree_start
First polynomial degree to start the loop. If diffusion, must be at least 1.
FlowSolverParam flow_solver_param
Contains the parameters for simulation cases (flow solver test)
Files for the baseline physics.
Definition: ADTypes.hpp:10
Performs grid convergence for various polynomial degrees.
unsigned int poly_degree
Polynomial order (P) of the basis functions for DG.
Main parameter class that contains the various other sub-parameter classes.
ManufacturedConvergenceStudyParam manufactured_convergence_study_param
Contains parameters for manufactured convergence study.
static std::unique_ptr< FlowSolver< dim, nstate > > select_flow_case(const Parameters::AllParameters *const parameters_input, const dealii::ParameterHandler &parameter_handler_input)
Factory to return the correct flow solver given input file.
const Parameters::AllParameters *const all_parameters
Pointer to all parameters.
Definition: tests.h:20
int number_of_mesh_refinements
Number of refinements for naca0012 and Gaussian bump based cases.
const dealii::ParameterHandler & parameter_handler
Parameter handler for storing the .prm file being ran.
unsigned int max_poly_degree_for_adaptation
Maximum polynomial order of the DG basis functions for adaptation.
Base class of all the tests.
Definition: tests.h:17
int run_test() const
Grid convergence on Euler Gaussian Bump.
EulerNACA0012()=delete
Constructor. Deleted the default constructor since it should not be used.