1 #ifndef VERTEXCFD_CLOSURE_INCOMPRESSIBLESHEARVARIABLES_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_INCOMPRESSIBLESHEARVARIABLES_IMPL_HPP
4 #include "utils/VertexCFD_Utils_SmoothMath.hpp"
5 #include "utils/VertexCFD_Utils_VectorField.hpp"
7 #include <Panzer_HierarchicParallelism.hpp>
13 namespace ClosureModel
16 template<
class EvalType,
class Traits,
int NumSpaceDim>
17 IncompressibleShearVariables<EvalType, Traits, NumSpaceDim>::
18 IncompressibleShearVariables(
const panzer::IntegrationRule& ir)
19 : _tau_w(
"wall_shear_stress", ir.dl_scalar)
20 , _u_tau(
"friction_velocity", ir.dl_scalar)
21 , _rho(
"density", ir.dl_scalar)
22 , _nu(
"kinematic_viscosity", ir.dl_scalar)
23 , _normals(
"Side Normal", ir.dl_vector)
26 this->addEvaluatedField(_tau_w);
27 this->addEvaluatedField(_u_tau);
30 this->addDependentField(_rho);
31 this->addDependentField(_nu);
32 this->addDependentField(_normals);
33 Utils::addDependentVectorField(
34 *
this, ir.dl_vector, _grad_velocity,
"GRAD_velocity_");
36 this->setName(
"Incompressible Shear/Friction Variables "
37 + std::to_string(num_space_dim) +
"D");
41 template<
class EvalType,
class Traits,
int NumSpaceDim>
42 void IncompressibleShearVariables<EvalType, Traits, NumSpaceDim>::evaluateFields(
43 typename Traits::EvalData workset)
45 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
47 Kokkos::parallel_for(this->getName(), policy, *
this);
51 template<
class EvalType,
class Traits,
int NumSpaceDim>
52 KOKKOS_INLINE_FUNCTION
void
53 IncompressibleShearVariables<EvalType, Traits, NumSpaceDim>::operator()(
54 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
56 const int cell = team.league_rank();
57 const int num_point = _u_tau.extent(1);
58 const int num_grad_dim = _normals.extent(2);
63 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
65 _tau_w(cell, point) = 0.0;
66 for (
int i = 0; i < num_space_dim; ++i)
68 scalar_type sum = 0.0;
69 for (
int dim = 0; dim < num_grad_dim; ++dim)
71 sum += _normals(cell, point, dim)
72 * _grad_velocity[i](cell, point, dim);
74 _tau_w(cell, point) += pow(sum, 2.0);
77 = sqrt(SmoothMath::max(_tau_w(cell, point), 1.0E-10, 0.0));
78 _tau_w(cell, point) *= _nu(cell, point);
81 _u_tau(cell, point) = sqrt(_tau_w(cell, point));
82 _tau_w(cell, point) *= _rho(cell, point);
91 #endif // end VERTEXCFD_CLOSURE_INCOMPRESSIBLESHEARVARIABLES_IMPL_HPP