VertexCFD  0.0-dev
VertexCFD_Closure_LorentzForce_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_LORENTZFORCE_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_LORENTZFORCE_IMPL_HPP
3 
4 #include "utils/VertexCFD_Utils_VectorField.hpp"
5 
6 #include <Panzer_HierarchicParallelism.hpp>
7 
8 #include <string>
9 
10 namespace VertexCFD
11 {
12 namespace ClosureModel
13 {
14 //---------------------------------------------------------------------------//
15 template<class EvalType, class Traits, int NumSpaceDim>
16 LorentzForce<EvalType, Traits, NumSpaceDim>::LorentzForce(
17  const panzer::IntegrationRule& ir)
18  : _sigma("electrical_conductivity", ir.dl_scalar)
19  , _grad_electric_potential("GRAD_electric_potential", ir.dl_vector)
20 {
21  // Evaluated fields
22  Utils::addEvaluatedVectorField(
23  *this, ir.dl_scalar, _lorentz_force, "VOLUMETRIC_SOURCE_momentum_");
24 
25  // Dependent fields
26  this->addDependentField(_sigma);
27  this->addDependentField(_grad_electric_potential);
28  Utils::addDependentVectorField(*this, ir.dl_scalar, _velocity, "velocity_");
29  Utils::addDependentVectorField(
30  *this, ir.dl_scalar, _ext_magn_field, "external_magnetic_field_");
31 
32  this->setName("Electric Potential Lorentz Force "
33  + std::to_string(num_space_dim) + "D");
34 }
35 
36 //---------------------------------------------------------------------------//
37 template<class EvalType, class Traits, int NumSpaceDim>
38 void LorentzForce<EvalType, Traits, NumSpaceDim>::evaluateFields(
39  typename Traits::EvalData workset)
40 {
41  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
42  workset.num_cells);
43  Kokkos::parallel_for(this->getName(), policy, *this);
44 }
45 
46 //---------------------------------------------------------------------------//
47 template<class EvalType, class Traits, int NumSpaceDim>
48 KOKKOS_INLINE_FUNCTION void
49 LorentzForce<EvalType, Traits, NumSpaceDim>::operator()(
50  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
51 {
52  const int cell = team.league_rank();
53  const int num_point = _grad_electric_potential.extent(1);
54  const int num_grad_dim = _grad_electric_potential.extent(2);
55 
56  Kokkos::parallel_for(
57  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
58  // Cross product term
59  _lorentz_force[0](cell, point)
60  = -_grad_electric_potential(cell, point, 1)
61  * _ext_magn_field[2](cell, point);
62 
63  _lorentz_force[1](cell, point)
64  = _grad_electric_potential(cell, point, 0)
65  * _ext_magn_field[2](cell, point);
66 
67  if (num_grad_dim == 3)
68  {
69  // x-component
70  _lorentz_force[0](cell, point)
71  += _grad_electric_potential(cell, point, 2)
72  * _ext_magn_field[1](cell, point);
73  // y-component
74  _lorentz_force[1](cell, point)
75  -= _grad_electric_potential(cell, point, 2)
76  * _ext_magn_field[0](cell, point);
77  // z-component
78  _lorentz_force[2](cell, point)
79  = -_grad_electric_potential(cell, point, 0)
80  * _ext_magn_field[1](cell, point)
81  + _grad_electric_potential(cell, point, 1)
82  * _ext_magn_field[0](cell, point);
83  }
84 
85  // Compute B norm and \vec{u} \cdot \vec{B}. NOTE: the external
86  // magnetic field is assumed of the same dimension as the mesh in
87  // this implementation.
88  scalar_type B2 = 0.0;
89  scalar_type B_dot_u = 0.0;
90  for (int dim = 0; dim < num_grad_dim; ++dim)
91  {
92  B2 += _ext_magn_field[dim](cell, point)
93  * _ext_magn_field[dim](cell, point);
94  B_dot_u += _ext_magn_field[dim](cell, point)
95  * _velocity[dim](cell, point);
96  }
97 
98  // Adding dot product terms
99  for (int dim = 0; dim < num_grad_dim; ++dim)
100  {
101  _lorentz_force[dim](cell, point)
102  += _ext_magn_field[dim](cell, point) * B_dot_u;
103  _lorentz_force[dim](cell, point)
104  -= B2 * _velocity[dim](cell, point);
105  _lorentz_force[dim](cell, point) *= _sigma(cell, point);
106  }
107  });
108 }
109 
110 //---------------------------------------------------------------------------//
111 
112 } // end namespace ClosureModel
113 } // end namespace VertexCFD
114 
115 #endif // end
116  // VERTEXCFD_CLOSURE_LORENTZFORCE_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23