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