1 #ifndef VERTEXCFD_CLOSURE_TOTALMAGNETICFIELDGRADIENT_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_TOTALMAGNETICFIELDGRADIENT_IMPL_HPP
4 #include "utils/VertexCFD_Utils_MagneticLayout.hpp"
5 #include "utils/VertexCFD_Utils_VectorField.hpp"
7 #include <Panzer_HierarchicParallelism.hpp>
11 namespace ClosureModel
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)
25 Utils::addDependentVectorField(
28 _grad_induced_magnetic_field,
29 gradient_prefix +
"GRAD_induced_magnetic_field_");
30 if (!_uniform_external_field)
32 Utils::addDependentVectorField(*
this,
34 _grad_external_magnetic_field,
35 "GRAD_external_magnetic_field_");
39 this->addEvaluatedField(_grad_total_magnetic_field);
40 this->addEvaluatedField(_divergence_total_magnetic_field);
43 this->setName(
"Total Magnetic Field Gradient"
44 + std::to_string(num_space_dim) +
"D");
48 template<
class EvalType,
class Traits,
int NumSpaceDim>
49 void TotalMagneticFieldGradient<EvalType, Traits, NumSpaceDim>::evaluateFields(
50 typename Traits::EvalData workset)
52 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
54 Kokkos::parallel_for(this->getName(), policy, *
this);
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
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);
68 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
69 for (
int field_dim = 0; field_dim < num_space_dim; ++field_dim)
71 for (
int grad_dim = 0; grad_dim < num_grad_dim; ++grad_dim)
73 _grad_total_magnetic_field(cell, point, grad_dim, field_dim)
74 = _grad_induced_magnetic_field[field_dim](
75 cell, point, grad_dim);
79 for (
int field_dim = num_space_dim;
80 field_dim < num_magnetic_field_dim;
83 for (
int grad_dim = 0; grad_dim < num_grad_dim; ++grad_dim)
85 _grad_total_magnetic_field(cell, point, grad_dim, field_dim)
90 if (!_uniform_external_field)
92 for (
int field_dim = 0; field_dim < num_magnetic_field_dim;
95 for (
int grad_dim = 0; grad_dim < num_grad_dim; ++grad_dim)
97 _grad_total_magnetic_field(
98 cell, point, grad_dim, field_dim)
99 += _grad_external_magnetic_field[field_dim](
100 cell, point, grad_dim);
105 _divergence_total_magnetic_field(cell, point) = 0.0;
106 for (
int grad_dim = 0; grad_dim < num_grad_dim; ++grad_dim)
108 _divergence_total_magnetic_field(cell, point)
109 += _grad_total_magnetic_field(
110 cell, point, grad_dim, grad_dim);