1 #ifndef VERTEXCFD_CLOSURE_JOULEHEATINGSOURCE_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_JOULEHEATINGSOURCE_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 JouleHeatingSource<EvalType, Traits, NumSpaceDim>::JouleHeatingSource(
17 const panzer::IntegrationRule& ir)
18 : _jh_source(
"VOLUMETRIC_SOURCE_energy", ir.dl_scalar)
19 , _sigma(
"electrical_conductivity", ir.dl_scalar)
22 this->addEvaluatedField(_jh_source);
25 this->addDependentField(_sigma);
26 Utils::addDependentVectorField(*
this,
28 _electric_current_density,
29 "electric_current_density_");
31 this->setName(
"Joule Heating Source Term " + std::to_string(num_space_dim)
36 template<
class EvalType,
class Traits,
int NumSpaceDim>
37 void JouleHeatingSource<EvalType, Traits, NumSpaceDim>::evaluateFields(
38 typename Traits::EvalData workset)
40 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
42 Kokkos::parallel_for(this->getName(), policy, *
this);
46 template<
class EvalType,
class Traits,
int NumSpaceDim>
47 KOKKOS_INLINE_FUNCTION
void
48 JouleHeatingSource<EvalType, Traits, NumSpaceDim>::operator()(
49 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
51 const int cell = team.league_rank();
52 const int num_point = _jh_source.extent(1);
55 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
57 _jh_source(cell, point) = 0.0;
58 for (
int dim = 0; dim < num_space_dim; ++dim)
60 _jh_source(cell, point)
61 += _electric_current_density[dim](cell, point)
62 * _electric_current_density[dim](cell, point)
63 / _sigma(cell, point);