1 #ifndef VERTEXCFD_CLOSURE_INCOMPRESSIBLESUPGFLUX_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_INCOMPRESSIBLESUPGFLUX_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 IncompressibleSUPGFlux<EvalType, Traits, NumSpaceDim>::IncompressibleSUPGFlux(
15 const panzer::IntegrationRule& ir,
16 const Teuchos::ParameterList& fluid_params,
17 const Teuchos::ParameterList& closure_params)
18 : _continuity_flux(
"VISCOUS_FLUX_continuity", ir.dl_vector)
19 , _energy_flux(
"VISCOUS_FLUX_energy", ir.dl_vector)
20 , _use_source(closure_params.isType<std::string>(
"Prefix for source "
22 , _source_prefix(_use_source
23 ? closure_params.get<std::string>(
"Prefix for source "
26 , _grad_lagrange_pressure(
"GRAD_lagrange_pressure", ir.dl_vector)
27 , _rho(
"density", ir.dl_scalar)
28 , _cp(
"specific_heat_capacity", ir.dl_scalar)
29 , _tau_supg_cont(
"tau_supg_continuity", ir.dl_scalar)
30 , _tau_supg_mom(
"tau_supg_momentum", ir.dl_scalar)
31 , _tau_supg_energy(
"tau_supg_energy", ir.dl_scalar)
32 , _dxdt_temperature(
"DXDT_temperature", ir.dl_scalar)
33 , _temperature(
"temperature", ir.dl_scalar)
34 , _energy_source(_source_prefix +
"_SOURCE_energy", ir.dl_scalar)
35 , _grad_temperature(
"GRAD_temperature", ir.dl_vector)
36 , _solve_temp(fluid_params.get<bool>(
"Build Temperature Equation"))
39 this->addContributedField(_continuity_flux);
40 this->addContributedField(_energy_flux);
41 Utils::addContributedVectorField(*
this,
48 Utils::addDependentVectorField(
49 *
this, ir.dl_scalar, _dxdt_velocity,
"DXDT_velocity_");
50 Utils::addDependentVectorField(*
this, ir.dl_scalar, _velocity,
"velocity_");
51 Utils::addDependentVectorField(
52 *
this, ir.dl_vector, _grad_velocity,
"GRAD_velocity_");
55 Utils::addDependentVectorField(*
this,
58 _source_prefix +
"_SOURCE_"
61 this->addDependentField(_rho);
62 this->addDependentField(_grad_lagrange_pressure);
63 this->addDependentField(_tau_supg_cont);
64 this->addDependentField(_tau_supg_mom);
67 this->addDependentField(_cp);
68 this->addDependentField(_dxdt_temperature);
69 this->addDependentField(_temperature);
70 this->addDependentField(_grad_temperature);
72 this->addDependentField(_energy_source);
73 this->addDependentField(_tau_supg_energy);
76 this->setName(
"Incompressible SUPG Flux " + std::to_string(num_space_dim)
81 template<
class EvalType,
class Traits,
int NumSpaceDim>
82 void IncompressibleSUPGFlux<EvalType, Traits, NumSpaceDim>::evaluateFields(
83 typename Traits::EvalData workset)
87 if constexpr (Sacado::IsADType<scalar_type>::value)
89 const int fad_size = Kokkos::dimension_scalar(_tau_supg_mom.get_view());
90 bytes = scratch_view::shmem_size(
91 _tau_supg_mom.extent(1), NUM_TMPS, fad_size);
95 bytes = scratch_view::shmem_size(_tau_supg_mom.extent(1), NUM_TMPS);
98 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
100 policy.set_scratch_size(0, Kokkos::PerTeam(bytes));
101 Kokkos::parallel_for(this->getName(), policy, *
this);
105 template<
class EvalType,
class Traits,
int NumSpaceDim>
106 KOKKOS_INLINE_FUNCTION
void
107 IncompressibleSUPGFlux<EvalType, Traits, NumSpaceDim>::operator()(
108 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
110 const int cell = team.league_rank();
111 const int num_point = _tau_supg_mom.extent(1);
114 scratch_view tmp_field;
115 if constexpr (Sacado::IsADType<scalar_type>::value)
117 const int fad_size = Kokkos::dimension_scalar(_tau_supg_mom.get_view());
119 = scratch_view(team.team_shmem(), num_point, NUM_TMPS, fad_size);
123 tmp_field = scratch_view(team.team_shmem(), num_point, NUM_TMPS);
126 Kokkos::parallel_for(
127 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
128 auto&& vel_res = tmp_field(point, VEL_RES);
129 auto&& temp_res = tmp_field(point, TEMP_RES);
135 temp_res = _rho(cell, point) * _cp(cell, point)
136 * _dxdt_temperature(cell, point);
138 temp_res -= _energy_source(cell, point);
140 for (
int i = 0; i < num_space_dim; ++i)
142 temp_res += _rho(cell, point) * _cp(cell, point)
143 * _velocity[i](cell, point)
144 * _grad_temperature(cell, point, i);
150 for (
int j = 0; j < num_space_dim; ++j)
154 vel_res = _dxdt_velocity[j](cell, point);
155 for (
int i = 0; i < num_space_dim; ++i)
157 vel_res += _velocity[i](cell, point)
158 * _grad_velocity[j](cell, point, i);
160 vel_res *= _rho(cell, point);
161 vel_res += _grad_lagrange_pressure(cell, point, j);
163 vel_res -= _momentum_source[j](cell, point);
166 _continuity_flux(cell, point, j) += _tau_supg_cont(cell, point)
170 for (
int i = 0; i < num_space_dim; ++i)
172 _momentum_flux[j](cell, point, i)
173 += _tau_supg_mom(cell, point) * vel_res
174 * _velocity[i](cell, point);
179 _energy_flux(cell, point, j)
180 += _tau_supg_energy(cell, point) * temp_res
181 * _velocity[j](cell, point);
192 #endif // end VERTEXCFD_CLOSURE_INCOMPRESSIBLESUPGFLUX_IMPL_HPP