13 #include "libmesh/nonlinear_solver.h" 22 " terms of turbulent kinetic energy dissipation (TKED).");
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>(
"k",
"Coupled turbulent kinetic energy.");
31 params.
addParam<std::vector<BoundaryName>>(
32 "walls", {},
"Boundaries that correspond to solid walls.");
36 "Boolean to determine if the problem should be used in a linear or nonlinear solve");
37 MooseEnum wall_treatment(
"eq_newton eq_incremental eq_linearized neq",
"neq");
40 "The method used for computing the wall functions " 41 "'eq_newton', 'eq_incremental', 'eq_linearized', 'neq'");
42 params.
addParam<
Real>(
"C1_eps", 1.44,
"First epsilon coefficient");
43 params.
addParam<
Real>(
"C2_eps", 1.92,
"Second epsilon coefficient");
44 params.
addParam<
Real>(
"C_mu", 0.09,
"Coupled turbulent kinetic energy closure.");
45 params.
addParam<
Real>(
"C_pl", 10.0,
"Production limiter constant multiplier.");
46 params.
set<
unsigned short>(
"ghost_layers") = 2;
47 params.
addParam<
bool>(
"newton_solve",
false,
"Whether a Newton nonlinear solve is being used");
54 _dim(_subproblem.
mesh().dimension()),
55 _u_var(getFunctor<
ADReal>(
"u")),
56 _v_var(params.isParamValid(
"v") ? &(getFunctor<
ADReal>(
"v")) : nullptr),
57 _w_var(params.isParamValid(
"w") ? &(getFunctor<
ADReal>(
"w")) : nullptr),
62 _wall_boundary_names(getParam<
std::vector<BoundaryName>>(
"walls")),
63 _linearized_model(getParam<bool>(
"linearized_model")),
65 _C1_eps(getParam<
Real>(
"C1_eps")),
66 _C2_eps(getParam<
Real>(
"C2_eps")),
67 _C_mu(getParam<
Real>(
"C_mu")),
68 _C_pl(getParam<
Real>(
"C_pl")),
69 _newton_solve(getParam<bool>(
"newton_solve"))
72 paramError(
"v",
"In two or more dimensions, the v velocity must be supplied!");
75 paramError(
"w",
"In three or more dimensions, the w velocity must be supplied!");
95 const auto old_state =
97 const auto mu =
_mu(elem_arg, state);
98 const auto rho =
_rho(elem_arg, state);
100 _newton_solve ? std::max(
_k(elem_arg, old_state), 1e-10) :
_k(elem_arg, old_state);
105 std::vector<ADReal> y_plus_vec;
107 Real tot_weight = 0.0;
111 velocity(1) = (*_v_var)(elem_arg, state);
113 velocity(2) = (*_w_var)(elem_arg, state);
117 mooseAssert(distance_vec.size(),
"Should have found a distance vector");
118 mooseAssert(distance_vec.size() == face_info_vec.size(),
119 "Should be as many distance vectors as face info vectors");
121 for (
unsigned int i = 0; i < distance_vec.size(); i++)
123 const auto distance = distance_vec[i];
124 mooseAssert(
distance > 0,
"Should be at a non-zero distance");
132 velocity -
velocity * face_info_vec[i]->normal() * face_info_vec[i]->normal());
137 y_plus_vec.push_back(y_plus);
144 const auto y_plus = y_plus_vec[i];
148 const auto fi = face_info_vec[i];
150 const Elem *
const loc_elem = defined_on_elem_side ? &fi->elem() : fi->neighborPtr();
153 destruction += 2.0 * TKE_old *
_mu(facearg, state) / rho /
154 Utility::pow<2>(distance_vec[i]) / tot_weight;
165 const auto & grad_u =
_u_var.gradient(elem_arg, state);
166 const auto Sij_xx = 2.0 * grad_u(0);
173 const auto grad_xx = grad_u(0);
183 auto trace = Sij_xx / 3.0;
187 const auto & grad_v = (*_v_var).gradient(elem_arg, state);
188 Sij_xy = grad_u(1) + grad_v(0);
189 Sij_yy = 2.0 * grad_v(1);
195 trace += Sij_yy / 3.0;
199 const auto & grad_w = (*_w_var).gradient(elem_arg, state);
201 Sij_xz = grad_u(2) + grad_w(0);
202 Sij_yz = grad_v(2) + grad_w(1);
203 Sij_zz = 2.0 * grad_w(2);
211 trace += Sij_zz / 3.0;
215 const auto symmetric_strain_tensor_sq_norm =
216 (Sij_xx - trace) * grad_xx + Sij_xy * grad_xy + Sij_xz * grad_xz + Sij_xy * grad_yx +
217 (Sij_yy - trace) * grad_yy + Sij_yz * grad_yz + Sij_xz * grad_zx + Sij_yz * grad_zy +
218 (Sij_zz - trace) * grad_zz;
220 ADReal production_k =
_mu_t(elem_arg, state) * symmetric_strain_tensor_sq_norm;
224 const ADReal production_limit =
_C_pl * rho * eps_old;
226 production_k = std::min(production_k, production_limit);
229 production =
_C1_eps * production_k / time_scale;
230 destruction =
_C2_eps * rho *
_var(elem_arg, state) / time_scale;
232 residual = destruction - production;
std::map< const Elem *, bool > _wall_bounded
Maps for wall treatment.
const Moose::Functor< ADReal > * _w_var
z-velocity
static constexpr Real von_karman_constant
static const std::string mu_t
const std::vector< BoundaryName > & _wall_boundary_names
Wall boundaries.
void getWallBoundedElements(const std::vector< BoundaryName > &wall_boundary_name, const FEProblemBase &fe_problem, const SubProblem &subproblem, const std::set< SubdomainID > &block_ids, std::map< const Elem *, bool > &wall_bounded_map)
Map marking wall bounded elements The map passed in wall_bounded_map gets cleared and re-populated...
static InputParameters validParams()
Moose::StateArg determineState() const
const Moose::Functor< ADReal > & _k
Turbulent kinetic energy.
NS::WallTreatmentEnum _wall_treatment
Method used for wall treatment.
static const std::string density
static const std::string TKE
WallTreatmentEnum
Wall treatment options.
const Moose::Functor< ADReal > * _v_var
y-velocity
virtual const std::set< SubdomainID > & blockIDs() const
Computes the source and sink terms for the turbulent kinetic energy dissipation rate.
DualNumber< Real, DNDerivativeType, true > ADReal
ADReal computeQpResidual() override
Real distance(const Point &p)
const Real _C2_eps
Value of the second epsilon closure coefficient.
Moose::ElemArg makeElemArg(const Elem *elem, bool correct_skewnewss=false) const
const Real _C1_eps
Value of the first epsilon closure coefficient.
const Moose::Functor< ADReal > & _rho
Density.
registerMooseObject("NavierStokesApp", INSFVTKEDSourceSink)
const Elem *const & _current_elem
static const std::string mu
INSFVTKEDSourceSink(const InputParameters ¶meters)
static InputParameters validParams()
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...
std::map< const Elem *, std::vector< Real > > _dist
FEProblemBase & _fe_problem
std::map< const Elem *, std::vector< const FaceInfo * > > _face_infos
void paramError(const std::string ¶m, Args... args) const
void getWallDistance(const std::vector< BoundaryName > &wall_boundary_name, const FEProblemBase &fe_problem, const SubProblem &subproblem, const std::set< SubdomainID > &block_ids, std::map< const Elem *, std::vector< Real >> &dist_map)
Map storing wall ditance for near-wall marked elements The map passed in dist_map gets cleared and re...
const bool _newton_solve
Whether a nonlinear Newton-like solver is being used (as opposed to a linearized solver) ...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const unsigned int _dim
The dimension of the simulation.
const Real _C_mu
C_mu constant.
const bool _linearized_model
If the user wants to use the linearized model.
const Moose::Functor< ADReal > & _mu
Dynamic viscosity.
void getElementFaceArgs(const std::vector< BoundaryName > &wall_boundary_name, const FEProblemBase &fe_problem, const SubProblem &subproblem, const std::set< SubdomainID > &block_ids, std::map< const Elem *, std::vector< const FaceInfo *>> &face_info_map)
Map storing face arguments to wall bounded faces The map passed in face_info_map gets cleared and re-...
const Moose::Functor< ADReal > & _mu_t
Turbulent dynamic viscosity.
static const std::string velocity
ADReal computeSpeed(const ADRealVectorValue &velocity)
Compute the speed (velocity norm) given the supplied velocity.
virtual void initialSetup() override
MooseUnits pow(const MooseUnits &, int)
auto index_range(const T &sizable)
const Moose::Functor< ADReal > & _u_var
x-velocity
MooseVariableFV< Real > & _var
virtual bool hasFaceSide(const FaceInfo &fi, const bool fi_elem_side) const override