ExaDG
Loading...
Searching...
No Matches
error_calculation.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2021 by the ExaDG authors
6 *
7 * This program is free software: you can redistribute it and/or modify
8 * it under the terms of the GNU General Public License as published by
9 * the Free Software Foundation, either version 3 of the License, or
10 * (at your option) any later version.
11 *
12 * This program is distributed in the hope that it will be useful,
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 * GNU General Public License for more details.
16 *
17 * You should have received a copy of the GNU General Public License
18 * along with this program. If not, see <https://www.gnu.org/licenses/>.
19 * ______________________________________________________________________
20 */
21
22#ifndef EXADG_POSTPROCESSOR_ERROR_CALCULATION_H_
23#define EXADG_POSTPROCESSOR_ERROR_CALCULATION_H_
24
25// deal.II
26#include <deal.II/base/function.h>
27#include <deal.II/dofs/dof_handler.h>
28#include <deal.II/fe/mapping_q.h>
29#include <deal.II/lac/la_parallel_vector.h>
30
31// ExaDG
32#include <exadg/postprocessor/time_control.h>
33#include <exadg/utilities/print_functions.h>
34
35namespace ExaDG
36{
37template<int dim>
38struct ErrorCalculationData
39{
40 ErrorCalculationData()
41 : calculate_relative_errors(true),
42 calculate_H1_seminorm_error(false),
43 write_errors_to_file(false),
44 compute_convergence_table(false),
45 spatially_weight_error(false),
46 weight(nullptr),
47 additional_quadrature_points(3),
48 directory("output/"),
49 name("all fields")
50 {
51 }
52
53 void
54 print(dealii::ConditionalOStream & pcout, bool unsteady)
55 {
56 print_parameter(pcout, "Error calculation", unsteady == true and analytical_solution);
57 if(unsteady == true and time_control_data.is_active)
58 {
59 print(pcout, unsteady, time_control_data);
60 print_parameter(pcout, "Calculate relative errors", calculate_relative_errors);
61 print_parameter(pcout, "Calculate H1-seminorm error", calculate_H1_seminorm_error);
62 print_parameter(pcout, "Write errors to file", write_errors_to_file);
63 if(write_errors_to_file)
64 print_parameter(pcout, "Directory", directory);
65 print_parameter(pcout, "Name", name);
66 }
67 }
68
69 std::shared_ptr<dealii::Function<dim>> analytical_solution;
70
71 // relative or absolute errors?
72 // If calculate_relative_errors == false, this implies that absolute errors are calculated
73 bool calculate_relative_errors;
74
75 // by default, only the L2-error is computed. Other norms have to be explicitly specified by the
76 // user.
77 bool calculate_H1_seminorm_error;
78
79 // data used to control the output
80 TimeControlData time_control_data;
81
82 // write errors to file?
83 bool write_errors_to_file;
84
85 // process computed files generating a convergence table
86 bool compute_convergence_table;
87
88 // If true, a spatially weighted norm is computed.
89 bool spatially_weight_error;
90 // Weight used to compute spatially weighted error.
91 std::shared_ptr<dealii::Function<dim>> weight;
92
93 // Number of quadrature points in 1D: fe_degree + additional_quadrature_points
94 unsigned int additional_quadrature_points;
95
96 // directory and name (used as filename and as identifier for screen output)
97 std::string directory;
98 std::string name;
99};
100
101template<int dim, typename Number>
102class ErrorCalculator
103{
104public:
105 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
106
107 ErrorCalculator(MPI_Comm const & comm);
108
109 void
110 setup(dealii::DoFHandler<dim> const & dof_handler,
111 dealii::Mapping<dim> const & mapping,
112 ErrorCalculationData<dim> const & error_data);
113
114 void
115 evaluate(VectorType const & solution, double const time, bool const unsteady);
116
117 TimeControl time_control;
118
119private:
120 void
121 do_evaluate(VectorType const & solution_vector, double const time);
122
123 std::string
124 filename_from_filename_base(std::string const & directory,
125 std::string const & filename_base,
126 bool const initial_call);
127
128 MPI_Comm const mpi_comm;
129
130 bool clear_files_L2, clear_files_H1_seminorm;
131
132 dealii::ObserverPointer<dealii::DoFHandler<dim> const> dof_handler;
133 dealii::ObserverPointer<dealii::Mapping<dim> const> mapping;
134
135 ErrorCalculationData<dim> error_data;
136};
137
138} // namespace ExaDG
139
140#endif /* EXADG_POSTPROCESSOR_ERROR_CALCULATION_H_ */
Definition time_control.h:64
Definition driver.cpp:33
Definition error_calculation.h:39
Definition time_control.h:40