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