1 #ifndef VERTEXCFD_CLOSURE_CONDUCTIONVOLUMETRICSOURCE_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_CONDUCTIONVOLUMETRICSOURCE_IMPL_HPP
4 #include <Panzer_HierarchicParallelism.hpp>
5 #include <Panzer_Workset_Utilities.hpp>
7 #include <Teuchos_StandardParameterEntryValidators.hpp>
11 namespace ClosureModel
14 template<
class EvalType,
class Traits>
16 const panzer::IntegrationRule& ir,
17 const Teuchos::ParameterList& closure_params)
18 : _source(
"SOURCE_energy", ir.dl_scalar)
19 , _num_space_dim(ir.spatial_dimension)
20 , _ir_degree(ir.cubature_degree)
21 , _q(closure_params.get<double>(
"Volumetric Heat Source Value"))
22 , _x_left(std::numeric_limits<double>::signaling_NaN())
23 , _x_right(std::numeric_limits<double>::signaling_NaN())
24 , _heat_source_type(HeatSourceType::constant)
27 this->addEvaluatedField(
_source);
30 if (closure_params.isType<std::string>(
"Heat Source Type"))
32 const auto type_validator = Teuchos::rcp(
33 new Teuchos::StringToIntegralParameterEntryValidator<HeatSourceType>(
34 Teuchos::tuple<std::string>(
"constant",
"xlinear"),
36 _heat_source_type = type_validator->getIntegralValue(
37 closure_params.get<std::string>(
"Heat Source Type"));
41 if (_heat_source_type == HeatSourceType::xlinear)
43 _x_left = closure_params.get<
double>(
"X-value of the left boundary");
44 _x_right = closure_params.get<
double>(
"X-value of the right boundary");
48 this->setName(
"Conduction Volumetric Heat Source "
49 + std::to_string(_num_space_dim) +
"D");
53 template<
class EvalType,
class Traits>
55 typename Traits::SetupData sd, PHX::FieldManager<Traits>&)
57 _ir_index = panzer::getIntegrationRuleIndex(_ir_degree, (*sd.worksets_)[0]);
61 template<
class EvalType,
class Traits>
63 typename Traits::EvalData workset)
65 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
68 _ip_coords = workset.int_rules[_ir_index]->ip_coordinates;
69 Kokkos::parallel_for(this->getName(), policy, *
this);
73 template<
class EvalType,
class Traits>
74 KOKKOS_INLINE_FUNCTION
void
76 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
78 const int cell = team.league_rank();
79 const int num_point = _source.extent(1);
81 Kokkos::parallel_for(Kokkos::TeamThreadRange(team, 0, num_point),
82 [&](
const int point) {
83 if (_heat_source_type == HeatSourceType::constant)
85 _source(cell, point) = _q;
89 const double x = _ip_coords(cell, point, 0);
90 _source(cell, point) = _q * (x - _x_left)
91 / (_x_right - _x_left);
101 #endif // end VERTEXCFD_CLOSURE_CONDUCTIONVOLUMETRICSOURCE_IMPL_HPP