VertexCFD  0.0-dev
VertexCFD_BoundaryState_ElectricPotentialInsulatingWall_impl.hpp
1 #ifndef VERTEXCFD_BOUNDARYSTATE_ELECTRICPOTENTIALINSULATINGWALL_IMPL_HPP
2 #define VERTEXCFD_BOUNDARYSTATE_ELECTRICPOTENTIALINSULATINGWALL_IMPL_HPP
3 
4 #include <Panzer_HierarchicParallelism.hpp>
5 
6 namespace VertexCFD
7 {
8 namespace BoundaryCondition
9 {
10 //---------------------------------------------------------------------------//
11 template<class EvalType, class Traits>
12 ElectricPotentialInsulatingWall<EvalType, Traits>::ElectricPotentialInsulatingWall(
13  const panzer::IntegrationRule& ir)
14  : _boundary_electric_potential("BOUNDARY_electric_potential", ir.dl_scalar)
15  , _boundary_grad_electric_potential("BOUNDARY_GRAD_electric_potential",
16  ir.dl_vector)
17  , _electric_potential("electric_potential", ir.dl_scalar)
18  , _grad_electric_potential("GRAD_electric_potential", ir.dl_vector)
19  , _normals("Side Normal", ir.dl_vector)
20 {
21  this->addDependentField(_electric_potential);
22  this->addEvaluatedField(_boundary_electric_potential);
23  this->addDependentField(_grad_electric_potential);
24  this->addEvaluatedField(_boundary_grad_electric_potential);
25 
26  this->addDependentField(_normals);
27 
28  this->setName("Boundary State Electric Potential Insulating Wall");
29 }
30 
31 //---------------------------------------------------------------------------//
32 template<class EvalType, class Traits>
33 void ElectricPotentialInsulatingWall<EvalType, Traits>::evaluateFields(
34  typename Traits::EvalData workset)
35 {
36  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
37  workset.num_cells);
38  Kokkos::parallel_for(this->getName(), policy, *this);
39 }
40 
41 //---------------------------------------------------------------------------//
42 template<class EvalType, class Traits>
43 KOKKOS_INLINE_FUNCTION void
44 ElectricPotentialInsulatingWall<EvalType, Traits>::operator()(
45  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
46 {
47  const int cell = team.league_rank();
48  const int num_point = _electric_potential.extent(1);
49  const int num_grad_dim = _normals.extent(2);
50 
51  Kokkos::parallel_for(
52  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
53  // Compute \grad(\phi) \cdot \vec{n}
54  scalar_type grad_phi_dot_n = 0.0;
55  for (int dim = 0; dim < num_grad_dim; ++dim)
56  {
57  grad_phi_dot_n += _grad_electric_potential(cell, point, dim)
58  * _normals(cell, point, dim);
59  }
60 
61  _boundary_electric_potential(cell, point)
62  = _electric_potential(cell, point);
63 
64  // Set and boundary gradients
65  for (int dim = 0; dim < num_grad_dim; ++dim)
66  {
67  _boundary_grad_electric_potential(cell, point, dim)
68  = _grad_electric_potential(cell, point, dim)
69  - grad_phi_dot_n * _normals(cell, point, dim);
70  }
71  });
72 }
73 
74 //---------------------------------------------------------------------------//
75 
76 } // end namespace BoundaryCondition
77 } // end namespace VertexCFD
78 
79 #endif // VERTEXCFD_BOUNDARYSTATE_ELECTRICPOTENTIALINSULATINGWALL_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23