Reference documentation for deal.II version 8.5.1
convergence_table.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2017 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/base/convergence_table.h>
17 #include <cmath>
18 
19 DEAL_II_NAMESPACE_OPEN
20 
22 {}
23 
24 
25 void ConvergenceTable::evaluate_convergence_rates(const std::string &data_column_key,
26  const std::string &reference_column_key,
27  const RateMode rate_mode,
28  const unsigned int dim)
29 {
30  Assert(columns.count(data_column_key),
31  ExcColumnNotExistent(data_column_key));
32  Assert(columns.count(reference_column_key),
33  ExcColumnNotExistent(reference_column_key));
34 
35  if (rate_mode == none)
36  return;
37 
38  // reset the auto fill mode flag since we are going to fill columns from
39  // the top that don't yet exist
40  set_auto_fill_mode (false);
41 
42  std::vector<internal::TableEntry> &entries=columns[data_column_key].entries;
43  std::vector<internal::TableEntry> &ref_entries=columns[reference_column_key].entries;
44  std::string rate_key = data_column_key+"...";
45 
46  const unsigned int n = entries.size();
47  const unsigned int n_ref = ref_entries.size();
48  Assert(n == n_ref, ExcDimensionMismatch(n, n_ref));
49 
50  std::vector<double> values(n);
51  std::vector<double> ref_values(n_ref);
52 
53  for (unsigned int i=0; i<n; ++i)
54  {
55  values[i] = entries[i].get_numeric_value();
56  ref_values[i] = ref_entries[i].get_numeric_value();
57  }
58 
59  unsigned int no_rate_entries = 0;
60 
61  switch (rate_mode)
62  {
63  // case none: already considered above
64  case reduction_rate:
65  rate_key += "red.rate";
66  no_rate_entries = columns[rate_key].entries.size();
67  // Calculate all missing rate values:
68  for (unsigned int i = no_rate_entries; i<n; ++i)
69  {
70  if (i == 0)
71  {
72  // no value available for the first row
73  add_value(rate_key, std::string("-"));
74  }
75  else
76  {
77  add_value(rate_key, values[i-1]/values[i] *
78  ref_values[i]/ref_values[i-1]);
79  }
80  }
81  break;
83  rate_key += "red.rate.log2";
84  no_rate_entries = columns[rate_key].entries.size();
85  // Calculate all missing rate values:
86  for (unsigned int i = no_rate_entries; i<n; ++i)
87  {
88  if (i == 0)
89  {
90  // no value available for the first row
91  add_value(rate_key, std::string("-"));
92  }
93  else
94  {
95  add_value(rate_key, dim*std::log(std::fabs(values[i-1]/values[i])) /
96  std::log(std::fabs(ref_values[i]/ref_values[i-1])));
97  }
98  }
99  break;
100  default:
101  AssertThrow(false, ExcNotImplemented());
102  }
103 
104  Assert(columns.count(rate_key), ExcInternalError());
105  columns[rate_key].flag = 1;
106  set_precision(rate_key, 2);
107 
108  std::string superkey = data_column_key;
109  if (!supercolumns.count(superkey))
110  {
111  add_column_to_supercolumn(data_column_key, superkey);
112  set_tex_supercaption(superkey, columns[data_column_key].tex_caption);
113  }
114 
115  // only add rate_key to the supercolumn once
116  if (no_rate_entries == 0)
117  {
118  add_column_to_supercolumn(rate_key, superkey);
119  }
120 
121 }
122 
123 
124 
125 void
126 ConvergenceTable::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 
169  case reduction_rate_log2:
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, std::log(std::fabs(values[i-1]/values[i]))/std::log(2.0));
183  }
184  }
185  break;
186 
187  default:
188  AssertThrow(false, ExcNotImplemented());
189  }
190 
191  Assert(columns.count(rate_key), ExcInternalError());
192  columns[rate_key].flag=1;
193  set_precision(rate_key, 2);
194 
195  // set the superkey equal to the key
196  std::string superkey=data_column_key;
197  // and set the tex caption of the supercolumn to the tex caption of the
198  // data_column.
199  if (!supercolumns.count(superkey))
200  {
201  add_column_to_supercolumn(data_column_key, superkey);
202  set_tex_supercaption(superkey, columns[data_column_key].tex_caption);
203  }
204 
205  // only add rate_key to the supercolumn once
206  if (no_rate_entries == 0)
207  {
208  add_column_to_supercolumn(rate_key, superkey);
209  }
210 }
211 
212 
213 
214 void
216 {
217  Assert(columns.count(key), ExcColumnNotExistent(key));
218 
219  const std::map<std::string, Column>::iterator col_iter=columns.find(key);
220  col_iter->second.flag=1;
221 }
222 
223 
224 
225 void
226 ConvergenceTable::evaluate_all_convergence_rates(const std::string &reference_column_key,
227  const RateMode rate_mode)
228 {
229  for (std::map<std::string, Column>::const_iterator col_iter=columns.begin();
230  col_iter!=columns.end(); ++col_iter)
231  if (!col_iter->second.flag)
232  evaluate_convergence_rates(col_iter->first, reference_column_key, rate_mode);
233 }
234 
235 
236 
237 void
239 {
240  for (std::map<std::string, Column>::const_iterator col_iter=columns.begin();
241  col_iter!=columns.end(); ++col_iter)
242  if (!col_iter->second.flag)
243  evaluate_convergence_rates(col_iter->first, rate_mode);
244 }
245 
246 DEAL_II_NAMESPACE_CLOSE
void set_precision(const std::string &key, const unsigned int precision)
void set_tex_supercaption(const std::string &superkey, const std::string &tex_supercaption)
void add_value(const std::string &key, const T value)
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
void add_column_to_supercolumn(const std::string &key, const std::string &superkey)
void omit_column_from_convergence_rate_evaluation(const std::string &key)
void evaluate_all_convergence_rates(const std::string &reference_column_key, const RateMode rate_mode)
#define Assert(cond, exc)
Definition: exceptions.h:313
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
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)
static ::ExceptionBase & ExcColumnNotExistent(std::string arg1)
static ::ExceptionBase & ExcNotImplemented()
std::map< std::string, Column > columns
std::map< std::string, std::vector< std::string > > supercolumns
static ::ExceptionBase & ExcInternalError()
void set_auto_fill_mode(const bool state)