VertexCFD  0.0-dev
VertexCFD_Closure_JouleHeatingSource_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_JOULEHEATINGSOURCE_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_JOULEHEATINGSOURCE_IMPL_HPP
3 
4 #include "utils/VertexCFD_Utils_VectorField.hpp"
5 
6 #include <Panzer_HierarchicParallelism.hpp>
7 
8 #include <string>
9 
10 namespace VertexCFD
11 {
12 namespace ClosureModel
13 {
14 //---------------------------------------------------------------------------//
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)
20 {
21  // Evaluated fields
22  this->addEvaluatedField(_jh_source);
23 
24  // Dependent fields
25  this->addDependentField(_sigma);
26  Utils::addDependentVectorField(*this,
27  ir.dl_scalar,
28  _electric_current_density,
29  "electric_current_density_");
30 
31  this->setName("Joule Heating Source Term " + std::to_string(num_space_dim)
32  + "D");
33 }
34 
35 //---------------------------------------------------------------------------//
36 template<class EvalType, class Traits, int NumSpaceDim>
37 void JouleHeatingSource<EvalType, Traits, NumSpaceDim>::evaluateFields(
38  typename Traits::EvalData workset)
39 {
40  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
41  workset.num_cells);
42  Kokkos::parallel_for(this->getName(), policy, *this);
43 }
44 
45 //---------------------------------------------------------------------------//
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
50 {
51  const int cell = team.league_rank();
52  const int num_point = _jh_source.extent(1);
53 
54  Kokkos::parallel_for(
55  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
56  // Joule heating source term
57  _jh_source(cell, point) = 0.0;
58  for (int dim = 0; dim < num_space_dim; ++dim)
59  {
60  _jh_source(cell, point)
61  += _electric_current_density[dim](cell, point)
62  * _electric_current_density[dim](cell, point)
63  / _sigma(cell, point);
64  }
65  });
66 }
67 
68 //---------------------------------------------------------------------------//
69 
70 } // end namespace ClosureModel
71 } // end namespace VertexCFD
72 
73 #endif // end
74  // VERTEXCFD_CLOSURE_JOULEHEATINGSOURCE_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23