1 #ifndef VERTEXCFD_CLOSURE_CONDUCTIONTIMESTEPSIZE_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_CONDUCTIONTIMESTEPSIZE_IMPL_HPP
4 #include <Panzer_HierarchicParallelism.hpp>
11 template<
class EvalType,
class Traits>
13 const panzer::IntegrationRule& ir)
14 : _local_dt(
"local_dt", ir.dl_scalar)
15 , _thermal_conductivity(
"thermal_conductivity", ir.dl_scalar)
16 , _specific_heat(
"solid_specific_heat_capacity", ir.dl_scalar)
17 , _solid_density(
"solid_density", ir.dl_scalar)
18 , _element_length(
"element_length", ir.dl_vector)
24 this->addDependentField(_element_length);
25 this->addDependentField(_solid_density);
26 this->addDependentField(_thermal_conductivity);
27 this->addDependentField(_specific_heat);
29 this->setName(
"Conduction Time-Step Size");
33 template<
class EvalType,
class Traits>
35 typename Traits::EvalData workset)
37 auto policy = panzer::HP::inst().teamPolicy<
scalar_type, PHX::Device>(
39 Kokkos::parallel_for(this->getName(), policy, *
this);
43 template<
class EvalType,
class Traits>
44 KOKKOS_INLINE_FUNCTION
void
46 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
48 const int cell = team.league_rank();
49 const int num_point = _local_dt.extent(1);
50 const int num_space_dim = _element_length.extent(2);
53 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
54 _local_dt(cell, point) = 0.0;
56 for (
int dim = 0; dim < num_space_dim; ++dim)
58 _local_dt(cell, point)
59 += _element_length(cell, point, dim)
60 * _element_length(cell, point, dim)
61 * _solid_density(cell, point)
62 * _specific_heat(cell, point)
63 / (2.0 * _thermal_conductivity(cell, point));
73 #endif // end VERTEXCFD_CLOSURE_CONDUCTIONTIMESTEPSIZE_IMPL_HPP