ExaDG
Loading...
Searching...
No Matches
boundary_descriptor.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_ACOUSTIC_CONSERVATION_EQUATIONS_USER_INTERFACE_BOUNDARY_DESCRIPTOR_H_
23#define EXADG_ACOUSTIC_CONSERVATION_EQUATIONS_USER_INTERFACE_BOUNDARY_DESCRIPTOR_H_
24
25#include <set>
26
27#include <deal.II/base/function.h>
28#include <deal.II/base/types.h>
29
30#include <exadg/functions_and_boundary_conditions/verify_boundary_conditions.h>
31
32namespace ExaDG
33{
34namespace Acoustics
35{
36/* Boundary conditions:
37 *
38 * +---------------------------+-------------------------+-------------------------+
39 * | example | pressure | velocity |
40 * +---------------------------+-------------------------+-------------------------+
41 * | prescribe pressure values | Dirichlet: | |
42 * | | prescribe g_p | no BCs to be prescribed |
43 * +-----------------------------------------------------+-------------------------+
44 * | prescribe velocity values | | Dirichlet: |
45 * | | no BCs to be prescribed | prescribe g_u |
46 * +---------------------------+-------------------------+-------------------------+
47 * | admittance (Y) BC | | Admittance: |
48 * | rho * u * n = Y / c * p | no BCs to be prescribed | prescribe Y |
49 * | Y = 0: sound hard | | |
50 * | Y = 1: first order ABC | | |
51 * +---------------------------+-------------------------+-------------------------+
52 */
53
54enum class BoundaryType
55{
56 Undefined,
57 PressureDirichlet,
58 VelocityDirichlet,
59 Admittance
60};
61
62template<int dim>
64{
65 using boundary_type = BoundaryType;
66
67 static constexpr int dimension = dim;
68
69 // Dirichlet: prescribe pressure
70 std::map<dealii::types::boundary_id, std::shared_ptr<dealii::Function<dim>>> pressure_dbc;
71
72 // Dirichlet: prescribe velocity
73 std::map<dealii::types::boundary_id, std::shared_ptr<dealii::Function<dim>>> velocity_dbc;
74
75 // BC for Admittance Y:
76 // Special cases are:
77 // Y = 0: sound hard (perfectly reflecting)
78 // Y = 1: first order ABC
79 std::map<dealii::types::boundary_id, std::shared_ptr<dealii::Function<dim>>> admittance_bc;
80
81 // return the boundary type
82 inline DEAL_II_ALWAYS_INLINE //
83 BoundaryType
84 get_boundary_type(dealii::types::boundary_id const & boundary_id) const
85 {
86 if(this->pressure_dbc.find(boundary_id) != this->pressure_dbc.end())
87 return BoundaryType::PressureDirichlet;
88
89 if(this->velocity_dbc.find(boundary_id) != this->velocity_dbc.end())
90 return BoundaryType::VelocityDirichlet;
91
92 if(this->admittance_bc.find(boundary_id) != this->admittance_bc.end())
93 return BoundaryType::Admittance;
94
95 AssertThrow(false, dealii::ExcMessage("Boundary type of face is invalid or not implemented."));
96
97 return BoundaryType::Undefined;
98 }
99
100 inline DEAL_II_ALWAYS_INLINE //
101 void
102 verify_boundary_conditions(
103 dealii::types::boundary_id const boundary_id,
104 std::set<dealii::types::boundary_id> const & periodic_boundary_ids) const
105 {
106 unsigned int counter = 0;
107 if(this->pressure_dbc.find(boundary_id) != this->pressure_dbc.end())
108 counter++;
109
110 if(this->velocity_dbc.find(boundary_id) != this->velocity_dbc.end())
111 counter++;
112
113 if(this->admittance_bc.find(boundary_id) != this->admittance_bc.end())
114 counter++;
115
116 if(periodic_boundary_ids.find(boundary_id) != periodic_boundary_ids.end())
117 counter++;
118
119 AssertThrow(counter == 1,
120 dealii::ExcMessage("Boundary face with non-unique boundary type found."));
121 }
122};
123
124} // namespace Acoustics
125} // namespace ExaDG
126
127#endif /* EXADG_ACOUSTIC_CONSERVATION_EQUATIONS_USER_INTERFACE_BOUNDARY_DESCRIPTOR_H_ */
Definition driver.cpp:33
Definition boundary_descriptor.h:64