Reference documentation for deal.II version 9.5.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
Public Member Functions | Private Member Functions | Private Attributes | List of all members
ParsedConvergenceTable Class Reference

The ParsedConvergenceTable class. More...

#include <deal.II/base/parsed_convergence_table.h>

Public Member Functions

 ParsedConvergenceTable (const std::vector< std::string > &component_names={"u"}, const std::vector< std::set< VectorTools::NormType > > &list_of_error_norms={ {VectorTools::H1_norm, VectorTools::L2_norm, VectorTools::Linfty_norm}})
 
 ParsedConvergenceTable (const std::vector< std::string > &component_names, const std::vector< std::set< VectorTools::NormType > > &list_of_error_norms, const double exponent, const std::set< std::string > &extra_columns, const std::string &rate_key, const std::string &rate_mode, const std::string &error_file_name, const unsigned int precision, const bool compute_error)
 
void add_parameters (ParameterHandler &prm)
 
template<int dim, int spacedim, typename VectorType >
void error_from_exact (const DoFHandler< dim, spacedim > &vspace, const VectorType &solution, const Function< spacedim > &exact, const Function< spacedim > *weight=nullptr)
 
template<int dim, int spacedim, typename VectorType >
void error_from_exact (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &vspace, const VectorType &solution, const Function< spacedim > &exact, const Function< spacedim > *weight=nullptr)
 
void add_extra_column (const std::string &column_name, const std::function< double()> &custom_function, const bool compute_rate=true)
 
template<int dim, int spacedim, typename VectorType >
void difference (const DoFHandler< dim, spacedim > &, const VectorType &, const VectorType &, const Function< spacedim > *weight=nullptr)
 
template<int dim, int spacedim, typename VectorType >
void difference (const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &, const VectorType &, const VectorType &, const Function< spacedim > *weight=nullptr)
 
void output_table (std::ostream &out)
 
void output_table ()
 

Private Member Functions

void prepare_table_for_output ()
 

Private Attributes

const std::vector< std::string > component_names
 
const std::vector< std::string > unique_component_names
 
const std::vector< ComponentMaskunique_component_masks
 
std::map< std::string, std::pair< std::function< double()>, bool > > extra_column_functions
 
std::vector< std::set< VectorTools::NormType > > norms_per_unique_component
 
double exponent
 
ConvergenceTable table
 
std::set< std::string > extra_columns
 
std::string rate_key
 
std::string rate_mode
 
unsigned int precision
 
std::string error_file_name
 
bool compute_error
 

Detailed Description

The ParsedConvergenceTable class.

This class simplifies the construction of convergence tables, reading the options for the generation of the table from a parameter file. It provides a series of methods that can be used to compute the error given a reference exact solution, or the difference between two numerical solutions, or any other custom computation of the error, given via std::function objects.

An example usage of this class is given by

table.add_parameters(prm);
for (unsigned int i = 0; i < n_cycles; ++i)
{
... // do some computations
table.error_from_exact(dof_handler, solution, exact_solution);
}
table.output_table(std::cout);
The ParsedConvergenceTable class.

The above code constructs a ParsedConvergenceTable that works for scalar problems, and will produce an error table with H1_norm, L2_norm, and Linfty_norm norms of the error.

Whenever a call to the methods error_from_exact() or difference() is made, the instance of this class inspects its parameters, computes all norms specified by the parameter given at construction time, possibly modified via a parameter file, computes all extra column entries specified using the method add_extra_column(), and writes one row of the convergence table.

Once you have finished with the computations, a call to output_table() will generate a formatted convergence table on the provided stream, and to the file (if any) specified in the parameter file.

With a small modification, the same code can be used to estimate the errors of mixed or multi-physics problems, e.g.:

using namespace VectorTools;
ParsedConvergenceTable table({"u,u,p"},{{H1_norm, L2_norm}, {L2_norm}});
table.add_parameters(prm);
for (unsigned int i = 0; i < n_cycles; ++i)
{
... // do some computations
table.error_from_exact(dof_handler, solution, exact_solution);
}
table.output_table(std::cout);

The above code assumes that you are solving a Stokes problem with three components. Two components for the vector velocity field u, and one component for the pressure field p, and will produce an error table with H1 and L2 norm of the error in the velocity field (first two, components) and L2 error in the pressure field.

