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.");
25 params.
addParam<MooseFunctorName>(
NS::TKE,
"The turbulent kinetic energy.");
26 params.
addParam<
Real>(
"C_mu", 0.09,
"Coupled turbulent kinetic energy closure.");
28 MooseEnum wall_treatment(
"eq_newton eq_incremental eq_linearized neq",
"neq");
30 "wall_treatment", wall_treatment,
"The method used for computing the wall functions");
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),
44 _C_mu(getParam<
Real>(
"C_mu")),
59 const auto mu =
_mu(re, old_state);
60 const auto rho =
_rho(re, old_state);
65 velocity(1) = (*_v_var)(re, old_state);
67 velocity(2) = (*_w_var)(re, old_state);
82 y_plus = wall_dist * u_tau * rho /
mu;
83 mu_wall = rho * Utility::pow<2>(u_tau) * wall_dist / parallel_speed;
98 const Real c_c = parallel_speed;
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;
107 y_plus =
std::pow(
_C_mu, 0.25) * wall_dist * std::sqrt(
_k(re, old_state)) * rho /
mu;
113 "For `INSFVTurbulentViscosityWallFunction` , wall treatment should not reach here");
115 const Real mut_log = mu_wall -
mu;
122 else if (y_plus >= 30.0)
128 const auto blending_function = (y_plus - 5.0) / 25.0;
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.");
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.");
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.");
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.");
static constexpr Real von_karman_constant
static const std::string mu_t
Real computeTurbulentViscosity() const
Real computeCellToFaceDistance() const
Moose::StateArg determineState() const
virtual Real computeBoundaryValueMatrixContribution() const override
static const std::string density
static const std::string TKE
WallTreatmentEnum
Wall treatment options.
LinearFVTurbulentViscosityWallFunctionBC(const InputParameters ¶meters)
Class constructor.
template Real findUStar< Real >(const Real &mu, const Real &rho, const Real &u, const Real dist)
Real distance(const Point &p)
Class implementing a Dirichlet boundary condition for the turbulent viscosity wall function in a RANS...
virtual Real computeBoundaryNormalGradient() const override
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)
virtual Real computeBoundaryGradientMatrixContribution() const override
const Elem * neighborPtr() const
const Moose::Functor< Real > * _v_var
y-velocity
static const std::string mu
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.
const Moose::Functor< Real > * _w_var
z-velocity
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual Real computeBoundaryGradientRHSContribution() const override
static constexpr Real E_turb_constant
const FaceInfo * _current_face_info
const NS::WallTreatmentEnum _wall_treatment
Method used for wall treatment.
void mooseError(Args &&... args) const
registerMooseObject("NavierStokesApp", LinearFVTurbulentViscosityWallFunctionBC)
static const std::string velocity
virtual Real computeBoundaryValue() const override
template Real computeSpeed< Real >(const libMesh::VectorValue< Real > &velocity)
const Real _C_mu
C_mu turbulent coefficient.
static constexpr Real mu_t_low_limit
static InputParameters validParams()
static InputParameters validParams()
const Moose::Functor< Real > & _u_var
x-velocity
MooseUnits pow(const MooseUnits &, int)
const Moose::Functor< Real > & _rho
Density.
virtual Real computeBoundaryValueRHSContribution() const override