[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
parameters_artificial_dissipation.cpp
1 #include "parameters/parameters_artificial_dissipation.h"
2 
3 namespace PHiLiP {
4 namespace Parameters {
5 
6 void ArtificialDissipationParam::declare_parameters (dealii::ParameterHandler &prm)
7 {
8  prm.enter_subsection("artificial dissipation");
9  {
10 
11  prm.declare_entry("add_artificial_dissipation", "false",
12  dealii::Patterns::Bool(),
13  "Persson's subscell shock capturing artificial dissipation.");
14 
15  prm.declare_entry("artificial_dissipation_type", "laplacian",
16  dealii::Patterns::Selection(
17  "laplacian |"
18  "physical |"
19  "enthalpy_conserving_laplacian |"),
20  "Type of artificial dissipation we want to implement. Choices are laplacian, physical and enthalpy_conserving_laplacian");
21 
22  prm.declare_entry("artificial_dissipation_test_type", "poly_order_convergence",
23  dealii::Patterns::Selection(
24  "residual_convergence |"
25  "discontinuity_sensor_activation |"
26  "enthalpy_conservation |"
27  "poly_order_convergence |"),
28  "Type of artificial dissipation test type we want to implement. Choices are residual_convergence, discontinuity_sensor_activation, poly_order_convergence");
29 
30  prm.declare_entry("mu_artificial_dissipation", "1.0",
31  dealii::Patterns::Double(-1e20,1e20),
32  "Mu (viscosity) from Persson's subcell shock capturing.");
33 
34  prm.declare_entry("kappa_artificial_dissipation", "1.0",
35  dealii::Patterns::Double(-1e20,1e20),
36  "Kappa from Persson's subcell shock capturing");
37 
38  prm.declare_entry("use_enthalpy_error", "false",
39  dealii::Patterns::Bool(),
40  "By default we calculate the entropy error from the conservative variables. Otherwise, compute the enthalpy error. An example is in Euler Gaussian bump.");
41 
42  }
43  prm.leave_subsection();
44 }
45 
46 void ArtificialDissipationParam::parse_parameters (dealii::ParameterHandler &prm)
47 {
48  prm.enter_subsection("artificial dissipation");
49  {
50 
51  add_artificial_dissipation = prm.get_bool("add_artificial_dissipation");
52  use_enthalpy_error = prm.get_bool("use_enthalpy_error");
53 
54  const std::string artificial_dissipation_string = prm.get("artificial_dissipation_type");
55  if (artificial_dissipation_string == "laplacian")
56  {
57  artificial_dissipation_type = laplacian;
58  }
59  else if (artificial_dissipation_string == "physical")
60  {
61  artificial_dissipation_type = physical;
62  }
63  else if (artificial_dissipation_string == "enthalpy_conserving_laplacian")
64  {
65  artificial_dissipation_type = enthalpy_conserving_laplacian;
66  }
67 
68  const std::string artificial_dissipation_test_string = prm.get("artificial_dissipation_test_type");
69  if (artificial_dissipation_test_string == "residual_convergence")
70  {
71  artificial_dissipation_test_type = residual_convergence;
72  }
73  else if (artificial_dissipation_test_string == "discontinuity_sensor_activation")
74  {
75  artificial_dissipation_test_type = discontinuity_sensor_activation;
76  }
77  else if (artificial_dissipation_test_string == "enthalpy_conservation")
78  {
79  artificial_dissipation_test_type = enthalpy_conservation;
80  }
81  else if (artificial_dissipation_test_string == "poly_order_convergence")
82  {
83  artificial_dissipation_test_type = poly_order_convergence;
84  }
85 
86  mu_artificial_dissipation = prm.get_double("mu_artificial_dissipation");
87  kappa_artificial_dissipation = prm.get_double("kappa_artificial_dissipation");
88  }
89  prm.leave_subsection();
90 }
91 
92 } // Parameters namespace
93 } // PHiLiP namespace
Files for the baseline physics.
Definition: ADTypes.hpp:10
static void declare_parameters(dealii::ParameterHandler &prm)
Function to declare parameters.
ArtificialDissipationTestType artificial_dissipation_test_type
Selected dissipation test type.
double kappa_artificial_dissipation
Parameter kappa from Persson and Peraire, 2008.
double mu_artificial_dissipation
Parameter mu from Persson & Peraire, 2008.
bool add_artificial_dissipation
Flag to add artificial dissipation from Persson's shock capturing paper.
ArtificialDissipationType artificial_dissipation_type
Selected artificial dissipation type specified in the input.
void parse_parameters(dealii::ParameterHandler &prm)
Function to parse parameters.