-
Notifications
You must be signed in to change notification settings - Fork 36
/
all_parameters.cpp
285 lines (242 loc) · 14.4 KB
/
all_parameters.cpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
#include <deal.II/base/mpi.h>
#include <deal.II/base/utilities.h>
#include "parameters/all_parameters.h"
namespace PHiLiP {
namespace Parameters {
AllParameters::AllParameters ()
: manufactured_convergence_study_param(ManufacturedConvergenceStudyParam())
, ode_solver_param(ODESolverParam())
, linear_solver_param(LinearSolverParam())
, euler_param(EulerParam())
, grid_refinement_study_param(GridRefinementStudyParam())
, pcout(std::cout, dealii::Utilities::MPI::this_mpi_process(MPI_COMM_WORLD)==0)
{ }
void AllParameters::declare_parameters (dealii::ParameterHandler &prm)
{
const int mpi_rank = dealii::Utilities::MPI::this_mpi_process(MPI_COMM_WORLD);
dealii::ConditionalOStream pcout(std::cout, mpi_rank==0);
pcout << "Declaring inputs." << std::endl;
prm.declare_entry("dimension", "-1",
dealii::Patterns::Integer(),
"Number of dimensions");
prm.declare_entry("mesh_type", "default_triangulation",
dealii::Patterns::Selection(
" default_triangulation | "
" triangulation | "
" parallel_shared_triangulation | "
" parallel_distributed_triangulation"),
"Type of triangulation to be used."
"Note: parralel_distributed_triangulation not availible int 1D."
" <default_triangulation | "
" triangulation | "
" parallel_shared_triangulation |"
" parallel_distributed_triangulation>.");
prm.declare_entry("overintegration", "0",
dealii::Patterns::Integer(),
"Number of extra quadrature points to use."
"If overintegration=0, then we use n_quad = soln_degree + 1.");
prm.declare_entry("use_weak_form", "true",
dealii::Patterns::Bool(),
"Use weak form by default. If false, use strong form.");
prm.declare_entry("use_collocated_nodes", "false",
dealii::Patterns::Bool(),
"Use Gauss-Legendre by default. Otherwise, use Gauss-Lobatto to collocate.");
prm.declare_entry("use_split_form", "false",
dealii::Patterns::Bool(),
"Use original form by defualt. Otherwise, split the fluxes.");
prm.declare_entry("use_periodic_bc", "false",
dealii::Patterns::Bool(),
"Use other boundary conditions by default. Otherwise use periodic (for 1d burgers only");
prm.declare_entry("use_energy", "false",
dealii::Patterns::Bool(),
"Not calculate energy by default. Otherwise, get energy per iteration.");
prm.declare_entry("use_L2_norm", "false",
dealii::Patterns::Bool(),
"Not calculate L2 norm by default (M+K). Otherwise, get L2 norm per iteration.");
prm.declare_entry("use_classical_Flux_Reconstruction", "false",
dealii::Patterns::Bool(),
"Not use Classical Flux Reconstruction by default. Otherwise, use Classical Flux Reconstruction.");
prm.declare_entry("flux_reconstruction", "cDG",
dealii::Patterns::Selection("cDG | cSD | cHU | cNegative | cNegative2 | cPlus | cPlus1D |c10Thousand | cHULumped"),
"Flux Reconstruction. "
"Choices are <cDG | cSD | cHU | cNegative | cNegative2 | cPlus | cPlus1D | c10Thousand | cHULumped>.");
prm.declare_entry("flux_reconstruction_aux", "kDG",
dealii::Patterns::Selection("kDG | kSD | kHU | kNegative | kNegative2 | kPlus | k10Thousand"),
"Flux Reconstruction for Auxiliary Equation. "
"Choices are <kDG | kSD | kHU | kNegative | kNegative2 | kPlus | k10Thousand>.");
prm.declare_entry("add_artificial_dissipation", "false",
dealii::Patterns::Bool(),
"Persson's subscell shock capturing artificial dissipation.");
prm.declare_entry("sipg_penalty_factor", "1.0",
dealii::Patterns::Double(1.0,1e200),
"Scaling of Symmetric Interior Penalty term to ensure coercivity.");
prm.declare_entry("test_type", "run_control",
dealii::Patterns::Selection(
" run_control | "
" grid_refinement_study | "
" burgers_energy_stability | "
" diffusion_exact_adjoint | "
" optimization_inverse_manufactured | "
" euler_gaussian_bump | "
" euler_gaussian_bump_adjoint | "
" euler_cylinder | "
" euler_cylinder_adjoint | "
" euler_vortex | "
" euler_entropy_waves | "
" euler_split_taylor_green | "
" euler_bump_optimization | "
" euler_naca_optimization | "
" shock_1d | "
" euler_naca0012 | "
" advection_periodicity"),
"The type of test we want to solve. "
"Choices are (only run control has been coded up for now)"
" <run_control | "
" grid_refinement_study | "
" burgers_energy_stability | "
" diffusion_exact_adjoint | "
" optimization_inverse_manufactured | "
" euler_gaussian_bump | "
" euler_gaussian_bump_adjoint | "
" euler_cylinder | "
" euler_cylinder_adjoint | "
" euler_vortex | "
" euler_entropy_waves | "
" euler_split_taylor_green |"
" euler_bump_optimization | "
" euler_naca_optimization | "
" shock_1d | "
" euler_naca0012 | "
" advection_periodicity >.");
prm.declare_entry("pde_type", "advection",
dealii::Patterns::Selection(
" advection | "
" diffusion | "
" convection_diffusion | "
" advection_vector | "
" burgers_inviscid | "
" euler |"
" mhd"),
"The PDE we want to solve. "
"Choices are "
" <advection | "
" diffusion | "
" convection_diffusion | "
" advection_vector | "
" burgers_inviscid | "
" euler | "
" mhd>.");
prm.declare_entry("conv_num_flux", "lax_friedrichs",
dealii::Patterns::Selection("lax_friedrichs | roe | split_form"),
"Convective numerical flux. "
"Choices are <lax_friedrichs | roe | split_form>.");
prm.declare_entry("diss_num_flux", "symm_internal_penalty",
dealii::Patterns::Selection("symm_internal_penalty | bassi_rebay_2"),
"Dissipative numerical flux. "
"Choices are <symm_internal_penalty | bassi_rebay_2>.");
Parameters::LinearSolverParam::declare_parameters (prm);
Parameters::ManufacturedConvergenceStudyParam::declare_parameters (prm);
Parameters::ODESolverParam::declare_parameters (prm);
Parameters::EulerParam::declare_parameters (prm);
Parameters::GridRefinementStudyParam::declare_parameters (prm);
pcout << "Done declaring inputs." << std::endl;
}
void AllParameters::parse_parameters (dealii::ParameterHandler &prm)
{
pcout << "Parsing main input..." << std::endl;
dimension = prm.get_integer("dimension");
const std::string mesh_type_string = prm.get("mesh_type");
if (mesh_type_string == "default_triangulation") { mesh_type = default_triangulation; }
else if (mesh_type_string == "triangulation") { mesh_type = triangulation; }
else if (mesh_type_string == "parallel_shared_triangulation") { mesh_type = parallel_shared_triangulation; }
else if (mesh_type_string == "parallel_distributed_triangulation") { mesh_type = parallel_distributed_triangulation; }
const std::string test_string = prm.get("test_type");
if (test_string == "run_control") { test_type = run_control; }
else if (test_string == "grid_refinement_study") { test_type = grid_refinement_study; }
else if (test_string == "burgers_energy_stability") { test_type = burgers_energy_stability; }
else if (test_string == "diffusion_exact_adjoint") { test_type = diffusion_exact_adjoint; }
else if (test_string == "euler_gaussian_bump") { test_type = euler_gaussian_bump; }
else if (test_string == "euler_gaussian_bump_adjoint") { test_type = euler_gaussian_bump_adjoint; }
else if (test_string == "euler_cylinder") { test_type = euler_cylinder; }
else if (test_string == "euler_cylinder_adjoint") { test_type = euler_cylinder_adjoint; }
else if (test_string == "euler_vortex") { test_type = euler_vortex; }
else if (test_string == "euler_entropy_waves") { test_type = euler_entropy_waves; }
else if (test_string == "advection_periodicity") { test_type = advection_periodicity; }
else if (test_string == "euler_split_taylor_green") { test_type = euler_split_taylor_green; }
else if (test_string == "euler_bump_optimization") { test_type = euler_bump_optimization; }
else if (test_string == "euler_naca_optimization") { test_type = euler_naca_optimization; }
else if (test_string == "shock_1d") { test_type = shock_1d; }
else if (test_string == "euler_naca0012") { test_type = euler_naca0012; }
else if (test_string == "optimization_inverse_manufactured") {test_type = optimization_inverse_manufactured; }
const std::string pde_string = prm.get("pde_type");
if (pde_string == "advection") {
pde_type = advection;
nstate = 1;
} else if (pde_string == "advection_vector") {
pde_type = advection_vector;
nstate = 2;
} else if (pde_string == "diffusion") {
pde_type = diffusion;
nstate = 1;
} else if (pde_string == "convection_diffusion") {
pde_type = convection_diffusion;
nstate = 1;
} else if (pde_string == "burgers_inviscid") {
pde_type = burgers_inviscid;
nstate = dimension;
} else if (pde_string == "euler") {
pde_type = euler;
nstate = dimension+2;
}
overintegration = prm.get_integer("overintegration");
use_weak_form = prm.get_bool("use_weak_form");
use_collocated_nodes = prm.get_bool("use_collocated_nodes");
use_split_form = prm.get_bool("use_split_form");
use_periodic_bc = prm.get_bool("use_periodic_bc");
use_energy = prm.get_bool("use_energy");
use_L2_norm = prm.get_bool("use_L2_norm");
use_classical_FR = prm.get_bool("use_classical_Flux_Reconstruction");
add_artificial_dissipation = prm.get_bool("add_artificial_dissipation");
sipg_penalty_factor = prm.get_double("sipg_penalty_factor");
const std::string conv_num_flux_string = prm.get("conv_num_flux");
if (conv_num_flux_string == "lax_friedrichs") conv_num_flux_type = lax_friedrichs;
if (conv_num_flux_string == "split_form") conv_num_flux_type = split_form;
if (conv_num_flux_string == "roe") conv_num_flux_type = roe;
const std::string diss_num_flux_string = prm.get("diss_num_flux");
if (diss_num_flux_string == "symm_internal_penalty") diss_num_flux_type = symm_internal_penalty;
if (diss_num_flux_string == "bassi_rebay_2") {
diss_num_flux_type = bassi_rebay_2;
sipg_penalty_factor = 0.0;
}
const std::string flux_reconstruction_string = prm.get("flux_reconstruction");
if (flux_reconstruction_string == "cDG") flux_reconstruction_type = cDG;
if (flux_reconstruction_string == "cSD") flux_reconstruction_type = cSD;
if (flux_reconstruction_string == "cHU") flux_reconstruction_type = cHU;
if (flux_reconstruction_string == "cNegative") flux_reconstruction_type = cNegative;
if (flux_reconstruction_string == "cNegative2") flux_reconstruction_type = cNegative2;
if (flux_reconstruction_string == "cPlus") flux_reconstruction_type = cPlus;
if (flux_reconstruction_string == "cPlus1D") flux_reconstruction_type = cPlus1D;
if (flux_reconstruction_string == "c10Thousand") flux_reconstruction_type = c10Thousand;
if (flux_reconstruction_string == "cHULumped") flux_reconstruction_type = cHULumped;
const std::string flux_reconstruction_aux_string = prm.get("flux_reconstruction_aux");
if (flux_reconstruction_aux_string == "kDG") flux_reconstruction_aux_type = kDG;
if (flux_reconstruction_aux_string == "kSD") flux_reconstruction_aux_type = kSD;
if (flux_reconstruction_aux_string == "kHU") flux_reconstruction_aux_type = kHU;
if (flux_reconstruction_aux_string == "kNegative") flux_reconstruction_aux_type = kNegative;
if (flux_reconstruction_aux_string == "kNegative2") flux_reconstruction_aux_type = kNegative2;
if (flux_reconstruction_aux_string == "kPlus") flux_reconstruction_aux_type = kPlus;
if (flux_reconstruction_aux_string == "k10Thousand") flux_reconstruction_aux_type = k10Thousand;
pcout << "Parsing linear solver subsection..." << std::endl;
linear_solver_param.parse_parameters (prm);
pcout << "Parsing ODE solver subsection..." << std::endl;
ode_solver_param.parse_parameters (prm);
pcout << "Parsing manufactured convergence study subsection..." << std::endl;
manufactured_convergence_study_param.parse_parameters (prm);
pcout << "Parsing euler subsection..." << std::endl;
euler_param.parse_parameters (prm);
pcout << "Parsing grid refinement study subsection..." << std::endl;
grid_refinement_study_param.parse_parameters (prm);
pcout << "Done parsing." << std::endl;
}
} // Parameters namespace
} // PHiLiP namespace