You may also call table.output_table() without arguments, to write the table only to the file specified in the parameter file.

By calling the method add_parameters() passing a ParameterHandler object, the following options will be defined in the given ParameterHandler object (in the current level of the ParameterHandler object, i.e., whatever level you have entered with the ParameterHandler::enter_subsection() method), and can be modified at run time through a parameter file :

set Enable computation of the errors = true
set Error file name =
set Error precision = 3
set Exponent for p-norms = 2
set Extra columns = dofs, cells
set List of error norms to compute = Linfty_norm, L2_norm, H1_norm
set Rate key = dofs
set Rate mode = reduction_rate_log2

When using this class, please cite

@article{SartoriGiulianiBardelloni-2018-a,
Author = {Sartori, Alberto and Giuliani, Nicola and
Bardelloni, Mauro and Heltai, Luca},
Journal = {SoftwareX},
Pages = {318--327},
Title = {{deal2lkit: A toolkit library for high performance
programming in deal.II}},
Doi = {10.1016/j.softx.2018.09.004},
Volume = {7},
Year = {2018}}

Definition at line 136 of file parsed_convergence_table.h.

Constructor & Destructor Documentation

◆ ParsedConvergenceTable() [1/2]

ParsedConvergenceTable::ParsedConvergenceTable ( const std::vector< std::string > &  component_names = {"u"},
const std::vector< std::set< VectorTools::NormType > > &  list_of_error_norms = { {VectorTools::H1_normVectorTools::L2_normVectorTools::Linfty_norm}} 
)

Minimal constructor for ParsedConvergenceTable objects.

The number of components must match the number of components of the finite element space that is used to compute the errors. If a component name is repeated, than it is interpreted as a vector field, and the errors of the repeated components are grouped together.

The size of the vector list_of_error_norms must match the number of unique component names, and may contain zero or more comma separated identifiers for the norm to compute for each component (see the documentation of VectorTools::NormType for the available options).

For example, the following constructor

using namespace VectorTools;
ParsedConvergenceTable table({"u", "v", "v"},
{{Linfty_norm}, {L2_norm, H1_norm}});

would produce (if the parameter file is left untouched) a table similar to

cells dofs u_Linfty_norm v_L2_norm v_H1_norm
4 9 1.183e-01 - 5.156e-02 - 2.615e-01 -
16 25 3.291e-02 2.50 1.333e-02 2.65 1.272e-01 1.41
64 81 8.449e-03 2.31 3.360e-03 2.34 6.313e-02 1.19
256 289 2.126e-03 2.17 8.418e-04 2.18 3.150e-02 1.09
1024 1089 5.325e-04 2.09 2.106e-04 2.09 1.574e-02 1.05

See the other constructor for a documentation of all the parameters you can change.

Parameters
component_namesSpecify the names of the components;
list_of_error_normsSpecify what error norms to compute for each unique component name.

Definition at line 63 of file parsed_convergence_table.cc.

◆ ParsedConvergenceTable() [2/2]

ParsedConvergenceTable::ParsedConvergenceTable ( const std::vector< std::string > &  component_names,
const std::vector< std::set< VectorTools::NormType > > &  list_of_error_norms,
const double  exponent,
const std::set< std::string > &  extra_columns,
const std::string &  rate_key,
const std::string &  rate_mode,
const std::string &  error_file_name,
const unsigned int  precision,
const bool  compute_error 
)

Full constructor for ParsedConvergenceTable.

Parameters
component_namesNames of the components. Repeated consecutive names are interpreted as components of a vector valued field;
list_of_error_normsSpecify what error norms to compute for each unique component name;
exponentThe exponent to use in p-norms;
extra_columnsExtra columns to add. These may be "cells" or "dofs";
rate_keySpecify the extra column by which we will compute the error rates. This key can either be one of "cells" or "dofs", or, if you add extra columns to the table via the method add_extra_column(), it may be one of the extra columns you added;
rate_modeSpecify the rate mode to use when computing error rates. This maybe either "reduction_rate", "reduction_rate_log2", or "none". See the documentation of ConvergenceTable::RateMode for an explanation of how each of this mode behaves;
error_file_nameName of error output file (with extension txt, gpl, tex, or org). If different from the empty string, than output_table() also writes in this file in the format deduced from its extension;
precisionHow many digits to use when writing the error;
compute_errorControl whether the filling of the table is enabled or not. This flag may be used to disable at run time any error computation;

