1 #ifndef VERTEXCFD_CLOSURE_MAGNETICPRESSURE_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_MAGNETICPRESSURE_IMPL_HPP
4 #include "utils/VertexCFD_Utils_MagneticLayout.hpp"
6 #include <Panzer_HierarchicParallelism.hpp>
10 namespace ClosureModel
13 template<
class EvalType,
class Traits>
14 MagneticPressure<EvalType, Traits>::MagneticPressure(
15 const panzer::IntegrationRule& ir,
16 const MHDProperties::FullInductionMHDProperties& mhd_props)
17 : _magnetic_permeability(mhd_props.vacuumMagneticPermeability())
18 , _total_magnetic_field(
19 "total_magnetic_field",
20 Utils::buildMagneticLayout(ir.dl_scalar, num_magnetic_field_dim))
21 , _magnetic_pressure(
"magnetic_pressure", ir.dl_scalar)
24 this->addDependentField(_total_magnetic_field);
27 this->addEvaluatedField(_magnetic_pressure);
30 this->setName(
"Magnetic Pressure");
34 template<
class EvalType,
class Traits>
35 void MagneticPressure<EvalType, Traits>::evaluateFields(
36 typename Traits::EvalData workset)
38 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
40 Kokkos::parallel_for(this->getName(), policy, *
this);
44 template<
class EvalType,
class Traits>
45 KOKKOS_INLINE_FUNCTION
void MagneticPressure<EvalType, Traits>::operator()(
46 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
48 const int cell = team.league_rank();
49 const int num_point = _magnetic_pressure.extent(1);
52 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
53 _magnetic_pressure(cell, point) = 0.0;
54 for (
int dim = 0; dim < num_magnetic_field_dim; ++dim)
56 _magnetic_pressure(cell, point)
57 += _total_magnetic_field(cell, point, dim)
58 * _total_magnetic_field(cell, point, dim);
60 _magnetic_pressure(cell, point) /= 2.0 * _magnetic_permeability;