1 #ifndef VERTEXCFD_BOUNDARYSTATE_TURBULENCEKEPSILONWALLFUNCTION_IMPL_HPP
2 #define VERTEXCFD_BOUNDARYSTATE_TURBULENCEKEPSILONWALLFUNCTION_IMPL_HPP
4 #include "utils/VertexCFD_Utils_SmoothMath.hpp"
5 #include "utils/VertexCFD_Utils_VectorField.hpp"
7 #include <Panzer_HierarchicParallelism.hpp>
11 namespace BoundaryCondition
14 template<
class EvalType,
class Traits,
int NumSpaceDim>
15 TurbulenceKEpsilonWallFunction<EvalType, Traits, NumSpaceDim>::
16 TurbulenceKEpsilonWallFunction(
const panzer::IntegrationRule& ir,
17 const Teuchos::ParameterList& bc_params)
18 : _boundary_k(
"BOUNDARY_turb_kinetic_energy", ir.dl_scalar)
19 , _boundary_e(
"BOUNDARY_turb_dissipation_rate", ir.dl_scalar)
20 , _boundary_grad_k(
"BOUNDARY_GRAD_turb_kinetic_energy", ir.dl_vector)
21 , _boundary_grad_e(
"BOUNDARY_GRAD_turb_dissipation_rate", ir.dl_vector)
22 , _boundary_u_tau(
"BOUNDARY_friction_velocity", ir.dl_scalar)
23 , _boundary_y_plus(
"BOUNDARY_y_plus", ir.dl_scalar)
24 , _wall_func_nu_t(
"wall_func_turbulent_eddy_viscosity", ir.dl_scalar)
25 , _k(
"turb_kinetic_energy", ir.dl_scalar)
26 , _e(
"turb_dissipation_rate", ir.dl_scalar)
27 , _grad_k(
"GRAD_turb_kinetic_energy", ir.dl_vector)
28 , _grad_e(
"GRAD_turb_dissipation_rate", ir.dl_vector)
29 , _normals(
"Side Normal", ir.dl_vector)
30 , _nu(
"kinematic_viscosity", ir.dl_scalar)
37 if (bc_params.isType<std::string>(
"Epsilon Condition Type"))
39 const std::string bc_type
40 = bc_params.get<std::string>(
"Epsilon Condition Type");
42 if (bc_type ==
"Neumann")
46 else if (bc_type ==
"Dirichlet")
52 std::string msg =
"Unknown Epsilon Condition Type " + bc_type;
53 msg +=
"\nPlease choose from Dirichlet (default) or Neumann.\n";
55 throw std::runtime_error(msg);
60 this->addEvaluatedField(_boundary_k);
61 this->addEvaluatedField(_boundary_e);
62 this->addEvaluatedField(_boundary_grad_k);
63 this->addEvaluatedField(_boundary_grad_e);
64 this->addEvaluatedField(_boundary_u_tau);
65 this->addEvaluatedField(_boundary_y_plus);
66 this->addEvaluatedField(_wall_func_nu_t);
69 this->addDependentField(_k);
70 this->addDependentField(_e);
71 this->addDependentField(_grad_k);
72 this->addDependentField(_grad_e);
73 this->addDependentField(_normals);
74 this->addDependentField(_nu);
75 Utils::addDependentVectorField(*
this, ir.dl_scalar, _velocity,
"velocity_");
77 this->setName(
"Boundary State Turbulence Model K-Epsilon Wall Function "
78 + std::to_string(_num_grad_dim) +
"D");
82 template<
class EvalType,
class Traits,
int NumSpaceDim>
83 void TurbulenceKEpsilonWallFunction<EvalType, Traits, NumSpaceDim>::evaluateFields(
84 typename Traits::EvalData workset)
86 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
88 Kokkos::parallel_for(this->getName(), policy, *
this);
92 template<
class EvalType,
class Traits,
int NumSpaceDim>
93 KOKKOS_INLINE_FUNCTION
void
94 TurbulenceKEpsilonWallFunction<EvalType, Traits, NumSpaceDim>::operator()(
95 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
97 const int cell = team.league_rank();
98 const int num_point = _grad_k.extent(1);
99 const double max_tol = 1e-10;
103 Kokkos::parallel_for(
104 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
106 scalar_type mag_u = 0.0;
108 for (
int d = 0; d < num_space_dim; d++)
110 mag_u += pow(_velocity[d](cell, point), 2.0);
113 mag_u = pow(SmoothMath::max(mag_u, max_tol, 0.0), 0.5);
116 _boundary_k(cell, point) = _k(cell, point);
118 for (
int d = 0; d < num_space_dim; ++d)
120 _boundary_grad_k(cell, point, d) = _grad_k(cell, point, d);
123 for (
int grad_dim = 0; grad_dim < num_space_dim; ++grad_dim)
125 _boundary_grad_k(cell, point, d)
126 -= _grad_k(cell, point, grad_dim)
127 * _normals(cell, point, grad_dim)
128 * _normals(cell, point, d);
134 const scalar_type u_tau = SmoothMath::max(
136 * pow(SmoothMath::max(_k(cell, point), max_tol, 0.0), 0.5),
140 const scalar_type nu_t_w = _kappa * _yp_tr * _nu(cell, point);
146 _boundary_e(cell, point) = _e(cell, point);
149 for (
int d = 0; d < num_space_dim; ++d)
151 _boundary_grad_e(cell, point, d)
152 = _grad_e(cell, point, d)
153 + (_kappa * pow(u_tau, 5.0) / pow(nu_t_w, 2.0))
154 * _normals(cell, point, d);
157 for (
int grad_dim = 0; grad_dim < num_space_dim; ++grad_dim)
159 _boundary_grad_e(cell, point, d)
160 -= _grad_e(cell, point, grad_dim)
161 * _normals(cell, point, grad_dim)
162 * _normals(cell, point, d);
167 _wall_func_nu_t(cell, point) = nu_t_w;
172 _boundary_e(cell, point) = pow(u_tau, 4.0) / nu_t_w;
175 for (
int d = 0; d < num_space_dim; ++d)
177 _boundary_grad_e(cell, point, d) = _grad_e(cell, point, d);
181 _wall_func_nu_t(cell, point)
182 = _C_mu * pow(_k(cell, point), 2.0)
183 / SmoothMath::max(_e(cell, point), max_tol, 0.0);
188 _boundary_u_tau(cell, point) = u_tau;
189 _boundary_y_plus(cell, point) = _yp_tr;
200 #endif // end VERTEXCFD_BOUNDARYSTATE_TURBULENCEKEPSILONWALLFUNCTION_IMPL_HPP