https://mooseframework.inl.gov
INSFVTurbulentTemperatureWallFunction.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 #include "NS.h"
13 
15 
18 {
20 
21  params.addClassDescription("Adds turbulent temperature wall function.");
22  params.addRequiredParam<MooseFunctorName>("T_w", "The wall temperature.");
23  params.addRequiredParam<MooseFunctorName>("u", "The velocity in the x direction.");
24  params.addParam<MooseFunctorName>("v", "The velocity in the y direction.");
25  params.addParam<MooseFunctorName>("w", "The velocity in the z direction.");
26  params.addRequiredParam<MooseFunctorName>(NS::density, "Density");
27  params.addRequiredParam<MooseFunctorName>(NS::mu, "Dynamic viscosity.");
28  params.addRequiredParam<MooseFunctorName>(NS::cp, "The specific heat at constant pressure.");
29  params.addParam<MooseFunctorName>(NS::kappa, "The thermal conductivity.");
30  params.addParam<MooseFunctorName>("Pr_t", 0.58, "The turbulent Prandtl number.");
31  params.addParam<MooseFunctorName>("k", "Turbulent kinetic energy functor.");
32  params.deprecateParam("k", NS::TKE, "01/01/2025");
33  params.addParam<Real>("C_mu", 0.09, "Coupled turbulent kinetic energy closure.");
34  MooseEnum wall_treatment("eq_newton eq_incremental eq_linearized neq", "neq");
35  params.addParam<MooseEnum>(
36  "wall_treatment", wall_treatment, "The method used for computing the wall functions");
37 
38  params.addParam<bool>("newton_solve", false, "Whether a Newton nonlinear solve is being used");
39  params.addParamNamesToGroup("newton_solve", "Advanced");
40  return params;
41 }
42 
44  const InputParameters & parameters)
45  : FVFluxBC(parameters),
46  _dim(_subproblem.mesh().dimension()),
47  _T_w(getFunctor<ADReal>("T_w")),
48  _u_var(getFunctor<ADReal>("u")),
49  _v_var(parameters.isParamValid("v") ? &(getFunctor<ADReal>("v")) : nullptr),
50  _w_var(parameters.isParamValid("w") ? &(getFunctor<ADReal>("w")) : nullptr),
51  _rho(getFunctor<ADReal>(NS::density)),
52  _mu(getFunctor<ADReal>(NS::mu)),
53  _cp(getFunctor<ADReal>(NS::cp)),
54  _kappa(getFunctor<ADReal>(NS::kappa)),
55  _Pr_t(getFunctor<ADReal>("Pr_t")),
56  _k(getFunctor<ADReal>(NS::TKE)),
57  _C_mu(getParam<Real>("C_mu")),
58  _wall_treatment(getParam<MooseEnum>("wall_treatment")),
59  _newton_solve(getParam<bool>("newton_solve"))
60 {
61 }
62 
63 ADReal
65 {
66  // Useful parameters
67  const FaceInfo & fi = *_face_info;
68  const Real wall_dist = std::abs((fi.elemCentroid() - fi.faceCentroid()) * fi.normal());
69  const Elem & _current_elem = fi.elem();
70  const auto current_argument = makeElemArg(&_current_elem);
71  const auto state = determineState();
73  const auto mu = _mu(current_argument, state);
74  const auto rho = _rho(current_argument, state);
75  const auto cp = _cp(current_argument, state);
76  const auto kappa = _kappa(current_argument, state);
77 
78  // Get the velocity vector
79  ADRealVectorValue velocity(_u_var(current_argument, state));
80  if (_v_var)
81  velocity(1) = (*_v_var)(current_argument, state);
82  if (_w_var)
83  velocity(2) = (*_w_var)(current_argument, state);
84 
85  // Compute the velocity and direction of the velocity component that is parallel to the wall
86  const ADReal parallel_speed =
87  NS::computeSpeed(velocity - velocity * (fi.normal()) * (fi.normal()));
88 
89  // Computing friction velocity and y+
90  ADReal u_tau, y_plus;
91 
92  if (_wall_treatment == "eq_newton")
93  {
94  // Full Newton-Raphson solve to find the wall quantities from the law of the wall
95  u_tau = NS::findUStar(mu, rho, parallel_speed, wall_dist);
96  y_plus = wall_dist * u_tau * rho / mu;
97  }
98  else if (_wall_treatment == "eq_incremental")
99  {
100  // Incremental solve on y_plus to get the near-wall quantities
101  y_plus = NS::findyPlus(mu, rho, std::max(parallel_speed, 1e-10), wall_dist);
102  u_tau = parallel_speed / (std::log(std::max(NS::E_turb_constant * y_plus, 1.0 + 1e-4)) /
104  }
105  else if (_wall_treatment == "eq_linearized")
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 = 1.0 / NS::von_karman_constant *
110  (std::log(NS::E_turb_constant * std::max(wall_dist, 1.0) / mu) + 1.0);
111  const ADReal c_c = parallel_speed;
112 
113  u_tau = (-b_c + std::sqrt(std::pow(b_c, 2) + 4.0 * a_c * c_c)) / (2.0 * a_c);
114  y_plus = wall_dist * u_tau * rho / mu;
115  }
116  else if (_wall_treatment == "neq")
117  {
118  // Assign non-equilibrium wall function value
119  y_plus = wall_dist * std::sqrt(std::sqrt(_C_mu) * _k(current_argument, old_state)) * rho / mu;
120  u_tau = parallel_speed / (std::log(std::max(NS::E_turb_constant * y_plus, 1.0 + 1e-4)) /
122  }
123 
124  ADReal alpha;
125  if (y_plus <= 5.0) // sub-laminar layer
126  {
127  alpha = kappa / (rho * cp);
128  }
129  else if (y_plus >= 30.0) // log-layer
130  {
131  const auto Pr = cp * mu / kappa;
132  const auto Pr_ratio = Pr / _Pr_t(current_argument, state);
133  const auto jayatilleke_P =
134  9.24 * (std::pow(Pr_ratio, 0.75) - 1.0) * (1.0 + 0.28 * std::exp(-0.007 * Pr_ratio));
135  const auto wall_scaling =
136  1.0 / NS::von_karman_constant * std::log(NS::E_turb_constant * y_plus) + jayatilleke_P;
137  alpha = u_tau * wall_dist / (_Pr_t(current_argument, state) * wall_scaling);
138  }
139  else // buffer layer
140  {
141  const auto alpha_lam = kappa / (rho * cp);
142  const auto Pr = cp * mu / kappa;
143  const auto Pr_t = _Pr_t(current_argument, state);
144  const auto Pr_ratio = Pr / Pr_t;
145  const auto jayatilleke_P =
146  9.24 * (std::pow(Pr_ratio, 0.75) - 1.0) * (1.0 + 0.28 * std::exp(-0.007 * Pr_ratio));
147  const auto wall_scaling =
148  1.0 / NS::von_karman_constant * std::log(NS::E_turb_constant * y_plus) + jayatilleke_P;
149  const auto alpha_turb = u_tau * wall_dist / (Pr_t * wall_scaling);
150  const auto blending_function = (y_plus - 5.0) / 25.0;
151  alpha = blending_function * alpha_turb + (1.0 - blending_function) * alpha_lam;
152  }
153 
154  // To make sure new derivatives are not introduced as the solve progresses
155  if (_newton_solve)
156  alpha += 0 * kappa * (rho * cp) + 0 * u_tau * _Pr_t(current_argument, state) * y_plus;
157 
158  const auto face_arg = singleSidedFaceArg();
159  return -rho * cp * alpha * (_T_w(face_arg, state) - _var(current_argument, state)) / wall_dist;
160 }
static constexpr Real von_karman_constant
Definition: NS.h:191
INSFVTurbulentTemperatureWallFunction(const InputParameters &parameters)
const FaceInfo * _face_info
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
static InputParameters validParams()
const Moose::Functor< ADReal > & _mu
Dynamic viscosity.
const Elem & elem() const
const bool _newton_solve
For Newton solves we want to add extra zero-valued terms regardless of y-plus to avoid sparsity patte...
Moose::StateArg determineState() const
const Point & faceCentroid() const
registerADMooseObject("NavierStokesApp", INSFVTurbulentTemperatureWallFunction)
MeshBase & mesh
static const std::string density
Definition: NS.h:33
static const std::string TKE
Definition: NS.h:176
Moose::FaceArg singleSidedFaceArg(const FaceInfo *fi=nullptr, Moose::FV::LimiterType limiter_type=Moose::FV::LimiterType::CentralDifference, bool correct_skewness=false, const Moose::StateArg *state_limiter=nullptr) const
MooseVariableFV< Real > & _var
DualNumber< Real, DNDerivativeType, true > ADReal
void addRequiredParam(const std::string &name, const std::string &doc_string)
Moose::ElemArg makeElemArg(const Elem *elem, bool correct_skewnewss=false) const
static const std::string cp
Definition: NS.h:121
const Moose::Functor< ADReal > & _rho
Density.
const Moose::Functor< ADReal > & _u_var
x-velocity
const Point & elemCentroid() const
void deprecateParam(const std::string &old_name, const std::string &new_name, const std::string &removal_date)
static const std::string mu
Definition: NS.h:123
const MooseEnum _wall_treatment
Method used for wall treatment.
ADReal findyPlus(const ADReal &mu, const ADReal &rho, const ADReal &u, Real dist)
Finds the non-dimensional wall distance normalized with the friction velocity Implements a fixed-poin...
const Moose::Functor< ADReal > & _T_w
Wall Temperature.
const Point & normal() const
const Moose::Functor< ADReal > & _kappa
Thermal conductivity.
const Moose::Functor< ADReal > & _Pr_t
Turbulent Prandtl number near the wall.
const Moose::Functor< ADReal > & _cp
The specific heat at constant pressure.
static const std::string kappa
Definition: NS.h:116
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Moose::Functor< ADReal > & _k
Turbulent kinetic energy.
static constexpr Real E_turb_constant
Definition: NS.h:192
static const std::string alpha
Definition: NS.h:134
const Moose::Functor< ADReal > * _w_var
z-velocity
void addClassDescription(const std::string &doc_string)
static const std::string velocity
Definition: NS.h:45
ADReal findUStar(const ADReal &mu, const ADReal &rho, const ADReal &u, Real dist)
Finds the friction velocity using standard velocity wall functions formulation.
const Moose::Functor< ADReal > * _v_var
y-velocity
ADReal computeSpeed(const ADRealVectorValue &velocity)
Compute the speed (velocity norm) given the supplied velocity.
This boundary condition applies a wall function for the energy equation for turbulent flows...
MooseUnits pow(const MooseUnits &, int)
void addParamNamesToGroup(const std::string &space_delim_names, const std::string group_name)