VertexCFD  0.0-dev
VertexCFD_Closure_IncompressibleLocalTimeStepSize_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_INCOMPRESSIBLELOCALTIMESTEPSIZE_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_INCOMPRESSIBLELOCALTIMESTEPSIZE_IMPL_HPP
3 
4 #include <utils/VertexCFD_Utils_VectorField.hpp>
5 
6 #include <utils/VertexCFD_Utils_SmoothMath.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 IncompressibleLocalTimeStepSize<EvalType, Traits, NumSpaceDim>::
19  IncompressibleLocalTimeStepSize(const panzer::IntegrationRule& ir)
20  : _local_dt("local_dt", ir.dl_scalar)
21  , _element_length("element_length", ir.dl_vector)
22 {
23  // Add evaluated field
24  this->addEvaluatedField(_local_dt);
25 
26  // Add dependent fields
27  this->addDependentField(_element_length);
28  Utils::addDependentVectorField(*this, ir.dl_scalar, _velocity, "velocity_");
29 
30  this->setName("Incompressible Local Time Step Size");
31 }
32 
33 //---------------------------------------------------------------------------//
34 template<class EvalType, class Traits, int NumSpaceDim>
35 void IncompressibleLocalTimeStepSize<EvalType, Traits, NumSpaceDim>::evaluateFields(
36  typename Traits::EvalData workset)
37 {
38  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
39  workset.num_cells);
40  Kokkos::parallel_for(this->getName(), policy, *this);
41 }
42 
43 //---------------------------------------------------------------------------//
44 template<class EvalType, class Traits, int NumSpaceDim>
45 KOKKOS_INLINE_FUNCTION void
46 IncompressibleLocalTimeStepSize<EvalType, Traits, NumSpaceDim>::operator()(
47  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
48 {
49  const int cell = team.league_rank();
50  const int num_point = _local_dt.extent(1);
51  const int num_grad_dim = _element_length.extent(2);
52 
53  const double tol = 1.0e-8;
54  Kokkos::parallel_for(
55  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
56  // Accumulate 1/dt
57  auto&& one_over_dt = _local_dt(cell, point);
58  one_over_dt = 0.0;
59  for (int dim = 0; dim < num_grad_dim; ++dim)
60  {
61  one_over_dt += SmoothMath::abs(_velocity[dim](cell, point), tol)
62  / _element_length(cell, point, dim);
63  }
64 
65  // Invert to compute dt
66  _local_dt(cell, point) = 1.0 / one_over_dt;
67  });
68 }
69 
70 //---------------------------------------------------------------------------//
71 
72 } // end namespace ClosureModel
73 } // end namespace VertexCFD
74 
75 #endif // end VERTEXCFD_CLOSURE_INCOMPRESSIBLELOCALTIMESTEPSIZE_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23