1 #ifndef VERTEXCFD_CLOSURE_INDUCTIONCONSTANTSOURCE_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_INDUCTIONCONSTANTSOURCE_IMPL_HPP
4 #include "utils/VertexCFD_Utils_VectorField.hpp"
6 #include <Panzer_HierarchicParallelism.hpp>
12 namespace ClosureModel
15 template<
class EvalType,
class Traits,
int NumSpaceDim>
16 InductionConstantSource<EvalType, Traits, NumSpaceDim>::InductionConstantSource(
17 const panzer::IntegrationRule& ir,
18 const Teuchos::ParameterList& closure_params)
20 const auto ind_input_source
21 = closure_params.get<Teuchos::Array<double>>(
"Induction Source");
23 for (
int dim = 0; dim < num_space_dim; ++dim)
24 _ind_input_source[dim] = ind_input_source[dim];
27 Utils::addEvaluatedVectorField(
28 *
this, ir.dl_scalar, _induction_source,
"CONSTANT_SOURCE_induction_");
30 this->setName(
"Induction Constant Source " + std::to_string(num_space_dim)
35 template<
class EvalType,
class Traits,
int NumSpaceDim>
36 void InductionConstantSource<EvalType, Traits, NumSpaceDim>::evaluateFields(
37 typename Traits::EvalData workset)
39 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
41 Kokkos::parallel_for(this->getName(), policy, *
this);
45 template<
class EvalType,
class Traits,
int NumSpaceDim>
46 KOKKOS_INLINE_FUNCTION
void
47 InductionConstantSource<EvalType, Traits, NumSpaceDim>::operator()(
48 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
50 const int cell = team.league_rank();
51 const int num_point = _induction_source[0].extent(1);
54 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
55 for (
int dim = 0; dim < num_space_dim; ++dim)
57 _induction_source[dim](cell, point) = _ind_input_source[dim];
67 #endif // end VERTEXCFD_CLOSURE_INDUCTIONCONSTANTSOURCE_IMPL_HPP