https://mooseframework.inl.gov
INSFVTKEDWallFunctionBC.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 #include "NS.h"
14 
16 
19 {
21  params.addClassDescription("Adds Reichardt extrapolated wall values to set up directly the"
22  "Dirichlet BC for the turbulent kinetic energy dissipation rate.");
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::mu_t, "The turbulent viscosity.");
29  params.addRequiredParam<MooseFunctorName>(NS::TKE, "The turbulent kinetic energy.");
30  params.addParam<MooseFunctorName>("C_mu", 0.09, "Coupled turbulent kinetic energy closure.");
31  params.addParam<bool>("newton_solve", false, "Whether a Newton nonlinear solve is being used");
32  params.addParamNamesToGroup("newton_solve", "Advanced");
33  return params;
34 }
35 
37  : FVDirichletBCBase(params),
38  _dim(_subproblem.mesh().dimension()),
39  _u_var(getFunctor<ADReal>("u")),
40  _v_var(params.isParamValid("v") ? &(getFunctor<ADReal>("v")) : nullptr),
41  _w_var(params.isParamValid("w") ? &(getFunctor<ADReal>("w")) : nullptr),
42  _rho(getFunctor<ADReal>(NS::density)),
43  _mu(getFunctor<ADReal>(NS::mu)),
44  _mu_t(getFunctor<ADReal>(NS::mu_t)),
45  _k(getFunctor<ADReal>(NS::TKE)),
46  _C_mu(getFunctor<ADReal>("C_mu")),
47  _newton_solve(getParam<bool>("newton_solve"))
48 {
49 }
50 
51 ADReal
53 {
54  using std::pow, std::abs;
55 
56  const bool use_elem = (fi.faceType(std::make_pair(_var.number(), _var.sys().number())) ==
58  const Real dist = abs(
59  ((use_elem ? fi.elemCentroid() : fi.neighborCentroid()) - fi.faceCentroid()) * fi.normal());
60  const Elem * const elem_ptr = use_elem ? fi.elemPtr() : fi.neighborPtr();
61  const auto elem_arg = makeElemArg(elem_ptr);
62  const auto mu = _mu(elem_arg, state);
63  const auto rho = _rho(elem_arg, state);
64 
65  // Assign boundary weights to element
66  // This is based on the theory of linear TKE development for each boundary
67  // This is, it assumes no interaction across turbulence production from boundaries
68  Real weight = 0.0;
69  for (unsigned int i_side = 0; i_side < elem_ptr->n_sides(); ++i_side)
70  weight += static_cast<Real>(_subproblem.mesh().getBoundaryIDs(elem_ptr, i_side).size());
71 
72  // Get the velocity vector
73  ADRealVectorValue velocity(_u_var(elem_arg, state));
74  if (_v_var)
75  velocity(1) = (*_v_var)(elem_arg, state);
76  if (_w_var)
77  velocity(2) = (*_w_var)(elem_arg, state);
78 
79  // Compute the velocity and direction of the velocity component that is parallel to the wall
80  const ADReal parallel_speed =
82 
83  // Get friction velocity
84  const ADReal u_star = NS::findUStar<ADReal>(mu, rho, parallel_speed, dist);
85 
86  // Get associated non-dimensional wall distance
87  const ADReal y_plus = dist * u_star * rho / mu;
88 
89  const auto TKE = _k(elem_arg, state);
90 
91  if (y_plus <= 5.0) // sub-laminar layer
92  {
93  const auto laminar_value = 2.0 * weight * TKE * mu / pow(dist, 2);
94  if (!_newton_solve)
95  return laminar_value;
96  else
97  // Additional zero term to make sure new derivatives are not introduced as y_plus
98  // changes
99  return laminar_value + 0 * _mu_t(elem_arg, state);
100  }
101  else if (y_plus >= 30.0) // log-layer
102  {
103  const auto turbulent_value =
104  weight * _C_mu(elem_arg, state) * pow(abs(TKE), 1.5) / (_mu_t(elem_arg, state) * dist);
105  if (!_newton_solve)
106  return turbulent_value;
107  else
108  // Additional zero term to make sure new derivatives are not introduced as y_plus changes
109  return turbulent_value + 0 * mu;
110  }
111  else // blending function
112  {
113  const auto laminar_value = 2.0 * weight * TKE * mu / pow(dist, 2);
114  const auto turbulent_value =
115  weight * _C_mu(elem_arg, state) * pow(abs(TKE), 1.5) / (_mu_t(elem_arg, state) * dist);
116  const auto interpolation_coef = (y_plus - 5.0) / 25.0;
117  return (interpolation_coef * (turbulent_value - laminar_value) + laminar_value);
118  }
119 }
const Moose::Functor< ADReal > & _rho
Density.
virtual MooseMesh & mesh()=0
MetaPhysicL::DualNumber< V, D, asd > abs(const MetaPhysicL::DualNumber< V, D, asd > &a)
ADReal boundaryValue(const FaceInfo &fi, const Moose::StateArg &state) const override
const Moose::Functor< ADReal > & _mu
Dynamic viscosity.
static const std::string mu_t
Definition: NS.h:125
const Moose::Functor< ADReal > & _mu_t
Turbulent dynamic viscosity.
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
unsigned int number() const
const Point & faceCentroid() const
MeshBase & mesh
static const std::string density
Definition: NS.h:33
static const std::string TKE
Definition: NS.h:176
const Moose::Functor< ADReal > & _C_mu
C_mu turbulent coefficient.
static InputParameters validParams()
MooseVariableFV< Real > & _var
const Moose::Functor< ADReal > * _v_var
y-velocity
const Point & neighborCentroid() const
DualNumber< Real, DNDerivativeType, true > ADReal
const Moose::Functor< ADReal > & _k
Turbulent kinetic energy.
void addRequiredParam(const std::string &name, const std::string &doc_string)
const Moose::Functor< ADReal > & _u_var
x-velocity
Moose::ElemArg makeElemArg(const Elem *elem, bool correct_skewnewss=false) const
const Elem * neighborPtr() const
const Point & elemCentroid() const
static InputParameters validParams()
dof_id_type weight(const MeshBase &mesh, const processor_id_type pid)
static const std::string mu
Definition: NS.h:123
SubProblem & _subproblem
template ADReal findUStar< ADReal >(const ADReal &mu, const ADReal &rho, const ADReal &u, const Real dist)
const Point & normal() const
unsigned int number() const
const Elem * elemPtr() const
ExpressionBuilder::EBTerm pow(const ExpressionBuilder::EBTerm &left, T exponent)
const Moose::Functor< ADReal > * _w_var
z-velocity
registerMooseObject("NavierStokesApp", INSFVTKEDWallFunctionBC)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
template ADReal computeSpeed< ADReal >(const libMesh::VectorValue< ADReal > &velocity)
const bool _newton_solve
Whether we are using a newton solve.
Applies a wall function to the turbulent kinetic energy dissipation rate.
void addClassDescription(const std::string &doc_string)
static const std::string velocity
Definition: NS.h:45
std::vector< BoundaryID > getBoundaryIDs(const Elem *const elem, const unsigned short int side) const
MooseUnits pow(const MooseUnits &, int)
INSFVTKEDWallFunctionBC(const InputParameters &parameters)
VarFaceNeighbors faceType(const std::pair< unsigned int, unsigned int > &var_sys) const
void addParamNamesToGroup(const std::string &space_delim_names, const std::string group_name)