1 #ifndef VERTEXCFD_BOUNDARYSTATE_FULLINDUCTIONFIXED_IMPL_HPP
2 #define VERTEXCFD_BOUNDARYSTATE_FULLINDUCTIONFIXED_IMPL_HPP
4 #include "utils/VertexCFD_Utils_VectorField.hpp"
6 #include <Panzer_HierarchicParallelism.hpp>
10 namespace BoundaryCondition
19 template<
class EvalType,
class Traits,
int NumSpaceDim>
20 FullInductionFixed<EvalType, Traits, NumSpaceDim>::FullInductionFixed(
21 const panzer::IntegrationRule& ir,
22 const Teuchos::ParameterList& bc_params,
23 const MHDProperties::FullInductionMHDProperties& mhd_props)
24 : _boundary_scalar_magnetic_potential(
"BOUNDARY_scalar_magnetic_potential",
26 , _build_magn_corr(mhd_props.buildMagnCorr())
27 , _dirichlet_scalar_magn_pot(bc_params.isType<double>(
"scalar_magnetic_"
29 , _bnd_scalar_magn_pot(std::numeric_limits<double>::signaling_NaN())
30 , _scalar_magnetic_potential(
"scalar_magnetic_potential", ir.dl_scalar)
32 for (
int dim = 0; dim < num_space_dim; ++dim)
34 const std::string magn_string =
"induced_magnetic_field_"
35 + std::to_string(dim);
36 _bnd_magn_field[dim] = bc_params.get<
double>(magn_string);
40 Utils::addEvaluatedVectorField(*
this,
42 _boundary_induced_magnetic_field,
43 "BOUNDARY_induced_magnetic_field_");
45 Utils::addEvaluatedVectorField(*
this,
47 _boundary_grad_induced_magnetic_field,
48 "BOUNDARY_GRAD_induced_magnetic_field_");
50 Utils::addDependentVectorField(*
this,
52 _grad_induced_magnetic_field,
53 "GRAD_induced_magnetic_field_");
57 this->addEvaluatedField(_boundary_scalar_magnetic_potential);
58 if (_dirichlet_scalar_magn_pot)
61 = bc_params.get<
double>(
"scalar_magnetic_potential");
65 this->addDependentField(_scalar_magnetic_potential);
69 this->setName(
"Boundary State Full Induction Fixed "
70 + std::to_string(num_space_dim) +
"D");
74 template<
class EvalType,
class Traits,
int NumSpaceDim>
75 void FullInductionFixed<EvalType, Traits, NumSpaceDim>::evaluateFields(
76 typename Traits::EvalData workset)
78 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
80 Kokkos::parallel_for(this->getName(), policy, *
this);
84 template<
class EvalType,
class Traits,
int NumSpaceDim>
85 KOKKOS_INLINE_FUNCTION
void
86 FullInductionFixed<EvalType, Traits, NumSpaceDim>::operator()(
87 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
89 const int cell = team.league_rank();
90 const int num_point = _boundary_induced_magnetic_field[0].extent(1);
91 const int num_grad_dim = _boundary_grad_induced_magnetic_field[0].extent(2);
94 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
96 for (
int dim = 0; dim < num_space_dim; ++dim)
98 _boundary_induced_magnetic_field[dim](cell, point)
99 = _bnd_magn_field[dim];
102 if (_build_magn_corr)
104 if (_dirichlet_scalar_magn_pot)
106 _boundary_scalar_magnetic_potential(cell, point)
107 = _bnd_scalar_magn_pot;
111 _boundary_scalar_magnetic_potential(cell, point)
112 = _scalar_magnetic_potential(cell, point);
117 for (
int d = 0; d < num_grad_dim; ++d)
119 for (
int field_dim = 0; field_dim < num_space_dim; ++field_dim)
121 _boundary_grad_induced_magnetic_field[field_dim](
123 = _grad_induced_magnetic_field[field_dim](
135 #endif // end VERTEXCFD_BOUNDARYSTATE_FULLINDUCTIONFIXED_IMPL_HPP