https://mooseframework.inl.gov
INSFVTurbulentViscosityWallFunction.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 "Function.h"
12 #include "NavierStokesMethods.h"
13 
15 
18 {
20  params.addClassDescription("Adds Dirichlet BC for wall values of the turbulent viscosity.");
21  params.addRequiredParam<MooseFunctorName>("u", "The velocity in the x direction.");
22  params.addParam<MooseFunctorName>("v", "The velocity in the y direction.");
23  params.addParam<MooseFunctorName>("w", "The velocity in the z direction.");
24  params.addRequiredParam<MooseFunctorName>(NS::density, "Density");
25  params.addRequiredParam<MooseFunctorName>(NS::mu, "Dynamic viscosity.");
26  params.addRequiredParam<MooseFunctorName>(NS::mu_t, "The turbulent viscosity.");
27  params.addParam<MooseFunctorName>(NS::TKE, "The turbulent kinetic energy.");
28  params.addParam<Real>("C_mu", 0.09, "Coupled turbulent kinetic energy closure.");
29 
30  MooseEnum wall_treatment("eq_newton eq_incremental eq_linearized neq", "neq");
31  params.addParam<MooseEnum>(
32  "wall_treatment", wall_treatment, "The method used for computing the wall functions");
33  return params;
34 }
35 
37  const InputParameters & params)
38  : FVDirichletBCBase(params),
39  _dim(_subproblem.mesh().dimension()),
40  _u_var(getFunctor<ADReal>("u")),
41  _v_var(params.isParamValid("v") ? &(getFunctor<ADReal>("v")) : nullptr),
42  _w_var(params.isParamValid("w") ? &(getFunctor<ADReal>("w")) : nullptr),
43  _rho(getFunctor<ADReal>(NS::density)),
44  _mu(getFunctor<ADReal>(NS::mu)),
45  _mu_t(getFunctor<ADReal>(NS::mu_t)),
46  _k(getFunctor<ADReal>(NS::TKE)),
47  _C_mu(getParam<Real>("C_mu")),
48  _wall_treatment(getParam<MooseEnum>("wall_treatment").getEnum<NS::WallTreatmentEnum>()),
49  _preserve_sparsity_pattern(_fv_problem.preserveMatrixSparsityPattern())
50 {
51 }
52 
53 ADReal
55  const Moose::StateArg & /* state */) const
56 {
57  using std::abs, std::max, std::sqrt, std::pow, std::log;
58 
59  // if on an internal face (internal to the mesh, but an external boundary of the flow area),
60  // we have to select the element on which the flow variables / viscosity is defined
62  const Real wall_dist = abs(
63  ((use_elem ? fi.elemCentroid() : fi.neighborCentroid()) - fi.faceCentroid()) * fi.normal());
64  const Elem * const elem_ptr = use_elem ? fi.elemPtr() : fi.neighborPtr();
65  const auto elem_arg = makeElemArg(elem_ptr);
66 
67  // Get the previous non linear iteration values
69  const auto mu = _mu(elem_arg, old_state);
70  const auto rho = _rho(elem_arg, old_state);
71 
72  // Get the velocity vector
73  ADRealVectorValue velocity(_u_var(elem_arg, old_state));
74  if (_v_var)
75  velocity(1) = (*_v_var)(elem_arg, old_state);
76  if (_w_var)
77  velocity(2) = (*_w_var)(elem_arg, old_state);
78 
79  // Compute the velocity and direction of the velocity component that is parallel to the wall
80  const auto parallel_speed =
82 
83  // Switch for determining the near wall quantities
84  // wall_treatment can be: "eq_newton eq_incremental eq_linearized neq"
85  ADReal y_plus;
86  ADReal mut_log; // turbulent log-layer viscosity
87  ADReal mu_wall; // total wall viscosity to obtain the shear stress at the wall
88 
90  {
91  // Full Newton-Raphson solve to find the wall quantities from the law of the wall
92  const auto u_tau = NS::findUStar<ADReal>(mu, rho, parallel_speed, wall_dist);
93  y_plus = wall_dist * u_tau * rho / mu;
94  mu_wall = rho * Utility::pow<2>(u_tau) * wall_dist / parallel_speed;
95  mut_log = mu_wall - mu;
96  }
98  {
99  // Incremental solve on y_plus to get the near-wall quantities
100  y_plus = NS::findyPlus<ADReal>(mu, rho, max(parallel_speed, 1e-10), wall_dist);
101  mu_wall =
102  mu * (NS::von_karman_constant * y_plus / log(max(NS::E_turb_constant * y_plus, 1 + 1e-4)));
103  mut_log = mu_wall - mu;
104  }
106  {
107  // Linearized approximation to the wall function to find the near-wall quantities faster
108  const ADReal a_c = 1 / NS::von_karman_constant;
109  const ADReal b_c =
110  1 / NS::von_karman_constant * (log(NS::E_turb_constant * max(wall_dist, 1.0) / mu) + 1.0);
111  const ADReal c_c = parallel_speed;
112 
113  const auto u_tau = (-b_c + sqrt(pow(b_c, 2) + 4.0 * a_c * c_c)) / (2.0 * a_c);
114  y_plus = wall_dist * u_tau * rho / mu;
115  mu_wall = rho * Utility::pow<2>(u_tau) * wall_dist / parallel_speed;
116  mut_log = mu_wall - mu;
117  }
119  {
120  // Assign non-equilibrium wall function value
121  y_plus = pow(_C_mu, 0.25) * wall_dist * sqrt(_k(elem_arg, old_state)) * rho / mu;
122  mu_wall = mu * (NS::von_karman_constant * y_plus /
123  log(max(NS::E_turb_constant * y_plus, 1.0 + 1e-4)));
124  mut_log = mu_wall - mu;
125  }
126  else
127  mooseAssert(false,
128  "For `INSFVTurbulentViscosityWallFunction` , wall treatment should not reach here");
129 
130  ADReal residual = 0;
131  // To keep the same sparsity pattern for all y_plus
133  residual = 0 * mut_log * y_plus;
134 
135  if (y_plus <= 5.0)
136  // sub-laminar layer
137  residual += 0.0;
138  else if (y_plus >= 30.0)
139  // log-layer
140  residual += max(mut_log, NS::mu_t_low_limit);
141  else
142  {
143  // buffer layer
144  const auto blending_function = (y_plus - 5.0) / 25.0;
145  // the blending depends on the mut_log at y+=30
146  const auto mut_log = mu * _mut_30;
147  residual += max(blending_function * mut_log, NS::mu_t_low_limit);
148  }
149  return residual;
150 }
MetaPhysicL::DualNumber< V, D, asd > abs(const MetaPhysicL::DualNumber< V, D, asd > &a)
static constexpr Real von_karman_constant
Definition: NS.h:201
static const std::string mu_t
Definition: NS.h:125
const std::pair< unsigned int, unsigned int > _var_sys_numbers_pair
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
const Moose::Functor< ADReal > * _v_var
y-velocity
Applies a wall function to the turbulent viscosity field.
const Point & faceCentroid() const
MeshBase & mesh
static const std::string density
Definition: NS.h:33
static const std::string TKE
Definition: NS.h:176
WallTreatmentEnum
Wall treatment options.
Definition: NS.h:182
const Moose::Functor< ADReal > * _w_var
z-velocity
const Point & neighborCentroid() const
DualNumber< Real, DNDerivativeType, true > ADReal
void addRequiredParam(const std::string &name, const std::string &doc_string)
auto max(const L &left, const R &right)
Moose::ElemArg makeElemArg(const Elem *elem, bool correct_skewnewss=false) const
INSFVTurbulentViscosityWallFunction(const InputParameters &parameters)
ADReal boundaryValue(const FaceInfo &fi, const Moose::StateArg &state) const override
const Elem * neighborPtr() const
const Point & elemCentroid() const
static InputParameters validParams()
static const std::string mu
Definition: NS.h:123
const Moose::Functor< ADReal > & _k
Turbulent kinetic energy.
registerMooseObject("NavierStokesApp", INSFVTurbulentViscosityWallFunction)
NS::WallTreatmentEnum _wall_treatment
Method used for wall treatment.
const Moose::Functor< ADReal > & _u_var
x-velocity
template ADReal findUStar< ADReal >(const ADReal &mu, const ADReal &rho, const ADReal &u, const Real dist)
const Point & normal() const
auto log(const T &)
const Elem * elemPtr() const
const Real _C_mu
C_mu turbulent coefficient.
ExpressionBuilder::EBTerm pow(const ExpressionBuilder::EBTerm &left, T exponent)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
template ADReal computeSpeed< ADReal >(const libMesh::VectorValue< ADReal > &velocity)
static constexpr Real E_turb_constant
Definition: NS.h:202
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template * sqrt(_arg)) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(tanh
template ADReal findyPlus< ADReal >(const ADReal &mu, const ADReal &rho, const ADReal &u, Real dist)
void addClassDescription(const std::string &doc_string)
const bool _preserve_sparsity_pattern
For Newton solves we want to add extra zero-valued terms regardless of y-plus to avoid sparsity patte...
static const std::string velocity
Definition: NS.h:45
const Moose::Functor< ADReal > & _rho
Density.
static constexpr Real mu_t_low_limit
Definition: NS.h:204
const Moose::Functor< ADReal > & _mu
Dynamic viscosity.
MooseUnits pow(const MooseUnits &, int)
VarFaceNeighbors faceType(const std::pair< unsigned int, unsigned int > &var_sys) const