1 #ifndef VERTEXCFD_CLOSURE_INCOMPRESSIBLEVISCOUSFLUX_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_INCOMPRESSIBLEVISCOUSFLUX_IMPL_HPP
4 #include "utils/VertexCFD_Utils_VectorField.hpp"
6 #include <Panzer_HierarchicParallelism.hpp>
10 namespace ClosureModel
13 template<
class EvalType,
class Traits,
int NumSpaceDim>
14 IncompressibleViscousFlux<EvalType, Traits, NumSpaceDim>::IncompressibleViscousFlux(
15 const panzer::IntegrationRule& ir,
16 const Teuchos::ParameterList& fluid_params,
17 const Teuchos::ParameterList& user_params,
18 const bool use_turbulence_model,
19 const std::string& flux_prefix,
20 const std::string& gradient_prefix)
21 : _continuity_flux(flux_prefix +
"VISCOUS_FLUX_continuity", ir.dl_vector)
22 , _energy_flux(flux_prefix +
"VISCOUS_FLUX_energy", ir.dl_vector)
23 , _rho(
"density", ir.dl_scalar)
24 , _nu(
"kinematic_viscosity", ir.dl_scalar)
25 , _cp(
"specific_heat_capacity", ir.dl_scalar)
26 , _k(
"thermal_conductivity", ir.dl_scalar)
27 , _grad_press(gradient_prefix +
"GRAD_lagrange_pressure", ir.dl_vector)
28 , _grad_temp(gradient_prefix +
"GRAD_temperature", ir.dl_vector)
29 , _nu_t(flux_prefix +
"turbulent_eddy_viscosity", ir.dl_scalar)
30 , _gamma(std::numeric_limits<double>::quiet_NaN())
31 , _beta(fluid_params.get<double>(
"Artificial compressibility"))
32 , _solve_temp(fluid_params.get<bool>(
"Build Temperature Equation"))
33 , _use_turbulence_model(use_turbulence_model)
34 , _continuity_model_name(fluid_params.isType<std::string>(
"Continuity "
36 ? fluid_params.get<std::string>(
"Continuity "
39 , _is_edac(_continuity_model_name ==
"EDAC" ? true : false)
41 ? (user_params.isType<double>(
"Turbulent Prandtl Number")
42 ? user_params.get<double>(
"Turbulent Prandtl Number")
44 : std::numeric_limits<double>::quiet_NaN())
47 this->addEvaluatedField(_continuity_flux);
48 Utils::addEvaluatedVectorField(*
this, ir.dl_vector, _momentum_flux,
49 flux_prefix +
"VISCOUS_FLUX_"
51 this->addEvaluatedField(_energy_flux);
54 Utils::addDependentVectorField(*
this,
57 gradient_prefix +
"GRAD_velocity_");
59 this->addDependentField(_rho);
60 this->addDependentField(_nu);
63 this->addDependentField(_grad_press);
65 _gamma = fluid_params.get<
double>(
"Heat capacity ratio");
70 this->addDependentField(_k);
71 this->addDependentField(_cp);
72 this->addDependentField(_grad_temp);
75 if (_use_turbulence_model)
77 this->addDependentField(_nu_t);
80 this->setName(
"Incompressible Viscous Flux "
81 + std::to_string(num_space_dim) +
"D");
85 template<
class EvalType,
class Traits,
int NumSpaceDim>
86 void IncompressibleViscousFlux<EvalType, Traits, NumSpaceDim>::evaluateFields(
87 typename Traits::EvalData workset)
89 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
91 Kokkos::parallel_for(this->getName(), policy, *
this);
95 template<
class EvalType,
class Traits,
int NumSpaceDim>
96 KOKKOS_INLINE_FUNCTION
void
97 IncompressibleViscousFlux<EvalType, Traits, NumSpaceDim>::operator()(
98 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
100 const int cell = team.league_rank();
101 const int num_point = _continuity_flux.extent(1);
103 Kokkos::parallel_for(
104 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
106 for (
int i = 0; i < num_space_dim; ++i)
113 _continuity_flux(cell, point, i)
114 = _gamma * _k(cell, point)
115 * _grad_press(cell, point, i)
116 / (_beta * _rho(cell, point) * _cp(cell, point));
120 _continuity_flux(cell, point, i)
121 = _rho(cell, point) * _nu(cell, point)
122 * _grad_press(cell, point, i) / _beta;
128 _continuity_flux(cell, point, i) = 0.0;
134 _energy_flux(cell, point, i)
135 = _k(cell, point) * _grad_temp(cell, point, i);
136 if (_use_turbulence_model)
138 _energy_flux(cell, point, i)
139 += _nu_t(cell, point) * _rho(cell, point)
140 * _cp(cell, point) * _grad_temp(cell, point, i)
146 for (
int j = 0; j < num_space_dim; ++j)
148 _momentum_flux[j](cell, point, i)
149 = _rho(cell, point) * _nu(cell, point)
150 * _grad_velocity[j](cell, point, i);
151 if (_use_turbulence_model)
153 _momentum_flux[j](cell, point, i)
154 += _rho(cell, point) * _nu_t(cell, point)
155 * _grad_velocity[j](cell, point, i);
167 #endif // end VERTEXCFD_CLOSURE_INCOMPRESSIBLEVISCOUSFLUX_IMPL_HPP