VertexCFD  0.0-dev
VertexCFD_Closure_IncompressibleSpalartAllmarasEddyViscosity_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_INCOMPRESSIBLESPALARTALLMARASEDDYVISCOSITY_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_INCOMPRESSIBLESPALARTALLMARASEDDYVISCOSITY_IMPL_HPP
3 
4 #include "utils/VertexCFD_Utils_SmoothMath.hpp"
5 
6 #include <Panzer_HierarchicParallelism.hpp>
7 
8 namespace VertexCFD
9 {
10 namespace ClosureModel
11 {
12 //---------------------------------------------------------------------------//
13 template<class EvalType, class Traits>
14 IncompressibleSpalartAllmarasEddyViscosity<EvalType, Traits>::
15  IncompressibleSpalartAllmarasEddyViscosity(const panzer::IntegrationRule& ir)
16  : _sa_var("spalart_allmaras_variable", ir.dl_scalar)
17  , _nu("kinematic_viscosity", ir.dl_scalar)
18  , _cv1(7.1)
19  , _num_grad_dim(ir.spatial_dimension)
20  , _nu_t("turbulent_eddy_viscosity", ir.dl_scalar)
21 {
22  // Add dependent fields
23  this->addDependentField(_sa_var);
24  this->addDependentField(_nu);
25 
26  // Add evaluated fields
27  this->addEvaluatedField(_nu_t);
28 
29  // Closure model name
30  this->setName("Spalart-Allmaras Incompressible Turbulent Eddy Viscosity "
31  + std::to_string(_num_grad_dim) + "D");
32 }
33 
34 //---------------------------------------------------------------------------//
35 template<class EvalType, class Traits>
36 void IncompressibleSpalartAllmarasEddyViscosity<EvalType, Traits>::evaluateFields(
37  typename Traits::EvalData workset)
38 {
39  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
40  workset.num_cells);
41  Kokkos::parallel_for(this->getName(), policy, *this);
42 }
43 
44 //---------------------------------------------------------------------------//
45 template<class EvalType, class Traits>
46 KOKKOS_INLINE_FUNCTION void
47 IncompressibleSpalartAllmarasEddyViscosity<EvalType, Traits>::operator()(
48  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
49 {
50  const int cell = team.league_rank();
51  const int num_point = _nu_t.extent(1);
52  using std::pow;
53 
54  Kokkos::parallel_for(
55  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
56  const auto sa_var = _sa_var(cell, point);
57  if (sa_var >= 0.0)
58  {
59  const scalar_type xi3 = pow(sa_var / _nu(cell, point), 3.0);
60  const scalar_type f_v1 = xi3 / (xi3 + _cv1 * _cv1 * _cv1);
61  _nu_t(cell, point) = sa_var * f_v1;
62  }
63  else
64  {
65  _nu_t(cell, point) = 0.0;
66  }
67  });
68 }
69 
70 //---------------------------------------------------------------------------//
71 
72 } // end namespace ClosureModel
73 } // end namespace VertexCFD
74 
75 #endif // end
76  // VERTEXCFD_CLOSURE_INCOMPRESSIBLESPALARTALLMARASEDDYVISCOSITY_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23