VertexCFD  0.0-dev
VertexCFD_Closure_ConductionTimeDerivative_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_CONDUCTIONTIMEDERIVATIVE_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_CONDUCTIONTIMEDERIVATIVE_IMPL_HPP
3 
4 #include <Panzer_HierarchicParallelism.hpp>
5 
6 namespace VertexCFD
7 {
8 namespace ClosureModel
9 {
10 //---------------------------------------------------------------------------//
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)
18 {
19  // Dependent fields
20  this->addDependentField(_dxdt_temperature);
21  this->addDependentField(_density);
22  this->addDependentField(_specific_heat);
23 
24  // Evaluated field
25  this->addEvaluatedField(_dqdt_energy);
26 
27  this->setName("Conduction Time Derivative");
28 }
29 
30 //---------------------------------------------------------------------------//
31 template<class EvalType, class Traits>
32 void ConductionTimeDerivative<EvalType, Traits>::evaluateFields(
33  typename Traits::EvalData workset)
34 {
35  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
36  workset.num_cells);
37  Kokkos::parallel_for(this->getName(), policy, *this);
38 }
39 
40 //---------------------------------------------------------------------------//
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
45 {
46  const int cell = team.league_rank();
47  const int num_point = _density.extent(1);
48 
49  Kokkos::parallel_for(
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);
54  });
55 }
56 
57 //---------------------------------------------------------------------------//
58 
59 } // end namespace ClosureModel
60 } // namespace VertexCFD
61 
62 #endif // end VERTEXCFD_CLOSURE_CONDUCTIONTIMEDERIVATIVE_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23