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_
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>
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)
62 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> value_m;
64 if(operator_type == OperatorType::full or operator_type == OperatorType::homogeneous)
66 value_m = integrator.get_value(q);
68 else if(operator_type == OperatorType::inhomogeneous)
74 AssertThrow(
false, dealii::ExcMessage(
"Specified OperatorType is not implemented!"));
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,
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,
93 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> value_p;
95 if(boundary_type == BoundaryTypeU::Dirichlet or boundary_type == BoundaryTypeU::DirichletCached)
97 if(operator_type == OperatorType::full or operator_type == OperatorType::inhomogeneous)
99 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> g;
101 if(boundary_type == BoundaryTypeU::Dirichlet)
103 auto bc = boundary_descriptor->dirichlet_bc.find(boundary_id)->second;
104 auto q_points = integrator.quadrature_point(q);
106 g = FunctionEvaluator<1, dim, Number>::value(*bc, q_points, time);
108 else if(boundary_type == BoundaryTypeU::DirichletCached)
110 auto bc = boundary_descriptor->get_dirichlet_cached_data();
111 g = FunctionEvaluator<1, dim, Number>::value(*bc,
112 integrator.get_current_cell_index(),
114 integrator.get_quadrature_index());
118 AssertThrow(
false, dealii::ExcMessage(
"Not implemented."));
121 value_p = -value_m + dealii::Tensor<1, dim, dealii::VectorizedArray<Number>>(2.0 * g);
123 else if(operator_type == OperatorType::homogeneous)
129 AssertThrow(
false, dealii::ExcMessage(
"Specified OperatorType is not implemented!"));
132 else if(boundary_type == BoundaryTypeU::Neumann)
136 else if(boundary_type == BoundaryTypeU::Symmetry)
138 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> normal_m =
139 integrator.get_normal_vector(q);
141 value_p = value_m - 2.0 * (value_m * normal_m) * normal_m;
145 AssertThrow(
false, dealii::ExcMessage(
"Boundary type of face is invalid or not implemented."));
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,
182 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> u_p;
184 if(boundary_type == BoundaryTypeU::Dirichlet or boundary_type == BoundaryTypeU::DirichletCached)
186 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> g;
188 if(operator_type == OperatorType::full or operator_type == OperatorType::inhomogeneous)
190 if(boundary_type == BoundaryTypeU::Dirichlet)
192 auto bc = boundary_descriptor->dirichlet_bc.find(boundary_id)->second;
193 auto q_points = integrator.quadrature_point(q);
195 g = FunctionEvaluator<1, dim, Number>::value(*bc, q_points, time);
197 else if(boundary_type == BoundaryTypeU::DirichletCached)
199 auto bc = boundary_descriptor->get_dirichlet_cached_data();
200 g = FunctionEvaluator<1, dim, Number>::value(*bc,
201 integrator.get_current_cell_index(),
203 integrator.get_quadrature_index());
207 AssertThrow(
false, dealii::ExcMessage(
"Not implemented."));
211 if(type_dirichlet_bc == TypeDirichletBCs::Mirror)
213 u_p = -u_m + dealii::make_vectorized_array<Number>(2.0) * g;
215 else if(type_dirichlet_bc == TypeDirichletBCs::Direct)
221 AssertThrow(
false, dealii::ExcMessage(
"TypeDirichletBCs is not implemented."));
224 else if(boundary_type == BoundaryTypeU::Neumann)
228 else if(boundary_type == BoundaryTypeU::Symmetry)
230 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> normal_m =
231 integrator.get_normal_vector(q);
233 u_p = u_m - 2. * (u_m * normal_m) * normal_m;
237 AssertThrow(
false, dealii::ExcMessage(
"Boundary type of face is invalid or not implemented."));
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)
254 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> value_p;
256 if(boundary_type == BoundaryTypeU::Dirichlet or boundary_type == BoundaryTypeU::DirichletCached)
258 if(operator_type == OperatorType::full or operator_type == OperatorType::inhomogeneous)
260 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> g = integrator_bc.get_value(q);
262 value_p = -value_m + dealii::make_vectorized_array<Number>(2.0) * g;
264 else if(operator_type == OperatorType::homogeneous)
270 AssertThrow(
false, dealii::ExcMessage(
"Specified OperatorType is not implemented!"));
273 else if(boundary_type == BoundaryTypeU::Neumann)
277 else if(boundary_type == BoundaryTypeU::Symmetry)
279 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> normal_m =
280 integrator_bc.get_normal_vector(q);
282 value_p = value_m - 2.0 * (value_m * normal_m) * normal_m;
286 AssertThrow(
false, dealii::ExcMessage(
"Boundary type of face is invalid or not implemented."));
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)
317 dealii::VectorizedArray<Number> value_m = dealii::make_vectorized_array<Number>(0.0);
319 if(operator_type == OperatorType::full or operator_type == OperatorType::homogeneous)
321 value_m = integrator.get_value(q);
323 else if(operator_type == OperatorType::inhomogeneous)
329 AssertThrow(
false, dealii::ExcMessage(
"Specified OperatorType is not implemented!"));
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,
347 double const & inverse_scaling_factor)
349 dealii::VectorizedArray<Number> value_p = dealii::make_vectorized_array<Number>(0.0);
351 if(boundary_type == BoundaryTypeP::Dirichlet)
353 if(operator_type == OperatorType::full or operator_type == OperatorType::inhomogeneous)
355 auto bc = boundary_descriptor->dirichlet_bc.find(boundary_id)->second;
356 auto q_points = integrator.quadrature_point(q);
358 dealii::VectorizedArray<Number> g =
359 FunctionEvaluator<0, dim, Number>::value(*bc, q_points, time);
361 value_p = -value_m + 2.0 * inverse_scaling_factor * g;
363 else if(operator_type == OperatorType::homogeneous)
369 AssertThrow(
false, dealii::ExcMessage(
"Specified OperatorType is not implemented!"));
372 else if(boundary_type == BoundaryTypeP::Neumann)
378 AssertThrow(
false, dealii::ExcMessage(
"Boundary type of face is invalid or not implemented."));
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)
394 dealii::VectorizedArray<Number> value_p = dealii::make_vectorized_array<Number>(0.0);
396 if(boundary_type == BoundaryTypeP::Dirichlet)
398 if(operator_type == OperatorType::full or operator_type == OperatorType::inhomogeneous)
400 dealii::VectorizedArray<Number> g = integrator_bc.get_value(q);
402 value_p = -value_m + 2.0 * inverse_scaling_factor * g;
404 else if(operator_type == OperatorType::homogeneous)
410 AssertThrow(
false, dealii::ExcMessage(
"Specified OperatorType is not implemented!"));
413 else if(boundary_type == BoundaryTypeP::Neumann)
419 AssertThrow(
false, dealii::ExcMessage(
"Boundary type of face is invalid or not implemented."));
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,
469 bool const variable_normal_vector)
471 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> normal_gradient_p;
473 if(boundary_type == BoundaryTypeU::Dirichlet or boundary_type == BoundaryTypeU::DirichletCached)
475 normal_gradient_p = normal_gradient_m;
477 else if(boundary_type == BoundaryTypeU::Neumann)
479 if(operator_type == OperatorType::full or operator_type == OperatorType::inhomogeneous)
481 auto bc = boundary_descriptor->neumann_bc.find(boundary_id)->second;
482 auto q_points = integrator.quadrature_point(q);
484 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> h;
485 if(variable_normal_vector ==
false)
487 h = FunctionEvaluator<1, dim, Number>::value(*bc, q_points, time);
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);
496 normal_gradient_p = -normal_gradient_m + 2.0 * h;
498 else if(operator_type == OperatorType::homogeneous)
500 normal_gradient_p = -normal_gradient_m;
504 AssertThrow(
false, dealii::ExcMessage(
"Specified OperatorType is not implemented!"));
507 else if(boundary_type == BoundaryTypeU::Symmetry)
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;
514 AssertThrow(
false, dealii::ExcMessage(
"Boundary type of face is invalid or not implemented."));
517 return normal_gradient_p;
Definition boundary_descriptor.h:191
Definition boundary_descriptor.h:81