VertexCFD  0.0-dev
VertexCFD_BoundaryState_FullInductionFixed_impl.hpp
1 #ifndef VERTEXCFD_BOUNDARYSTATE_FULLINDUCTIONFIXED_IMPL_HPP
2 #define VERTEXCFD_BOUNDARYSTATE_FULLINDUCTIONFIXED_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 // This function sets fixed boundary values for the induced magnetic field
14 // components from input. When solving the magnetic correction potential
15 // equation for divergence cleaning, the boundary value for the scalar
16 // magnetic potential is extrapolated (free) unless a value is set in the
17 // input BC parameters.
18 //---------------------------------------------------------------------------//
19 template<class EvalType, class Traits, int NumSpaceDim>
20 FullInductionFixed<EvalType, Traits, NumSpaceDim>::FullInductionFixed(
21  const panzer::IntegrationRule& ir,
22  const Teuchos::ParameterList& bc_params,
23  const MHDProperties::FullInductionMHDProperties& mhd_props)
24  : _boundary_scalar_magnetic_potential("BOUNDARY_scalar_magnetic_potential",
25  ir.dl_scalar)
26  , _build_magn_corr(mhd_props.buildMagnCorr())
27  , _dirichlet_scalar_magn_pot(bc_params.isType<double>("scalar_magnetic_"
28  "potential"))
29  , _bnd_scalar_magn_pot(std::numeric_limits<double>::signaling_NaN())
30  , _scalar_magnetic_potential("scalar_magnetic_potential", ir.dl_scalar)
31 {
32  for (int dim = 0; dim < num_space_dim; ++dim)
33  {
34  const std::string magn_string = "induced_magnetic_field_"
35  + std::to_string(dim);
36  _bnd_magn_field[dim] = bc_params.get<double>(magn_string);
37  }
38 
39  // Add evaluated/dependent fields
40  Utils::addEvaluatedVectorField(*this,
41  ir.dl_scalar,
42  _boundary_induced_magnetic_field,
43  "BOUNDARY_induced_magnetic_field_");
44 
45  Utils::addEvaluatedVectorField(*this,
46  ir.dl_vector,
47  _boundary_grad_induced_magnetic_field,
48  "BOUNDARY_GRAD_induced_magnetic_field_");
49 
50  Utils::addDependentVectorField(*this,
51  ir.dl_vector,
52  _grad_induced_magnetic_field,
53  "GRAD_induced_magnetic_field_");
54 
55  if (_build_magn_corr)
56  {
57  this->addEvaluatedField(_boundary_scalar_magnetic_potential);
58  if (_dirichlet_scalar_magn_pot)
59  {
60  _bnd_scalar_magn_pot
61  = bc_params.get<double>("scalar_magnetic_potential");
62  }
63  else
64  {
65  this->addDependentField(_scalar_magnetic_potential);
66  }
67  }
68 
69  this->setName("Boundary State Full Induction Fixed "
70  + std::to_string(num_space_dim) + "D");
71 }
72 
73 //---------------------------------------------------------------------------//
74 template<class EvalType, class Traits, int NumSpaceDim>
75 void FullInductionFixed<EvalType, Traits, NumSpaceDim>::evaluateFields(
76  typename Traits::EvalData workset)
77 {
78  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
79  workset.num_cells);
80  Kokkos::parallel_for(this->getName(), policy, *this);
81 }
82 
83 //---------------------------------------------------------------------------//
84 template<class EvalType, class Traits, int NumSpaceDim>
85 KOKKOS_INLINE_FUNCTION void
86 FullInductionFixed<EvalType, Traits, NumSpaceDim>::operator()(
87  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
88 {
89  const int cell = team.league_rank();
90  const int num_point = _boundary_induced_magnetic_field[0].extent(1);
91  const int num_grad_dim = _boundary_grad_induced_magnetic_field[0].extent(2);
92 
93  Kokkos::parallel_for(
94  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
95  // Set the fixed boundary values
96  for (int dim = 0; dim < num_space_dim; ++dim)
97  {
98  _boundary_induced_magnetic_field[dim](cell, point)
99  = _bnd_magn_field[dim];
100  }
101 
102  if (_build_magn_corr)
103  {
104  if (_dirichlet_scalar_magn_pot)
105  {
106  _boundary_scalar_magnetic_potential(cell, point)
107  = _bnd_scalar_magn_pot;
108  }
109  else
110  {
111  _boundary_scalar_magnetic_potential(cell, point)
112  = _scalar_magnetic_potential(cell, point);
113  }
114  }
115 
116  // Set gradients
117  for (int d = 0; d < num_grad_dim; ++d)
118  {
119  for (int field_dim = 0; field_dim < num_space_dim; ++field_dim)
120  {
121  _boundary_grad_induced_magnetic_field[field_dim](
122  cell, point, d)
123  = _grad_induced_magnetic_field[field_dim](
124  cell, point, d);
125  }
126  }
127  });
128 }
129 
130 //---------------------------------------------------------------------------//
131 
132 } // end namespace BoundaryCondition
133 } // end namespace VertexCFD
134 
135 #endif // end VERTEXCFD_BOUNDARYSTATE_FULLINDUCTIONFIXED_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23