1 #include <deal.II/base/mpi.h>     2 #include <deal.II/base/utilities.h>     3 #include <deal.II/base/patterns.h>     5 #include "parameters/all_parameters.h"    12 namespace Parameters {
    31     , pcout(
std::cout, dealii::Utilities::MPI::this_mpi_process(MPI_COMM_WORLD)==0)
    36     const int mpi_rank = dealii::Utilities::MPI::this_mpi_process(MPI_COMM_WORLD);
    37     dealii::ConditionalOStream 
pcout(std::cout, mpi_rank==0);
    38     pcout << 
"Declaring inputs." << std::endl;
    39     prm.declare_entry(
"dimension", 
"-1",
    40                       dealii::Patterns::Integer(),
    41                       "Number of dimensions");
    43     prm.declare_entry(
"run_type", 
"integration_test",
    44                       dealii::Patterns::Selection(
    45                       " integration_test | "    47                       "Type of run (default is integration_test). "    48                       "Choices are  <integration_test | flow_simulation>.");
    50     prm.declare_entry(
"mesh_type", 
"default_triangulation",
    51                       dealii::Patterns::Selection(
    52                       " default_triangulation | "    54                       " parallel_shared_triangulation | "    55                       " parallel_distributed_triangulation"),
    56                       "Type of triangulation to be used."    57                       "Note: parralel_distributed_triangulation not availible int 1D."    58                       " <default_triangulation | "    60                       "  parallel_shared_triangulation |"    61                       "  parallel_distributed_triangulation>.");
    63     prm.declare_entry(
"overintegration", 
"0",
    64                       dealii::Patterns::Integer(),
    65                       "Number of extra quadrature points to use."    66                       "If overintegration=0, then we use n_quad = soln_degree + 1.");
    68     prm.declare_entry(
"use_weak_form", 
"true",
    69                       dealii::Patterns::Bool(),
    70                       "Use weak form by default. If false, use strong form.");
    72     prm.declare_entry(
"flux_nodes_type", 
"GL",
    73                       dealii::Patterns::Selection(
    75                       "Flux nodes type, default is GL for uncollocated. NOTE: Solution nodes are type GLL."    76                       "Choices are <GL | GLL>.");
    78     prm.declare_entry(
"use_split_form", 
"false",
    79                       dealii::Patterns::Bool(),
    80                       "Use original form by defualt. Otherwise, split the fluxes.");
    82     prm.declare_entry(
"two_point_num_flux_type", 
"KG",
    83                       dealii::Patterns::Selection(
    85                       "Two point flux type. "    86                       "Choices are <KG | IR | CH | Ra>.");
    88     prm.declare_entry(
"use_curvilinear_split_form", 
"false",
    89                       dealii::Patterns::Bool(),
    90                       "Use original form by defualt. Otherwise, split the curvilinear fluxes.");
    92     prm.declare_entry(
"store_residual_cpu_time", 
"false",
    93                       dealii::Patterns::Bool(),
    94                       "Do not store the residual local processor cpu time by default. Store the residual cpu time if true.");
    96     prm.declare_entry(
"use_weight_adjusted_mass", 
"false",
    97                       dealii::Patterns::Bool(),
    98                       "Use original form by defualt. Otherwise, use the weight adjusted low storage mass matrix for curvilinear.");
   100     prm.declare_entry(
"use_periodic_bc", 
"false",
   101                       dealii::Patterns::Bool(),
   102                       "Use other boundary conditions by default. Otherwise use periodic (for 1d burgers only");
   104     prm.declare_entry(
"use_curvilinear_grid", 
"false",
   105                       dealii::Patterns::Bool(),
   106                       "Use straight grid by default. Curvilinear is true. Only used in taylor_green_scaling test.");
   108     prm.declare_entry(
"use_energy", 
"false",
   109                       dealii::Patterns::Bool(),
   110                       "Not calculate energy by default. Otherwise, get energy per iteration.");
   112     prm.declare_entry(
"use_L2_norm", 
"false",
   113                       dealii::Patterns::Bool(),
   114                       "Not calculate L2 norm by default (M+K). Otherwise, get L2 norm per iteration.");
   116     prm.declare_entry(
"use_classical_FR", 
"false",
   117                       dealii::Patterns::Bool(),
   118                       "Not use Classical Flux Reconstruction by default. Otherwise, use Classical Flux Reconstruction.");
   120     prm.declare_entry(
"flux_reconstruction", 
"cDG",
   121                       dealii::Patterns::Selection(
   122                       "cDG | cSD | cHU | cNegative | cNegative2 | cPlus | c10Thousand | cHULumped | user_specified_value"),
   123                       "Flux Reconstruction. "   125                       " <cDG | cSD | cHU | cNegative | cNegative2 | cPlus | c10Thousand | cHULumped | user_specified_value>.");
   127     prm.declare_entry(
"FR_user_specified_correction_parameter_value", 
"0.0",
   128                       dealii::Patterns::Double(-dealii::Patterns::Double::max_double_value, dealii::Patterns::Double::max_double_value),
   129                       "User specified flux recontruction correction parameter value. "   130                       "Enter a 1D correction parameter. "   131                       "Internally, the input c value is divided by 2 to account for the basis and adjusted for the deal.ii reference element. "   132                       "Default value is 0.0. ");
   134     prm.declare_entry(
"flux_reconstruction_aux", 
"kDG",
   135                       dealii::Patterns::Selection(
   136                       "kDG | kSD | kHU | kNegative | kNegative2 | kPlus | k10Thousand"),
   137                       "Flux Reconstruction for Auxiliary Equation. "   138                       "Choices are <kDG | kSD | kHU | kNegative | kNegative2 | kPlus | k10Thousand>.");
   140     prm.declare_entry(
"sipg_penalty_factor", 
"1.0",
   141                       dealii::Patterns::Double(1.0,1e200),
   142                       "Scaling of Symmetric Interior Penalty term to ensure coercivity.");
   144     prm.declare_entry(
"use_invariant_curl_form", 
"false",
   145                       dealii::Patterns::Bool(),
   146                       "Use conservative curl form for metric cofactor by default. If true, then use invariant curl form.");
   148     prm.declare_entry(
"use_inverse_mass_on_the_fly", 
"false",
   149                       dealii::Patterns::Bool(),
   150                       "Build global mass inverse matrix and apply it. Otherwise, use inverse mass on-the-fly by default for explicit timestepping.");
   152     prm.declare_entry(
"check_valid_metric_Jacobian", 
"true",
   153                       dealii::Patterns::Bool(),
   154                       "Check validty of metric Jacobian when high-order grid is constructed by default. Do not check if false. Not checking is useful if the metric terms are built on the fly with operators, it reduces the memory cost for high polynomial grids. The metric Jacobian is never checked for strong form, regardless of the user input.");
   156     prm.declare_entry(
"energy_file", 
"energy_file",
   157                       dealii::Patterns::FileName(dealii::Patterns::FileName::FileType::input),
   158                       "Input file for energy test.");
   160     prm.declare_entry(
"test_type", 
"run_control",
   161                       dealii::Patterns::Selection(
   163                       " grid_refinement_study | "   164                       " stability_fr_parameter_range | "   165                       " advection_limiter | "   166                       " burgers_limiter | "   167                       " burgers_energy_stability | "   168                       " diffusion_exact_adjoint | "   169                       " optimization_inverse_manufactured | "   170                       " euler_gaussian_bump | "   171                       " euler_gaussian_bump_enthalpy | "   172                       " euler_gaussian_bump_adjoint | "   174                       " euler_cylinder_adjoint | "   176                       " euler_entropy_waves | "   177                       " euler_split_taylor_green | "   178                       " taylor_green_scaling | "   179                       " euler_bump_optimization | "   180                       " euler_naca_optimization | "   184                       " unsteady_reduced_order |"   185                       " convection_diffusion_periodicity |"   187                       " POD_adaptive_sampling_run | "   188                       " adaptive_sampling_testing | "   189                       " finite_difference_sensitivity | "   190                       " advection_periodicity | "   191                       " dual_weighted_residual_mesh_adaptation | "   192                       " anisotropic_mesh_adaptation | "   193                       " taylor_green_vortex_energy_check | "   194                       " taylor_green_vortex_restart_check | "   195                       " homogeneous_isotropic_turbulence_initialization_check | "   196                       " time_refinement_study | "   197                       " time_refinement_study_reference | "   198                       " rrk_numerical_entropy_conservation_check | "   199                       " euler_entropy_conserving_split_forms_check | "   200                       " h_refinement_study_isentropic_vortex | "   201                       " build_NNLS_problem |"   202                       " hyper_reduction_comparison |"   203                       " hyper_adaptive_sampling_run |"   204                       " hyper_reduction_post_sampling |"   205                       " ROM_error_post_sampling |"   206                       " HROM_error_post_sampling | "   207                       " hyper_adaptive_sampling_new_error |"   208                       " halton_sampling_run |"   209                       " naca0012_unsteady_check_quick | "   212                       "The type of test we want to solve. "   215                       "  grid_refinement_study | "   216                       "  stability_fr_parameter_range | "   217                       "  advection_limiter | "   218                       "  burgers_limiter | "   219                       "  burgers_energy_stability | "   220                       "  diffusion_exact_adjoint | "   221                       "  optimization_inverse_manufactured | "   222                       "  euler_gaussian_bump | "   223                       "  euler_gaussian_bump_enthalpy | "   224                       "  euler_gaussian_bump_adjoint | "   226                       "  euler_cylinder_adjoint | "   228                       "  euler_entropy_waves | "   229                       "  euler_split_taylor_green |"   230                       " taylor_green_scaling | "   231                       "  euler_bump_optimization | "   232                       "  euler_naca_optimization | "   235                       "  convection_diffusion_periodicity |"   237                       "  unsteady_reduced_order | "   239                       "  POD_adaptive_sampling_run | "   240                       "  adaptive_sampling_testing | "   241                       "  finite_difference_sensitivity | "   242                       "  advection_periodicity | "   243                       "  dual_weighted_residual_mesh_adaptation | "   244                       "  anisotropic_mesh_adaptation | "   245                       "  taylor_green_vortex_energy_check | "   246                       "  taylor_green_vortex_restart_check | "   247                       "  homogeneous_isotropic_turbulence_initialization_check | "   248                       "  time_refinement_study | "   249                       "  time_refinement_study_reference | "   250                       "  rrk_numerical_entropy_conservation_check | "   251                       "  euler_entropy_conserving_split_forms_check | "   252                       "  h_refinement_study_isentropic_vortex | "   253                       "  build_NNLS_problem |"   254                       "  hyper_reduction_comparison |"   255                       "  hyper_adaptive_sampling_run |"   256                       "  hyper_reduction_post_sampling |"   257                       "  ROM_error_post_sampling |"   258                       "  HROM_error_post_sampling | "   259                       "  hyper_adaptive_sampling_new_error |"   260                       "  halton_sampling_run |"   261                       "  naca0012_unsteady_check_quick | "   265     prm.declare_entry(
"pde_type", 
"advection",
   266                       dealii::Patterns::Selection(
   269                       " convection_diffusion | "   270                       " advection_vector | "   271                       " burgers_inviscid | "   272                       " burgers_viscous | "   273                       " burgers_rewienski | "   278                       "The PDE we want to solve. "   282                       "  convection_diffusion | "   283                       "  advection_vector | "   284                       "  burgers_inviscid | "   285                       "  burgers_viscous | "   286                       "  burgers_rewienski | "   292     prm.declare_entry(
"model_type", 
"large_eddy_simulation",
   293                       dealii::Patterns::Selection(
   294                       "large_eddy_simulation | reynolds_averaged_navier_stokes"),
   295                       "Enum of physics models "   296                       "(i.e. model equations and/or terms additional to Navier-Stokes or a chosen underlying baseline physics)."   298                       " <large_eddy_simulation | reynolds_averaged_navier_stokes>.");
   300     prm.declare_entry(
"conv_num_flux", 
"lax_friedrichs",
   301                       dealii::Patterns::Selection(
   307                       " two_point_flux_with_lax_friedrichs_dissipation | "   308                       " two_point_flux_with_roe_dissipation | "   309                       " two_point_flux_with_l2roe_dissipation"),
   310                       "Convective numerical flux. "   312                       " <lax_friedrichs | "   317                       " two_point_flux_with_lax_friedrichs_dissipation | "   318                       " two_point_flux_with_roe_dissipation | "   319                       " two_point_flux_with_l2roe_dissipation>.");
   321     prm.declare_entry(
"diss_num_flux", 
"symm_internal_penalty",
   322                       dealii::Patterns::Selection(
"symm_internal_penalty | bassi_rebay_2 | central_visc_flux"),
   323                       "Dissipative numerical flux. "   324                       "Choices are <symm_internal_penalty | bassi_rebay_2 | central_visc_flux>.");
   326     prm.declare_entry(
"non_physical_behavior", 
"return_big_number",
   327                       dealii::Patterns::Selection(
"return_big_number | abort_run | print_warning"),
   328                       "Behavior when a nonphysical result is detected in physics, "   329                       "For example negative density or NaN. "   330                       "return_big_number will set the quantity to BIG_NUMBER without any warnings "   331                       "abort_run will std::abort() "   332                       "print_warning will return BIG_NUMBER and print a warning to console. "   333                       "Choices are <return_big_number | abort_run | print_warning>.");
   335     prm.declare_entry(
"solution_vtk_files_directory_name", 
".",
   336                       dealii::Patterns::FileName(dealii::Patterns::FileName::FileType::input),
   337                       "Name of directory for writing solution vtk files. Current directory by default.");
   339     prm.declare_entry(
"output_high_order_grid", 
"false",
   340                       dealii::Patterns::Bool(),
   341                       "Outputs the high-order mesh vtu files. False by default");
   343     prm.declare_entry(
"enable_higher_order_vtk_output", 
"true",
   344                       dealii::Patterns::Bool(),
   345                       "Enable writing of higher-order vtk files. True by default; "   346                       "number of subdivisions is chosen according to the max of grid_degree and poly_degree.");
   348     prm.declare_entry(
"output_face_results_vtk", 
"false",
   349                       dealii::Patterns::Bool(),
   350                       "Outputs the surface solution vtk files. False by default");
   352     prm.declare_entry(
"do_renumber_dofs", 
"true",
   353                       dealii::Patterns::Bool(),
   354                       "Flag for renumbering DOFs using Cuthill-McKee renumbering. True by default. Set to false if doing 3D unsteady flow simulations.");
   356     prm.declare_entry(
"renumber_dofs_type", 
"CuthillMckee",
   357                       dealii::Patterns::Selection(
   359                       "Renumber the dof handler type. Currently the only choice is Cuthill-Mckee.");
   361     prm.declare_entry(
"matching_surface_jac_det_tolerance", 
"1.3e-11",
   362                       dealii::Patterns::Double(0, dealii::Patterns::Double::max_double_value),
   363                       "Tolerance for checking that the determinant of surface jacobians at element faces matches. "   364                       "Note: Currently only used in weak dg.");
   383     pcout << 
"Done declaring inputs." << std::endl;
   388     pcout << 
"Parsing main input..." << std::endl;
   390     dimension = prm.get_integer(
"dimension");
   392     const std::string run_type_string = prm.get(
"run_type");
   393     if      (run_type_string == 
"integration_test") { 
run_type = integration_test; }
   394     else if (run_type_string == 
"flow_simulation")  { 
run_type = flow_simulation; }
   396     const std::string mesh_type_string = prm.get(
"mesh_type");
   397     if      (mesh_type_string == 
"default_triangulation")              { 
mesh_type = default_triangulation; }
   398     else if (mesh_type_string == 
"triangulation")                      { 
mesh_type = triangulation; }
   399     else if (mesh_type_string == 
"parallel_shared_triangulation")      { 
mesh_type = parallel_shared_triangulation; }
   400     else if (mesh_type_string == 
"parallel_distributed_triangulation") { 
mesh_type = parallel_distributed_triangulation; }
   402 const std::string test_string = prm.get(
"test_type");
   403     if      (test_string == 
"run_control")                              { 
test_type = run_control; }
   404     else if (test_string == 
"grid_refinement_study")                    { 
test_type = grid_refinement_study; }
   405     else if (test_string == 
"stability_fr_parameter_range")             { 
test_type = stability_fr_parameter_range; }
   406     else if (test_string == 
"advection_limiter")                        { 
test_type = advection_limiter; }
   407     else if (test_string == 
"burgers_limiter")                          { 
test_type = burgers_limiter; }
   408     else if (test_string == 
"burgers_energy_stability")                 { 
test_type = burgers_energy_stability; }
   409     else if (test_string == 
"diffusion_exact_adjoint")                  { 
test_type = diffusion_exact_adjoint; }
   410     else if (test_string == 
"euler_gaussian_bump")                      { 
test_type = euler_gaussian_bump; }
   411     else if (test_string == 
"euler_gaussian_bump_enthalpy")             { 
test_type = euler_gaussian_bump_enthalpy; }
   412     else if (test_string == 
"euler_gaussian_bump_adjoint")              { 
test_type = euler_gaussian_bump_adjoint; }
   413     else if (test_string == 
"euler_cylinder")                           { 
test_type = euler_cylinder; }
   414     else if (test_string == 
"euler_cylinder_adjoint")                   { 
test_type = euler_cylinder_adjoint; }
   415     else if (test_string == 
"euler_vortex")                             { 
test_type = euler_vortex; }
   416     else if (test_string == 
"euler_entropy_waves")                      { 
test_type = euler_entropy_waves; }
   417     else if (test_string == 
"advection_periodicity")                    { 
test_type = advection_periodicity; }
   418     else if (test_string == 
"convection_diffusion_periodicity")         { 
test_type = convection_diffusion_periodicity; }
   419     else if (test_string == 
"euler_split_taylor_green")                 { 
test_type = euler_split_taylor_green; }
   420     else if (test_string == 
"taylor_green_scaling")                     { 
test_type = taylor_green_scaling; }
   421     else if (test_string == 
"euler_bump_optimization")                  { 
test_type = euler_bump_optimization; }
   422     else if (test_string == 
"euler_naca_optimization")                  { 
test_type = euler_naca_optimization; }
   423     else if (test_string == 
"shock_1d")                                 { 
test_type = shock_1d; }
   424     else if (test_string == 
"reduced_order")                            { 
test_type = reduced_order; }
   425     else if (test_string == 
"unsteady_reduced_order")                   { 
test_type = unsteady_reduced_order; }
   426     else if (test_string == 
"POD_adaptation")                           { 
test_type = POD_adaptation; }
   427     else if (test_string == 
"POD_adaptive_sampling_run")                { 
test_type = POD_adaptive_sampling_run; }
   428     else if (test_string == 
"adaptive_sampling_testing")                { 
test_type = adaptive_sampling_testing; }
   429     else if (test_string == 
"finite_difference_sensitivity")            { 
test_type = finite_difference_sensitivity; }
   430     else if (test_string == 
"euler_naca0012")                           { 
test_type = euler_naca0012; }
   431     else if (test_string == 
"optimization_inverse_manufactured")        { 
test_type = optimization_inverse_manufactured; }
   432     else if (test_string == 
"dual_weighted_residual_mesh_adaptation")   { 
test_type = dual_weighted_residual_mesh_adaptation; }
   433     else if (test_string == 
"anisotropic_mesh_adaptation")              { 
test_type = anisotropic_mesh_adaptation; }
   434     else if (test_string == 
"taylor_green_vortex_energy_check")         { 
test_type = taylor_green_vortex_energy_check; }
   435     else if (test_string == 
"taylor_green_vortex_restart_check")        { 
test_type = taylor_green_vortex_restart_check; }
   436     else if (test_string == 
"homogeneous_isotropic_turbulence_initialization_check")
   437                                                                         { 
test_type = homogeneous_isotropic_turbulence_initialization_check; }
   438     else if (test_string == 
"time_refinement_study")                    { 
test_type = time_refinement_study; }
   439     else if (test_string == 
"time_refinement_study_reference")          { 
test_type = time_refinement_study_reference; }
   440     else if (test_string == 
"h_refinement_study_isentropic_vortex")     { 
test_type = h_refinement_study_isentropic_vortex; }
   441     else if (test_string == 
"rrk_numerical_entropy_conservation_check") { 
test_type = rrk_numerical_entropy_conservation_check; }
   442     else if (test_string == 
"euler_entropy_conserving_split_forms_check") 
   443                                                                         { 
test_type = euler_entropy_conserving_split_forms_check; }
   444     else if (test_string == 
"h_refinement_study_isentropic_vortex")     { 
test_type = h_refinement_study_isentropic_vortex; }
   445     else if (test_string == 
"khi_robustness")                           { 
test_type = khi_robustness; }
   446     else if (test_string == 
"build_NNLS_problem")                       { 
test_type = build_NNLS_problem; }
   447     else if (test_string == 
"hyper_reduction_comparison")               { 
test_type = hyper_reduction_comparison; }
   448     else if (test_string == 
"hyper_adaptive_sampling_run")              { 
test_type = hyper_adaptive_sampling_run; }
   449     else if (test_string == 
"hyper_reduction_post_sampling")            { 
test_type = hyper_reduction_post_sampling; }
   450     else if (test_string == 
"ROM_error_post_sampling")                  { 
test_type = ROM_error_post_sampling; }
   451     else if (test_string == 
"HROM_error_post_sampling")                 { 
test_type = HROM_error_post_sampling; }
   452     else if (test_string == 
"hyper_adaptive_sampling_new_error")        { 
test_type = hyper_adaptive_sampling_new_error; }
   453     else if (test_string == 
"halton_sampling_run")                      { 
test_type = halton_sampling_run; }
   454     else if (test_string == 
"low_density")                              { 
test_type = low_density; }
   455     else if (test_string == 
"naca0012_unsteady_check_quick")            { 
test_type = naca0012_unsteady_check_quick; }
   461     const std::string flux_nodes_string = prm.get(
"flux_nodes_type");
   469     const std::string two_point_num_flux_string = prm.get(
"two_point_num_flux_type");
   493     const std::string conv_num_flux_string = prm.get(
"conv_num_flux");
   494     if (conv_num_flux_string == 
"lax_friedrichs")                                          { 
conv_num_flux_type = ConvectiveNumericalFlux::lax_friedrichs; }
   495     if (conv_num_flux_string == 
"roe")                                                     { 
conv_num_flux_type = ConvectiveNumericalFlux::roe; }
   496     if (conv_num_flux_string == 
"l2roe")                                                   { 
conv_num_flux_type = ConvectiveNumericalFlux::l2roe; }
   497     if (conv_num_flux_string == 
"central_flux")                                            { 
conv_num_flux_type = ConvectiveNumericalFlux::central_flux; }
   498     if (conv_num_flux_string == 
"two_point_flux")                                 { 
conv_num_flux_type = ConvectiveNumericalFlux::two_point_flux; }
   499     if (conv_num_flux_string == 
"two_point_flux_with_lax_friedrichs_dissipation") { 
conv_num_flux_type = ConvectiveNumericalFlux::two_point_flux_with_lax_friedrichs_dissipation; }
   500     if (conv_num_flux_string == 
"two_point_flux_with_roe_dissipation")            { 
conv_num_flux_type = ConvectiveNumericalFlux::two_point_flux_with_roe_dissipation; }
   501     if (conv_num_flux_string == 
"two_point_flux_with_l2roe_dissipation")          { 
conv_num_flux_type = ConvectiveNumericalFlux::two_point_flux_with_l2roe_dissipation; }
   503     const std::string diss_num_flux_string = prm.get(
"diss_num_flux");
   504     if (diss_num_flux_string == 
"symm_internal_penalty") { 
diss_num_flux_type = symm_internal_penalty; }
   505     if (diss_num_flux_string == 
"bassi_rebay_2") {
   509     if (diss_num_flux_string == 
"central_visc_flux") 
diss_num_flux_type = central_visc_flux;
   511     const std::string flux_reconstruction_string = prm.get(
"flux_reconstruction");
   520     if (flux_reconstruction_string == 
"user_specified_value") 
   525         pcout << 
"Warning: User-specified FR parameter has been set, but flux_reconstruction_type is " << std::endl
   526               << 
"not chosen as user_specified_value. This may be unintended." << std::endl;
   529     const std::string flux_reconstruction_aux_string = prm.get(
"flux_reconstruction_aux");
   538     const std::string non_physical_behavior_string = prm.get(
"non_physical_behavior");
   539     if (non_physical_behavior_string == 
"return_big_number") { 
non_physical_behavior_type = NonPhysicalBehaviorEnum::return_big_number;}
   545     struct stat info_vtk;
   548                   << 
"Please create the directory and restart. Aborting..." << std::endl;
   556     const std::string renumber_dofs_type_string = prm.get(
"renumber_dofs_type");
   557     if (renumber_dofs_type_string == 
"CuthillMckee") { 
renumber_dofs_type = RenumberDofsType::CuthillMckee; }
   561     pcout << 
"Parsing linear solver subsection..." << std::endl;
   564     pcout << 
"Parsing ODE solver subsection..." << std::endl;
   567     pcout << 
"Parsing manufactured convergence study subsection..." << std::endl;
   570     pcout << 
"Parsing euler subsection..." << std::endl;
   573     pcout << 
"Parsing navier stokes subsection..." << std::endl;
   576     pcout << 
"Parsing reduced order subsection..." << std::endl;
   579     pcout << 
"Parsing hyperreduction subsection..." << std::endl;
   582     pcout << 
"Parsing Burgers subsection..." << std::endl;
   585     pcout << 
"Parsing physics model subsection..." << std::endl;
   588     pcout << 
"Parsing grid refinement study subsection..." << std::endl;
   591     pcout << 
"Parsing artificial dissipation subsection..." << std::endl;
   594     pcout << 
"Parsing limiter subsection..." << std::endl;
   597     pcout << 
"Parsing flow solver subsection..." << std::endl;
   600     pcout << 
"Parsing mesh adaptation subsection..." << std::endl;
   603     pcout << 
"Parsing functional subsection..." << std::endl;
   607     const std::string model_string = prm.get(
"model_type");
   608     if (model_string == 
"large_eddy_simulation") { 
model_type = large_eddy_simulation; }
   609     else if (model_string == 
"reynolds_averaged_navier_stokes") { 
model_type = reynolds_averaged_navier_stokes; }
   611     const std::string pde_string = prm.get(
"pde_type");
   612     if (pde_string == 
"advection") {
   615     } 
else if (pde_string == 
"advection_vector") {
   618     } 
else if (pde_string == 
"diffusion") {
   621     } 
else if (pde_string == 
"convection_diffusion") {
   624     } 
else if (pde_string == 
"burgers_inviscid") {
   627     } 
else if (pde_string == 
"burgers_viscous") {
   630     } 
else if (pde_string == 
"burgers_rewienski") {
   633     } 
else if (pde_string == 
"euler") {
   637     else if (pde_string == 
"navier_stokes") {
   641     else if (pde_string == 
"physics_model") {
   647         else if (
model_type == reynolds_averaged_navier_stokes)
   654     pcout << 
"Parsing time refinement study subsection..." << std::endl;
   657     pcout << 
"Done parsing." << std::endl;
 void parse_parameters(dealii::ParameterHandler &prm)
Parses input file and sets the variables. 
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. 
void parse_parameters(dealii::ParameterHandler &prm)
Parses input file and sets the variables. 
static void declare_parameters(dealii::ParameterHandler &prm)
Declares the possible variables and sets the defaults. 
TimeRefinementStudyParam time_refinement_study_param
Contains the parameters for time refinement study. 
NonPhysicalBehaviorEnum non_physical_behavior_type
Specify behavior on nonphysical results. 
AllParameters()
Constructor. 
FlowSolverParam flow_solver_param
Contains the parameters for simulation cases (flow solver test) 
static void declare_parameters(dealii::ParameterHandler &prm)
Declares the possible variables and sets the defaults. 
Parameters related to the manufactured convergence study. 
bool output_face_results_vtk
Flag for outputting the surface solution vtk files. 
double matching_surface_jac_det_tolerance
static void declare_parameters(dealii::ParameterHandler &prm)
Declares the possible variables and sets the defaults. 
void parse_parameters(dealii::ParameterHandler &prm)
Parses input file and sets the variables. 
Parameters related to time refinement studies. 
void parse_parameters(dealii::ParameterHandler &prm)
Parses input file and sets the variables. 
void parse_parameters(dealii::ParameterHandler &prm)
Parses input file and sets the variables. 
bool use_energy
Flag to use an energy monotonicity test. 
static void declare_parameters(dealii::ParameterHandler &prm)
Declares the possible variables and sets the defaults. 
void parse_parameters(dealii::ParameterHandler &prm)
Function to parse parameters. 
Files for the baseline physics. 
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. 
void parse_parameters(dealii::ParameterHandler &prm)
Parses input file and sets the variables. 
void parse_parameters(dealii::ParameterHandler &prm)
Parses input file and sets the variables. 
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_Aux flux_reconstruction_aux_type
Store flux reconstruction type for the auxiliary variables. 
EulerParam euler_param
Contains parameters for the Euler equations non-dimensionalization. 
void parse_parameters(dealii::ParameterHandler &prm)
Parse parameters. 
Parameters related to the limiter. 
unsigned int dimension
Number of dimensions. Note that it has to match the executable PHiLiP_xD. 
RunType run_type
Selected RunType from the input file. 
MeshAdaptationParam mesh_adaptation_param
Constains parameters for mesh adaptation. 
std::string energy_file
Energy file. 
static void declare_parameters(dealii::ParameterHandler &prm)
Declares the possible variables and sets the defaults. 
static void declare_parameters(dealii::ParameterHandler &prm)
Function to declare parameters. 
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. 
static void declare_parameters(dealii::ParameterHandler &prm)
Declares the possible variables and sets the defaults. 
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. 
NavierStokesParam navier_stokes_param
Contains parameters for the Navier-Stokes equations non-dimensionalization. 
void parse_parameters(dealii::ParameterHandler &prm)
Parses input file and sets the variables. 
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. 
static void declare_parameters(dealii::ParameterHandler &prm)
Declares the possible variables and sets the defaults. 
bool use_L2_norm
Flag to use an L2 energy monotonicity test (for FR) 
double FR_user_specified_correction_parameter_value
User specified flux recontruction correction parameter value. 
Parameters related to the Navier Stokes equations. 
bool use_split_form
Flag to use split form. 
Parameterse related to the functional object. 
void parse_parameters(dealii::ParameterHandler &prm)
Parses input file and sets the variables. 
void parse_parameters(dealii::ParameterHandler &prm)
Parses input file and sets the variables. 
static void declare_parameters(dealii::ParameterHandler &prm)
Declares the possible variables and sets the defaults. 
Parameters for Mesh Adaptation. 
static void declare_parameters(dealii::ParameterHandler &prm)
Declares the possible variables and sets the defaults. 
Parameters related to Physics Models. 
void parse_parameters(dealii::ParameterHandler &prm)
Parses input file and sets the variables. 
void parse_parameters(dealii::ParameterHandler &prm)
Parses input file and sets the variables. 
ReynoldsAveragedNavierStokesModel RANS_model_type
Store the Reynolds-averaged Navier-Stokes (RANS) model type. 
ConvectiveNumericalFlux conv_num_flux_type
Store convective flux type. 
static void declare_parameters(dealii::ParameterHandler &prm)
Declares the possible variables and sets the defaults. 
bool use_inverse_mass_on_the_fly
Flag to use inverse mass matrix on-the-fly for explicit solves. 
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. 
void parse_parameters(dealii::ParameterHandler &prm)
Parses input file and sets the variables. 
FunctionalParam functional_param
Contains parameters for functional. 
static void declare_parameters(dealii::ParameterHandler &prm)
Declares the possible variables and sets the defaults. 
DissipativeNumericalFlux diss_num_flux_type
Store diffusive flux type. 
static void declare_parameters(dealii::ParameterHandler &prm)
Declares the possible variables and sets the defaults. 
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 use_periodic_bc
Flag to use periodic BC. 
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. 
static void declare_parameters(dealii::ParameterHandler &prm)
Declare parameters. 
void parse_parameters(dealii::ParameterHandler &prm)
Function to parse parameters. 
static void declare_parameters(dealii::ParameterHandler &prm)
Function to declare parameters. 
bool store_residual_cpu_time
Flag to store the residual local processor cput time. 
Parameters related to reduced-order model. 
static void declare_parameters(dealii::ParameterHandler &prm)
Declares the possible variables and sets the defaults.