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_INCOMPRESSIBLE_NAVIER_STOKES_SPATIAL_DISCRETIZATION_OPERATORS_WEAK_BOUNDARY_CONDITIONS_H_
23#define INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_SPATIAL_DISCRETIZATION_OPERATORS_WEAK_BOUNDARY_CONDITIONS_H_
24
25#include <exadg/functions_and_boundary_conditions/evaluate_functions.h>
26#include <exadg/incompressible_navier_stokes/user_interface/boundary_descriptor.h>
27#include <exadg/incompressible_navier_stokes/user_interface/enum_types.h>
28#include <exadg/matrix_free/integrators.h>
29#include <exadg/operators/operator_type.h>
30
31namespace ExaDG
32{
33namespace IncNS
34{
35// clang-format off
36/*
37 * Velocity:
38 *
39 * The following two functions calculate the interior/exterior value for boundary faces depending on the
40 * operator type, the type of the boundary face and the given boundary conditions.
41 *
42 * +-------------------------+--------------------+------------------------------+
43 * | Dirichlet boundaries | Neumann boundaries | symmetry boundaries |
44 * +-------------------------+-------------------------+--------------------+------------------------------+
45 * | full operator | u⁺ = -u⁻ + 2g | u⁺ = u⁻ | u⁺ = u⁻ - 2 (u⁻*n)n |
46 * +-------------------------+-------------------------+--------------------+------------------------------+
47 * | homogeneous operator | u⁺ = -u⁻ | u⁺ = u⁻ | u⁺ = u⁻ - 2 (u⁻*n)n |
48 * +-------------------------+-------------------------+--------------------+------------------------------+
49 * | inhomogeneous operator | u⁺ = -u⁻ + 2g , u⁻ = 0 | u⁺ = u⁻ , u⁻ = 0 | u⁺ = u⁻ - 2 (u⁻*n)n , u⁻ = 0 |
50 * +-------------------------+-------------------------+--------------------+------------------------------+
51 *
52 */
53// clang-format on
54template<int dim, typename Number>
55inline DEAL_II_ALWAYS_INLINE //
56 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>>
57 calculate_interior_value(unsigned int const q,
58 FaceIntegrator<dim, dim, Number> const & integrator,
59 OperatorType const & operator_type)
60{
61 // element e⁻
62 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> value_m;
63
64 if(operator_type == OperatorType::full or operator_type == OperatorType::homogeneous)
65 {
66 value_m = integrator.get_value(q);
67 }
68 else if(operator_type == OperatorType::inhomogeneous)
69 {
70 // do nothing, value_m is already initialized with zeros
71 }
72 else
73 {
74 AssertThrow(false, dealii::ExcMessage("Specified OperatorType is not implemented!"));
75 }
76
77 return value_m;
78}
79
80template<int dim, typename Number>
81inline DEAL_II_ALWAYS_INLINE //
82 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>>
83 calculate_exterior_value(dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> const & value_m,
84 unsigned int const q,
85 FaceIntegrator<dim, dim, Number> const & integrator,
86 OperatorType const & operator_type,
87 BoundaryTypeU const & boundary_type,
88 dealii::types::boundary_id const boundary_id,
89 std::shared_ptr<BoundaryDescriptorU<dim> const> boundary_descriptor,
90 double const & time)
91{
92 // element e⁺
93 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> value_p;
94
95 if(boundary_type == BoundaryTypeU::Dirichlet or boundary_type == BoundaryTypeU::DirichletCached)
96 {
97 if(operator_type == OperatorType::full or operator_type == OperatorType::inhomogeneous)
98 {
99 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> g;
100
101 if(boundary_type == BoundaryTypeU::Dirichlet)
102 {
103 auto bc = boundary_descriptor->dirichlet_bc.find(boundary_id)->second;
104 auto q_points = integrator.quadrature_point(q);
105
106 g = FunctionEvaluator<1, dim, Number>::value(*bc, q_points, time);
107 }
108 else if(boundary_type == BoundaryTypeU::DirichletCached)
109 {
110 auto bc = boundary_descriptor->get_dirichlet_cached_data();
111 g = FunctionEvaluator<1, dim, Number>::value(*bc,
112 integrator.get_current_cell_index(),
113 q,
114 integrator.get_quadrature_index());
115 }
116 else
117 {
118 AssertThrow(false, dealii::ExcMessage("Not implemented."));
119 }
120
121 value_p = -value_m + dealii::Tensor<1, dim, dealii::VectorizedArray<Number>>(2.0 * g);
122 }
123 else if(operator_type == OperatorType::homogeneous)
124 {
125 value_p = -value_m;
126 }
127 else
128 {
129 AssertThrow(false, dealii::ExcMessage("Specified OperatorType is not implemented!"));
130 }
131 }
132 else if(boundary_type == BoundaryTypeU::Neumann)
133 {
134 value_p = value_m;
135 }
136 else if(boundary_type == BoundaryTypeU::Symmetry)
137 {
138 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> normal_m =
139 integrator.get_normal_vector(q);
140
141 value_p = value_m - 2.0 * (value_m * normal_m) * normal_m;
142 }
143 else
144 {
145 AssertThrow(false, dealii::ExcMessage("Boundary type of face is invalid or not implemented."));
146 }
147
148 return value_p;
149}
150
151/*
152 * This function calculates the exterior velocity on boundary faces
153 * according to:
154 *
155 * Dirichlet boundary: u⁺ = -u⁻ + 2g (mirror) or g (direct)
156 * Neumann boundary: u⁺ = u⁻
157 * symmetry boundary: u⁺ = u⁻ -(u⁻*n)n - (u⁻*n)n = u⁻ - 2 (u⁻*n)n
158 *
159 * The name "convective" indicates that this function is used when
160 * evaluating the convective operator.
161 *
162 * OperatorType:
163 * - To evaluate the nonlinear convective term (residual evaluation), use OperatorType::full.
164 * - To evaluate the linearization of the nonlinear convective term, use OperatorType::homogeneous.
165 * - To evaluate the linearly implicit convective term, use OperatorType::homogeneous for the
166 * "matrix-vector" product and OperatorType::inhomogeneous for the "rhs" contribution.
167 */
168template<int dim, typename Number>
169inline DEAL_II_ALWAYS_INLINE //
170 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>>
171 calculate_exterior_value_convective(
172 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> const & u_m,
173 unsigned int const q,
174 FaceIntegrator<dim, dim, Number> & integrator,
175 OperatorType const & operator_type,
176 BoundaryTypeU const & boundary_type,
177 TypeDirichletBCs const & type_dirichlet_bc,
178 dealii::types::boundary_id const boundary_id,
179 std::shared_ptr<BoundaryDescriptorU<dim> const> boundary_descriptor,
180 double const & time)
181{
182 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> u_p;
183
184 if(boundary_type == BoundaryTypeU::Dirichlet or boundary_type == BoundaryTypeU::DirichletCached)
185 {
186 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> g;
187
188 if(operator_type == OperatorType::full or operator_type == OperatorType::inhomogeneous)
189 {
190 if(boundary_type == BoundaryTypeU::Dirichlet)
191 {
192 auto bc = boundary_descriptor->dirichlet_bc.find(boundary_id)->second;
193 auto q_points = integrator.quadrature_point(q);
194
195 g = FunctionEvaluator<1, dim, Number>::value(*bc, q_points, time);
196 }
197 else if(boundary_type == BoundaryTypeU::DirichletCached)
198 {
199 auto bc = boundary_descriptor->get_dirichlet_cached_data();
200 g = FunctionEvaluator<1, dim, Number>::value(*bc,
201 integrator.get_current_cell_index(),
202 q,
203 integrator.get_quadrature_index());
204 }
205 else
206 {
207 AssertThrow(false, dealii::ExcMessage("Not implemented."));
208 }
209 }
210
211 if(type_dirichlet_bc == TypeDirichletBCs::Mirror)
212 {
213 u_p = -u_m + dealii::make_vectorized_array<Number>(2.0) * g;
214 }
215 else if(type_dirichlet_bc == TypeDirichletBCs::Direct)
216 {
217 u_p = g;
218 }
219 else
220 {
221 AssertThrow(false, dealii::ExcMessage("TypeDirichletBCs is not implemented."));
222 }
223 }
224 else if(boundary_type == BoundaryTypeU::Neumann)
225 {
226 u_p = u_m;
227 }
228 else if(boundary_type == BoundaryTypeU::Symmetry)
229 {
230 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> normal_m =
231 integrator.get_normal_vector(q);
232
233 u_p = u_m - 2. * (u_m * normal_m) * normal_m;
234 }
235 else
236 {
237 AssertThrow(false, dealii::ExcMessage("Boundary type of face is invalid or not implemented."));
238 }
239
240 return u_p;
241}
242
243template<int dim, typename Number>
244inline DEAL_II_ALWAYS_INLINE //
245 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>>
246 calculate_exterior_value_from_dof_vector(
247 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> const & value_m,
248 unsigned int const q,
249 FaceIntegrator<dim, dim, Number> const & integrator_bc,
250 OperatorType const & operator_type,
251 BoundaryTypeU const & boundary_type)
252{
253 // element e⁺
254 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> value_p;
255
256 if(boundary_type == BoundaryTypeU::Dirichlet or boundary_type == BoundaryTypeU::DirichletCached)
257 {
258 if(operator_type == OperatorType::full or operator_type == OperatorType::inhomogeneous)
259 {
260 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> g = integrator_bc.get_value(q);
261
262 value_p = -value_m + dealii::make_vectorized_array<Number>(2.0) * g;
263 }
264 else if(operator_type == OperatorType::homogeneous)
265 {
266 value_p = -value_m;
267 }
268 else
269 {
270 AssertThrow(false, dealii::ExcMessage("Specified OperatorType is not implemented!"));
271 }
272 }
273 else if(boundary_type == BoundaryTypeU::Neumann)
274 {
275 value_p = value_m;
276 }
277 else if(boundary_type == BoundaryTypeU::Symmetry)
278 {
279 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> normal_m =
280 integrator_bc.get_normal_vector(q);
281
282 value_p = value_m - 2.0 * (value_m * normal_m) * normal_m;
283 }
284 else
285 {
286 AssertThrow(false, dealii::ExcMessage("Boundary type of face is invalid or not implemented."));
287 }
288
289 return value_p;
290}
291
292/*
293 * Pressure:
294 *
295 * These two function calculate the interior/exterior value on boundary faces depending on the
296 * operator type, the type of the boundary face and the given boundary conditions.
297 *
298 * +--------------------+----------------------+
299 * | Neumann boundaries | Dirichlet boundaries |
300 * +-------------------------+--------------------+----------------------+
301 * | full operator | p⁺ = p⁻ | p⁺ = - p⁻ + 2g |
302 * +-------------------------+--------------------+----------------------+
303 * | homogeneous operator | p⁺ = p⁻ | p⁺ = - p⁻ |
304 * +-------------------------+--------------------+----------------------+
305 * | inhomogeneous operator | p⁺ = 0 , p⁻ = 0 | p⁺ = 2g , p⁻ = 0 |
306 * +-------------------------+--------------------+----------------------+
307 *
308 */
309template<int dim, typename Number>
310inline DEAL_II_ALWAYS_INLINE //
311 dealii::VectorizedArray<Number>
312 calculate_interior_value(unsigned int const q,
313 FaceIntegrator<dim, 1, Number> const & integrator,
314 OperatorType const & operator_type)
315{
316 // element e⁻
317 dealii::VectorizedArray<Number> value_m = dealii::make_vectorized_array<Number>(0.0);
318
319 if(operator_type == OperatorType::full or operator_type == OperatorType::homogeneous)
320 {
321 value_m = integrator.get_value(q);
322 }
323 else if(operator_type == OperatorType::inhomogeneous)
324 {
325 // do nothing, value_m is already initialized with zeros
326 }
327 else
328 {
329 AssertThrow(false, dealii::ExcMessage("Specified OperatorType is not implemented!"));
330 }
331
332 return value_m;
333}
334
335
336template<int dim, typename Number>
337inline DEAL_II_ALWAYS_INLINE //
338 dealii::VectorizedArray<Number>
339 calculate_exterior_value(dealii::VectorizedArray<Number> const & value_m,
340 unsigned int const q,
341 FaceIntegrator<dim, 1, Number> const & integrator,
342 OperatorType const & operator_type,
343 BoundaryTypeP const & boundary_type,
344 dealii::types::boundary_id const boundary_id,
345 std::shared_ptr<BoundaryDescriptorP<dim> const> boundary_descriptor,
346 double const & time,
347 double const & inverse_scaling_factor)
348{
349 dealii::VectorizedArray<Number> value_p = dealii::make_vectorized_array<Number>(0.0);
350
351 if(boundary_type == BoundaryTypeP::Dirichlet)
352 {
353 if(operator_type == OperatorType::full or operator_type == OperatorType::inhomogeneous)
354 {
355 auto bc = boundary_descriptor->dirichlet_bc.find(boundary_id)->second;
356 auto q_points = integrator.quadrature_point(q);
357
358 dealii::VectorizedArray<Number> g =
359 FunctionEvaluator<0, dim, Number>::value(*bc, q_points, time);
360
361 value_p = -value_m + 2.0 * inverse_scaling_factor * g;
362 }
363 else if(operator_type == OperatorType::homogeneous)
364 {
365 value_p = -value_m;
366 }
367 else
368 {
369 AssertThrow(false, dealii::ExcMessage("Specified OperatorType is not implemented!"));
370 }
371 }
372 else if(boundary_type == BoundaryTypeP::Neumann)
373 {
374 value_p = value_m;
375 }
376 else
377 {
378 AssertThrow(false, dealii::ExcMessage("Boundary type of face is invalid or not implemented."));
379 }
380
381 return value_p;
382}
383
384template<int dim, typename Number>
385inline DEAL_II_ALWAYS_INLINE //
386 dealii::VectorizedArray<Number>
387 calculate_exterior_value_from_dof_vector(dealii::VectorizedArray<Number> const & value_m,
388 unsigned int const q,
389 FaceIntegrator<dim, 1, Number> const & integrator_bc,
390 OperatorType const & operator_type,
391 BoundaryTypeP const & boundary_type,
392 double const & inverse_scaling_factor)
393{
394 dealii::VectorizedArray<Number> value_p = dealii::make_vectorized_array<Number>(0.0);
395
396 if(boundary_type == BoundaryTypeP::Dirichlet)
397 {
398 if(operator_type == OperatorType::full or operator_type == OperatorType::inhomogeneous)
399 {
400 dealii::VectorizedArray<Number> g = integrator_bc.get_value(q);
401
402 value_p = -value_m + 2.0 * inverse_scaling_factor * g;
403 }
404 else if(operator_type == OperatorType::homogeneous)
405 {
406 value_p = -value_m;
407 }
408 else
409 {
410 AssertThrow(false, dealii::ExcMessage("Specified OperatorType is not implemented!"));
411 }
412 }
413 else if(boundary_type == BoundaryTypeP::Neumann)
414 {
415 value_p = value_m;
416 }
417 else
418 {
419 AssertThrow(false, dealii::ExcMessage("Boundary type of face is invalid or not implemented."));
420 }
421
422 return value_p;
423}
424
425// clang-format off
426/*
427 * This function calculates the exterior velocity gradient in normal direction depending
428 * on the operator type, the type of the boundary face and the given boundary conditions.
429 *
430 * The function argument normal_gradient_m has the meaning of F(u⁻)*n with
431 *
432 * Divergence formulation: F(u) = F_nu(u) / nu = ( grad(u) + grad(u)^T )
433 * Laplace formulation: F(u) = F_nu(u) / nu = grad(u)
434 *
435 * +---------------------------------+---------------------------------------+----------------------------------------------------+
436 * | Dirichlet boundaries | Neumann boundaries | symmetry boundaries |
437 * +-------------------------+---------------------------------+---------------------------------------+----------------------------------------------------+
438 * | full operator | F(u⁺)*n = F(u⁻)*n | F(u⁺)*n = -F(u⁻)*n + 2h | F(u⁺)*n = -F(u⁻)*n + 2*[(F(u⁻)*n)*n]n |
439 * +-------------------------+---------------------------------+---------------------------------------+----------------------------------------------------+
440 * | homogeneous operator | F(u⁺)*n = F(u⁻)*n | F(u⁺)*n = -F(u⁻)*n | F(u⁺)*n = -F(u⁻)*n + 2*[(F(u⁻)*n)*n]n |
441 * +-------------------------+---------------------------------+---------------------------------------+----------------------------------------------------+
442 * | inhomogeneous operator | F(u⁺)*n = F(u⁻)*n, F(u⁻)*n = 0 | F(u⁺)*n = -F(u⁻)*n + 2h , F(u⁻)*n = 0 | F(u⁺)*n = -F(u⁻)*n + 2*[(F(u⁻)*n)*n]n, F(u⁻)*n = 0 |
443 * +-------------------------+---------------------------------+---------------------------------------+----------------------------------------------------+
444 *
445 * +---------------------------------+---------------------------------------+----------------------------------------------------+
446 * | Dirichlet boundaries | Neumann boundaries | symmetry boundaries |
447 * +-------------------------+---------------------------------+---------------------------------------+----------------------------------------------------+
448 * | full operator | {{F(u)}}*n = F(u⁻)*n | {{F(u)}}*n = h | {{F(u)}}*n = [(F(u⁻)*n)*n]n |
449 * +-------------------------+---------------------------------+---------------------------------------+----------------------------------------------------+
450 * | homogeneous operator | {{F(u)}}*n = F(u⁻)*n | {{F(u)}}*n = 0 | {{F(u)}}*n = [(F(u⁻)*n)*n]n |
451 * +-------------------------+---------------------------------+---------------------------------------+----------------------------------------------------+
452 * | inhomogeneous operator | {{F(u)}}*n = 0 | {{F(u)}}*n = h | {{F(u)}}*n = 0 |
453 * +-------------------------+---------------------------------+---------------------------------------+----------------------------------------------------+
454 */
455// clang-format on
456
457template<int dim, typename Number>
458inline DEAL_II_ALWAYS_INLINE //
459 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>>
460 calculate_exterior_normal_gradient(
461 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> const & normal_gradient_m,
462 unsigned int const q,
463 FaceIntegrator<dim, dim, Number> const & integrator,
464 OperatorType const & operator_type,
465 BoundaryTypeU const & boundary_type,
466 dealii::types::boundary_id const boundary_id,
467 std::shared_ptr<BoundaryDescriptorU<dim> const> boundary_descriptor,
468 double const & time,
469 bool const variable_normal_vector)
470{
471 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> normal_gradient_p;
472
473 if(boundary_type == BoundaryTypeU::Dirichlet or boundary_type == BoundaryTypeU::DirichletCached)
474 {
475 normal_gradient_p = normal_gradient_m;
476 }
477 else if(boundary_type == BoundaryTypeU::Neumann)
478 {
479 if(operator_type == OperatorType::full or operator_type == OperatorType::inhomogeneous)
480 {
481 auto bc = boundary_descriptor->neumann_bc.find(boundary_id)->second;
482 auto q_points = integrator.quadrature_point(q);
483
484 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> h;
485 if(variable_normal_vector == false)
486 {
487 h = FunctionEvaluator<1, dim, Number>::value(*bc, q_points, time);
488 }
489 else
490 {
491 auto normals_m = integrator.get_normal_vector(q);
492 h = FunctionEvaluator<1, dim, Number>::value(
493 *(std::dynamic_pointer_cast<FunctionWithNormal<dim>>(bc)), q_points, normals_m, time);
494 }
495
496 normal_gradient_p = -normal_gradient_m + 2.0 * h;
497 }
498 else if(operator_type == OperatorType::homogeneous)
499 {
500 normal_gradient_p = -normal_gradient_m;
501 }
502 else
503 {
504 AssertThrow(false, dealii::ExcMessage("Specified OperatorType is not implemented!"));
505 }
506 }
507 else if(boundary_type == BoundaryTypeU::Symmetry)
508 {
509 auto normal_m = integrator.get_normal_vector(q);
510 normal_gradient_p = -normal_gradient_m + 2.0 * (normal_gradient_m * normal_m) * normal_m;
511 }
512 else
513 {
514 AssertThrow(false, dealii::ExcMessage("Boundary type of face is invalid or not implemented."));
515 }
516
517 return normal_gradient_p;
518}
519
520} // namespace IncNS
521} // namespace ExaDG
522
523#endif /* INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_SPATIAL_DISCRETIZATION_OPERATORS_WEAK_BOUNDARY_CONDITIONS_H_ \
524 */
Definition driver.cpp:33
Definition boundary_descriptor.h:191
Definition boundary_descriptor.h:81