Reference documentation for deal.II version GIT relicensing-233-g802318d791 2024-03-28 20:20:02+00:00
\(\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
convergence_table.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 1999 - 2022 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
16
17#include <cmath>
18
20
21void
23 const std::string &data_column_key,
24 const std::string &reference_column_key,
25 const RateMode rate_mode,
26 const unsigned int dim)
27{
28 Assert(columns.count(data_column_key), ExcColumnNotExistent(data_column_key));
29 Assert(columns.count(reference_column_key),
30 ExcColumnNotExistent(reference_column_key));
31
32 if (rate_mode == none)
33 return;
34
35 // reset the auto fill mode flag since we are going to fill columns from
36 // the top that don't yet exist
37 set_auto_fill_mode(false);
38
39 std::vector<internal::TableEntry> &entries = columns[data_column_key].entries;
40 std::vector<internal::TableEntry> &ref_entries =
41 columns[reference_column_key].entries;
42 std::string rate_key = data_column_key + "...";
43
44 const unsigned int n = entries.size();
45 const unsigned int n_ref = ref_entries.size();
46 Assert(n == n_ref, ExcDimensionMismatch(n, n_ref));
47
48 std::vector<double> values(n);
49 std::vector<double> ref_values(n_ref);
50
51 for (unsigned int i = 0; i < n; ++i)
52 {
53 values[i] = entries[i].get_numeric_value();
54 ref_values[i] = ref_entries[i].get_numeric_value();
55 }
56
57 unsigned int no_rate_entries = 0;
58
59 switch (rate_mode)
60 {
61 // case none: already considered above
62 case reduction_rate:
63 rate_key += "red.rate";
64 no_rate_entries = columns[rate_key].entries.size();
65 // Calculate all missing rate values:
66 for (unsigned int i = no_rate_entries; i < n; ++i)
67 {
68 if (i == 0)
69 {
70 // no value available for the first row
71 add_value(rate_key, std::string("-"));
72 }
73 else
74 {
75 add_value(rate_key,
76 values[i - 1] / values[i] * ref_values[i] /
77 ref_values[i - 1]);
78 }
79 }
80 break;
82 rate_key += "red.rate.log2";
83 no_rate_entries = columns[rate_key].entries.size();
84 // Calculate all missing rate values:
85 for (unsigned int i = no_rate_entries; i < n; ++i)
86 {
87 if (i == 0)
88 {
89 // no value available for the first row
90 add_value(rate_key, std::string("-"));
91 }
92 else
93 {
94 add_value(rate_key,
95 dim * std::log(std::fabs(values[i - 1] / values[i])) /
97 std::fabs(ref_values[i] / ref_values[i - 1])));
98 }
99 }
100 break;
101 default:
103 }
104
105 Assert(columns.count(rate_key), ExcInternalError());
106 columns[rate_key].flag = 1;
107 set_precision(rate_key, 2);
108
109 const std::string &superkey = data_column_key;
110 if (supercolumns.count(superkey) == 0u)
111 {
112 add_column_to_supercolumn(data_column_key, superkey);
113 set_tex_supercaption(superkey, columns[data_column_key].tex_caption);
114 }
115
116 // only add rate_key to the supercolumn once
117 if (no_rate_entries == 0)
118 {
119 add_column_to_supercolumn(rate_key, superkey);
120 }
121}
122
123
124
125void
126ConvergenceTable::evaluate_convergence_rates(const std::string &data_column_key,
127 const RateMode rate_mode)
128{
129 Assert(columns.count(data_column_key), ExcColumnNotExistent(data_column_key));
130
131 // reset the auto fill mode flag since we are going to fill columns from
132 // the top that don't yet exist
133 set_auto_fill_mode(false);
134
135 std::vector<internal::TableEntry> &entries = columns[data_column_key].entries;
136 std::string rate_key = data_column_key + "...";
137
138 const unsigned int n = entries.size();
139
140 std::vector<double> values(n);
141 for (unsigned int i = 0; i < n; ++i)
142 values[i] = entries[i].get_numeric_value();
143
144 unsigned int no_rate_entries = 0;
145
146 switch (rate_mode)
147 {
148 case none:
149 break;
150
151 case reduction_rate:
152 rate_key += "red.rate";
153 no_rate_entries = columns[rate_key].entries.size();
154 // Calculate all missing rate values:
155 for (unsigned int i = no_rate_entries; i < n; ++i)
156 {
157 if (i == 0)
158 {
159 // no value available for the first row
160 add_value(rate_key, std::string("-"));
161 }
162 else
163 {
164 add_value(rate_key, values[i - 1] / values[i]);
165 }
166 }
167 break;
168
170 rate_key += "red.rate.log2";
171 no_rate_entries = columns[rate_key].entries.size();
172 // Calculate all missing rate values:
173 for (unsigned int i = no_rate_entries; i < n; ++i)
174 {
175 if (i == 0)
176 {
177 // no value available for the first row
178 add_value(rate_key, std::string("-"));
179 }
180 else
181 {
182 add_value(rate_key,
183 std::log(std::fabs(values[i - 1] / values[i])) /
184 std::log(2.0));
185 }
186 }
187 break;
188
189 default:
191 }
192
193 Assert(columns.count(rate_key), ExcInternalError());
194 columns[rate_key].flag = 1;
195 set_precision(rate_key, 2);
196
197 // set the superkey equal to the key
198 const std::string &superkey = data_column_key;
199 // and set the tex caption of the supercolumn to the tex caption of the
200 // data_column.
201 if (supercolumns.count(superkey) == 0u)
202 {
203 add_column_to_supercolumn(data_column_key, superkey);
204 set_tex_supercaption(superkey, columns[data_column_key].tex_caption);
205 }
206
207 // only add rate_key to the supercolumn once
208 if (no_rate_entries == 0)
209 {
210 add_column_to_supercolumn(rate_key, superkey);
211 }
212}
213
214
215
216void
218 const std::string &key)
219{
220 Assert(columns.count(key), ExcColumnNotExistent(key));
221
222 const std::map<std::string, Column>::iterator col_iter = columns.find(key);
223 col_iter->second.flag = 1;
224}
225
226
227
228void
230 const std::string &reference_column_key,
231 const RateMode rate_mode)
232{
233 for (std::map<std::string, Column>::const_iterator col_iter = columns.begin();
234 col_iter != columns.end();
235 ++col_iter)
236 if (col_iter->second.flag == 0u)
237 evaluate_convergence_rates(col_iter->first,
238 reference_column_key,
239 rate_mode);
240}
241
242
243
244void
246{
247 for (std::map<std::string, Column>::const_iterator col_iter = columns.begin();
248 col_iter != columns.end();
249 ++col_iter)
250 if (col_iter->second.flag == 0u)
251 evaluate_convergence_rates(col_iter->first, rate_mode);
252}
253
void evaluate_all_convergence_rates(const std::string &reference_column_key, const RateMode rate_mode)
void evaluate_convergence_rates(const std::string &data_column_key, const std::string &reference_column_key, const RateMode rate_mode, const unsigned int dim=2)
void omit_column_from_convergence_rate_evaluation(const std::string &key)
void set_tex_supercaption(const std::string &superkey, const std::string &tex_supercaption)
void set_auto_fill_mode(const bool state)
void add_value(const std::string &key, const T value)
std::map< std::string, std::vector< std::string > > supercolumns
void add_column_to_supercolumn(const std::string &key, const std::string &superkey)
void set_precision(const std::string &key, const unsigned int precision)
std::map< std::string, Column > columns
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
static ::ExceptionBase & ExcColumnNotExistent(std::string arg1)
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
#define AssertThrow(cond, exc)
::VectorizedArray< Number, width > log(const ::VectorizedArray< Number, width > &)