VertexCFD  0.0-dev
VertexCFD_Closure_SolidFullInductionLocalTimeStepSize_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_SOLIDFULLINDUCTIONLOCALTIMESTEPSIZE_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_SOLIDFULLINDUCTIONLOCALTIMESTEPSIZE_IMPL_HPP
3 
4 #include "utils/VertexCFD_Utils_MagneticLayout.hpp"
5 #include "utils/VertexCFD_Utils_SmoothMath.hpp"
6 
7 #include <Panzer_HierarchicParallelism.hpp>
8 
9 #include <cmath>
10 
11 namespace VertexCFD
12 {
13 namespace ClosureModel
14 {
15 //---------------------------------------------------------------------------//
16 template<class EvalType, class Traits>
17 SolidFullInductionLocalTimeStepSize<EvalType, Traits>::
18  SolidFullInductionLocalTimeStepSize(
19  const panzer::IntegrationRule& ir,
20  const MHDProperties::FullInductionMHDProperties& mhd_props)
21  : _local_dt("local_dt", ir.dl_scalar)
22  , _magnetic_permeability(mhd_props.vacuumMagneticPermeability())
23  , _c_h(mhd_props.hyperbolicDivergenceCleaningSpeed())
24  , _element_length("element_length", ir.dl_vector)
25  , _solid_density("solid_density", ir.dl_scalar)
26  , _total_magnetic_field(
27  "total_magnetic_field",
28  Utils::buildMagneticLayout(ir.dl_scalar, num_magnetic_field_dim))
29 {
30  // Add evaluated field
31  this->addEvaluatedField(_local_dt);
32 
33  // Add dependent fields
34  this->addDependentField(_element_length);
35  this->addDependentField(_solid_density);
36  this->addDependentField(_total_magnetic_field);
37 
38  this->setName("Solid Full Induction Local Time Step Size");
39 }
40 
41 //---------------------------------------------------------------------------//
42 template<class EvalType, class Traits>
43 void SolidFullInductionLocalTimeStepSize<EvalType, Traits>::evaluateFields(
44  typename Traits::EvalData workset)
45 {
46  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
47  workset.num_cells);
48  Kokkos::parallel_for(this->getName(), policy, *this);
49 }
50 
51 //---------------------------------------------------------------------------//
52 template<class EvalType, class Traits>
53 KOKKOS_INLINE_FUNCTION void
54 SolidFullInductionLocalTimeStepSize<EvalType, Traits>::operator()(
55  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
56 {
57  const int cell = team.league_rank();
58  const int num_point = _local_dt.extent(1);
59  const int num_grad_dim = _element_length.extent(2);
60 
61  const double tol = 1.0e-8;
62  Kokkos::parallel_for(
63  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
64  _local_dt(cell, point) = 0.0;
65  const auto& rho = _solid_density(cell, point);
66  for (int dim = 0; dim < num_grad_dim; ++dim)
67  {
68  using Kokkos::sqrt;
69  using SmoothMath::abs;
70  using SmoothMath::max;
71  _local_dt(cell, point)
72  += max(abs(_total_magnetic_field(cell, point, dim), tol)
73  / sqrt(_magnetic_permeability * rho),
74  _c_h,
75  tol)
76  / _element_length(cell, point, dim);
77  }
78  _local_dt(cell, point) = 1.0 / _local_dt(cell, point);
79  });
80 }
81 
82 //---------------------------------------------------------------------------//
83 
84 } // end namespace ClosureModel
85 } // end namespace VertexCFD
86 
87 #endif // end VERTEXCFD_CLOSURE_SOLIDFULLINDUCTIONLOCALTIMESTEPSIZE_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23