The parameters you specify with this constructor can be written to a ParameterHandler object by calling the add_parameters() method. Once you call the add_parameters() method, the following options will be defined in the given ParameterHandler object, and the parameters of the instance of this class will follow the modification you make to the ParameterHandler object at run time:

# Listing of Parameters
# ---------------------
# When set to false, no computations are performed.
set Enable computation of the errors = true
# Set this to a filename with extension .txt, .gpl, .org, or .tex to enable
# writing the convergence table to a file.
set Error file name =
# Number of digits to use when printing the error.
set Error precision = 3
# Extra columns to add to the table. Available options are dofs and cells.
set Extra columns = dofs, cells
# The exponent to use when computing p-norms.
set Exponent for p-norms = 2
# Each component is separated by a semicolon and each norm by a comma. See
# the documentation of VectorTools::NormType for a list of implemented
# norms. If you want to skip a component, leave its entry empty.
set List of error norms to compute = Linfty_norm, L2_norm, H1_norm
# Key to use when computing convergence rates. If this is set to a
# column that is not present, or to none, then no error rates are computed.
set Rate key = dofs
# What type of error rate to compute. Available options are
# reduction_rate_log2, reduction_rate, and none.
set Rate mode = reduction_rate_log2

Definition at line 79 of file parsed_convergence_table.cc.

Member Function Documentation

◆ add_parameters()

void ParsedConvergenceTable::add_parameters ( ParameterHandler prm)

Attach all the parameters in this class to entries of the parameter handler prm. Whenever the content of prm changes, the parameters of this class will be updated.

Definition at line 108 of file parsed_convergence_table.cc.

◆ error_from_exact() [1/2]

template<int dim, int spacedim, typename VectorType >
void ParsedConvergenceTable::error_from_exact ( const DoFHandler< dim, spacedim > &  vspace,
const VectorType &  solution,
const Function< spacedim > &  exact,
const Function< spacedim > *  weight = nullptr 
)

Add a row to the error table, containing the error between solution and the exact function, in the norm(s) specified in the parameter file.

If you specify a weight function during this call, then this is used to compute weighted errors. The weight function can be either a scalar function (which will be used for all components), or a vector function. When it is a vector function, an assertion is triggered if the number of components does not coincide with the number of components of the underlying finite element space.

◆ error_from_exact() [2/2]

template<int dim, int spacedim, typename VectorType >
void ParsedConvergenceTable::error_from_exact ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  vspace,
const VectorType &  solution,
const Function< spacedim > &  exact,
const Function< spacedim > *  weight = nullptr 
)

Same as above, with a different mapping.

◆ add_extra_column()

void ParsedConvergenceTable::add_extra_column ( const std::string &  column_name,
const std::function< double()> &  custom_function,
const bool  compute_rate = true 
)

Add an additional column (with name column_name) to the table, by invoking the function custom_function, when calling error_from_exact() or difference().

You can call this method as many times as you want. If column_name was already used in a previous call, then calling this method with the same name will overwrite whatever function you had previously specified. If you use a lambda function for this call, make sure that the variables used internally in the lambda function remain valid until the call to error_from_exact() or difference().

Make sure you add all extra columns before the first call to error_from_exact() or difference(). Adding additional columns to the convergence table after you already started filling the table will trigger an exception.

This method may be used, for example, to compute the error w.r.t. to time step increments in time, for example:

using namespace VectorTools;
ParsedConvergenceTable table({"u"}, {{L2_norm}});
double dt = .5;
auto dt_function = [&]() {
return dt;
};
table.add_extra_column("dt", dt_function, false);
for (unsigned int i = 0; i < n_cycles; ++i)
{
// ... compute solution at the current dt
table.error_from_exact(dof_handler, solution, exact_solution);
dt /= 2.0;
}
table.output_table(std::cout);

will produce a table similar to

dt u_L2_norm
5.000e-1 5.156e-02 -
2.500e-2 1.333e-02 2.65
1.250e-2 3.360e-03 2.34
6.250e-3 8.418e-04 2.18

provided that you use the following parameter file (only non default entries are shown here):

