VertexCFD  0.0-dev
VertexCFD_Closure_TotalMagneticField_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_TOTALMAGNETICFIELD_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_TOTALMAGNETICFIELD_IMPL_HPP
3 
4 #include "utils/VertexCFD_Utils_MagneticLayout.hpp"
5 #include "utils/VertexCFD_Utils_VectorField.hpp"
6 
7 #include <Panzer_HierarchicParallelism.hpp>
8 
9 namespace VertexCFD
10 {
11 namespace ClosureModel
12 {
13 //---------------------------------------------------------------------------//
14 template<class EvalType, class Traits, int NumSpaceDim>
15 TotalMagneticField<EvalType, Traits, NumSpaceDim>::TotalMagneticField(
16  const panzer::IntegrationRule& ir, const std::string& field_prefix)
17  : _total_magnetic_field(
18  "total_magnetic_field",
19  Utils::buildMagneticLayout(ir.dl_scalar, num_magnetic_field_dim))
20 {
21  // Add dependent fields
22  Utils::addDependentVectorField(*this,
23  ir.dl_scalar,
24  _induced_magnetic_field,
25  field_prefix + "induced_magnetic_field_");
26  Utils::addDependentVectorField(*this,
27  ir.dl_scalar,
28  _external_magnetic_field,
29  "external_magnetic_field_");
30  // Add evaluated fields
31  this->addEvaluatedField(_total_magnetic_field);
32 
33  // Closure model name
34  this->setName("Total Magnetic Field " + std::to_string(num_space_dim)
35  + "D");
36 }
37 
38 //---------------------------------------------------------------------------//
39 template<class EvalType, class Traits, int NumSpaceDim>
40 void TotalMagneticField<EvalType, Traits, NumSpaceDim>::evaluateFields(
41  typename Traits::EvalData workset)
42 {
43  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
44  workset.num_cells);
45  Kokkos::parallel_for(this->getName(), policy, *this);
46 }
47 
48 //---------------------------------------------------------------------------//
49 template<class EvalType, class Traits, int NumSpaceDim>
50 KOKKOS_INLINE_FUNCTION void
51 TotalMagneticField<EvalType, Traits, NumSpaceDim>::operator()(
52  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
53 {
54  const int cell = team.league_rank();
55  const int num_point = _total_magnetic_field.extent(1);
56 
57  Kokkos::parallel_for(
58  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
59  for (int dim = 0; dim < num_space_dim; ++dim)
60  {
61  _total_magnetic_field(cell, point, dim)
62  = _induced_magnetic_field[dim](cell, point)
63  + _external_magnetic_field[dim](cell, point);
64  }
65 
66  for (int dim = num_space_dim; dim < num_magnetic_field_dim; ++dim)
67  {
68  _total_magnetic_field(cell, point, dim)
69  = _external_magnetic_field[dim](cell, point);
70  }
71  });
72 }
73 
74 //---------------------------------------------------------------------------//
75 
76 } // end namespace ClosureModel
77 } // end namespace VertexCFD
78 
79 #endif // end
80  // VERTEXCFD_CLOSURE_TOTALMAGNETICFIELD_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23