[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
all_parameters.h
1 #ifndef __ALL_PARAMETERS_H__
2 #define __ALL_PARAMETERS_H__
3 
4 #include <deal.II/base/conditional_ostream.h>
5 #include <deal.II/base/parameter_handler.h>
6 #include <deal.II/base/tensor.h>
7 
8 #include "parameters.h"
9 #include "parameters/parameters_ode_solver.h"
10 #include "parameters/parameters_linear_solver.h"
11 #include "parameters/parameters_manufactured_convergence_study.h"
12 
13 #include "parameters/parameters_euler.h"
14 #include "parameters/parameters_navier_stokes.h"
15 #include "parameters/parameters_physics_model.h"
16 
17 #include "parameters/parameters_reduced_order.h"
18 #include "parameters/parameters_hyper_reduction.h"
19 #include "parameters/parameters_burgers.h"
20 #include "parameters/parameters_grid_refinement_study.h"
21 #include "parameters/parameters_grid_refinement.h"
22 #include "parameters/parameters_artificial_dissipation.h"
23 #include "parameters/parameters_limiter.h"
24 #include "parameters/parameters_flow_solver.h"
25 #include "parameters/parameters_mesh_adaptation.h"
26 #include "parameters/parameters_functional.h"
27 #include "parameters/parameters_time_refinement_study.h"
28 
29 namespace PHiLiP {
30 namespace Parameters {
31 
34 {
35 public:
37  AllParameters();
38 
71 
73  unsigned int dimension;
74 
76  enum RunType {
77  integration_test,
78  flow_simulation
79  };
81 
83  enum MeshType {
84  default_triangulation,
85  triangulation,
86  parallel_shared_triangulation,
87  parallel_distributed_triangulation,
88  };
91 
93 
96 
99 
101  enum FluxNodes { GL, GLL };
104 
107 
110 
112  enum TwoPointNumericalFlux { KG, IR, CH, Ra };
115 
118 
121 
124 
126 
129 
132 
135 
138 
140  //The default ESFR scheme is the Nonlinearly Stable FR where the volume is also reconstructed
142 
145 
148 
151 
154 
157 
159  std::string energy_file;
160 
162  int nstate;
163 
165  enum TestType {
166  run_control,
167  grid_refinement_study,
168  advection_limiter,
169  burgers_limiter,
170  burgers_energy_stability,
171  diffusion_exact_adjoint,
172  euler_gaussian_bump,
173  euler_gaussian_bump_enthalpy,
174  euler_gaussian_bump_adjoint,
175  euler_cylinder,
176  euler_cylinder_adjoint,
177  euler_vortex,
178  euler_entropy_waves,
179  euler_split_taylor_green,
180  taylor_green_scaling,
181  burgers_split_form,
182  optimization_inverse_manufactured,
183  euler_bump_optimization,
184  euler_naca_optimization,
185  shock_1d,
186  euler_naca0012,
187  reduced_order,
188  unsteady_reduced_order,
189  convection_diffusion_periodicity,
190  POD_adaptation,
191  POD_adaptive_sampling_run,
192  adaptive_sampling_testing,
193  finite_difference_sensitivity,
194  advection_periodicity,
195  dual_weighted_residual_mesh_adaptation,
196  anisotropic_mesh_adaptation,
197  taylor_green_vortex_energy_check,
198  taylor_green_vortex_restart_check,
199  time_refinement_study,
200  time_refinement_study_reference,
201  rrk_numerical_entropy_conservation_check,
202  euler_entropy_conserving_split_forms_check,
203  h_refinement_study_isentropic_vortex,
204  khi_robustness,
205  naca0012_unsteady_check_quick,
206  homogeneous_isotropic_turbulence_initialization_check,
207  build_NNLS_problem,
208  hyper_reduction_comparison,
209  hyper_adaptive_sampling_run,
210  hyper_reduction_post_sampling,
211  ROM_error_post_sampling,
212  HROM_error_post_sampling,
213  hyper_adaptive_sampling_new_error,
214  halton_sampling_run,
215  low_density
216  };
219 
222  advection,
223  diffusion,
224  convection_diffusion,
225  advection_vector,
226  burgers_inviscid,
227  burgers_viscous,
228  burgers_rewienski,
229  euler,
230  mhd,
231  navier_stokes,
232  physics_model,
233  };
236 
238  enum ModelType {
239  large_eddy_simulation,
240  reynolds_averaged_navier_stokes,
241  };
244 
247  manufactured_dirichlet,
248  manufactured_neumann,
249  manufactured_inout_flow,
250  };
251 
253  enum SourceTerm {
254  zero,
255  manufactured,
256  };
257 
260  lax_friedrichs,
261  roe,
262  l2roe,
263  central_flux,
264  two_point_flux,
265  two_point_flux_with_lax_friedrichs_dissipation,
266  two_point_flux_with_roe_dissipation,
267  two_point_flux_with_l2roe_dissipation
268  };
271 
273  enum DissipativeNumericalFlux { symm_internal_penalty, bassi_rebay_2, central_visc_flux };
276 
278  enum Flux_Reconstruction {cDG, cSD, cHU, cNegative, cNegative2, cPlus, c10Thousand, cHULumped};
281 
283  enum Flux_Reconstruction_Aux {kDG, kSD, kHU, kNegative, kNegative2, kPlus, k10Thousand};
286 
288  enum NonPhysicalBehaviorEnum {return_big_number, abort_run, print_warning};
291 
294 
297 
300 
303 
306 
308  enum RenumberDofsType { CuthillMckee };
311 
315 
317 
320  static void declare_parameters (dealii::ParameterHandler &prm);
321 
323 
326  void parse_parameters (dealii::ParameterHandler &prm);
327 
328  //FunctionParser<dim> initial_conditions;
329  //BoundaryConditions boundary_conditions[max_n_boundaries];
330 protected:
331  dealii::ConditionalOStream pcout;
332 
333 };
334 
335 } // Parameters namespace
336 } // PHiLiP namespace
337 
338 #endif
339 
TwoPointNumericalFlux two_point_num_flux_type
Store selected TwoPointNumericalFlux from the input file.
PartialDifferentialEquation pde_type
Store the PDE type to be solved.
bool use_curvilinear_grid
Flag to use curvilinear grid.
LimiterParam limiter_param
Contains parameters for limiter.
Parameters related to the linear solver.
TimeRefinementStudyParam time_refinement_study_param
Contains the parameters for time refinement study.
NonPhysicalBehaviorEnum non_physical_behavior_type
Specify behavior on nonphysical results.
FlowSolverParam flow_solver_param
Contains the parameters for simulation cases (flow solver test)
TestType
Possible integration tests to run.
Parameters related to the manufactured convergence study.
bool output_face_results_vtk
Flag for outputting the surface solution vtk files.
Parameters related to time refinement studies.
bool use_energy
Flag to use an energy monotonicity test.
PartialDifferentialEquation
Possible Partial Differential Equations to solve.
Files for the baseline physics.
Definition: ADTypes.hpp:10
GridRefinementStudyParam grid_refinement_study_param
Contains the parameters for grid refinement study.
bool use_weak_form
Flag to use weak or strong form of DG.
dealii::ConditionalOStream pcout
Parallel std::cout that only outputs on mpi_rank==0.
bool use_invariant_curl_form
Flag to use invariant curl form for metric cofactor operator.
Parameters related to reduced-order model.
ReducedOrderModelParam reduced_order_param
Contains parameters for the Reduced-Order model.
Flux_Reconstruction
Type of correction in Flux Reconstruction.
Flux_Reconstruction_Aux flux_reconstruction_aux_type
Store flux reconstruction type for the auxiliary variables.
EulerParam euler_param
Contains parameters for the Euler equations non-dimensionalization.
ModelType
Types of models available.
Parameters related to the limiter.
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.
RunType run_type
Selected RunType from the input file.
MeshAdaptationParam mesh_adaptation_param
Constains parameters for mesh adaptation.
std::string energy_file
Energy file.
NonPhysicalBehaviorEnum
Enum of nonphysical behavior.
ManufacturedConvergenceStudyParam manufactured_convergence_study_param
Contains parameters for manufactured convergence study.
bool enable_higher_order_vtk_output
Enable writing of higher-order vtk results.
Flux_Reconstruction flux_reconstruction_type
Store flux reconstruction type.
ODESolverParam ode_solver_param
Contains parameters for ODE solver.
Parameters related to reduced-order model.
bool use_classical_FR
Flag to use a Classical ESFR scheme where only the surface is reconstructed.
TwoPointNumericalFlux
Two point numerical flux type for split form.
NavierStokesParam navier_stokes_param
Contains parameters for the Navier-Stokes equations non-dimensionalization.
DissipativeNumericalFlux
Possible dissipative numerical flux types.
LinearSolverParam linear_solver_param
Contains parameters for linear solver.
bool do_renumber_dofs
Flag for renumbering DOFs.
Parameters related to the flow solver.
TestType test_type
Store selected TestType from the input file.
bool check_valid_metric_Jacobian
Flag to check if the metric Jacobian is valid when high-order grid is constructed.
bool use_L2_norm
Flag to use an L2 energy monotonicity test (for FR)
Parameters related to the Navier Stokes equations.
bool use_split_form
Flag to use split form.
Parameterse related to the functional object.
ConvectiveNumericalFlux
Possible convective numerical flux types.
MeshType
Mesh type to be used in defining the triangulation.
Parameters related to Physics Models.
Flux_Reconstruction_Aux
Type of correction in Flux Reconstruction for the auxiliary variables.
ConvectiveNumericalFlux conv_num_flux_type
Store convective flux type.
bool use_inverse_mass_on_the_fly
Flag to use inverse mass matrix on-the-fly for explicit solves.
RenumberDofsType
Renumber dofs type.
int nstate
Number of state variables. Will depend on PDE.
void parse_parameters(dealii::ParameterHandler &prm)
Retrieve parameters from dealii::ParameterHandler.
std::string solution_vtk_files_directory_name
Name of directory for writing solution vtk files.
Parameters related to collection of grid refinement runs.
FluxNodes flux_nodes_type
Store selected FluxNodes from the input file.
Parameters related to the ODE solver.
BoundaryType
Possible boundary types, NOT IMPLEMENTED YET.
FunctionalParam functional_param
Contains parameters for functional.
DissipativeNumericalFlux diss_num_flux_type
Store diffusive flux type.
Parameters related to the artificial dissipation.
static void declare_parameters(dealii::ParameterHandler &prm)
Declare parameters that can be set as inputs and set up the default options.
double sipg_penalty_factor
Scaling of Symmetric Interior Penalty term to ensure coercivity.
bool use_weight_adjusted_mass
Flag to use weight-adjusted Mass Matrix for curvilinear elements.
bool store_global_mass_matrix
Flag to store global mass matrix.
bool use_periodic_bc
Flag to use periodic BC.
SourceTerm
Possible source terms, NOT IMPLEMENTED YET.
ModelType model_type
Store the model type.
HyperReductionParam hyper_reduction_param
Contains parameters for Hyperreduction.
ArtificialDissipationParam artificial_dissipation_param
Contains parameters for artificial dissipation.
bool output_high_order_grid
Flag for outputting the high-order grid vtk files.
RenumberDofsType renumber_dofs_type
Store selected RenumberDofsType from the input file.
bool use_curvilinear_split_form
Flag to use curvilinear metric split form.
PhysicsModelParam physics_model_param
Contains parameters for Physics Model.
MeshType mesh_type
Store selected MeshType from the input file.
BurgersParam burgers_param
Contains parameters for Burgers equation.
bool use_collocated_nodes
Flag for using collocated nodes; determined based on flux_nodes_type and overintegration input parame...
Parameters related to the linear solver.
int overintegration
Number of additional quadrature points to use.
bool store_residual_cpu_time
Flag to store the residual local processor cput time.
Parameters related to reduced-order model.