ExaDG
Loading...
Searching...
No Matches
weak_boundary_conditions.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 INCLUDE_EXADG_CONVECTION_DIFFUSION_SPATIAL_DISCRETIZATION_OPERATORS_WEAK_BOUNDARY_CONDITIONS_H_
23#define INCLUDE_EXADG_CONVECTION_DIFFUSION_SPATIAL_DISCRETIZATION_OPERATORS_WEAK_BOUNDARY_CONDITIONS_H_
24
25#include <exadg/convection_diffusion/user_interface/boundary_descriptor.h>
26#include <exadg/functions_and_boundary_conditions/evaluate_functions.h>
27#include <exadg/matrix_free/integrators.h>
28#include <exadg/operators/operator_type.h>
29
30namespace ExaDG
31{
32namespace ConvDiff
33{
34/*
35 * The following two functions calculate the interior_value/exterior_value
36 * depending on the operator type, the type of the boundary face
37 * and the given boundary conditions.
38 *
39 * +----------------------+--------------------+
40 * | Dirichlet boundaries | Neumann boundaries |
41 * +-------------------------+----------------------+--------------------+
42 * | full operator | phi⁺ = -phi⁻ + 2g | phi⁺ = phi⁻ |
43 * +-------------------------+----------------------+--------------------+
44 * | homogeneous operator | phi⁺ = -phi⁻ | phi⁺ = phi⁻ |
45 * +-------------------------+----------------------+--------------------+
46 * | inhomogeneous operator | phi⁻ = 0, phi⁺ = 2g | phi⁻ = 0, phi⁺ = 0 |
47 * +-------------------------+----------------------+--------------------+
48 */
49template<int dim, typename Number>
50inline DEAL_II_ALWAYS_INLINE //
51 dealii::VectorizedArray<Number>
52 calculate_interior_value(unsigned int const q,
53 FaceIntegrator<dim, 1, Number> const & integrator,
54 OperatorType const & operator_type)
55{
56 dealii::VectorizedArray<Number> value_m = dealii::make_vectorized_array<Number>(0.0);
57
58 if(operator_type == OperatorType::full or operator_type == OperatorType::homogeneous)
59 {
60 value_m = integrator.get_value(q);
61 }
62 else if(operator_type == OperatorType::inhomogeneous)
63 {
64 // do nothing (value_m already initialized with 0.0)
65 }
66 else
67 {
68 AssertThrow(false, dealii::ExcMessage("Specified OperatorType is not implemented!"));
69 }
70
71 return value_m;
72}
73
74template<int dim, typename Number>
75inline DEAL_II_ALWAYS_INLINE //
76 dealii::VectorizedArray<Number>
77 calculate_exterior_value(dealii::VectorizedArray<Number> const & value_m,
78 unsigned int const q,
79 FaceIntegrator<dim, 1, Number> const & integrator,
80 OperatorType const & operator_type,
81 BoundaryType const & boundary_type,
82 dealii::types::boundary_id const boundary_id,
83 std::shared_ptr<BoundaryDescriptor<dim> const> boundary_descriptor,
84 double const & time)
85{
86 dealii::VectorizedArray<Number> value_p = dealii::make_vectorized_array<Number>(0.0);
87
88 if(boundary_type == BoundaryType::Dirichlet)
89 {
90 if(operator_type == OperatorType::full or operator_type == OperatorType::inhomogeneous)
91 {
92 dealii::VectorizedArray<Number> g;
93
94 auto bc = boundary_descriptor->dirichlet_bc.find(boundary_id)->second;
95 auto q_points = integrator.quadrature_point(q);
96
97 g = FunctionEvaluator<0, dim, Number>::value(*bc, q_points, time);
98
99 value_p = -value_m + 2.0 * g;
100 }
101 else if(operator_type == OperatorType::homogeneous)
102 {
103 value_p = -value_m;
104 }
105 else
106 {
107 AssertThrow(false, dealii::ExcMessage("Specified OperatorType is not implemented!"));
108 }
109 }
110 else if(boundary_type == BoundaryType::Neumann)
111 {
112 value_p = value_m;
113 }
114 else
115 {
116 AssertThrow(false, dealii::ExcMessage("Boundary type of face is invalid or not implemented."));
117 }
118
119 return value_p;
120}
121
122// clang-format off
123 /*
124 * The following two functions calculate the interior/exterior gradient
125 * in normal direction depending on the operator type, the type of the boundary face
126 * and the given boundary conditions.
127 *
128 * +-----------------------------------------------+------------------------------------------------------+
129 * | Dirichlet boundaries | Neumann boundaries |
130 * +-------------------------+-----------------------------------------------+------------------------------------------------------+
131 * | full operator | grad(phi⁺)*n = grad(phi⁻)*n | grad(phi⁺)*n = -grad(phi⁻)*n + 2h |
132 * +-------------------------+-----------------------------------------------+------------------------------------------------------+
133 * | homogeneous operator | grad(phi⁺)*n = grad(phi⁻)*n | grad(phi⁺)*n = -grad(phi⁻)*n |
134 * +-------------------------+-----------------------------------------------+------------------------------------------------------+
135 * | inhomogeneous operator | grad(phi⁺)*n = grad(phi⁻)*n, grad(phi⁻)*n = 0 | grad(phi⁺)*n = -grad(phi⁻)*n + 2h, grad(phi⁻)*n = 0 |
136 * +-------------------------+-----------------------------------------------+------------------------------------------------------+
137 *
138 * +-----------------------------------------------+------------------------------------------------------+
139 * | Dirichlet boundaries | Neumann boundaries |
140 * +-------------------------+-----------------------------------------------+------------------------------------------------------+
141 * | full operator | {{grad(phi)}}*n = grad(phi⁻)*n | {{grad(phi)}}*n = h |
142 * +-------------------------+-----------------------------------------------+------------------------------------------------------+
143 * | homogeneous operator | {{grad(phi)}}*n = grad(phi⁻)*n | {{grad(phi)}}*n = 0 |
144 * +-------------------------+-----------------------------------------------+------------------------------------------------------+
145 * | inhomogeneous operator | {{grad(phi)}}*n = 0 | {{grad(phi)}}*n = h |
146 * +-------------------------+-----------------------------------------------+------------------------------------------------------+
147 */
148// clang-format on
149template<int dim, typename Number>
150inline DEAL_II_ALWAYS_INLINE //
151 dealii::VectorizedArray<Number>
152 calculate_interior_normal_gradient(unsigned int const q,
153 FaceIntegrator<dim, 1, Number> const & integrator,
154 OperatorType const & operator_type)
155{
156 dealii::VectorizedArray<Number> normal_gradient_m = dealii::make_vectorized_array<Number>(0.0);
157
158 if(operator_type == OperatorType::full or operator_type == OperatorType::homogeneous)
159 {
160 normal_gradient_m = integrator.get_normal_derivative(q);
161 }
162 else if(operator_type == OperatorType::inhomogeneous)
163 {
164 // do nothing (normal_gradient_m already initialized with 0.0)
165 }
166 else
167 {
168 AssertThrow(false, dealii::ExcMessage("Specified OperatorType is not implemented!"));
169 }
170
171 return normal_gradient_m;
172}
173
174template<int dim, typename Number>
175inline DEAL_II_ALWAYS_INLINE //
176 dealii::VectorizedArray<Number>
177 calculate_exterior_normal_gradient(
178 dealii::VectorizedArray<Number> const & normal_gradient_m,
179 unsigned int const q,
180 FaceIntegrator<dim, 1, Number> const & integrator,
181 OperatorType const & operator_type,
182 BoundaryType const & boundary_type,
183 dealii::types::boundary_id const boundary_id,
184 std::shared_ptr<BoundaryDescriptor<dim> const> boundary_descriptor,
185 double const & time)
186{
187 dealii::VectorizedArray<Number> normal_gradient_p = dealii::make_vectorized_array<Number>(0.0);
188
189 if(boundary_type == BoundaryType::Dirichlet)
190 {
191 normal_gradient_p = normal_gradient_m;
192 }
193 else if(boundary_type == BoundaryType::Neumann)
194 {
195 if(operator_type == OperatorType::full or operator_type == OperatorType::inhomogeneous)
196 {
197 auto bc = boundary_descriptor->neumann_bc.find(boundary_id)->second;
198 auto q_points = integrator.quadrature_point(q);
199
200 auto h = FunctionEvaluator<0, dim, Number>::value(*bc, q_points, time);
201
202 normal_gradient_p = -normal_gradient_m + 2.0 * h;
203 }
204 else if(operator_type == OperatorType::homogeneous)
205 {
206 normal_gradient_p = -normal_gradient_m;
207 }
208 else
209 {
210 AssertThrow(false, dealii::ExcMessage("Specified OperatorType is not implemented!"));
211 }
212 }
213 else
214 {
215 AssertThrow(false, dealii::ExcMessage("Boundary type of face is invalid or not implemented."));
216 }
217
218 return normal_gradient_p;
219}
220
221} // namespace ConvDiff
222} // namespace ExaDG
223
224#endif /* INCLUDE_EXADG_CONVECTION_DIFFUSION_SPATIAL_DISCRETIZATION_OPERATORS_WEAK_BOUNDARY_CONDITIONS_H_ \
225 */
Definition driver.cpp:33