VertexCFD  0.0-dev
VertexCFD_Closure_MagneticPressure_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_MAGNETICPRESSURE_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_MAGNETICPRESSURE_IMPL_HPP
3 
4 #include "utils/VertexCFD_Utils_MagneticLayout.hpp"
5 
6 #include <Panzer_HierarchicParallelism.hpp>
7 
8 namespace VertexCFD
9 {
10 namespace ClosureModel
11 {
12 //---------------------------------------------------------------------------//
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)
22 {
23  // Add dependent fields
24  this->addDependentField(_total_magnetic_field);
25 
26  // Add evaluated fields
27  this->addEvaluatedField(_magnetic_pressure);
28 
29  // Closure model name
30  this->setName("Magnetic Pressure");
31 }
32 
33 //---------------------------------------------------------------------------//
34 template<class EvalType, class Traits>
35 void MagneticPressure<EvalType, Traits>::evaluateFields(
36  typename Traits::EvalData workset)
37 {
38  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
39  workset.num_cells);
40  Kokkos::parallel_for(this->getName(), policy, *this);
41 }
42 
43 //---------------------------------------------------------------------------//
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
47 {
48  const int cell = team.league_rank();
49  const int num_point = _magnetic_pressure.extent(1);
50 
51  Kokkos::parallel_for(
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)
55  {
56  _magnetic_pressure(cell, point)
57  += _total_magnetic_field(cell, point, dim)
58  * _total_magnetic_field(cell, point, dim);
59  }
60  _magnetic_pressure(cell, point) /= 2.0 * _magnetic_permeability;
61  });
62 }
63 
64 //---------------------------------------------------------------------------//
65 
66 } // end namespace ClosureModel
67 } // end namespace VertexCFD
68 
69 #endif // end
70  // VERTEXCFD_CLOSURE_MAGNETICPRESSURE_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23