Reference documentation for deal.II version 8.4.2
convergence_table.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2014 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:
64  break;
65  case reduction_rate:
66  rate_key += "red.rate";
67  no_rate_entries = columns[rate_key].entries.size();
68  // Calculate all missing rate values:
69  for (unsigned int i = no_rate_entries; i<n; ++i)
70  {
71  if (i == 0)
72  {
73  // no value available for the first row
74  add_value(rate_key, std::string("-"));
75  }
76  else
77  {
78  add_value(rate_key, values[i-1]/values[i] *
79  ref_values[i]/ref_values[i-1]);
80  }
81  }
82  break;
84  rate_key += "red.rate.log2";
85  no_rate_entries = columns[rate_key].entries.size();
86  // Calculate all missing rate values:
87  for (unsigned int i = no_rate_entries; i<n; ++i)
88  {
89  if (i == 0)
90  {
91  // no value available for the first row
92  add_value(rate_key, std::string("-"));
93  }
94  else
95  {
96  add_value(rate_key, dim*std::log(std::fabs(values[i-1]/values[i])) /
97  std::log(std::fabs(ref_values[i]/ref_values[i-1])));
98  }
99  }
100  break;
101  default:
102  AssertThrow(false, ExcNotImplemented());
103  }
104 
105  Assert(columns.count(rate_key), ExcInternalError());
106  columns[rate_key].flag = 1;
107  set_precision(rate_key, 2);
108 
109  std::string superkey = data_column_key;
110  if (!supercolumns.count(superkey))
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 
125 
126 void
127 ConvergenceTable::evaluate_convergence_rates(const std::string &data_column_key,
128  const RateMode rate_mode)
129 {
130  Assert(columns.count(data_column_key), ExcColumnNotExistent(data_column_key));
131 
132  // reset the auto fill mode flag since we are going to fill columns from
133  // the top that don't yet exist
134  set_auto_fill_mode (false);
135 
136  std::vector<internal::TableEntry> &entries = columns[data_column_key].entries;
137  std::string rate_key=data_column_key+"...";
138 
139  const unsigned int n=entries.size();
140 
141  std::vector<double> values(n);
142  for (unsigned int i=0; i<n; ++i)
143  values[i] = entries[i].get_numeric_value();
144 
145  unsigned int no_rate_entries = 0;
146 
147  switch (rate_mode)
148  {
149  case none:
150  break;
151 
152  case reduction_rate:
153  rate_key+="red.rate";
154  no_rate_entries = columns[rate_key].entries.size();
155  // Calculate all missing rate values:
156  for (unsigned int i = no_rate_entries; i<n; ++i)
157  {
158  if (i == 0)
159  {
160  // no value available for the first row
161  add_value(rate_key, std::string("-"));
162  }
163  else
164  {
165  add_value(rate_key, values[i-1]/values[i]);
166  }
167  }
168  break;
169 
170  case reduction_rate_log2:
171  rate_key+="red.rate.log2";
172  no_rate_entries = columns[rate_key].entries.size();
173  // Calculate all missing rate values:
174  for (unsigned int i = no_rate_entries; i<n; ++i)
175  {
176  if (i == 0)
177  {
178  // no value available for the first row
179  add_value(rate_key, std::string("-"));
180  }
181  else
182  {
183  add_value(rate_key, std::log(std::fabs(values[i-1]/values[i]))/std::log(2.0));
184  }
185  }
186  break;
187 
188  default:
189  AssertThrow(false, ExcNotImplemented());
190  }
191 
192  Assert(columns.count(rate_key), ExcInternalError());
193  columns[rate_key].flag=1;
194  set_precision(rate_key, 2);
195 
196  // set the superkey equal to the key
197  std::string superkey=data_column_key;
198  // and set the tex caption of the supercolumn to the tex caption of the
199  // data_column.
200  if (!supercolumns.count(superkey))
201  {
202  add_column_to_supercolumn(data_column_key, superkey);
203  set_tex_supercaption(superkey, columns[data_column_key].tex_caption);
204  }
205 
206  // only add rate_key to the supercolumn once
207  if (no_rate_entries == 0)
208  {
209  add_column_to_supercolumn(rate_key, superkey);
210  }
211 }
212 
213 
214 
215 void
217 {
218  Assert(columns.count(key), ExcColumnNotExistent(key));
219 
220  const std::map<std::string, Column>::iterator col_iter=columns.find(key);
221  col_iter->second.flag=1;
222 }
223 
224 
225 
226 void
227 ConvergenceTable::evaluate_all_convergence_rates(const std::string &reference_column_key,
228  const RateMode rate_mode)
229 {
230  for (std::map<std::string, Column>::const_iterator col_iter=columns.begin();
231  col_iter!=columns.end(); ++col_iter)
232  if (!col_iter->second.flag)
233  evaluate_convergence_rates(col_iter->first, reference_column_key, rate_mode);
234 }
235 
236 
237 
238 void
240 {
241  for (std::map<std::string, Column>::const_iterator col_iter=columns.begin();
242  col_iter!=columns.end(); ++col_iter)
243  if (!col_iter->second.flag)
244  evaluate_convergence_rates(col_iter->first, rate_mode);
245 }
246 
247 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:358
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:294
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)
std::map< std::string, Column > columns
std::map< std::string, std::vector< std::string > > supercolumns
void set_auto_fill_mode(const bool state)