Reference documentation for deal.II version 9.1.1
|
#include <deal.II/base/convergence_table.h>
Public Types | |
enum | RateMode { none, reduction_rate, reduction_rate_log2 } |
Public Types inherited from TableHandler | |
enum | TextOutputFormat { table_with_headers, table_with_separate_column_description, simple_table_with_separate_column_description, org_mode_table } |
Public Member Functions | |
ConvergenceTable ()=default | |
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 | evaluate_convergence_rates (const std::string &data_column_key, const RateMode rate_mode) |
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) |
void | evaluate_all_convergence_rates (const RateMode rate_mode) |
Public Member Functions inherited from TableHandler | |
TableHandler () | |
void | declare_column (const std::string &key) |
template<typename T > | |
void | add_value (const std::string &key, const T value) |
void | start_new_row () |
void | set_auto_fill_mode (const bool state) |
void | add_column_to_supercolumn (const std::string &key, const std::string &superkey) |
void | set_column_order (const std::vector< std::string > &new_order) |
void | set_precision (const std::string &key, const unsigned int precision) |
void | set_scientific (const std::string &key, const bool scientific) |
void | set_tex_caption (const std::string &key, const std::string &tex_caption) |
void | set_tex_table_caption (const std::string &table_caption) |
void | set_tex_table_label (const std::string &table_label) |
void | set_tex_supercaption (const std::string &superkey, const std::string &tex_supercaption) |
void | set_tex_format (const std::string &key, const std::string &format="c") |
void | write_text (std::ostream &out, const TextOutputFormat format=table_with_headers) const |
void | write_tex (std::ostream &file, const bool with_header=true) const |
void | clear () |
void | clear_current_row () |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Static Public Member Functions | |
static ::ExceptionBase & | ExcRateColumnAlreadyExists (std::string arg1) |
Static Public Member Functions inherited from TableHandler | |
static ::ExceptionBase & | ExcColumnNotExistent (std::string arg1) |
static ::ExceptionBase & | ExcSuperColumnNotExistent (std::string arg1) |
static ::ExceptionBase & | ExcColumnOrSuperColumnNotExistent (std::string arg1) |
static ::ExceptionBase & | ExcWrongNumberOfDataEntries (std::string arg1, int arg2, std::string arg3, int arg4) |
static ::ExceptionBase & | ExcUndefinedTexFormat (std::string arg1) |
Additional Inherited Members | |
Protected Member Functions inherited from TableHandler | |
void | get_selected_columns (std::vector< std::string > &sel_columns) const |
unsigned int | n_rows () const |
Protected Attributes inherited from TableHandler | |
std::vector< std::string > | column_order |
std::map< std::string, Column > | columns |
std::map< std::string, std::vector< std::string > > | supercolumns |
std::map< std::string, std::string > | tex_supercaptions |
std::string | tex_table_caption |
std::string | tex_table_label |
bool | auto_fill_mode |
The ConvergenceTable class is an application to the TableHandler class and stores some convergence data, such as residuals of the cg-method, or some evaluated L2-errors of discrete solutions, etc, and evaluates convergence rates or orders.
The already implemented RateMode's are reduction_rate, where the convergence rate is the quotient of two following rows, and reduction_rate_log2, that evaluates the order of convergence. These standard evaluations are useful for global refinement, for local refinement this may not be an appropriate method, as the convergence rates should be set in relation to the number of cells or the number of DoFs. The implementations of these non-standard methods is left to a user.
The number of cells and the number of DoFs may be added to the table by calling e.g. add_value("n cells", n_cells)
. The table data is also added by calling add_value(). Before the output of the table the functions evaluate_convergence_rates() and evaluate_all_convergence_rates() may be called.
There are two possibilities of how to evaluate the convergence rates of multiple columns in the same RateMode.
A detailed discussion of this class can also be found in the step-7 and step-13 example programs.
Definition at line 63 of file convergence_table.h.
Rate in relation to the rows.
Enumerator | |
---|---|
none | Do not do anything. |
reduction_rate | Quotient of values in the previous row and in this row. |
reduction_rate_log2 | Logarithm of reduction_rate to the base 2 representing the order of convergence when halving the grid size, e.g. from h to h/2. |
Definition at line 74 of file convergence_table.h.
|
default |
Constructor.
void ConvergenceTable::evaluate_convergence_rates | ( | const std::string & | data_column_key, |
const std::string & | reference_column_key, | ||
const RateMode | rate_mode, | ||
const unsigned int | dim = 2 |
||
) |
Evaluate the convergence rates of the data column data_column_key
due to the RateMode in relation to the reference column reference_column_key
. Be sure that the value types of the table entries of the data column and the reference data column is a number, i.e. double, float, (unsigned) int, and so on.
As this class has no information on the space dimension upon which the reference column vs. the value column is based upon, it needs to be passed as last argument to this method. The default dimension for the reference column is 2, which is appropriate for the number of cells in 2D. If you work in 3D, set the number to 3. If the reference column is \(1/h\), remember to set the dimension to 1 also when working in 3D to get correct rates.
The new rate column and the data column will be merged to a supercolumn. The tex caption of the supercolumn will be (by default) the same as the one of the data column. This may be changed by using the set_tex_supercaption (...)
function of the base class TableHandler.
This method behaves in the following way:
If RateMode is reduction_rate, then the computed output is \( \frac{e_{n-1}/k_{n-1}}{e_n/k_n}, \) where \(k\) is the reference column (no dimension dependence!).
If RateMode is reduction_rate_log2, then the computed output is \( dim \frac{\log |e_{n-1}/e_{n}|}{\log |k_n/k_{n-1}|} \).
This is useful, for example, if we use as reference key the number of degrees of freedom or better, the number of cells. Assuming that the error is proportional to \( C (1/\sqrt{k})^r \) in 2D, then this method will produce the rate \(r\) as a result. For general dimension, as described by the last parameter of this function, the formula needs to be \( C (1/\sqrt[dim]{k})^r \).
Definition at line 23 of file convergence_table.cc.
void ConvergenceTable::evaluate_convergence_rates | ( | const std::string & | data_column_key, |
const RateMode | rate_mode | ||
) |
Evaluate the convergence rates of the data column data_column_key
due to the RateMode. Be sure that the value types of the table entries of the data column is a number, i.e. double, float, (unsigned) int, and so on.
The new rate column and the data column will be merged to a supercolumn. The tex caption of the supercolumn will be (by default) the same as the one of the data column. This may be changed by using the set_tex_supercaption() function of the base class TableHandler.
Definition at line 127 of file convergence_table.cc.
void ConvergenceTable::omit_column_from_convergence_rate_evaluation | ( | const std::string & | key | ) |
Omit this column key
(not supercolumn!) from the evaluation of the convergence rates of ‘all’ columns (see the following two functions).
The Column::flag==1 is reserved for omitting the column from convergence rate evaluation.
Definition at line 218 of file convergence_table.cc.
void ConvergenceTable::evaluate_all_convergence_rates | ( | const std::string & | reference_column_key, |
const RateMode | rate_mode | ||
) |
Evaluate convergence rates due to the rate_mode
in relation to the reference column reference_column_key
. This function evaluates the rates of ALL columns except of the columns that are to be omitted (see previous function) and except of the columns that are previously evaluated rate columns. This function allows to evaluate the convergence rate for almost all columns of a table without calling evaluate_convergence_rates() for each column separately.
Example: Columns like n cells
or n dofs
columns may be wanted to be omitted in the evaluation of the convergence rates. Hence they should omitted by calling the omit_column_from_convergence_rate_evaluation().
Definition at line 230 of file convergence_table.cc.
void ConvergenceTable::evaluate_all_convergence_rates | ( | const RateMode | rate_mode | ) |
Evaluate convergence rates due to the rate_mode
. This function evaluates the rates of ALL columns except of the columns that are to be omitted (see previous function) and except of the columns that are previously evaluated rate columns. This function allows to evaluate the convergence rate for almost all columns of a table without calling evaluate_convergence_rates() for each column separately.
Example: Columns like n cells
or n dofs
columns may be wanted to be omitted in the evaluation of the convergence rates. Hence they should omitted by calling the omit_column_from_convergence_rate_evaluation().
Definition at line 246 of file convergence_table.cc.