[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
parameters_grid_refinement_study.cpp
1 #include <array>
2 
3 #include "parameters_grid_refinement_study.h"
4 
5 #include "parameters/parameters_manufactured_convergence_study.h"
6 
7 namespace PHiLiP {
8 
9 namespace Parameters {
10 
11 void GridRefinementStudyParam::declare_parameters(dealii::ParameterHandler &prm)
12 {
13  prm.enter_subsection("grid refinement study");
14  {
17 
18  prm.enter_subsection("grid refinement");
19  {
21  }
22  prm.leave_subsection();
23 
24  // looping over the elements of the array to declare the different refinement
25  for(unsigned int i = 0; i < MAX_REFINEMENTS; ++i)
26  {
27  prm.enter_subsection("grid refinement [" + dealii::Utilities::int_to_string(i,1) + "]");
28  {
30  }
31  prm.leave_subsection();
32  }
33 
34  prm.declare_entry("poly_degree", "1",
35  dealii::Patterns::Integer(),
36  "Polynomial order of starting mesh.");
37 
38  prm.declare_entry("poly_degree_max", "5",
39  dealii::Patterns::Integer(),
40  "Maximum polynomial order.");
41 
42  prm.declare_entry("poly_degree_grid", "2",
43  dealii::Patterns::Integer(),
44  "Polynomial degree of the grid.");
45 
46  prm.declare_entry("num_refinements", "0",
47  dealii::Patterns::Integer(0, MAX_REFINEMENTS),
48  "Number of different refinements to be performed.");
49 
50  prm.declare_entry("grid_type", "hypercube",
51  dealii::Patterns::Selection("hypercube|sinehypercube|read_grid"),
52  "Enum of generated grid. "
53  "If read_grid, must have grids xxxx.msh."
54  "Choices are <hypercube|read_grid>.");
55 
56  prm.declare_entry("input_grid", "xxxx",
57  dealii::Patterns::FileName(dealii::Patterns::FileName::FileType::input),
58  "Name of Gmsh grid xxxx.msh used in the grid refinement study if read_grid is chosen as the grid_type.");
59 
60  prm.declare_entry("grid_left", "0.0",
61  dealii::Patterns::Double(),
62  "for grid_type hypercube, left bound of domain.");
63 
64  prm.declare_entry("grid_right", "1.0",
65  dealii::Patterns::Double(),
66  "for grid_type hypercube, right bound of domain.");
67 
68  prm.declare_entry("grid_size", "4",
69  dealii::Patterns::Integer(),
70  "Initial grid size (number of elements per side).");
71 
72  prm.declare_entry("use_interpolation", "false",
73  dealii::Patterns::Bool(),
74  "Indicates whether to interpolate the problem instead of solving with DG.");
75 
76  prm.declare_entry("approximate_functional", "false",
77  dealii::Patterns::Bool(),
78  "Indicates whether function is to be approximated from manufactured solution"
79  "or exact value read from functional_value parameter.");
80 
81  prm.declare_entry("functional_value", "0.0",
82  dealii::Patterns::Double(),
83  "Exact value of functional for goal-oriented convergence.");
84 
85  prm.declare_entry("output_vtk", "true",
86  dealii::Patterns::Bool(),
87  "Output flag for grid_refinement vtk files.");
88 
89  prm.declare_entry("output_adjoint_vtk", "false",
90  dealii::Patterns::Bool(),
91  "output flag for adjoint vtk files.");
92 
93  prm.declare_entry("output_solution_error", "true",
94  dealii::Patterns::Bool(),
95  "ouput the convergence table for the solution error.");
96 
97  prm.declare_entry("output_functional_error", "false",
98  dealii::Patterns::Bool(),
99  "ouput the convergence table for the functional error.");
100 
101  prm.declare_entry("output_gnuplot_solution", "true",
102  dealii::Patterns::Bool(),
103  "Output flag for gnuplot solution error figure.");
104 
105  prm.declare_entry("output_gnuplot_functional", "false",
106  dealii::Patterns::Bool(),
107  "Output flag for gnuplot functional error figure.");
108 
109  prm.declare_entry("refresh_gnuplot", "true",
110  dealii::Patterns::Bool(),
111  "Indicates whetherto output a new gnuplot figure at every iteration."
112  "Requires output_gnuplot == true.");
113 
114  prm.declare_entry("output_solution_time", "false",
115  dealii::Patterns::Bool(),
116  "Output flag for wall clock solution timing.");
117 
118  prm.declare_entry("output_adjoint_time", "false",
119  dealii::Patterns::Bool(),
120  "Output flag for wall clock adjoint timing.");
121  }
122  prm.leave_subsection();
123 }
124 
125 void GridRefinementStudyParam::parse_parameters(dealii::ParameterHandler &prm)
126 {
127  prm.enter_subsection("grid refinement study");
128  {
131 
132  // default section gets written into vector pos 1 by default
133  prm.enter_subsection("grid refinement");
134  {
135  grid_refinement_param_vector[0].parse_parameters(prm);
136  }
137  prm.leave_subsection();
138 
139  num_refinements = prm.get_integer("num_refinements");
140 
141  // overwrite if num_refinements > 0
142  for(unsigned int i = 0; i < num_refinements; ++i)
143  {
144  prm.enter_subsection("grid refinement [" + dealii::Utilities::int_to_string(i,1) + "]");
145  {
146  grid_refinement_param_vector[i].parse_parameters(prm);
147  }
148  prm.leave_subsection();
149  }
150 
151  poly_degree = prm.get_integer("poly_degree");
152  poly_degree_max = prm.get_integer("poly_degree_max");
153  poly_degree_grid = prm.get_integer("poly_degree_grid");
154 
155  const std::string grid_string = prm.get("grid_type");
157  if(grid_string == "hypercube") {grid_type = GridEnum::hypercube;}
158  else if(grid_string == "read_grid") {grid_type = GridEnum::read_grid;
159  input_grid = prm.get("input_grid");}
160 
161  grid_left = prm.get_double("grid_left");
162  grid_right = prm.get_double("grid_right");
163 
164  grid_size = prm.get_integer("grid_size");
165 
166  use_interpolation = prm.get_bool("use_interpolation");
167 
168  approximate_functional = prm.get_bool("approximate_functional");
169  functional_value = prm.get_double("functional_value");
170 
171  output_vtk = prm.get_bool("output_vtk");
172  output_adjoint_vtk = prm.get_bool("output_adjoint_vtk");
173 
174  output_solution_error = prm.get_bool("output_solution_error");
175  output_functional_error = prm.get_bool("output_functional_error");
176 
177  output_gnuplot_solution = prm.get_bool("output_gnuplot_solution");
178  output_gnuplot_functional = prm.get_bool("output_gnuplot_functional");
179  refresh_gnuplot = prm.get_bool("refresh_gnuplot");
180 
181  output_solution_time = prm.get_bool("output_solution_time");
182  output_adjoint_time = prm.get_bool("output_adjoint_time");
183  }
184  prm.leave_subsection();
185 }
186 
187 } // namespace Parameters
188 
189 } // namespace PHiLiP
bool output_gnuplot_functional
Flag to enable output of gnuplot graph of functional error convergence.
bool refresh_gnuplot
Flag to enable gnuplot refresh between iteration runs.
unsigned int poly_degree_grid
Initial grid polynomial degree.
bool use_interpolation
Flag to enable interpolation operation.
void parse_parameters(dealii::ParameterHandler &prm)
Parses input file and sets the variables.
static const unsigned int MAX_REFINEMENTS
Maximum number of different refinement procedures.
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.
unsigned int num_refinements
Number of different refinement procedures stored, 0 indicates to use the default pathway.
Files for the baseline physics.
Definition: ADTypes.hpp:10
bool output_solution_time
Flag to enable output of grid refinement wall-clock solution time.
ManufacturedSolutionParam manufactured_solution_param
Manufactured solution parameterse to be used with grid refinement study.
std::array< GridRefinementParam, MAX_REFINEMENTS > grid_refinement_param_vector
Array of grid refinement parameters to be run as part of grid refinement study.
unsigned int poly_degree_max
Maximimum allocated solution polynomial degree.
bool approximate_functional
Flag to enable approximation of the functional value on a fine grid before refinement run...
static void declare_parameters(dealii::ParameterHandler &prm)
Declares the possible variables and sets the defaults.
std::string input_grid
Input pathway for GridEnum::read_grid type.
bool output_solution_error
Flag to enable output of grid refinement solution error convergence.
bool output_gnuplot_solution
Flag to enable output of gnuplot graph of solution error convergence.
static void declare_parameters(dealii::ParameterHandler &prm)
Declares the possible variables and sets the defaults.
bool output_vtk
Flag to enable output of grid refinement .vtk file.
FunctionalParam functional_param
Functional parameters to be used with grid refinement study.
double functional_value
Specified exact functional value for comparison of error convergence.
bool output_functional_error
Flag to enable output of grid refinement functional error convergence.
unsigned int grid_size
Number of initial elements in each axis for GridEnum::hypercube type.
static void declare_parameters(dealii::ParameterHandler &prm)
Declares the possible variables and sets the defaults.
double grid_left
Lower coordinate bound for GridEnum::hypercube type.
bool output_adjoint_time
Flag to enable output of grid refinement wall-clock adjoint time.
unsigned int poly_degree
Initial solution polynomial degree.
double grid_right
Upper coordinate bound for GridEnum::hypercube type.
static void declare_parameters(dealii::ParameterHandler &prm)
Declares the possible variables and sets the defaults.
GridEnum
Types of grids that can be used for convergence study.
bool output_adjoint_vtk
Flag to enable output of adjoint .vtk file.