VertexCFD  0.0-dev
VertexCFD_BoundaryState_TurbulenceKEpsilonWallFunction_impl.hpp
1 #ifndef VERTEXCFD_BOUNDARYSTATE_TURBULENCEKEPSILONWALLFUNCTION_IMPL_HPP
2 #define VERTEXCFD_BOUNDARYSTATE_TURBULENCEKEPSILONWALLFUNCTION_IMPL_HPP
3 
4 #include "utils/VertexCFD_Utils_SmoothMath.hpp"
5 #include "utils/VertexCFD_Utils_VectorField.hpp"
6 
7 #include <Panzer_HierarchicParallelism.hpp>
8 
9 namespace VertexCFD
10 {
11 namespace BoundaryCondition
12 {
13 //---------------------------------------------------------------------------//
14 template<class EvalType, class Traits, int NumSpaceDim>
15 TurbulenceKEpsilonWallFunction<EvalType, Traits, NumSpaceDim>::
16  TurbulenceKEpsilonWallFunction(const panzer::IntegrationRule& ir,
17  const Teuchos::ParameterList& bc_params)
18  : _boundary_k("BOUNDARY_turb_kinetic_energy", ir.dl_scalar)
19  , _boundary_e("BOUNDARY_turb_dissipation_rate", ir.dl_scalar)
20  , _boundary_grad_k("BOUNDARY_GRAD_turb_kinetic_energy", ir.dl_vector)
21  , _boundary_grad_e("BOUNDARY_GRAD_turb_dissipation_rate", ir.dl_vector)
22  , _boundary_u_tau("BOUNDARY_friction_velocity", ir.dl_scalar)
23  , _boundary_y_plus("BOUNDARY_y_plus", ir.dl_scalar)
24  , _wall_func_nu_t("wall_func_turbulent_eddy_viscosity", ir.dl_scalar)
25  , _k("turb_kinetic_energy", ir.dl_scalar)
26  , _e("turb_dissipation_rate", ir.dl_scalar)
27  , _grad_k("GRAD_turb_kinetic_energy", ir.dl_vector)
28  , _grad_e("GRAD_turb_dissipation_rate", ir.dl_vector)
29  , _normals("Side Normal", ir.dl_vector)
30  , _nu("kinematic_viscosity", ir.dl_scalar)
31  , _C_mu(0.09)
32  , _kappa(0.41)
33  , _yp_tr(11.06)
34  , _neumann(false)
35 {
36  // Check for epsilon boundary specification
37  if (bc_params.isType<std::string>("Epsilon Condition Type"))
38  {
39  const std::string bc_type
40  = bc_params.get<std::string>("Epsilon Condition Type");
41 
42  if (bc_type == "Neumann")
43  {
44  _neumann = true;
45  }
46  else if (bc_type == "Dirichlet")
47  {
48  _neumann = false;
49  }
50  else
51  {
52  std::string msg = "Unknown Epsilon Condition Type " + bc_type;
53  msg += "\nPlease choose from Dirichlet (default) or Neumann.\n";
54 
55  throw std::runtime_error(msg);
56  }
57  }
58 
59  // Add evaluated fields
60  this->addEvaluatedField(_boundary_k);
61  this->addEvaluatedField(_boundary_e);
62  this->addEvaluatedField(_boundary_grad_k);
63  this->addEvaluatedField(_boundary_grad_e);
64  this->addEvaluatedField(_boundary_u_tau);
65  this->addEvaluatedField(_boundary_y_plus);
66  this->addEvaluatedField(_wall_func_nu_t);
67 
68  // Add dependent fields
69  this->addDependentField(_k);
70  this->addDependentField(_e);
71  this->addDependentField(_grad_k);
72  this->addDependentField(_grad_e);
73  this->addDependentField(_normals);
74  this->addDependentField(_nu);
75  Utils::addDependentVectorField(*this, ir.dl_scalar, _velocity, "velocity_");
76 
77  this->setName("Boundary State Turbulence Model K-Epsilon Wall Function "
78  + std::to_string(_num_grad_dim) + "D");
79 }
80 
81 //---------------------------------------------------------------------------//
82 template<class EvalType, class Traits, int NumSpaceDim>
83 void TurbulenceKEpsilonWallFunction<EvalType, Traits, NumSpaceDim>::evaluateFields(
84  typename Traits::EvalData workset)
85 {
86  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
87  workset.num_cells);
88  Kokkos::parallel_for(this->getName(), policy, *this);
89 }
90 
91 //---------------------------------------------------------------------------//
92 template<class EvalType, class Traits, int NumSpaceDim>
93 KOKKOS_INLINE_FUNCTION void
94 TurbulenceKEpsilonWallFunction<EvalType, Traits, NumSpaceDim>::operator()(
95  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
96 {
97  const int cell = team.league_rank();
98  const int num_point = _grad_k.extent(1);
99  const double max_tol = 1e-10;
100 
101  using std::pow;
102 
103  Kokkos::parallel_for(
104  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
105  // Calculate velocity magnitude
106  scalar_type mag_u = 0.0;
107 
108  for (int d = 0; d < num_space_dim; d++)
109  {
110  mag_u += pow(_velocity[d](cell, point), 2.0);
111  }
112 
113  mag_u = pow(SmoothMath::max(mag_u, max_tol, 0.0), 0.5);
114 
115  // Set turbulent kinetic energy boundary using Neumann condition
116  _boundary_k(cell, point) = _k(cell, point);
117 
118  for (int d = 0; d < num_space_dim; ++d)
119  {
120  _boundary_grad_k(cell, point, d) = _grad_k(cell, point, d);
121 
122  // Subtract boundary normal component
123  for (int grad_dim = 0; grad_dim < num_space_dim; ++grad_dim)
124  {
125  _boundary_grad_k(cell, point, d)
126  -= _grad_k(cell, point, grad_dim)
127  * _normals(cell, point, grad_dim)
128  * _normals(cell, point, d);
129  }
130  }
131 
132  // Calculate friction velocity and nu_t at the wall as suggested
133  // by Kuzmin et al. (2007)
134  const scalar_type u_tau = SmoothMath::max(
135  pow(_C_mu, 0.25)
136  * pow(SmoothMath::max(_k(cell, point), max_tol, 0.0), 0.5),
137  mag_u / _yp_tr,
138  0.0);
139 
140  const scalar_type nu_t_w = _kappa * _yp_tr * _nu(cell, point);
141 
142  // Set epsilon and nu_t at boundary according to user setting
143  if (_neumann)
144  {
145  // Set epsilon boundary equal to interior value
146  _boundary_e(cell, point) = _e(cell, point);
147 
148  // Calculate boundary gradient
149  for (int d = 0; d < num_space_dim; ++d)
150  {
151  _boundary_grad_e(cell, point, d)
152  = _grad_e(cell, point, d)
153  + (_kappa * pow(u_tau, 5.0) / pow(nu_t_w, 2.0))
154  * _normals(cell, point, d);
155 
156  // Subtract interior boundary normal component
157  for (int grad_dim = 0; grad_dim < num_space_dim; ++grad_dim)
158  {
159  _boundary_grad_e(cell, point, d)
160  -= _grad_e(cell, point, grad_dim)
161  * _normals(cell, point, grad_dim)
162  * _normals(cell, point, d);
163  }
164  }
165 
166  // Set boundary nu_t to analytical solution
167  _wall_func_nu_t(cell, point) = nu_t_w;
168  }
169  else
170  {
171  // Calculate epsilon boundary value
172  _boundary_e(cell, point) = pow(u_tau, 4.0) / nu_t_w;
173 
174  // Set boundary gradient equal to interior gradient
175  for (int d = 0; d < num_space_dim; ++d)
176  {
177  _boundary_grad_e(cell, point, d) = _grad_e(cell, point, d);
178  }
179 
180  // Calculate boundary nu_t with standard k-epsilon form
181  _wall_func_nu_t(cell, point)
182  = _C_mu * pow(_k(cell, point), 2.0)
183  / SmoothMath::max(_e(cell, point), max_tol, 0.0);
184  }
185 
186  // Fill in boundary friction velocity and y+ fields for
187  // post-processing
188  _boundary_u_tau(cell, point) = u_tau;
189  _boundary_y_plus(cell, point) = _yp_tr;
190 
191  // TODO: define TKE production term at boundary to enforce P = D?
192  });
193 }
194 
195 //---------------------------------------------------------------------------//
196 
197 } // end namespace BoundaryCondition
198 } // end namespace VertexCFD
199 
200 #endif // end VERTEXCFD_BOUNDARYSTATE_TURBULENCEKEPSILONWALLFUNCTION_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23