VertexCFD  0.0-dev
VertexCFD_Closure_ExternalMagneticField_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_EXTERNALMAGNETICFIELD_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_EXTERNALMAGNETICFIELD_IMPL_HPP
3 
4 #include "utils/VertexCFD_Utils_VectorField.hpp"
5 
6 #include <Panzer_HierarchicParallelism.hpp>
7 #include <Panzer_Workset_Utilities.hpp>
8 #include <Teuchos_StandardParameterEntryValidators.hpp>
9 
10 #include <string>
11 
12 namespace VertexCFD
13 {
14 namespace ClosureModel
15 {
16 //---------------------------------------------------------------------------//
17 template<class EvalType, class Traits>
18 ExternalMagneticField<EvalType, Traits>::ExternalMagneticField(
19  const panzer::IntegrationRule& ir,
20  const Teuchos::ParameterList& user_params)
21  : _ir_degree(ir.cubature_degree)
22  , _ext_magn_type(ExtMagnType::constant)
23 {
24  // Get external magnetic vector type
25  if (user_params.isType<std::string>("External Magnetic Field Type"))
26  {
27  const auto type_validator = Teuchos::rcp(
28  new Teuchos::StringToIntegralParameterEntryValidator<ExtMagnType>(
29  Teuchos::tuple<std::string>("constant", "toroidal"),
30  "constant"));
31  _ext_magn_type = type_validator->getIntegralValue(
32  user_params.get<std::string>("External Magnetic Field Type"));
33  }
34 
35  // Get external magnetic vector value
36  if (_ext_magn_type == ExtMagnType::constant)
37  {
38  const auto ext_magn_vct = user_params.get<Teuchos::Array<double>>(
39  "External Magnetic Field Value");
40 
41  const auto zero_array = Teuchos::Array<double>(field_size, 0.0);
42  const auto d_ext_magn_vct_dt
43  = user_params.isType<Teuchos::Array<double>>(
44  "External Magnetic Field Time Rate of "
45  "Change")
46  ? user_params.get<Teuchos::Array<double>>(
47  "External Magnetic Field Time Rate "
48  "of Change")
49  : zero_array;
50 
51  for (int dim = 0; dim < field_size; ++dim)
52  {
53  _ext_magn_vct[dim] = ext_magn_vct[dim];
54  _d_ext_magn_vct_dt[dim] = d_ext_magn_vct_dt[dim];
55  }
56  }
57  else if (_ext_magn_type == ExtMagnType::toroidal)
58  {
59  _toroidal_field_magn
60  = user_params.get<double>("Toroidal Field Magnitude");
61  }
62 
63  // Evaluated fields
64  Utils::addEvaluatedVectorField(
65  *this, ir.dl_scalar, _ext_magn_field, "external_magnetic_field_");
66 
67  this->setName("External Magnetic Field");
68 }
69 
70 //---------------------------------------------------------------------------//
71 template<class EvalType, class Traits>
72 void ExternalMagneticField<EvalType, Traits>::postRegistrationSetup(
73  typename Traits::SetupData sd, PHX::FieldManager<Traits>&)
74 {
75  _ir_index = panzer::getIntegrationRuleIndex(_ir_degree, (*sd.worksets_)[0]);
76 }
77 
78 //---------------------------------------------------------------------------//
79 template<class EvalType, class Traits>
80 void ExternalMagneticField<EvalType, Traits>::evaluateFields(
81  typename Traits::EvalData workset)
82 {
83  _time = workset.time;
84 
85  _ip_coords = workset.int_rules[_ir_index]->ip_coordinates;
86  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
87  workset.num_cells);
88  Kokkos::parallel_for(this->getName(), policy, *this);
89 }
90 
91 //---------------------------------------------------------------------------//
92 template<class EvalType, class Traits>
93 KOKKOS_INLINE_FUNCTION void ExternalMagneticField<EvalType, Traits>::operator()(
94  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
95 {
96  const int cell = team.league_rank();
97  const int num_point = _ext_magn_field[0].extent(1);
98  using std::sqrt;
99 
100  Kokkos::parallel_for(
101  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
102  if (_ext_magn_type == ExtMagnType::constant)
103  {
104  for (int dim = 0; dim < field_size; ++dim)
105  {
106  _ext_magn_field[dim](cell, point)
107  = _ext_magn_vct[dim] + _d_ext_magn_vct_dt[dim] * _time;
108  }
109  }
110  else if (_ext_magn_type == ExtMagnType::toroidal)
111  {
112  _ext_magn_field[0](cell, point)
113  = -_ip_coords(cell, point, 1) * _toroidal_field_magn
114  / (_ip_coords(cell, point, 0) * _ip_coords(cell, point, 0)
115  + _ip_coords(cell, point, 1)
116  * _ip_coords(cell, point, 1));
117 
118  _ext_magn_field[1](cell, point)
119  = _ip_coords(cell, point, 0) * _toroidal_field_magn
120  / (_ip_coords(cell, point, 0) * _ip_coords(cell, point, 0)
121  + _ip_coords(cell, point, 1)
122  * _ip_coords(cell, point, 1));
123 
124  _ext_magn_field[2](cell, point) = 0.0;
125  }
126  });
127 }
128 
129 //---------------------------------------------------------------------------//
130 
131 } // end namespace ClosureModel
132 } // end namespace VertexCFD
133 
134 #endif // end VERTEXCFD_CLOSURE_EXTERNALMAGNETICFIELD_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23