1 #ifndef VERTEXCFD_CLOSURE_INCOMPRESSIBLEREALIZABLEKEPSILONEDDYVISCOSITY_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_INCOMPRESSIBLEREALIZABLEKEPSILONEDDYVISCOSITY_IMPL_HPP
4 #include "utils/VertexCFD_Utils_SmoothMath.hpp"
5 #include "utils/VertexCFD_Utils_VectorField.hpp"
7 #include <Panzer_HierarchicParallelism.hpp>
11 namespace ClosureModel
14 template<
class EvalType,
class Traits,
int NumSpaceDim>
15 IncompressibleRealizableKEpsilonEddyViscosity<EvalType, Traits, NumSpaceDim>::
16 IncompressibleRealizableKEpsilonEddyViscosity(
17 const panzer::IntegrationRule& ir)
18 : _turb_kinetic_energy(
"turb_kinetic_energy", ir.dl_scalar)
19 , _turb_dissipation_rate(
"turb_dissipation_rate", ir.dl_scalar)
21 , _nu_t(
"turbulent_eddy_viscosity", ir.dl_scalar)
24 this->addDependentField(_turb_kinetic_energy);
25 this->addDependentField(_turb_dissipation_rate);
27 Utils::addDependentVectorField(
28 *
this, ir.dl_vector, _grad_velocity,
"GRAD_velocity_");
31 this->addEvaluatedField(_nu_t);
35 "Realizable K-Epsilon Incompressible Turbulent Eddy Viscosity "
36 + std::to_string(num_space_dim) +
"D");
40 template<
class EvalType,
class Traits,
int NumSpaceDim>
41 void IncompressibleRealizableKEpsilonEddyViscosity<EvalType, Traits, NumSpaceDim>::
42 evaluateFields(
typename Traits::EvalData workset)
44 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
46 Kokkos::parallel_for(this->getName(), policy, *
this);
50 template<
class EvalType,
class Traits,
int NumSpaceDim>
51 KOKKOS_INLINE_FUNCTION
void
52 IncompressibleRealizableKEpsilonEddyViscosity<EvalType, Traits, NumSpaceDim>::
53 operator()(
const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
60 const int cell = team.league_rank();
61 const int num_point = _nu_t.extent(1);
62 const auto tol = 1.0e-10;
63 const double sqrt_six = sqrt(6.0);
66 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
72 scalar_type Omega2 = 0.0;
75 for (
int i = 0; i < num_space_dim; ++i)
77 for (
int j = 0; j < num_space_dim; ++j)
80 * (_grad_velocity[i](cell, point, j)
81 + _grad_velocity[j](cell, point, i)),
84 * (_grad_velocity[i](cell, point, j)
85 - _grad_velocity[j](cell, point, i)),
87 for (
int k = 0; k < num_space_dim; ++k)
90 * (_grad_velocity[i](cell, point, j)
91 + _grad_velocity[j](cell, point, i))
92 * (_grad_velocity[j](cell, point, k)
93 + _grad_velocity[k](cell, point, j))
94 * (_grad_velocity[k](cell, point, i)
95 + _grad_velocity[i](cell, point, k));
100 S2 = SmoothMath::max(S2, tol, 0.0);
101 Omega2 = SmoothMath::max(Omega2, tol, 0.0);
105 const scalar_type Us = sqrt(S2 + Omega2);
106 const scalar_type phi
108 * acos(SmoothMath::max(
109 SmoothMath::min(sqrt_six * W, 1.0, 0.0), -1.0, 0.0));
110 const scalar_type As = sqrt_six * cos(phi);
111 const scalar_type C_nu
114 + (As * Us * _turb_kinetic_energy(cell, point)
116 _turb_dissipation_rate(cell, point), tol, 0.0)));
120 = C_nu * pow(_turb_kinetic_energy(cell, point), 2.0)
122 _turb_dissipation_rate(cell, point), tol, 0.0);