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_
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>
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)
56 dealii::VectorizedArray<Number> value_m = dealii::make_vectorized_array<Number>(0.0);
58 if(operator_type == OperatorType::full or operator_type == OperatorType::homogeneous)
60 value_m = integrator.get_value(q);
62 else if(operator_type == OperatorType::inhomogeneous)
68 AssertThrow(
false, dealii::ExcMessage(
"Specified OperatorType is not implemented!"));
74template<
int dim,
typename Number>
75inline DEAL_II_ALWAYS_INLINE
76 dealii::VectorizedArray<Number>
77 calculate_exterior_value(dealii::VectorizedArray<Number>
const & value_m,
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,
86 dealii::VectorizedArray<Number> value_p = dealii::make_vectorized_array<Number>(0.0);
88 if(boundary_type == BoundaryType::Dirichlet)
90 if(operator_type == OperatorType::full or operator_type == OperatorType::inhomogeneous)
92 dealii::VectorizedArray<Number> g;
94 auto bc = boundary_descriptor->dirichlet_bc.find(boundary_id)->second;
95 auto q_points = integrator.quadrature_point(q);
97 g = FunctionEvaluator<0, dim, Number>::value(*bc, q_points, time);
99 value_p = -value_m + 2.0 * g;
101 else if(operator_type == OperatorType::homogeneous)
107 AssertThrow(
false, dealii::ExcMessage(
"Specified OperatorType is not implemented!"));
110 else if(boundary_type == BoundaryType::Neumann)
116 AssertThrow(
false, dealii::ExcMessage(
"Boundary type of face is invalid or not implemented."));
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)
156 dealii::VectorizedArray<Number> normal_gradient_m = dealii::make_vectorized_array<Number>(0.0);
158 if(operator_type == OperatorType::full or operator_type == OperatorType::homogeneous)
160 normal_gradient_m = integrator.get_normal_derivative(q);
162 else if(operator_type == OperatorType::inhomogeneous)
168 AssertThrow(
false, dealii::ExcMessage(
"Specified OperatorType is not implemented!"));
171 return normal_gradient_m;
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,
187 dealii::VectorizedArray<Number> normal_gradient_p = dealii::make_vectorized_array<Number>(0.0);
189 if(boundary_type == BoundaryType::Dirichlet)
191 normal_gradient_p = normal_gradient_m;
193 else if(boundary_type == BoundaryType::Neumann)
195 if(operator_type == OperatorType::full or operator_type == OperatorType::inhomogeneous)
197 auto bc = boundary_descriptor->neumann_bc.find(boundary_id)->second;
198 auto q_points = integrator.quadrature_point(q);
200 auto h = FunctionEvaluator<0, dim, Number>::value(*bc, q_points, time);
202 normal_gradient_p = -normal_gradient_m + 2.0 * h;
204 else if(operator_type == OperatorType::homogeneous)
206 normal_gradient_p = -normal_gradient_m;
210 AssertThrow(
false, dealii::ExcMessage(
"Specified OperatorType is not implemented!"));
215 AssertThrow(
false, dealii::ExcMessage(
"Boundary type of face is invalid or not implemented."));
218 return normal_gradient_p;
Definition boundary_descriptor.h:46