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 spatially_weight_error(false),
45 weight(nullptr),
46 additional_quadrature_points(3),
47 directory("output/"),
48 name("all fields")
49 {
50 }
51
52 void
53 print(dealii::ConditionalOStream & pcout, bool unsteady)
54 {
55 print_parameter(pcout, "Error calculation", unsteady == true and analytical_solution);
56 if(unsteady == true and time_control_data.is_active)
57 {
58 print(pcout, unsteady, time_control_data);
59 print_parameter(pcout, "Calculate relative errors", calculate_relative_errors);
60 print_parameter(pcout, "Calculate H1-seminorm error", calculate_H1_seminorm_error);
61 print_parameter(pcout, "Write errors to file", write_errors_to_file);
62 if(write_errors_to_file)
63 print_parameter(pcout, "Directory", directory);
64 print_parameter(pcout, "Name", name);
65 }
66 }
67
68 std::shared_ptr<dealii::Function<dim>> analytical_solution;
69
70 // relative or absolute errors?
71 // If calculate_relative_errors == false, this implies that absolute errors are calculated
72 bool calculate_relative_errors;
73
74 // by default, only the L2-error is computed. Other norms have to be explicitly specified by the
75 // user.
76 bool calculate_H1_seminorm_error;
77
78 // data used to control the output
79 TimeControlData time_control_data;
80
81 // write errors to file?
82 bool write_errors_to_file;
83
84 // If true, a spatially weighted norm is computed.
85 bool spatially_weight_error;
86 // Weight used to compute spatially weighted error.
87 std::shared_ptr<dealii::Function<dim>> weight;
88 // Number of quadrature points in 1D: fe_degree + additional_quadrature_points
89 unsigned int additional_quadrature_points;
90
91 // directory and name (used as filename and as identifier for screen output)
92 std::string directory;
93 std::string name;
94};
95
96template<int dim, typename Number>
97class ErrorCalculator
98{
99public:
100 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
101
102 ErrorCalculator(MPI_Comm const & comm);
103
104 void
105 setup(dealii::DoFHandler<dim> const & dof_handler,
106 dealii::Mapping<dim> const & mapping,
107 ErrorCalculationData<dim> const & error_data);
108
109 void
110 evaluate(VectorType const & solution, double const time, bool const unsteady);
111
112 TimeControl time_control;
113
114private:
115 void
116 do_evaluate(VectorType const & solution_vector, double const time);
117
118 MPI_Comm const mpi_comm;
119
120 bool clear_files_L2, clear_files_H1_seminorm;
121
122 dealii::ObserverPointer<dealii::DoFHandler<dim> const> dof_handler;
123 dealii::ObserverPointer<dealii::Mapping<dim> const> mapping;
124
125 ErrorCalculationData<dim> error_data;
126};
127
128} // namespace ExaDG
129
130#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