Reference documentation for deal.II version GIT 00e6fe71c9 2023-06-04 19:35:01+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\}}\)
convergence_table.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2021 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.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_convergence_table_h
17 #define dealii_convergence_table_h
18 
19 
20 #include <deal.II/base/config.h>
21 
23 
25 
26 
65 {
66 public:
70  ConvergenceTable() = default;
71 
75  enum RateMode
76  {
90  };
91 
134  void
135  evaluate_convergence_rates(const std::string &data_column_key,
136  const std::string &reference_column_key,
137  const RateMode rate_mode,
138  const unsigned int dim = 2);
139 
140 
157  void
158  evaluate_convergence_rates(const std::string &data_column_key,
159  const RateMode rate_mode);
160 
168  void
169  omit_column_from_convergence_rate_evaluation(const std::string &key);
170 
185  void
186  evaluate_all_convergence_rates(const std::string &reference_column_key,
187  const RateMode rate_mode);
188 
202  void
203  evaluate_all_convergence_rates(const RateMode rate_mode);
204 
214  std::string,
215  << "Rate column <" << arg1 << "> does already exist.");
217 };
218 
219 
221 
222 #endif
ConvergenceTable()=default
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)
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:474
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:475
static ::ExceptionBase & ExcRateColumnAlreadyExists(std::string arg1)
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:510