VertexCFD  0.0-dev
VertexCFD_Closure_TotalMagneticFieldGradient_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_TOTALMAGNETICFIELDGRADIENT_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_TOTALMAGNETICFIELDGRADIENT_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 TotalMagneticFieldGradient<EvalType, Traits, NumSpaceDim>::TotalMagneticFieldGradient(
16  const panzer::IntegrationRule& ir, const std::string& gradient_prefix)
17  : _uniform_external_field(true)
18  , _grad_total_magnetic_field(
19  gradient_prefix + "GRAD_total_magnetic_field",
20  Utils::buildMagneticGradLayout(ir.dl_vector, num_magnetic_field_dim))
21  , _divergence_total_magnetic_field(
22  gradient_prefix + "divergence_total_magnetic_field", ir.dl_scalar)
23 {
24  // Add dependent fields
25  Utils::addDependentVectorField(
26  *this,
27  ir.dl_vector,
28  _grad_induced_magnetic_field,
29  gradient_prefix + "GRAD_induced_magnetic_field_");
30  if (!_uniform_external_field)
31  {
32  Utils::addDependentVectorField(*this,
33  ir.dl_vector,
34  _grad_external_magnetic_field,
35  "GRAD_external_magnetic_field_");
36  }
37 
38  // Add evaluated fields
39  this->addEvaluatedField(_grad_total_magnetic_field);
40  this->addEvaluatedField(_divergence_total_magnetic_field);
41 
42  // Closure model name
43  this->setName("Total Magnetic Field Gradient"
44  + std::to_string(num_space_dim) + "D");
45 }
46 
47 //---------------------------------------------------------------------------//
48 template<class EvalType, class Traits, int NumSpaceDim>
49 void TotalMagneticFieldGradient<EvalType, Traits, NumSpaceDim>::evaluateFields(
50  typename Traits::EvalData workset)
51 {
52  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
53  workset.num_cells);
54  Kokkos::parallel_for(this->getName(), policy, *this);
55 }
56 
57 //---------------------------------------------------------------------------//
58 template<class EvalType, class Traits, int NumSpaceDim>
59 KOKKOS_INLINE_FUNCTION void
60 TotalMagneticFieldGradient<EvalType, Traits, NumSpaceDim>::operator()(
61  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
62 {
63  const int cell = team.league_rank();
64  const int num_point = _grad_total_magnetic_field.extent(1);
65  const int num_grad_dim = _grad_total_magnetic_field.extent(2);
66 
67  Kokkos::parallel_for(
68  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
69  for (int field_dim = 0; field_dim < num_space_dim; ++field_dim)
70  {
71  for (int grad_dim = 0; grad_dim < num_grad_dim; ++grad_dim)
72  {
73  _grad_total_magnetic_field(cell, point, grad_dim, field_dim)
74  = _grad_induced_magnetic_field[field_dim](
75  cell, point, grad_dim);
76  }
77  }
78 
79  for (int field_dim = num_space_dim;
80  field_dim < num_magnetic_field_dim;
81  ++field_dim)
82  {
83  for (int grad_dim = 0; grad_dim < num_grad_dim; ++grad_dim)
84  {
85  _grad_total_magnetic_field(cell, point, grad_dim, field_dim)
86  = 0.0;
87  }
88  }
89 
90  if (!_uniform_external_field)
91  {
92  for (int field_dim = 0; field_dim < num_magnetic_field_dim;
93  ++field_dim)
94  {
95  for (int grad_dim = 0; grad_dim < num_grad_dim; ++grad_dim)
96  {
97  _grad_total_magnetic_field(
98  cell, point, grad_dim, field_dim)
99  += _grad_external_magnetic_field[field_dim](
100  cell, point, grad_dim);
101  }
102  }
103  }
104 
105  _divergence_total_magnetic_field(cell, point) = 0.0;
106  for (int grad_dim = 0; grad_dim < num_grad_dim; ++grad_dim)
107  {
108  _divergence_total_magnetic_field(cell, point)
109  += _grad_total_magnetic_field(
110  cell, point, grad_dim, grad_dim);
111  }
112  });
113 }
114 
115 //---------------------------------------------------------------------------//
116 
117 } // end namespace ClosureModel
118 } // end namespace VertexCFD
119 
120 #endif // end
121  // VERTEXCFD_CLOSURE_TOTALMAGNETICFIELDGRADIENT_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23