Reference documentation for deal.II version 9.6.0
\(\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.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 1999 - 2023 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
15#ifndef dealii_convergence_table_h
16#define dealii_convergence_table_h
17
18
19#include <deal.II/base/config.h>
20
22
24
25
64{
65public:
69 ConvergenceTable() = default;
70
90
133 void
134 evaluate_convergence_rates(const std::string &data_column_key,
135 const std::string &reference_column_key,
136 const RateMode rate_mode,
137 const unsigned int dim = 2);
138
139
156 void
157 evaluate_convergence_rates(const std::string &data_column_key,
158 const RateMode rate_mode);
159
167 void
168 omit_column_from_convergence_rate_evaluation(const std::string &key);
169
184 void
185 evaluate_all_convergence_rates(const std::string &reference_column_key,
186 const RateMode rate_mode);
187
201 void
203
213 std::string,
214 << "Rate column <" << arg1 << "> does already exist.");
216};
217
218
220
221#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:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
static ::ExceptionBase & ExcRateColumnAlreadyExists(std::string arg1)
#define DeclException1(Exception1, type1, outsequence)
Definition exceptions.h:516