Reference documentation for deal.II version 9.6.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
parsed_convergence_table.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2019 - 2023 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
19
21
22namespace
23{
24 std::vector<std::string>
25 get_unique_component_names(const std::vector<std::string> &component_names)
26 {
27 auto elements = component_names;
28 elements.erase(std::unique(elements.begin(), elements.end()),
29 elements.end());
30 return elements;
31 }
32
33
34
35 std::vector<ComponentMask>
36 get_unique_component_masks(const std::vector<std::string> &component_names)
37 {
38 const auto unique_component_names =
39 get_unique_component_names(component_names);
40
41 std::vector<ComponentMask> masks;
42
43 std::vector<std::vector<bool>> bools(
44 unique_component_names.size(),
45 std::vector<bool>(component_names.size(), false));
46
47 unsigned int j = 0;
48 for (unsigned int i = 0; i < component_names.size(); ++i)
49 {
50 if (unique_component_names[j] != component_names[i])
51 masks.emplace_back(bools[j++]);
52 bools[j][i] = true;
53 }
54 masks.emplace_back(bools[j++]);
55 AssertDimension(j, unique_component_names.size());
56 return masks;
57 }
58} // namespace
59
60
61
63 const std::vector<std::string> &solution_names,
64 const std::vector<std::set<VectorTools::NormType>> &list_of_error_norms)
65 : ParsedConvergenceTable(solution_names,
66 list_of_error_norms,
67 2.0,
68 {"cells", "dofs"},
69 "dofs",
70 "reduction_rate_log2",
71 "",
72 3,
73 true)
74{}
75
76
77
79 const std::vector<std::string> &component_names,
80 const std::vector<std::set<VectorTools::NormType>> &list_of_error_norms,
81 const double exponent,
82 const std::set<std::string> &extra_columns,
83 const std::string &rate_key,
84 const std::string &rate_mode,
85 const std::string &error_file_name,
86 const unsigned int precision,
87 const bool compute_error)
88 : component_names(component_names)
89 , unique_component_names(get_unique_component_names(component_names))
90 , unique_component_masks(get_unique_component_masks(component_names))
91 , norms_per_unique_component(list_of_error_norms)
92 , exponent(exponent)
93 , extra_columns(extra_columns)
94 , rate_key(rate_key)
95 , rate_mode(rate_mode)
96 , precision(precision)
97 , error_file_name(error_file_name)
98 , compute_error(compute_error)
99{
102}
103
104
105
106void
108{
109 prm.add_parameter("Enable computation of the errors",
111 "When set to false, no computations are performed.");
112
113 prm.add_parameter("Error precision",
114 precision,
115 "Number of digits to use when printing the error.",
117
118 prm.add_parameter("Error file name",
120 "Set this to a filename with extension .txt, .gpl, .org, "
121 "or .tex to enable writing the convergence table to a "
122 "file.");
123
124 prm.add_parameter(
125 "List of error norms to compute",
127 "Each component is separated by a semicolon "
128 "and each norm by a comma. See the documentation of VectorTools::NormType "
129 "for a list of implemented norms. If you want to skip a component, leave "
130 "its entry empty.");
131
132
133 prm.add_parameter("Exponent for p-norms",
134 exponent,
135 "The exponent to use when computing p-norms.",
137
138 prm.add_parameter("Extra columns",
140 "Extra columns to add to the table. Available options "
141 "are dofs and cells.",
142 Patterns::List(Patterns::Selection("dofs|cells")));
143
144 prm.add_parameter("Rate key",
145 rate_key,
146 "Key to use when computing convergence rates. If "
147 "this is set to a column that is not present, or to the "
148 "empty string, then no error rates are computed.");
149
150 prm.add_parameter("Rate mode",
151 rate_mode,
152 "What type of error rate to compute. Available options are "
153 "reduction_rate_log2, reduction_rate, and none.",
155 "reduction_rate|reduction_rate_log2|none"));
156}
157
158
159
160void
162{
163 if (compute_error)
164 {
165 // Add convergence rates if the rate_key is not empty
166 if (rate_key != "")
167 {
168 bool has_key = false;
169 for (const auto &col : extra_columns)
170 {
171 if (rate_key == col)
172 has_key = true;
173
174 if (col != "")
176 }
177
178 for (const auto &extra_col : extra_column_functions)
179 if (extra_col.second.second == false)
180 {
181 if (rate_key == extra_col.first)
182 has_key = true;
184 extra_col.first);
185 }
186
187 if (has_key)
188 {
189 if (rate_mode == "reduction_rate_log2")
192 else if (rate_mode == "reduction_rate")
195 else
196 {
197 Assert(rate_mode != "none", ExcInternalError());
198 }
199 }
200 else
201 {
202 AssertThrow(rate_key != "",
204 "You specified the key <" + rate_key +
205 "> to compute convergence rates, but that key does "
206 "not exist in the current table."));
207 }
208 }
209 }
210}
211
212
213
214void
216{
217 if (compute_error)
218 {
220 table.write_text(out);
221 output_table();
222 }
223}
224
225
226
227void
229{
230 if (compute_error && error_file_name != "")
231 {
233
234 const std::string error_file_format =
235 error_file_name.substr(error_file_name.find_last_of('.') + 1);
236
237 std::ofstream table_file(error_file_name);
238
239 if (error_file_format == "tex")
240 table.write_tex(table_file);
241 else if (error_file_format == "txt")
242 table.write_text(table_file);
243 else if (error_file_format == "gpl")
244 table.write_text(table_file,
246 else if (error_file_format == "org")
248 else
249 {
250 AssertThrow(false,
252 std::string("Unrecognized file format: ") +
253 error_file_format));
254 }
255 table_file.close();
256 }
257}
258
259
260
261void
263 const std::string &column_name,
264 const std::function<double()> &custom_function,
265 const bool compute_rate)
266{
267 extra_column_functions[column_name] = {custom_function, compute_rate};
268}
269
void evaluate_all_convergence_rates(const std::string &reference_column_key, const RateMode rate_mode)
void omit_column_from_convergence_rate_evaluation(const std::string &key)
void add_parameter(const std::string &entry, ParameterType &parameter, const std::string &documentation="", const Patterns::PatternBase &pattern= *Patterns::Tools::Convert< ParameterType >::to_pattern(), const bool has_to_be_set=false)
The ParsedConvergenceTable class.
void add_extra_column(const std::string &column_name, const std::function< double()> &custom_function, const bool compute_rate=true)
void add_parameters(ParameterHandler &prm)
std::set< std::string > extra_columns
ParsedConvergenceTable(const std::vector< std::string > &component_names={"u"}, const std::vector< std::set< VectorTools::NormType > > &list_of_error_norms={ {VectorTools::H1_norm, VectorTools::L2_norm, VectorTools::Linfty_norm}})
const std::vector< std::string > unique_component_names
std::vector< std::set< VectorTools::NormType > > norms_per_unique_component
std::map< std::string, std::pair< std::function< double()>, bool > > extra_column_functions
void write_text(std::ostream &out, const TextOutputFormat format=table_with_headers) const
@ table_with_separate_column_description
void write_tex(std::ostream &file, const bool with_header=true) const
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)