1 #ifndef VERTEXCFD_CLOSURE_CONDUCTIONTIMEDERIVATIVE_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_CONDUCTIONTIMEDERIVATIVE_IMPL_HPP
4 #include <Panzer_HierarchicParallelism.hpp>
11 template<
class EvalType,
class Traits>
12 ConductionTimeDerivative<EvalType, Traits>::ConductionTimeDerivative(
13 const panzer::IntegrationRule& ir)
14 : _dqdt_energy(
"DQDT_energy", ir.dl_scalar)
15 , _density(
"solid_density", ir.dl_scalar)
16 , _specific_heat(
"solid_specific_heat_capacity", ir.dl_scalar)
17 , _dxdt_temperature(
"DXDT_temperature", ir.dl_scalar)
20 this->addDependentField(_dxdt_temperature);
21 this->addDependentField(_density);
22 this->addDependentField(_specific_heat);
25 this->addEvaluatedField(_dqdt_energy);
27 this->setName(
"Conduction Time Derivative");
31 template<
class EvalType,
class Traits>
32 void ConductionTimeDerivative<EvalType, Traits>::evaluateFields(
33 typename Traits::EvalData workset)
35 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
37 Kokkos::parallel_for(this->getName(), policy, *
this);
41 template<
class EvalType,
class Traits>
42 KOKKOS_INLINE_FUNCTION
void
43 ConductionTimeDerivative<EvalType, Traits>::operator()(
44 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
46 const int cell = team.league_rank();
47 const int num_point = _density.extent(1);
50 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
51 _dqdt_energy(cell, point) = _density(cell, point)
52 * _specific_heat(cell, point)
53 * _dxdt_temperature(cell, point);
62 #endif // end VERTEXCFD_CLOSURE_CONDUCTIONTIMEDERIVATIVE_IMPL_HPP