https://mooseframework.inl.gov
LinearFVTKEDSourceSink.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
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 
10 #include "LinearFVTKEDSourceSink.h"
11 #include "Assembly.h"
12 #include "SubProblem.h"
13 #include "NavierStokesMethods.h"
14 
16 
19 {
21  params.addClassDescription("Elemental kernel to compute the production and destruction "
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>(NS::TKE, "Coupled turbulent kinetic energy.");
27  params.addRequiredParam<MooseFunctorName>(NS::density, "fluid density");
28  params.addRequiredParam<MooseFunctorName>(NS::mu, "Dynamic viscosity.");
29  params.addRequiredParam<MooseFunctorName>(NS::mu_t, "Turbulent viscosity.");
30 
31  params.addParam<std::vector<BoundaryName>>(
32  "walls", {}, "Boundaries that correspond to solid walls.");
33  params.addParam<bool>(
34  "linearized_model",
35  true,
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");
38  params.addParam<MooseEnum>("wall_treatment",
39  wall_treatment,
40  "The method used for computing the wall functions "
41  "'eq_newton', 'eq_incremental', 'eq_linearized', 'neq'");
42 
43  params.addParam<Real>("C1_eps", 1.44, "First epsilon coefficient");
44  params.addParam<Real>("C2_eps", 1.92, "Second epsilon coefficient");
45  params.addParam<Real>("C_mu", 0.09, "Coupled turbulent kinetic energy closure.");
46  params.addParam<Real>("C_pl", 10.0, "Production limiter constant multiplier.");
47 
48  return params;
49 }
50 
52  : LinearFVElementalKernel(params),
53  _dim(_subproblem.mesh().dimension()),
54  _u_var(getFunctor<Real>("u")),
55  _v_var(params.isParamValid("v") ? &(getFunctor<Real>("v")) : nullptr),
56  _w_var(params.isParamValid("w") ? &(getFunctor<Real>("w")) : nullptr),
57  _k(getFunctor<Real>(NS::TKE)),
58  _rho(getFunctor<Real>(NS::density)),
59  _mu(getFunctor<Real>(NS::mu)),
60  _mu_t(getFunctor<Real>(NS::mu_t)),
61  _wall_boundary_names(getParam<std::vector<BoundaryName>>("walls")),
62  _linearized_model(getParam<bool>("linearized_model")),
63  _wall_treatment(getParam<MooseEnum>("wall_treatment").getEnum<NS::WallTreatmentEnum>()),
64  _C1_eps(getParam<Real>("C1_eps")),
65  _C2_eps(getParam<Real>("C2_eps")),
66  _C_mu(getParam<Real>("C_mu")),
67  _C_pl(getParam<Real>("C_pl"))
68 {
69  if (_dim >= 2 && !_v_var)
70  paramError("v", "In two or more dimensions, the v velocity must be supplied!");
71 
72  if (_dim >= 3 && !_w_var)
73  paramError("w", "In three or more dimensions, the w velocity must be supplied!");
74 
75  // Strain tensor term requires velocity gradients;
76  if (dynamic_cast<const MooseLinearVariableFV<Real> *>(&_u_var))
77  requestVariableCellGradient(getParam<MooseFunctorName>("u"));
78  if (dynamic_cast<const MooseLinearVariableFV<Real> *>(_v_var))
79  requestVariableCellGradient(getParam<MooseFunctorName>("v"));
80  if (dynamic_cast<const MooseLinearVariableFV<Real> *>(_w_var))
81  requestVariableCellGradient(getParam<MooseFunctorName>("w"));
82 }
83 
84 void
86 {
92 }
93 
94 Real
96 {
97  if (_wall_bounded.find(_current_elem_info->elem()) != _wall_bounded.end())
98  // TKED value for near wall element will be directly assigned for this cell
99  return 1.0;
100  else
101  {
102  // Convenient definitions
103  const auto state = determineState();
104  const auto elem_arg = makeElemArg(_current_elem_info->elem());
105  const Real rho = _rho(elem_arg, state);
106  const Real TKE = _k(elem_arg, state);
107  const auto epsilon = _var.getElemValue(*_current_elem_info, state);
108 
109  // Compute destruction
110  const auto destruction = _C2_eps * rho * epsilon / TKE;
111 
112  // Assign to matrix (term gets multiplied by TKED)
113  return destruction * _current_elem_volume;
114  }
115 }
116 
117 Real
119 {
120  if (_wall_bounded.find(_current_elem_info->elem()) != _wall_bounded.end())
121  {
122  // Convenient definitions
123  const auto state = determineState();
124  const auto elem_arg = makeElemArg(_current_elem_info->elem());
125  const Real rho = _rho(elem_arg, state);
126  const Real mu = _mu(elem_arg, state);
127  const Real TKE = _k(elem_arg, state);
128 
129  // Convenient variables
130  Real destruction = 0.0;
131  std::vector<Real> y_plus_vec, velocity_grad_norm_vec;
132  Real tot_weight = 0.0;
133 
134  // Get velocity vector
135  RealVectorValue velocity(_u_var(elem_arg, state));
136  if (_v_var)
137  velocity(1) = (*_v_var)(elem_arg, state);
138  if (_w_var)
139  velocity(2) = (*_w_var)(elem_arg, state);
140 
141  // Get near wall faceInfo and distances from cell center to every wall
142  const auto & face_info_vec = libmesh_map_find(_face_infos, _current_elem_info->elem());
143  const auto & distance_vec = libmesh_map_find(_dist, _current_elem_info->elem());
144  mooseAssert(distance_vec.size(), "Should have found a distance vector");
145  mooseAssert(distance_vec.size() == face_info_vec.size(),
146  "Should be as many distance vectors as face info vectors");
147 
148  // Update y+ and wall face cell
149  for (unsigned int i = 0; i < distance_vec.size(); i++)
150  {
151  const auto distance = distance_vec[i];
152  mooseAssert(distance > 0, "Should be at a non-zero distance");
153 
154  Real y_plus;
155  if (_wall_treatment == NS::WallTreatmentEnum::NEQ) // Non-equilibrium / Non-iterative
156  y_plus = distance * std::sqrt(std::sqrt(_C_mu) * TKE) * rho / mu;
157  else // Equilibrium / Iterative
158  {
159  const auto parallel_speed = NS::computeSpeed<Real>(
160  velocity - velocity * face_info_vec[i]->normal() * face_info_vec[i]->normal());
161  y_plus = NS::findyPlus<Real>(mu, rho, std::max(parallel_speed, 1e-10), distance);
162  }
163 
164  y_plus_vec.push_back(y_plus);
165  tot_weight += 1.0;
166  }
167 
168  // Compute near wall epsilon value
169  for (const auto i : index_range(y_plus_vec))
170  {
171  const auto y_plus = y_plus_vec[i];
172 
173  if (y_plus < 11.25)
174  destruction += 2.0 * TKE * mu / rho / Utility::pow<2>(distance_vec[i]) / tot_weight;
175  else
176  destruction += std::pow(_C_mu, 0.75) * std::pow(TKE, 1.5) /
177  (NS::von_karman_constant * distance_vec[i]) / tot_weight;
178  }
179 
180  // Assign the computed value of TKED for element near the wall
181  return destruction;
182  }
183  else
184  {
185  // Convenient definitions
186  const auto state = determineState();
187  const auto elem_arg = makeElemArg(_current_elem_info->elem());
188  const Real rho = _rho(elem_arg, state);
189  const Real TKE = _k(elem_arg, state);
190  const Real TKED = _var.getElemValue(*_current_elem_info, state);
191 
192  // Compute production of TKE
193  const auto symmetric_strain_tensor_sq_norm =
195  Real production = _mu_t(elem_arg, state) * symmetric_strain_tensor_sq_norm;
196 
197  // Limit TKE production (needed for flows with stagnation zones)
198  const Real production_limit = _C_pl * rho * TKED;
199  production = std::min(production, production_limit);
200 
201  // Compute production - recasted with mu_t definition to avoid division by epsilon
202  const auto production_epsilon = _C1_eps * production * TKED / TKE;
203 
204  // Assign to matrix (term gets multiplied by TKED)
205  return production_epsilon * _current_elem_volume;
206  }
207 }
const ElemInfo * _current_elem_info
virtual Real computeMatrixContribution() override
virtual Real computeRightHandSideContribution() override
static constexpr Real von_karman_constant
Definition: NS.h:191
static const std::string mu_t
Definition: NS.h:125
SubProblem & _subproblem
const Moose::Functor< Real > * _v_var
y-velocity
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
void requestVariableCellGradient(const std::string &variable_name)
registerMooseObject("NavierStokesApp", LinearFVTKEDSourceSink)
Moose::StateArg determineState() const
const Elem * elem() const
const Moose::Functor< Real > & _u_var
x-velocity
static InputParameters validParams()
std::map< const Elem *, std::vector< Real > > _dist
MooseLinearVariableFV< Real > & _var
MeshBase & mesh
static const std::string density
Definition: NS.h:33
static const std::string TKE
Definition: NS.h:176
WallTreatmentEnum
Wall treatment options.
Definition: NS.h:182
NS::WallTreatmentEnum _wall_treatment
Method used for wall treatment.
void getWallBoundedElements(const std::vector< BoundaryName > &wall_boundary_name, const FEProblemBase &fe_problem, const SubProblem &subproblem, const std::set< SubdomainID > &block_ids, std::unordered_set< const Elem *> &wall_bounded)
Map marking wall bounded elements The map passed in wall_bounded_map gets cleared and re-populated...
const std::vector< BoundaryName > & _wall_boundary_names
Wall boundaries.
virtual const std::set< SubdomainID > & blockIDs() const
const Real _C2_eps
Value of the second epsilon closure coefficient.
const Moose::Functor< Real > & _rho
Density.
Real distance(const Point &p)
void addRequiredParam(const std::string &name, const std::string &doc_string)
const Real _C_pl
Production Limiter Constant.
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)
LinearFVTKEDSourceSink(const InputParameters &params)
Class constructor.
const unsigned int _dim
The dimension of the simulation.
static const std::string mu
Definition: NS.h:123
static InputParameters validParams()
std::map< const Elem *, std::vector< const FaceInfo * > > _face_infos
Kernel that adds contributions to the source and the sink of the turbulent kinetic energy dissipation...
virtual void initialSetup() override
void paramError(const std::string &param, Args... args) const
const Real _C1_eps
Value of the first epsilon closure coefficient.
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...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Real getElemValue(const ElemInfo &elem_info, const StateArg &state) const
const Moose::Functor< Real > & _mu_t
Turbulent dynamic viscosity.
static const std::string TKED
Definition: NS.h:177
void addClassDescription(const std::string &doc_string)
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-...
static const std::string velocity
Definition: NS.h:45
FEProblemBase & _fe_problem
const Moose::Functor< Real > & _mu
Dynamic viscosity.
template Real computeSpeed< Real >(const libMesh::VectorValue< Real > &velocity)
const Moose::Functor< Real > & _k
Turbulent kinetic energy.
virtual void initialSetup()
const Real _C_mu
C_mu constant.
MooseUnits pow(const MooseUnits &, int)
std::unordered_set< const Elem * > _wall_bounded
Maps for wall treatment.
auto index_range(const T &sizable)
const Moose::Functor< Real > * _w_var
z-velocity
template Real computeShearStrainRateNormSquared< Real >(const Moose::Functor< Real > &u, const Moose::Functor< Real > *v, const Moose::Functor< Real > *w, const Moose::ElemArg &elem_arg, const Moose::StateArg &state)