set Extra columns =
set List of error norms to compute = L2_norm
set Rate key = dt
Parameters
column_nameName of the column to add;
custom_functionFunction that will be called to fill the given entry. You need to make sure that the scope of this function is valid up to the call to error_from_exact() or difference();
compute_rateIf set to true, then this column will be included in the list of columns for which error rates are computed. You may want to set this to false if you want to compute error rates with respect to this column. In this case, you should also specify column_name as the rate key in the parameter file.

Definition at line 263 of file parsed_convergence_table.cc.

◆ difference() [1/2]

template<int dim, int spacedim, typename VectorType >
void ParsedConvergenceTable::difference ( const DoFHandler< dim, spacedim > &  ,
const VectorType &  ,
const VectorType &  ,
const Function< spacedim > *  weight = nullptr 
)

Difference between two solutions in the same vector space.

◆ difference() [2/2]

template<int dim, int spacedim, typename VectorType >
void ParsedConvergenceTable::difference ( const Mapping< dim, spacedim > &  mapping,
const DoFHandler< dim, spacedim > &  ,
const VectorType &  ,
const VectorType &  ,
const Function< spacedim > *  weight = nullptr 
)

Same as above, with a non default mapping.

◆ output_table() [1/2]

void ParsedConvergenceTable::output_table ( std::ostream &  out)

Write the error table to the out stream (in text format), and (possibly) to the file stream specified in the parameters (with the format deduced from the file name extension).

Definition at line 216 of file parsed_convergence_table.cc.

◆ output_table() [2/2]

void ParsedConvergenceTable::output_table ( )

Write the error table to the file stream specified in the parameters.

If the "Error file name" option in the parameter file is set to the empty string, no output is written.

Definition at line 229 of file parsed_convergence_table.cc.

◆ prepare_table_for_output()

void ParsedConvergenceTable::prepare_table_for_output ( )
private

Add rates to the output table.

Definition at line 162 of file parsed_convergence_table.cc.

Member Data Documentation

◆ component_names

const std::vector<std::string> ParsedConvergenceTable::component_names
private

Names of the solution components.

Definition at line 411 of file parsed_convergence_table.h.

◆ unique_component_names

const std::vector<std::string> ParsedConvergenceTable::unique_component_names
private

Same as above, but containing repeated component names only once.

Definition at line 416 of file parsed_convergence_table.h.

◆ unique_component_masks

const std::vector<ComponentMask> ParsedConvergenceTable::unique_component_masks
private

Masks for each unique component name.

Definition at line 421 of file parsed_convergence_table.h.

◆ extra_column_functions

std::map<std::string, std::pair<std::function<double()>, bool> > ParsedConvergenceTable::extra_column_functions
private

Additional methods to call when adding rows to the table.

Definition at line 427 of file parsed_convergence_table.h.

◆ norms_per_unique_component

std::vector<std::set<VectorTools::NormType> > ParsedConvergenceTable::norms_per_unique_component
private

Type of error to compute per components.

Definition at line 432 of file parsed_convergence_table.h.

◆ exponent

double ParsedConvergenceTable::exponent
private

Exponent to use in p-norm types.

Definition at line 437 of file parsed_convergence_table.h.

◆ table

ConvergenceTable ParsedConvergenceTable::table
private

The actual table

Definition at line 442 of file parsed_convergence_table.h.

◆ extra_columns

std::set<std::string> ParsedConvergenceTable::extra_columns
private

Extra columns to add to the table.

Definition at line 447 of file parsed_convergence_table.h.

◆ rate_key

std::string ParsedConvergenceTable::rate_key
private

The name of column with respect to which we compute convergence rates.

Definition at line 452 of file parsed_convergence_table.h.

◆ rate_mode

std::string ParsedConvergenceTable::rate_mode
private

Reduction rate mode. See ConvergenceTable::RateMode for a documentation.

Definition at line 457 of file parsed_convergence_table.h.

◆ precision

unsigned int ParsedConvergenceTable::precision
private

The precision used to output the table.

Definition at line 462 of file parsed_convergence_table.h.

◆ error_file_name

std::string ParsedConvergenceTable::error_file_name
private

Filename to use when writing to file.

Definition at line 467 of file parsed_convergence_table.h.

◆ compute_error

bool ParsedConvergenceTable::compute_error
private

Compute the error. If this is false, all methods that perform the computation of the error are disabled and don't do anything.

Definition at line 473 of file parsed_convergence_table.h.


The documentation for this class was generated from the following files: