22#ifndef EXADG_ACOUSTIC_CONSERVATION_EQUATIONS_SPATIAL_DISCRETIZATION_OPERATORS_WEAK_BOUNDARY_CONDITIONS_H_
23#define EXADG_ACOUSTIC_CONSERVATION_EQUATIONS_SPATIAL_DISCRETIZATION_OPERATORS_WEAK_BOUNDARY_CONDITIONS_H_
25#include <exadg/acoustic_conservation_equations/user_interface/boundary_descriptor.h>
26#include <exadg/functions_and_boundary_conditions/boundary_face_integrator_base.h>
27#include <exadg/functions_and_boundary_conditions/evaluate_functions.h>
28#include <exadg/matrix_free/integrators.h>
35template<
int dim,
typename Number>
36inline DEAL_II_ALWAYS_INLINE
37 dealii::VectorizedArray<Number>
38 calculate_exterior_value_pressure(
unsigned int const q,
39 FaceIntegrator<dim, 1, Number>
const & pressure_m,
40 BoundaryType
const & boundary_type,
41 dealii::types::boundary_id
const boundary_id,
45 if(boundary_type == BoundaryType::PressureDirichlet)
47 auto const g = FunctionEvaluator<0, dim, Number>::value(
48 *boundary_descriptor.pressure_dbc.find(boundary_id)->second,
49 pressure_m.quadrature_point(q),
51 return -pressure_m.get_value(q) + Number{2.0} * g;
53 else if(boundary_type == BoundaryType::VelocityDirichlet or
54 boundary_type == BoundaryType::Admittance)
56 return pressure_m.get_value(q);
60 AssertThrow(
false, dealii::ExcMessage(
"Boundary type of face is invalid or not implemented."));
65template<
int dim,
typename Number>
66inline DEAL_II_ALWAYS_INLINE
67 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>>
68 calculate_exterior_value_velocity(
unsigned int const q,
69 FaceIntegrator<dim, dim, Number>
const & velocity_m,
70 FaceIntegrator<dim, 1, Number>
const & pressure_m,
71 double const speed_of_sound,
72 BoundaryType
const & boundary_type,
73 dealii::types::boundary_id
const boundary_id,
77 if(boundary_type == BoundaryType::VelocityDirichlet)
79 auto const g = FunctionEvaluator<1, dim, Number>::value(
80 *boundary_descriptor.velocity_dbc.find(boundary_id)->second,
81 velocity_m.quadrature_point(q),
83 return -velocity_m.get_value(q) + Number{2.0} * g;
85 else if(boundary_type == BoundaryType::PressureDirichlet)
87 return velocity_m.get_value(q);
89 else if(boundary_type == BoundaryType::Admittance)
91 auto const Y = FunctionEvaluator<0, dim, Number>::value(
92 *boundary_descriptor.admittance_bc.find(boundary_id)->second,
93 velocity_m.quadrature_point(q),
96 auto const n = velocity_m.normal_vector(q);
97 auto const rho_um = velocity_m.get_value(q);
98 auto const pm = pressure_m.get_value(q);
100 auto const rho_um_normal = (rho_um * n) * n;
101 auto const rho_um_tangential = rho_um - rho_um_normal;
104 auto const rho_up_normal = (Number{2.0} * Y * pm / speed_of_sound) * n - rho_um_normal;
107 return rho_up_normal + rho_um_tangential;
111 AssertThrow(
false, dealii::ExcMessage(
"Boundary type of face is invalid or not implemented."));
119template<
int dim,
typename Number>
120class BoundaryFaceIntegratorP :
public BoundaryFaceIntegratorBase<BoundaryDescriptor<dim>, Number>
122 using FaceIntegratorP = FaceIntegrator<dim, 1, Number>;
125 BoundaryFaceIntegratorP(FaceIntegratorP
const & integrator_m_in,
127 : BoundaryFaceIntegratorBase<BoundaryDescriptor<dim>, Number>(integrator_m_in.get_matrix_free(),
128 boundary_descriptor_in),
129 integrator_m(integrator_m_in)
133 inline DEAL_II_ALWAYS_INLINE
134 typename FaceIntegratorP::value_type
135 get_value(
unsigned int const q)
const
137 return calculate_exterior_value_pressure(q,
141 this->boundary_descriptor,
142 this->evaluation_time);
146 FaceIntegratorP
const & integrator_m;
152template<
int dim,
typename Number>
153class BoundaryFaceIntegratorU :
public BoundaryFaceIntegratorBase<BoundaryDescriptor<dim>, Number>
155 using FaceIntegratorU = FaceIntegrator<dim, dim, Number>;
156 using FaceIntegratorP = FaceIntegrator<dim, 1, Number>;
159 BoundaryFaceIntegratorU(FaceIntegratorU
const & integrator_velocity_m_in,
160 FaceIntegratorP
const & integrator_pressure_m_in,
161 double const speed_of_sound,
163 : BoundaryFaceIntegratorBase<BoundaryDescriptor<dim>, Number>(
164 integrator_velocity_m_in.get_matrix_free(),
165 boundary_descriptor_in),
166 integrator_velocity_m(integrator_velocity_m_in),
167 integrator_pressure_m(integrator_pressure_m_in),
172 inline DEAL_II_ALWAYS_INLINE
173 typename FaceIntegratorU::value_type
174 get_value(
unsigned int const q)
const
176 return calculate_exterior_value_velocity(q,
177 integrator_velocity_m,
178 integrator_pressure_m,
182 this->boundary_descriptor,
183 this->evaluation_time);
187 FaceIntegratorU
const & integrator_velocity_m;
188 FaceIntegratorP
const & integrator_pressure_m;
Definition boundary_descriptor.h:64