VertexCFD  0.0-dev
VertexCFD_Closure_InductionConvectiveFlux_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_INDUCTIONCONVECTIVEFLUX_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_INDUCTIONCONVECTIVEFLUX_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 InductionConvectiveFlux<EvalType, Traits, NumSpaceDim>::InductionConvectiveFlux(
16  const panzer::IntegrationRule& ir,
17  const MHDProperties::FullInductionMHDProperties& mhd_props,
18  const std::string& flux_prefix,
19  const std::string& field_prefix)
20  : _magnetic_correction_potential_flux(
21  flux_prefix + "CONVECTIVE_FLUX_magnetic_correction_potential",
22  ir.dl_vector)
23  , _total_magnetic_field(
24  "total_magnetic_field",
25  Utils::buildMagneticLayout(ir.dl_scalar, num_magnetic_field_dim))
26  , _scalar_magnetic_potential(field_prefix + "scalar_magnetic_potential",
27  ir.dl_scalar)
28  , _magnetic_pressure("magnetic_pressure", ir.dl_scalar)
29  , _solve_magn_corr(mhd_props.buildMagnCorr())
30  , _magnetic_permeability(mhd_props.vacuumMagneticPermeability())
31  , _c_h(mhd_props.hyperbolicDivergenceCleaningSpeed())
32 {
33  // Evaluated fields
34  Utils::addEvaluatedVectorField(*this, ir.dl_vector, _induction_flux,
35  flux_prefix + "CONVECTIVE_FLUX_"
36  "induction_");
37 
38  if (_solve_magn_corr)
39  {
40  this->addEvaluatedField(_magnetic_correction_potential_flux);
41  }
42 
43  // Dependent fields
44  Utils::addDependentVectorField(
45  *this, ir.dl_scalar, _velocity, field_prefix + "velocity_");
46  this->addDependentField(_total_magnetic_field);
47  this->addDependentField(_magnetic_pressure);
48 
49  if (_solve_magn_corr)
50  {
51  this->addDependentField(_scalar_magnetic_potential);
52  }
53 
54  this->setName("Induction Convective Flux " + std::to_string(num_space_dim)
55  + "D");
56 }
57 
58 //---------------------------------------------------------------------------//
59 template<class EvalType, class Traits, int NumSpaceDim>
60 void InductionConvectiveFlux<EvalType, Traits, NumSpaceDim>::evaluateFields(
61  typename Traits::EvalData workset)
62 {
63  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
64  workset.num_cells);
65  Kokkos::parallel_for(this->getName(), policy, *this);
66 }
67 
68 //---------------------------------------------------------------------------//
69 template<class EvalType, class Traits, int NumSpaceDim>
70 KOKKOS_INLINE_FUNCTION void
71 InductionConvectiveFlux<EvalType, Traits, NumSpaceDim>::operator()(
72  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
73 {
74  const int cell = team.league_rank();
75  const int num_point = _induction_flux[0].extent(1);
76 
77  Kokkos::parallel_for(
78  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
79  for (int flux_dim = 0; flux_dim < num_space_dim; ++flux_dim)
80  {
81  for (int vec_dim = 0; vec_dim < num_space_dim; ++vec_dim)
82  {
83  if (vec_dim != flux_dim)
84  {
85  // Set the off-diagonal flux terms for the induction
86  // equation.
87  _induction_flux[vec_dim](cell, point, flux_dim)
88  = _velocity[flux_dim](cell, point)
89  * _total_magnetic_field(cell, point, vec_dim)
90  - _total_magnetic_field(cell, point, flux_dim)
91  * _velocity[vec_dim](cell, point);
92  }
93  }
94 
95  // Set diagonal flux terms for the induction equation, which
96  // are nonzero only when using divergence cleaning.
97  if (_solve_magn_corr)
98  {
99  _induction_flux[flux_dim](cell, point, flux_dim)
100  = _c_h * _scalar_magnetic_potential(cell, point);
101  _magnetic_correction_potential_flux(cell, point, flux_dim)
102  = _c_h * _total_magnetic_field(cell, point, flux_dim);
103  }
104  else
105  {
106  _induction_flux[flux_dim](cell, point, flux_dim) = 0.0;
107  }
108  }
109  });
110 }
111 
112 //---------------------------------------------------------------------------//
113 
114 } // end namespace ClosureModel
115 } // end namespace VertexCFD
116 
117 #endif // end VERTEXCFD_CLOSURE_INDUCTIONCONVECTIVEFLUX_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23