1 #ifndef VERTEXCFD_BOUNDARYSTATE_FULLINDUCTIONCONDUCTING_IMPL_HPP
2 #define VERTEXCFD_BOUNDARYSTATE_FULLINDUCTIONCONDUCTING_IMPL_HPP
4 #include "utils/VertexCFD_Utils_VectorField.hpp"
6 #include <Panzer_HierarchicParallelism.hpp>
10 namespace BoundaryCondition
19 template<
class EvalType,
class Traits,
int NumSpaceDim>
20 FullInductionConducting<EvalType, Traits, NumSpaceDim>::FullInductionConducting(
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",
26 , _build_magn_corr(mhd_props.buildMagnCorr())
27 , _build_resistive_flux(mhd_props.buildResistiveFlux())
28 , _dirichlet_scalar_magn_pot(bc_params.isType<double>(
"scalar_magnetic_"
30 , _bnd_scalar_magn_pot(std::numeric_limits<double>::signaling_NaN())
31 , _magnetic_permeability(mhd_props.vacuumMagneticPermeability())
32 , _normals(
"Side Normal", ir.dl_vector)
33 , _scalar_magnetic_potential(
"scalar_magnetic_potential", ir.dl_scalar)
34 , _resistivity(
"resistivity", ir.dl_scalar)
36 for (
int dim = 0; dim < num_space_dim; ++dim)
38 const std::string magn_string =
"induced_magnetic_field_"
39 + std::to_string(dim);
40 _bnd_magn_field[dim] = bc_params.get<
double>(magn_string);
44 Utils::addEvaluatedVectorField(*
this,
46 _boundary_induced_magnetic_field,
47 "BOUNDARY_induced_magnetic_field_");
49 Utils::addEvaluatedVectorField(*
this,
51 _boundary_grad_induced_magnetic_field,
52 "BOUNDARY_GRAD_induced_magnetic_field_");
54 Utils::addDependentVectorField(*
this,
56 _induced_magnetic_field,
57 "induced_magnetic_field_");
59 Utils::addDependentVectorField(*
this,
61 _grad_induced_magnetic_field,
62 "GRAD_induced_magnetic_field_");
66 this->addEvaluatedField(_boundary_scalar_magnetic_potential);
67 if (_dirichlet_scalar_magn_pot)
70 = bc_params.get<
double>(
"scalar_magnetic_potential");
74 this->addDependentField(_scalar_magnetic_potential);
77 this->addDependentField(_normals);
79 if (_build_resistive_flux)
81 this->addDependentField(_resistivity);
82 Utils::addDependentVectorField(
83 *
this, ir.dl_scalar, _boundary_velocity,
"BOUNDARY_velocity_");
84 Utils::addDependentVectorField(*
this,
86 _external_magnetic_field,
87 "external_magnetic_field_");
90 this->setName(
"Boundary State Full Induction Conducting "
91 + std::to_string(num_space_dim) +
"D");
95 template<
class EvalType,
class Traits,
int NumSpaceDim>
96 void FullInductionConducting<EvalType, Traits, NumSpaceDim>::evaluateFields(
97 typename Traits::EvalData workset)
101 if constexpr (Sacado::IsADType<scalar_type>::value)
103 const int fad_size = Kokkos::dimension_scalar(
104 _boundary_induced_magnetic_field[0].get_view());
105 bytes = scratch_view_B::shmem_size(
106 _boundary_induced_magnetic_field[0].extent(1), fad_size);
107 bytes += scratch_view_gradB::shmem_size(
108 _boundary_induced_magnetic_field[0].extent(1),
114 bytes = scratch_view_B::shmem_size(
115 _boundary_induced_magnetic_field[0].extent(1));
116 bytes += scratch_view_gradB::shmem_size(
117 _boundary_induced_magnetic_field[0].extent(1), num_space_dim);
120 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
122 policy.set_scratch_size(0, Kokkos::PerTeam(bytes));
123 Kokkos::parallel_for(this->getName(), policy, *
this);
127 template<
class EvalType,
class Traits,
int NumSpaceDim>
128 KOKKOS_INLINE_FUNCTION
void
129 FullInductionConducting<EvalType, Traits, NumSpaceDim>::operator()(
130 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
132 const int cell = team.league_rank();
133 const int num_point = _boundary_induced_magnetic_field[0].extent(1);
134 const int num_grad_dim = _boundary_grad_induced_magnetic_field[0].extent(2);
135 scratch_view_B scratch_data_B;
136 scratch_view_gradB scratch_data_gradB;
137 if constexpr (Sacado::IsADType<scalar_type>::value)
139 const int fad_size = Kokkos::dimension_scalar(
140 _boundary_induced_magnetic_field[0].get_view());
141 scratch_data_B = scratch_view_B(team.team_shmem(), num_point, fad_size);
142 scratch_data_gradB = scratch_view_gradB(
143 team.team_shmem(), num_point, num_space_dim, fad_size);
147 scratch_data_B = scratch_view_B(team.team_shmem(), num_point);
149 = scratch_view_gradB(team.team_shmem(), num_point, num_space_dim);
152 Kokkos::parallel_for(
153 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
155 auto&& B_dot_n = scratch_data_B(point);
158 = Kokkos::subview(scratch_data_gradB, point, Kokkos::ALL);
159 for (
int d = 0; d < num_space_dim; ++d)
161 for (
int grad_dim = 0; grad_dim < num_grad_dim; ++grad_dim)
163 B_dot_n += (_induced_magnetic_field[grad_dim](cell, point)
164 - _bnd_magn_field[grad_dim])
165 * _normals(cell, point, grad_dim);
166 for (
int field_dim = 0; field_dim < num_grad_dim; ++field_dim)
168 gradB_dot_n[field_dim]
169 += _grad_induced_magnetic_field[field_dim](
170 cell, point, grad_dim)
171 * _normals(cell, point, grad_dim);
176 for (
int dim = 0; dim < num_space_dim; ++dim)
178 _boundary_induced_magnetic_field[dim](cell, point)
179 = _induced_magnetic_field[dim](cell, point);
180 if (dim < num_grad_dim)
182 _boundary_induced_magnetic_field[dim](cell, point)
183 -= B_dot_n * _normals(cell, point, dim);
187 if (_build_magn_corr)
189 if (_dirichlet_scalar_magn_pot)
191 _boundary_scalar_magnetic_potential(cell, point)
192 = _bnd_scalar_magn_pot;
196 _boundary_scalar_magnetic_potential(cell, point)
197 = _scalar_magnetic_potential(cell, point);
203 if (_build_resistive_flux)
205 const scalar_type inv_eta = _magnetic_permeability
206 / _resistivity(cell, point);
207 for (
int d = 0; d < num_grad_dim; ++d)
209 for (
int fdim = 0; fdim < num_grad_dim; ++fdim)
214 += _normals(cell, point, d) * inv_eta
215 * (_boundary_velocity[fdim](cell, point)
216 * (_external_magnetic_field[d](cell, point)
217 + _bnd_magn_field[d])
218 - _boundary_velocity[d](cell, point)
219 * (_external_magnetic_field[fdim](
221 + _bnd_magn_field[fdim]));
227 for (
int d = 0; d < num_grad_dim; ++d)
229 for (
int fdim = 0; fdim < num_space_dim; ++fdim)
231 _boundary_grad_induced_magnetic_field[fdim](cell, point, d)
232 = _grad_induced_magnetic_field[fdim](cell, point, d);
233 if (fdim < num_grad_dim)
235 _boundary_grad_induced_magnetic_field[fdim](
237 -= gradB_dot_n[fdim] * _normals(cell, point, d);
249 #endif // end VERTEXCFD_BOUNDARYSTATE_FULLINDUCTIONCONDUCTING_IMPL_HPP