1 #ifndef VERTEXCFD_CLOSURE_CONDUCTIONTHERMALCONDUCTIVITY_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_CONDUCTIONTHERMALCONDUCTIVITY_IMPL_HPP
4 #include <Panzer_HierarchicParallelism.hpp>
6 #include <Teuchos_StandardParameterEntryValidators.hpp>
10 namespace ClosureModel
13 template<
class EvalType,
class Traits>
15 const panzer::IntegrationRule& ir,
16 const Teuchos::ParameterList& closure_params)
17 : _thermal_conductivity(
"thermal_conductivity", ir.dl_scalar)
18 , _temperature(
"temperature", ir.dl_scalar)
19 , _num_space_dim(ir.spatial_dimension)
20 , _k0(std::numeric_limits<double>::quiet_NaN())
21 , _therm_cond_type(ThermCondType::constant)
27 if (closure_params.isType<std::string>(
"Thermal Conductivity Type"))
29 const auto type_validator = Teuchos::rcp(
30 new Teuchos::StringToIntegralParameterEntryValidator<ThermCondType>(
31 Teuchos::tuple<std::string>(
"constant",
"inverse_proportional"),
33 _therm_cond_type = type_validator->getIntegralValue(
34 closure_params.get<std::string>(
"Thermal Conductivity Type"));
38 if (_therm_cond_type == ThermCondType::constant)
40 _k0 = closure_params.get<
double>(
"Thermal Conductivity Value");
44 _k0 = closure_params.get<
double>(
"Thermal Conductivity Coefficient");
48 this->addDependentField(_temperature);
50 this->setName(
"Conduction Thermal Conductivity "
51 + std::to_string(_num_space_dim) +
"D");
55 template<
class EvalType,
class Traits>
57 typename Traits::EvalData workset)
59 auto policy = panzer::HP::inst().teamPolicy<
scalar_type, PHX::Device>(
61 Kokkos::parallel_for(this->getName(), policy, *
this);
65 template<
class EvalType,
class Traits>
66 KOKKOS_INLINE_FUNCTION
void
68 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
70 const int cell = team.league_rank();
71 const int num_point = _thermal_conductivity.extent(1);
73 Kokkos::parallel_for(Kokkos::TeamThreadRange(team, 0, num_point),
74 [&](
const int point) {
75 if (_therm_cond_type == ThermCondType::constant)
77 _thermal_conductivity(cell, point) = _k0;
81 _thermal_conductivity(cell, point)
82 = _k0 / _temperature(cell, point);
92 #endif // end VERTEXCFD_CLOSURE_CONDUCTIONTHERMALCONDUCTIVITY_IMPL_HPP