VertexCFD  0.0-dev
VertexCFD_Closure_ConductionTimeStepSize_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_CONDUCTIONTIMESTEPSIZE_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_CONDUCTIONTIMESTEPSIZE_IMPL_HPP
3 
4 #include <Panzer_HierarchicParallelism.hpp>
5 
6 namespace VertexCFD
7 {
8 namespace ClosureModel
9 {
10 //---------------------------------------------------------------------------//
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)
19 {
20  // Add evaluated field
21  this->addEvaluatedField(_local_dt);
22 
23  // Add dependent fields
24  this->addDependentField(_element_length);
25  this->addDependentField(_solid_density);
26  this->addDependentField(_thermal_conductivity);
27  this->addDependentField(_specific_heat);
28 
29  this->setName("Conduction Time-Step Size");
30 }
31 
32 //---------------------------------------------------------------------------//
33 template<class EvalType, class Traits>
35  typename Traits::EvalData workset)
36 {
37  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
38  workset.num_cells);
39  Kokkos::parallel_for(this->getName(), policy, *this);
40 }
41 
42 //---------------------------------------------------------------------------//
43 template<class EvalType, class Traits>
44 KOKKOS_INLINE_FUNCTION void
46  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
47 {
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);
51 
52  Kokkos::parallel_for(
53  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
54  _local_dt(cell, point) = 0.0;
55 
56  for (int dim = 0; dim < num_space_dim; ++dim)
57  {
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));
64  }
65  });
66 }
67 
68 //---------------------------------------------------------------------------//
69 
70 } // end namespace ClosureModel
71 } // end namespace VertexCFD
72 
73 #endif // end VERTEXCFD_CLOSURE_CONDUCTIONTIMESTEPSIZE_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23
VertexCFD::ClosureModel::ConductionTimeStepSize::operator()
KOKKOS_INLINE_FUNCTION void operator()(const Kokkos::TeamPolicy< PHX::exec_space >::member_type &team) const
Functor invoked by Kokkos for each team of cells.
Definition: VertexCFD_Closure_ConductionTimeStepSize_impl.hpp:45
VertexCFD::ClosureModel::ConductionTimeStepSize::ConductionTimeStepSize
ConductionTimeStepSize(const panzer::IntegrationRule &ir)
Constructor.
Definition: VertexCFD_Closure_ConductionTimeStepSize_impl.hpp:12
VertexCFD::ClosureModel::ConductionTimeStepSize::evaluateFields
void evaluateFields(typename Traits::EvalData workset) override
Evaluate the time‑step size on the supplied workset.
Definition: VertexCFD_Closure_ConductionTimeStepSize_impl.hpp:34
VertexCFD::ClosureModel::ConductionTimeStepSize::_local_dt
PHX::MDField< scalar_type, panzer::Cell, panzer::Point > _local_dt
Output field: local time‑step size.
Definition: VertexCFD_Closure_ConductionTimeStepSize.hpp:75
VertexCFD::ClosureModel::ConductionTimeStepSize::scalar_type
typename EvalType::ScalarT scalar_type
Scalar type used by the evaluation type.
Definition: VertexCFD_Closure_ConductionTimeStepSize.hpp:36