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.");
27 params.
addParam<MooseFunctorName>(
NS::TKE,
"The turbulent kinetic energy.");
28 params.
addParam<
Real>(
"C_mu", 0.09,
"Coupled turbulent kinetic energy closure.");
30 MooseEnum wall_treatment(
"eq_newton eq_incremental eq_linearized neq",
"neq");
32 "wall_treatment", wall_treatment,
"The method used for computing the wall functions");
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),
47 _C_mu(getParam<
Real>(
"C_mu")),
49 _preserve_sparsity_pattern(_fv_problem.preserveMatrixSparsityPattern())
57 using std::abs, std::max, std::sqrt,
std::pow, std::log;
69 const auto mu =
_mu(elem_arg, old_state);
70 const auto rho =
_rho(elem_arg, old_state);
75 velocity(1) = (*_v_var)(elem_arg, old_state);
77 velocity(2) = (*_w_var)(elem_arg, old_state);
80 const auto parallel_speed =
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;
103 mut_log = mu_wall -
mu;
111 const ADReal c_c = parallel_speed;
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;
121 y_plus =
pow(
_C_mu, 0.25) * wall_dist *
sqrt(
_k(elem_arg, old_state)) * rho /
mu;
124 mut_log = mu_wall -
mu;
128 "For `INSFVTurbulentViscosityWallFunction` , wall treatment should not reach here");
133 residual = 0 * mut_log * y_plus;
138 else if (y_plus >= 30.0)
144 const auto blending_function = (y_plus - 5.0) / 25.0;
MetaPhysicL::DualNumber< V, D, asd > abs(const MetaPhysicL::DualNumber< V, D, asd > &a)
static constexpr Real von_karman_constant
static const std::string mu_t
const std::pair< unsigned int, unsigned int > _var_sys_numbers_pair
const Moose::Functor< ADReal > * _v_var
y-velocity
Applies a wall function to the turbulent viscosity field.
static InputParameters validParams()
const Point & faceCentroid() const
static const std::string density
static const std::string TKE
WallTreatmentEnum
Wall treatment options.
const Moose::Functor< ADReal > * _w_var
z-velocity
const Point & neighborCentroid() const
DualNumber< Real, DNDerivativeType, true > ADReal
auto max(const L &left, const R &right)
Moose::ElemArg makeElemArg(const Elem *elem, bool correct_skewnewss=false) const
INSFVTurbulentViscosityWallFunction(const InputParameters ¶meters)
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
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
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
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)
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
const Moose::Functor< ADReal > & _rho
Density.
static constexpr Real mu_t_low_limit
const Moose::Functor< ADReal > & _mu
Dynamic viscosity.
MooseUnits pow(const MooseUnits &, int)
VarFaceNeighbors faceType(const std::pair< unsigned int, unsigned int > &var_sys) const