1 #ifndef VERTEXCFD_CLOSURE_INCOMPRESSIBLEWALEEDDYVISCOSITY_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_INCOMPRESSIBLEWALEEDDYVISCOSITY_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 IncompressibleWALEEddyViscosity<EvalType, Traits, NumSpaceDim>::
16 IncompressibleWALEEddyViscosity(
const panzer::IntegrationRule& ir,
17 const Teuchos::ParameterList& user_params)
18 : _element_length(
"les_element_length", ir.dl_vector)
21 , _k_sgs(
"sgs_kinetic_energy", ir.dl_scalar)
22 , _nu_t(
"turbulent_eddy_viscosity", ir.dl_scalar)
25 if (user_params.isSublist(
"Turbulence Parameters"))
27 Teuchos::ParameterList turb_list
28 = user_params.sublist(
"Turbulence Parameters");
30 if (turb_list.isType<
double>(
"C_w"))
32 _C_w = turb_list.get<
double>(
"C_w");
35 if (turb_list.isType<
double>(
"C_k"))
37 _C_k = turb_list.get<
double>(
"C_k");
42 this->addDependentField(_element_length);
44 Utils::addDependentVectorField(
45 *
this, ir.dl_vector, _grad_velocity,
"GRAD_velocity_");
48 this->addEvaluatedField(_k_sgs);
49 this->addEvaluatedField(_nu_t);
52 this->setName(
"Incompressible WALE Turbulent Eddy Viscosity "
53 + std::to_string(num_space_dim) +
"D");
57 template<
class EvalType,
class Traits,
int NumSpaceDim>
58 void IncompressibleWALEEddyViscosity<EvalType, Traits, NumSpaceDim>::evaluateFields(
59 typename Traits::EvalData workset)
63 if constexpr (Sacado::IsADType<scalar_type>::value)
65 const int fad_size = Kokkos::dimension_scalar(_nu_t.get_view());
66 bytes = scratch_view::shmem_size(_nu_t.extent(1), NUM_TMPS, fad_size);
70 bytes = scratch_view::shmem_size(_nu_t.extent(1), NUM_TMPS);
72 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
74 policy.set_scratch_size(0, Kokkos::PerTeam(bytes));
75 Kokkos::parallel_for(this->getName(), policy, *
this);
79 template<
class EvalType,
class Traits,
int NumSpaceDim>
80 KOKKOS_INLINE_FUNCTION
void
81 IncompressibleWALEEddyViscosity<EvalType, Traits, NumSpaceDim>::operator()(
82 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
87 const int cell = team.league_rank();
88 const int num_point = _nu_t.extent(1);
89 const auto tol = 1.0e-10;
90 const double one_third = 1.0 / 3.0;
92 scratch_view tmp_field;
93 if constexpr (Sacado::IsADType<scalar_type>::value)
95 const int fad_size = Kokkos::dimension_scalar(_nu_t.get_view());
97 = scratch_view(team.team_shmem(), num_point, NUM_TMPS, fad_size);
101 tmp_field = scratch_view(team.team_shmem(), num_point, NUM_TMPS);
104 Kokkos::parallel_for(
105 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
108 auto&& mag_sqr_s = tmp_field(point, MAG_SQR_S);
109 auto&& mag_sqr_w = tmp_field(point, MAG_SQR_W);
114 for (
int i = 0; i < num_space_dim; ++i)
116 for (
int j = 0; j < num_space_dim; ++j)
120 * (_grad_velocity[i](cell, point, j)
121 + _grad_velocity[j](cell, point, i)),
126 * (_grad_velocity[i](cell, point, j)
127 - _grad_velocity[j](cell, point, i)),
131 delta += pow(_element_length(cell, point, i), 2.0);
134 mag_sqr_s = SmoothMath::max(mag_sqr_s, tol, 0.0);
135 mag_sqr_w = SmoothMath::max(mag_sqr_w, tol, 0.0);
136 delta = sqrt(SmoothMath::max(delta, tol, 0.0));
140 auto&& mag_sqr_sd = tmp_field(point, MAG_SQR_SD);
141 auto&& Sd_ij = tmp_field(point, SD_IJ);
144 for (
int i = 0; i < num_space_dim; ++i)
146 for (
int j = 0; j < num_space_dim; ++j)
150 for (
int k = 0; k < num_space_dim; ++k)
154 * (_grad_velocity[i](cell, point, k)
155 + _grad_velocity[k](cell, point, i))
156 * (_grad_velocity[k](cell, point, j)
157 + _grad_velocity[j](cell, point, k));
161 * (_grad_velocity[i](cell, point, k)
162 - _grad_velocity[k](cell, point, i))
163 * (_grad_velocity[k](cell, point, j)
164 - _grad_velocity[j](cell, point, k));
170 Sd_ij -= one_third * (mag_sqr_s - mag_sqr_w);
173 mag_sqr_sd += pow(Sd_ij, 2.0);
177 mag_sqr_sd = SmoothMath::max(mag_sqr_sd, tol, 0.0);
181 = pow(_C_w * delta, 2.0) * pow(mag_sqr_sd, 3.0 / 2.0)
182 / (pow(mag_sqr_s, 5.0 / 2.0) + pow(mag_sqr_sd, 5.0 / 4.0));
185 _k_sgs(cell, point) = pow(_nu_t(cell, point) / (_C_k * delta), 2.0);