VertexCFD  0.0-dev
VertexCFD_BoundaryState_ElectricCurrentDensityInsulating_impl.hpp
1 #ifndef VERTEXCFD_BOUNDARYSTATE_ELECTRICCURRENTDENSITYINSULATING_IMPL_HPP
2 #define VERTEXCFD_BOUNDARYSTATE_ELECTRICCURRENTDENSITYINSULATING_IMPL_HPP
3 
4 #include "utils/VertexCFD_Utils_VectorField.hpp"
5 
6 #include <Panzer_HierarchicParallelism.hpp>
7 
8 namespace VertexCFD
9 {
10 namespace BoundaryCondition
11 {
12 //---------------------------------------------------------------------------//
13 template<class EvalType, class Traits, int NumSpaceDim>
14 ElectricCurrentDensityInsulating<EvalType, Traits, NumSpaceDim>::
15  ElectricCurrentDensityInsulating(const panzer::IntegrationRule& ir)
16  : _boundary_electric_potential("BOUNDARY_electric_potential", ir.dl_scalar)
17  , _boundary_grad_electric_potential("BOUNDARY_GRAD_electric_potential",
18  ir.dl_vector)
19  , _electric_potential("electric_potential", ir.dl_scalar)
20  , _grad_electric_potential("GRAD_electric_potential", ir.dl_vector)
21  , _normals("Side Normal", ir.dl_vector)
22 {
23  this->addDependentField(_electric_potential);
24  this->addEvaluatedField(_boundary_electric_potential);
25 
26  this->addDependentField(_grad_electric_potential);
27  this->addEvaluatedField(_boundary_grad_electric_potential);
28 
29  this->addDependentField(_normals);
30 
31  Utils::addDependentVectorField(
32  *this, ir.dl_scalar, _boundary_velocity, "BOUNDARY_velocity_");
33  Utils::addDependentVectorField(*this, ir.dl_scalar, _velocity, "velocity_");
34  Utils::addDependentVectorField(
35  *this, ir.dl_scalar, _ext_magn_field, "external_magnetic_field_");
36 
37  this->setName("Boundary State Electric Current Density Insulating");
38 }
39 
40 //---------------------------------------------------------------------------//
41 template<class EvalType, class Traits, int NumSpaceDim>
42 void ElectricCurrentDensityInsulating<EvalType, Traits, NumSpaceDim>::evaluateFields(
43  typename Traits::EvalData workset)
44 {
45  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
46  workset.num_cells);
47  Kokkos::parallel_for(this->getName(), policy, *this);
48 }
49 
50 //---------------------------------------------------------------------------//
51 template<class EvalType, class Traits, int NumSpaceDim>
52 KOKKOS_INLINE_FUNCTION void
53 ElectricCurrentDensityInsulating<EvalType, Traits, NumSpaceDim>::operator()(
54  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
55 {
56  const int cell = team.league_rank();
57  const int num_point = _electric_potential.extent(1);
58  const int num_grad_dim = _normals.extent(2);
59 
60  Kokkos::parallel_for(
61  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
62  // using '_boundary_grad_electric_potential' as temporary variable
63  // to calculate the U \cross B
64 
65  // x-component
66  _boundary_grad_electric_potential(cell, point, 0)
67  = _velocity[1](cell, point) * _ext_magn_field[2](cell, point);
68  // y-component
69  _boundary_grad_electric_potential(cell, point, 1)
70  = -_velocity[0](cell, point) * _ext_magn_field[2](cell, point);
71  if (num_grad_dim == 3)
72  {
73  // contribution to x-y from z-component
74  _boundary_grad_electric_potential(cell, point, 0)
75  -= _velocity[2](cell, point)
76  * _ext_magn_field[1](cell, point);
77  _boundary_grad_electric_potential(cell, point, 1)
78  += _velocity[2](cell, point)
79  * _ext_magn_field[0](cell, point);
80  // z-component
81  _boundary_grad_electric_potential(cell, point, 2)
82  = _velocity[0](cell, point)
83  * _ext_magn_field[1](cell, point);
84  _boundary_grad_electric_potential(cell, point, 2)
85  -= _velocity[1](cell, point)
86  * _ext_magn_field[0](cell, point);
87  }
88 
89  // Compute (\grad(\phi) - U \cross B) \cdot \vec{n}
90  auto&& grad_phi_dot_n = _boundary_electric_potential(cell, point);
91  grad_phi_dot_n = 0.0;
92  for (int dim = 0; dim < num_grad_dim; ++dim)
93  {
94  grad_phi_dot_n
95  += (_grad_electric_potential(cell, point, dim)
96  - _boundary_grad_electric_potential(cell, point, dim))
97  * _normals(cell, point, dim);
98  }
99 
100  // using '_boundary_grad_electric_potential' as temporary variable
101  // to calculate the (U \cross B) - (U_{bc} \cross B)
102 
103  // x-component
104  _boundary_grad_electric_potential(cell, point, 0)
105  -= _boundary_velocity[1](cell, point)
106  * _ext_magn_field[2](cell, point);
107  // y-component
108  _boundary_grad_electric_potential(cell, point, 1)
109  -= -_boundary_velocity[0](cell, point)
110  * _ext_magn_field[2](cell, point);
111  if (num_grad_dim == 3)
112  {
113  // contribution to x-y from z-component
114  _boundary_grad_electric_potential(cell, point, 0)
115  += _boundary_velocity[2](cell, point)
116  * _ext_magn_field[1](cell, point);
117  _boundary_grad_electric_potential(cell, point, 1)
118  -= _boundary_velocity[2](cell, point)
119  * _ext_magn_field[0](cell, point);
120  // z-component
121  _boundary_grad_electric_potential(cell, point, 2)
122  -= _boundary_velocity[0](cell, point)
123  * _ext_magn_field[1](cell, point);
124  _boundary_grad_electric_potential(cell, point, 2)
125  += _boundary_velocity[1](cell, point)
126  * _ext_magn_field[0](cell, point);
127  }
128 
129  // Set and boundary gradients
130  for (int dim = 0; dim < num_grad_dim; ++dim)
131  {
132  _boundary_grad_electric_potential(cell, point, dim)
133  += _grad_electric_potential(cell, point, dim)
134  - grad_phi_dot_n * _normals(cell, point, dim);
135  }
136 
137  _boundary_electric_potential(cell, point)
138  = _electric_potential(cell, point);
139  });
140 }
141 
142 //---------------------------------------------------------------------------//
143 
144 } // end namespace BoundaryCondition
145 } // end namespace VertexCFD
146 
147 #endif // VERTEXCFD_BOUNDARYSTATE_ELECTRICCURRENTDENSITYINSULATING_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23