22#ifndef INCLUDE_EXADG_POISSON_SPATIAL_DISCRETIZATION_WEAK_BOUNDARY_CONDITIONS_H_
23#define INCLUDE_EXADG_POISSON_SPATIAL_DISCRETIZATION_WEAK_BOUNDARY_CONDITIONS_H_
25#include <exadg/functions_and_boundary_conditions/evaluate_functions.h>
26#include <exadg/matrix_free/integrators.h>
27#include <exadg/operators/operator_type.h>
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)
55 if(operator_type == OperatorType::full or operator_type == OperatorType::homogeneous)
57 return integrator.get_value(q);
59 else if(operator_type == OperatorType::inhomogeneous)
61 return dealii::Tensor<rank, dim, dealii::VectorizedArray<Number>>();
65 AssertThrow(
false, dealii::ExcMessage(
"Specified OperatorType is not implemented!"));
68 return dealii::Tensor<rank, dim, dealii::VectorizedArray<Number>>();
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,
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,
84 dealii::Tensor<rank, dim, dealii::VectorizedArray<Number>> value_p;
86 if(boundary_type == BoundaryType::Dirichlet or boundary_type == BoundaryType::DirichletCached)
88 if(operator_type == OperatorType::full or operator_type == OperatorType::inhomogeneous)
90 dealii::Tensor<rank, dim, dealii::VectorizedArray<Number>> g;
92 if(boundary_type == BoundaryType::Dirichlet)
94 auto bc = boundary_descriptor->dirichlet_bc.find(boundary_id)->second;
95 auto q_points = integrator.quadrature_point(q);
97 g = FunctionEvaluator<rank, dim, Number>::value(*bc, q_points, time);
99 else if(boundary_type == BoundaryType::DirichletCached)
101 auto bc = boundary_descriptor->get_dirichlet_cached_data();
102 g = FunctionEvaluator<rank, dim, Number>::value(*bc,
103 integrator.get_current_cell_index(),
105 integrator.get_quadrature_index());
109 AssertThrow(
false, dealii::ExcMessage(
"Not implemented."));
112 value_p = -value_m + dealii::Tensor<rank, dim, dealii::VectorizedArray<Number>>(2.0 * g);
114 else if(operator_type == OperatorType::homogeneous)
120 AssertThrow(
false, dealii::ExcMessage(
"Specified OperatorType is not implemented!"));
123 else if(boundary_type == BoundaryType::Neumann)
129 AssertThrow(
false, dealii::ExcMessage(
"Boundary type of face is invalid or not implemented."));
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)
169 dealii::Tensor<rank, dim, dealii::VectorizedArray<Number>> normal_gradient_m;
171 if(operator_type == OperatorType::full or operator_type == OperatorType::homogeneous)
173 normal_gradient_m = integrator.get_normal_derivative(q);
175 else if(operator_type == OperatorType::inhomogeneous)
181 AssertThrow(
false, dealii::ExcMessage(
"Specified OperatorType is not implemented!"));
184 return normal_gradient_m;
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,
200 dealii::Tensor<rank, dim, dealii::VectorizedArray<Number>> normal_gradient_p;
202 if(boundary_type == BoundaryType::Dirichlet or boundary_type == BoundaryType::DirichletCached)
204 normal_gradient_p = normal_gradient_m;
206 else if(boundary_type == BoundaryType::Neumann)
208 if(operator_type == OperatorType::full or operator_type == OperatorType::inhomogeneous)
210 auto bc = boundary_descriptor->neumann_bc.find(boundary_id)->second;
211 auto q_points = integrator.quadrature_point(q);
213 auto h = FunctionEvaluator<rank, dim, Number>::value(*bc, q_points, time);
216 -normal_gradient_m + dealii::Tensor<rank, dim, dealii::VectorizedArray<Number>>(2.0 * h);
218 else if(operator_type == OperatorType::homogeneous)
220 normal_gradient_p = -normal_gradient_m;
224 AssertThrow(
false, dealii::ExcMessage(
"Specified OperatorType is not implemented!"));
229 AssertThrow(
false, dealii::ExcMessage(
"Boundary type of face is invalid or not implemented."));
232 return normal_gradient_p;
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,
249 dealii::Tensor<rank, dim, dealii::VectorizedArray<Number>> normal_gradient;
251 if(boundary_type == BoundaryType::Neumann)
253 auto bc = boundary_descriptor->neumann_bc.find(boundary_id)->second;
254 auto q_points = integrator.quadrature_point(q);
256 normal_gradient = FunctionEvaluator<rank, dim, Number>::value(*bc, q_points, time);
262 Assert(boundary_type == BoundaryType::Dirichlet or
263 boundary_type == BoundaryType::DirichletCached,
264 dealii::ExcMessage(
"Boundary type of face is invalid or not implemented."));
267 return normal_gradient;
Definition boundary_descriptor.h:46