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