[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
parameters_reduced_order.cpp
1 #include "parameters/parameters_reduced_order.h"
2 namespace PHiLiP {
3 namespace Parameters {
4 
5 void ReducedOrderModelParam::declare_parameters (dealii::ParameterHandler &prm)
6 {
7  prm.enter_subsection("reduced order");
8  {
9  prm.declare_entry("adaptation_tolerance", "1",
10  dealii::Patterns::Double(0, dealii::Patterns::Double::max_double_value),
11  "Tolerance for POD adaptation");
12  prm.declare_entry("path_to_search", ".",
13  dealii::Patterns::FileName(dealii::Patterns::FileName::FileType::input),
14  "Path to search for saved snapshots or POD basis.");
15  prm.declare_entry("reduced_residual_tolerance", "1E-13",
16  dealii::Patterns::Double(0, dealii::Patterns::Double::max_double_value),
17  "Tolerance for nonlinear reduced residual");
18  prm.declare_entry("num_halton", "0",
19  dealii::Patterns::Integer(0, dealii::Patterns::Integer::max_int_value),
20  "Number of Halton sequence points to add to initial snapshot set");
21  prm.declare_entry("file_path_for_snapshot_locations", "",
22  dealii::Patterns::FileName(dealii::Patterns::FileName::FileType::input),
23  "Path to search for lhs snapshots (should contain snapshot_table)");
24  prm.declare_entry("recomputation_coefficient", "5",
25  dealii::Patterns::Integer(0, dealii::Patterns::Integer::max_int_value),
26  "Number of Halton sequence points to add to initial snapshot set");
27  prm.declare_entry("parameter_names", "mach, alpha",
28  dealii::Patterns::List(dealii::Patterns::Anything(), 0, 10, ","),
29  "Names of parameters for adaptive sampling");
30  prm.declare_entry("parameter_min_values", "0.4, 0",
31  dealii::Patterns::List(dealii::Patterns::Double(), 0, 10, ","),
32  "Minimum values for parameters");
33  prm.declare_entry("parameter_max_values", "0.7, 4",
34  dealii::Patterns::List(dealii::Patterns::Double(), 0, 10, ","),
35  "Maximum values for parameters");
36  prm.declare_entry("number_modes", "0",
37  dealii::Patterns::Integer(0, dealii::Patterns::Integer::max_int_value),
38  "Number of modes used in the POD Basis");
39  prm.declare_entry("singular_value_threshold", "1",
40  dealii::Patterns::Double(0, dealii::Patterns::Double::max_double_value),
41  "Threshold for the Singular Value cutoff");
42  prm.declare_entry("output_snapshot_every_x_timesteps","0",
43  dealii::Patterns::Integer(0,dealii::Patterns::Integer::max_int_value),
44  "Number of Timesteps before snapshot is added");
45  prm.declare_entry("FOM_error_linear_solver_type", "direct",
46  dealii::Patterns::Selection("direct|gmres"),
47  "Type of linear solver used for first adjoint problem (DWR between FOM and ROM)"
48  "Choices are <direct|gmres>.");
49  prm.declare_entry("residual_error_bool", "false",
50  dealii::Patterns::Bool(),
51  "Use residual/reduced residual for error indicator instead of DWR. False by default.");
52  }
53  prm.leave_subsection();
54 }
55 
56 void ReducedOrderModelParam::parse_parameters (dealii::ParameterHandler &prm)
57 {
58  prm.enter_subsection("reduced order");
59  {
60  adaptation_tolerance = prm.get_double("adaptation_tolerance");
61  reduced_residual_tolerance = prm.get_double("reduced_residual_tolerance");
62  singular_value_threshold = prm.get_double("singular_value_threshold");
63  num_halton = prm.get_integer("num_halton");
64  file_path_for_snapshot_locations = prm.get("file_path_for_snapshot_locations");
65  recomputation_coefficient = prm.get_integer("recomputation_coefficient");
66  number_modes = prm.get_integer("number_modes");
67  output_snapshot_every_x_timesteps = prm.get_integer("output_snapshot_every_x_timesteps");
68  path_to_search = prm.get("path_to_search");
69 
70  std::string parameter_names_string = prm.get("parameter_names");
71  std::unique_ptr<dealii::Patterns::PatternBase> ListPatternNames(new dealii::Patterns::List(dealii::Patterns::Anything(), 0, 10, ",")); //Note, in a future version of dealii, this may change from a unique_ptr to simply the object. Will need to use std::move(ListPattern) in next line.
72  parameter_names = dealii::Patterns::Tools::Convert<decltype(parameter_names)>::to_value(parameter_names_string, ListPatternNames);
73 
74  std::string parameter_min_string = prm.get("parameter_min_values");
75  std::unique_ptr<dealii::Patterns::PatternBase> ListPatternMin(new dealii::Patterns::List(dealii::Patterns::Double(), 0, 10, ",")); //Note, in a future version of dealii, this may change from a unique_ptr to simply the object. Will need to use std::move(ListPattern) in next line.
76  parameter_min_values = dealii::Patterns::Tools::Convert<decltype(parameter_min_values)>::to_value(parameter_min_string, ListPatternMin);
77 
78  std::string parameter_max_string = prm.get("parameter_max_values");
79  std::unique_ptr<dealii::Patterns::PatternBase> ListPatternMax(new dealii::Patterns::List(dealii::Patterns::Double(), 0, 10, ",")); //Note, in a future version of dealii, this may change from a unique_ptr to simply the object. Will need to use std::move(ListPattern) in next line.
80  parameter_max_values = dealii::Patterns::Tools::Convert<decltype(parameter_max_values)>::to_value(parameter_max_string, ListPatternMax);
81 
82 
83  const std::string solver_string = prm.get("FOM_error_linear_solver_type");
84  if (solver_string == "direct") FOM_error_linear_solver_type = LinearSolverEnum::direct;
85  if (solver_string == "gmres") FOM_error_linear_solver_type = LinearSolverEnum::gmres;
86  residual_error_bool = prm.get_bool("residual_error_bool");
87  }
88  prm.leave_subsection();
89 }
90 
91 } // Parameters namespace
92 } // PHiLiP namespace
bool residual_error_bool
Use residual/reduced residual for error indicator instead of DWR. False by default.
std::vector< std::string > parameter_names
Names of parameters.
static void declare_parameters(dealii::ParameterHandler &prm)
Declares the possible variables and sets the defaults.
int number_modes
Number of modes to include in the POD Basis.
int num_halton
Number of Halton sequence points to add to initial snapshot set.
void parse_parameters(dealii::ParameterHandler &prm)
Parses input file and sets the variables.
std::vector< double > parameter_min_values
Minimum value of parameters.
Files for the baseline physics.
Definition: ADTypes.hpp:10
std::string path_to_search
Path to search for snapshots or saved POD basis.
int recomputation_coefficient
Recomputation parameter for adaptive sampling algorithm.
std::string file_path_for_snapshot_locations
Path to search for file with pre-determined snapshot locations used to build POD (actual FOM snapshot...
std::vector< double > parameter_max_values
Maximum value of parameters.
double singular_value_threshold
Singular Value Threshold in the POD Basis.
LinearSolverEnum FOM_error_linear_solver_type
Type of linear solver used for first adjoint problem (DWR between FOM and ROM) (direct or gmres) ...
int output_snapshot_every_x_timesteps
Number of timesteps before putting solution in snapshot matrix.
double reduced_residual_tolerance
Tolerance of the reduced-order nonlinear residual.
double adaptation_tolerance
Tolerance for POD adaptation.