ExaDG
Loading...
Searching...
No Matches
boundary_face_integrator_base.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2023 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_FUNCTIONS_AND_BOUNDARY_CONDITIONS_BOUNDARY_FACE_INTEGRATOR_BASE_H_
23#define EXADG_FUNCTIONS_AND_BOUNDARY_CONDITIONS_BOUNDARY_FACE_INTEGRATOR_BASE_H_
24
25#include <deal.II/matrix_free/matrix_free.h>
26
27#include <exadg/utilities/enum_utilities.h>
28
29namespace ExaDG
30{
36template<typename BoundaryDescriptorType, typename Number>
38{
39 using BoundaryType = typename BoundaryDescriptorType::boundary_type;
40 static constexpr int dim = BoundaryDescriptorType::dimension;
41
42public:
43 void
44 reinit(unsigned int const face, Number const time)
45 {
46 evaluation_time = time;
47
48 auto const boundary_id_new = matrix_free.get_boundary_id(face);
49
50 // only update boundary_type if needed to avoid an unnecessary search in boundary_descriptor
51 if(boundary_id_new != boundary_id)
52 {
53 boundary_id = boundary_id_new;
54 boundary_type = boundary_descriptor.get_boundary_type(boundary_id);
55 }
56 }
57
58 // A corresponding function has to be implemented in the deriving class.
59 // The return type varies dependent on the value type.
60 // auto
61 // get_value();
62
63protected:
64 BoundaryFaceIntegratorBase(dealii::MatrixFree<dim, Number> const & matrix_free_in,
65 BoundaryDescriptorType const & boundary_descriptor_in)
66 : matrix_free(matrix_free_in),
67 boundary_descriptor(boundary_descriptor_in),
68 evaluation_time(Number{0.0}),
69 boundary_id(dealii::numbers::invalid_boundary_id),
70 boundary_type(Utilities::default_constructor<BoundaryType>())
71 {
72 }
73
74 dealii::MatrixFree<dim, Number> const & matrix_free;
75 BoundaryDescriptorType const & boundary_descriptor;
76
77 Number evaluation_time;
78 dealii::types::boundary_id boundary_id;
79 BoundaryType boundary_type;
80};
81
82} // namespace ExaDG
83
84#endif /*EXADG_FUNCTIONS_AND_BOUNDARY_CONDITIONS_BOUNDARY_FACE_INTEGRATOR_BASE_H_*/
Definition boundary_face_integrator_base.h:38
Definition driver.cpp:33