1 #ifndef VERTEXCFD_BOUNDARYSTATE_TURBULENCEINLETOUTLET_IMPL_HPP
2 #define VERTEXCFD_BOUNDARYSTATE_TURBULENCEINLETOUTLET_IMPL_HPP
4 #include "utils/VertexCFD_Utils_SmoothMath.hpp"
5 #include "utils/VertexCFD_Utils_VectorField.hpp"
7 #include <Panzer_HierarchicParallelism.hpp>
11 namespace BoundaryCondition
16 template<
class EvalType,
class Traits,
int NumSpaceDim>
17 TurbulenceInletOutlet<EvalType, Traits, NumSpaceDim>::TurbulenceInletOutlet(
18 const panzer::IntegrationRule& ir,
19 const Teuchos::ParameterList& bc_params,
20 const std::string variable_name)
21 : _boundary_variable(
"BOUNDARY_" + variable_name, ir.dl_scalar)
22 , _boundary_grad_variable(
"BOUNDARY_GRAD_" + variable_name, ir.dl_vector)
23 , _variable(variable_name, ir.dl_scalar)
24 , _grad_variable(
"GRAD_" + variable_name, ir.dl_vector)
25 , _normals(
"Side Normal", ir.dl_vector)
26 , _inlet_value(bc_params.get<double>(variable_name +
" Inlet Value"))
29 this->addEvaluatedField(_boundary_variable);
30 this->addEvaluatedField(_boundary_grad_variable);
33 Utils::addDependentVectorField(*
this, ir.dl_scalar, _velocity,
"velocity_");
34 this->addDependentField(_variable);
35 this->addDependentField(_grad_variable);
36 this->addDependentField(_normals);
38 this->setName(variable_name
39 +
" Boundary State Turbulence Model Inlet/Outlet "
40 + std::to_string(num_space_dim) +
"D");
44 template<
class EvalType,
class Traits,
int NumSpaceDim>
45 void TurbulenceInletOutlet<EvalType, Traits, NumSpaceDim>::evaluateFields(
46 typename Traits::EvalData workset)
48 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
50 Kokkos::parallel_for(this->getName(), policy, *
this);
54 template<
class EvalType,
class Traits,
int NumSpaceDim>
55 KOKKOS_INLINE_FUNCTION
void
56 TurbulenceInletOutlet<EvalType, Traits, NumSpaceDim>::operator()(
57 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
59 const int cell = team.league_rank();
60 const int num_point = _grad_variable.extent(1);
61 const int num_grad_dim = _grad_variable.extent(2);
62 const double smooth_ramp = 1.0e-8;
65 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
67 scalar_type vel_dot_n = 0.0;
68 for (
int dim = 0; dim < num_grad_dim; ++dim)
70 vel_dot_n += _velocity[dim](cell, point)
71 * _normals(cell, point, dim);
75 const scalar_type outlet
76 = SmoothMath::ramp(vel_dot_n, -smooth_ramp, smooth_ramp);
79 _boundary_variable(cell, point) = (1.0 - outlet) * _inlet_value
80 + outlet * _variable(cell, point);
83 for (
int d = 0; d < num_grad_dim; ++d)
85 _boundary_grad_variable(cell, point, d)
86 = _grad_variable(cell, point, d);
96 #endif // end VERTEXCFD_BOUNDARYSTATE_TURBULENCEINLETOUTLET_IMPL_HPP