1 #ifndef VERTEXCFD_CLOSURE_INCOMPRESSIBLELOCALTIMESTEPSIZE_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_INCOMPRESSIBLELOCALTIMESTEPSIZE_IMPL_HPP
4 #include <utils/VertexCFD_Utils_VectorField.hpp>
6 #include <utils/VertexCFD_Utils_SmoothMath.hpp>
8 #include <Panzer_HierarchicParallelism.hpp>
14 namespace ClosureModel
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)
24 this->addEvaluatedField(_local_dt);
27 this->addDependentField(_element_length);
28 Utils::addDependentVectorField(*
this, ir.dl_scalar, _velocity,
"velocity_");
30 this->setName(
"Incompressible Local Time Step Size");
34 template<
class EvalType,
class Traits,
int NumSpaceDim>
35 void IncompressibleLocalTimeStepSize<EvalType, Traits, NumSpaceDim>::evaluateFields(
36 typename Traits::EvalData workset)
38 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
40 Kokkos::parallel_for(this->getName(), policy, *
this);
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
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);
53 const double tol = 1.0e-8;
55 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
57 auto&& one_over_dt = _local_dt(cell, point);
59 for (
int dim = 0; dim < num_grad_dim; ++dim)
61 one_over_dt += SmoothMath::abs(_velocity[dim](cell, point), tol)
62 / _element_length(cell, point, dim);
66 _local_dt(cell, point) = 1.0 / one_over_dt;
75 #endif // end VERTEXCFD_CLOSURE_INCOMPRESSIBLELOCALTIMESTEPSIZE_IMPL_HPP