1 #ifndef VERTEXCFD_CLOSURE_VISCOUSHEAT_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_VISCOUSHEAT_IMPL_HPP
4 #include "utils/VertexCFD_Utils_VectorField.hpp"
6 #include <Panzer_HierarchicParallelism.hpp>
12 namespace ClosureModel
15 template<
class EvalType,
class Traits,
int NumSpaceDim>
16 IncompressibleViscousHeat<EvalType, Traits, NumSpaceDim>::IncompressibleViscousHeat(
17 const panzer::IntegrationRule& ir)
18 : _viscous_heat_continuity_source(
"VISCOUS_HEAT_SOURCE_continuity",
20 , _viscous_heat_energy_source(
"VISCOUS_HEAT_SOURCE_energy", ir.dl_scalar)
21 , _rho(
"density", ir.dl_scalar)
22 , _nu(
"kinematic_viscosity", ir.dl_scalar)
25 this->addEvaluatedField(_viscous_heat_continuity_source);
26 this->addEvaluatedField(_viscous_heat_energy_source);
28 Utils::addEvaluatedVectorField(*
this,
30 _viscous_heat_momentum_source,
31 "VISCOUS_HEAT_SOURCE_"
35 this->addDependentField(_rho);
36 this->addDependentField(_nu);
37 Utils::addDependentVectorField(
38 *
this, ir.dl_vector, _grad_velocity,
"GRAD_velocity_");
40 this->setName(
"Incompressible Viscous Heat "
41 + std::to_string(num_space_dim) +
"D");
45 template<
class EvalType,
class Traits,
int NumSpaceDim>
46 void IncompressibleViscousHeat<EvalType, Traits, NumSpaceDim>::evaluateFields(
47 typename Traits::EvalData workset)
49 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
51 Kokkos::parallel_for(this->getName(), policy, *
this);
55 template<
class EvalType,
class Traits,
int NumSpaceDim>
56 KOKKOS_INLINE_FUNCTION
void
57 IncompressibleViscousHeat<EvalType, Traits, NumSpaceDim>::operator()(
58 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
60 const int cell = team.league_rank();
61 const int num_point = _viscous_heat_continuity_source.extent(1);
64 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
65 _viscous_heat_continuity_source(cell, point) = 0.0;
68 _viscous_heat_energy_source(cell, point) = 0.0;
70 for (
int i = 0; i < num_space_dim; ++i)
72 for (
int j = 0; j < num_space_dim; ++j)
75 const scalar_type e_ij
77 * (_grad_velocity[j](cell, point, i)
78 + _grad_velocity[i](cell, point, j));
81 _viscous_heat_energy_source(cell, point)
82 += 2.0 * _rho(cell, point) * _nu(cell, point) * e_ij
87 _viscous_heat_momentum_source[i](cell, point) = 0.0;
97 #endif // end VERTEXCFD_CLOSURE_VISCOUSHEAT_IMPL_HPP