https://mooseframework.inl.gov
LinearFVTurbulentViscosityWallFunctionBC.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
11 #include "NavierStokesMethods.h"
12 
14 
17 {
19  params.addClassDescription("Adds Dirichlet BC for wall values of the turbulent viscosity.");
20  params.addRequiredParam<MooseFunctorName>("u", "The velocity in the x direction.");
21  params.addParam<MooseFunctorName>("v", "The velocity in the y direction.");
22  params.addParam<MooseFunctorName>("w", "The velocity in the z direction.");
23  params.addRequiredParam<MooseFunctorName>(NS::density, "Density");
24  params.addRequiredParam<MooseFunctorName>(NS::mu, "Dynamic viscosity.");
25  params.addParam<MooseFunctorName>(NS::TKE, "The turbulent kinetic energy.");
26  params.addParam<Real>("C_mu", 0.09, "Coupled turbulent kinetic energy closure.");
27 
28  MooseEnum wall_treatment("eq_newton eq_incremental eq_linearized neq", "neq");
29  params.addParam<MooseEnum>(
30  "wall_treatment", wall_treatment, "The method used for computing the wall functions");
31  return params;
32 }
33 
35  const InputParameters & params)
37  _dim(_subproblem.mesh().dimension()),
38  _u_var(getFunctor<Real>("u")),
39  _v_var(params.isParamValid("v") ? &(getFunctor<Real>("v")) : nullptr),
40  _w_var(params.isParamValid("w") ? &(getFunctor<Real>("w")) : nullptr),
41  _rho(getFunctor<Real>(NS::density)),
42  _mu(getFunctor<Real>(NS::mu)),
43  _k(getFunctor<Real>(NS::TKE)),
44  _C_mu(getParam<Real>("C_mu")),
45  _wall_treatment(getParam<MooseEnum>("wall_treatment").getEnum<NS::WallTreatmentEnum>())
46 {
47 }
48 
49 Real
51 {
52  // Utility functions
53  const auto wall_dist = computeCellToFaceDistance();
57  const auto re = makeElemArg(elem);
59  const auto mu = _mu(re, old_state);
60  const auto rho = _rho(re, old_state);
61 
62  // Get the velocity vector
63  RealVectorValue velocity(_u_var(re, old_state));
64  if (_v_var)
65  velocity(1) = (*_v_var)(re, old_state);
66  if (_w_var)
67  velocity(2) = (*_w_var)(re, old_state);
68 
69  // Compute the velocity and direction of the velocity component that is parallel to the wall
70  const auto parallel_speed = NS::computeSpeed<Real>(
72 
73  // Switch for determining the near wall quantities
74  // wall_treatment can be: "eq_newton eq_incremental eq_linearized neq"
75  Real y_plus = 0.0;
76  Real mu_wall = 0.0; // total wall viscosity to obtain the shear stress at the wall
77 
79  {
80  // Full Newton-Raphson solve to find the wall quantities from the law of the wall
81  const auto u_tau = NS::findUStar<Real>(mu, rho, parallel_speed, wall_dist);
82  y_plus = wall_dist * u_tau * rho / mu;
83  mu_wall = rho * Utility::pow<2>(u_tau) * wall_dist / parallel_speed;
84  }
86  {
87  // Incremental solve on y_plus to get the near-wall quantities
88  y_plus = NS::findyPlus<Real>(mu, rho, std::max(parallel_speed, 1e-10), wall_dist);
89  mu_wall = mu * (NS::von_karman_constant * y_plus /
90  std::log(std::max(NS::E_turb_constant * y_plus, 1 + 1e-4)));
91  }
93  {
94  // Linearized approximation to the wall function to find the near-wall quantities faster
95  const Real a_c = 1 / NS::von_karman_constant;
96  const Real b_c = 1 / NS::von_karman_constant *
97  (std::log(NS::E_turb_constant * std::max(wall_dist, 1.0) / mu) + 1.0);
98  const Real c_c = parallel_speed;
99 
100  const auto u_tau = (-b_c + std::sqrt(std::pow(b_c, 2) + 4.0 * a_c * c_c)) / (2.0 * a_c);
101  y_plus = wall_dist * u_tau * rho / mu;
102  mu_wall = rho * Utility::pow<2>(u_tau) * wall_dist / parallel_speed;
103  }
105  {
106  // Assign non-equilibrium wall function value
107  y_plus = std::pow(_C_mu, 0.25) * wall_dist * std::sqrt(_k(re, old_state)) * rho / mu;
108  mu_wall = mu * (NS::von_karman_constant * y_plus /
109  std::log(std::max(NS::E_turb_constant * y_plus, 1.0 + 1e-4)));
110  }
111  else
112  mooseAssert(false,
113  "For `INSFVTurbulentViscosityWallFunction` , wall treatment should not reach here");
114 
115  const Real mut_log = mu_wall - mu; // turbulent log-layer viscosity
116 
117  Real mu_t = 0;
118 
119  if (y_plus <= 5.0)
120  // sub-laminar layer
121  mu_t += 0.0;
122  else if (y_plus >= 30.0)
123  // log-layer
124  mu_t += std::max(mut_log, NS::mu_t_low_limit);
125  else
126  {
127  // buffer layer
128  const auto blending_function = (y_plus - 5.0) / 25.0;
129  // the blending depends on the mut_log at y+=30
130  const auto mut_log = mu * _mut_30;
131  mu_t += std::max(blending_function * mut_log, NS::mu_t_low_limit);
132  }
133  return mu_t;
134 }
135 
136 Real
138 {
139  return this->computeTurbulentViscosity();
140 }
141 
142 Real
144 {
149  return (this->computeTurbulentViscosity() - raw_value(_var(elem_arg, determineState()))) /
150  distance;
151 }
152 
153 Real
155 {
156  mooseError("We should not solve for the turbulent viscosity directly meaning that this should "
157  "contribute to neither vector nor a right hand side.");
158 }
159 
160 Real
162 {
163  mooseError("We should not solve for the turbulent viscosity directly meaning that this should "
164  "contribute to neither vector nor a right hand side.");
165 }
166 
167 Real
169 {
170  mooseError("We should not solve for the turbulent viscosity directly meaning that this should "
171  "contribute to neither vector nor a right hand side.");
172 }
173 
174 Real
176 {
177  mooseError("We should not solve for the turbulent viscosity directly meaning that this should "
178  "contribute to neither vector nor a right hand side.");
179 }
static constexpr Real von_karman_constant
Definition: NS.h:201
static const std::string mu_t
Definition: NS.h:125
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
Real computeCellToFaceDistance() const
Moose::StateArg determineState() const
MeshBase & mesh
static const std::string density
Definition: NS.h:33
auto raw_value(const Eigen::Map< T > &in)
static const std::string TKE
Definition: NS.h:176
WallTreatmentEnum
Wall treatment options.
Definition: NS.h:182
LinearFVTurbulentViscosityWallFunctionBC(const InputParameters &parameters)
Class constructor.
template Real findUStar< Real >(const Real &mu, const Real &rho, const Real &u, const Real dist)
Real distance(const Point &p)
void addRequiredParam(const std::string &name, const std::string &doc_string)
Class implementing a Dirichlet boundary condition for the turbulent viscosity wall function in a RANS...
Moose::ElemArg makeElemArg(const Elem *elem, bool correct_skewnewss=false) const
template Real findyPlus< Real >(const Real &mu, const Real &rho, const Real &u, Real dist)
const Elem * neighborPtr() const
static const std::string mu
Definition: NS.h:123
FaceInfo::VarFaceNeighbors _current_face_type
const Point & normal() const
MooseLinearVariableFV< Real > & _var
const Moose::Functor< Real > & _k
Turbulent kinetic energy.
const Elem * elemPtr() const
const Moose::Functor< Real > & _mu
Dynamic viscosity.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static constexpr Real E_turb_constant
Definition: NS.h:202
const FaceInfo * _current_face_info
const NS::WallTreatmentEnum _wall_treatment
Method used for wall treatment.
void mooseError(Args &&... args) const
void addClassDescription(const std::string &doc_string)
registerMooseObject("NavierStokesApp", LinearFVTurbulentViscosityWallFunctionBC)
static const std::string velocity
Definition: NS.h:45
template Real computeSpeed< Real >(const libMesh::VectorValue< Real > &velocity)
static constexpr Real mu_t_low_limit
Definition: NS.h:204
static InputParameters validParams()
MooseUnits pow(const MooseUnits &, int)