https://mooseframework.inl.gov
LinearFVTKESourceSink.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 "LinearFVTKESourceSink.h"
11 #include "Assembly.h"
12 #include "SubProblem.h"
13 #include "NavierStokesMethods.h"
14 
15 registerMooseObject("NavierStokesApp", LinearFVTKESourceSink);
16 
19 {
21  params.addClassDescription("Elemental kernel to compute the production and destruction "
22  " terms of turbulent kinetic energy (TKE).");
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::TKED,
27  "Coupled turbulent kinetic energy dissipation rate.");
28  params.addRequiredParam<MooseFunctorName>(NS::density, "Fluid density");
29  params.addRequiredParam<MooseFunctorName>(NS::mu, "Dynamic viscosity.");
30  params.addRequiredParam<MooseFunctorName>(NS::mu_t, "Turbulent viscosity.");
31 
32  params.addParam<std::vector<BoundaryName>>(
33  "walls", {}, "Boundaries that correspond to solid walls.");
34  params.addParam<bool>(
35  "linearized_model",
36  true,
37  "Boolean to determine if the problem should be use in a linear or nonlinear solve.");
38  MooseEnum wall_treatment("eq_newton eq_incremental eq_linearized neq", "neq");
39 
40  params.addParam<MooseEnum>("wall_treatment",
41  wall_treatment,
42  "The method used for computing the wall functions "
43  "'eq_newton', 'eq_incremental', 'eq_linearized', 'neq'");
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 
47  return params;
48 }
49 
51  : LinearFVElementalKernel(params),
52  _dim(_subproblem.mesh().dimension()),
53  _u_var(getFunctor<Real>("u")),
54  _v_var(params.isParamValid("v") ? &(getFunctor<Real>("v")) : nullptr),
55  _w_var(params.isParamValid("w") ? &(getFunctor<Real>("w")) : nullptr),
56  _epsilon(getFunctor<Real>(NS::TKED)),
57  _rho(getFunctor<Real>(NS::density)),
58  _mu(getFunctor<Real>(NS::mu)),
59  _mu_t(getFunctor<Real>(NS::mu_t)),
60  _wall_boundary_names(getParam<std::vector<BoundaryName>>("walls")),
61  _linearized_model(getParam<bool>("linearized_model")),
62  _wall_treatment(getParam<MooseEnum>("wall_treatment").getEnum<NS::WallTreatmentEnum>()),
63  _C_mu(getParam<Real>("C_mu")),
64  _C_pl(getParam<Real>("C_pl"))
65 {
66  if (_dim >= 2 && !_v_var)
67  paramError("v", "In two or more dimensions, the v velocity must be supplied!");
68 
69  if (_dim >= 3 && !_w_var)
70  paramError("w", "In three or more dimensions, the w velocity must be supplied!");
71 
72  // // Strain tensor term requires velocity gradients;
73  if (dynamic_cast<const MooseLinearVariableFV<Real> *>(&_u_var))
74  requestVariableCellGradient(getParam<MooseFunctorName>("u"));
75  if (dynamic_cast<const MooseLinearVariableFV<Real> *>(_v_var))
76  requestVariableCellGradient(getParam<MooseFunctorName>("v"));
77  if (dynamic_cast<const MooseLinearVariableFV<Real> *>(_w_var))
78  requestVariableCellGradient(getParam<MooseFunctorName>("w"));
79 }
80 
81 void
83 {
89 }
90 
91 Real
93 {
94  /*
95  Matrix contribution:
96  - Computes near-wall TKE destruction
97  - Computes bulk TKE destruction
98  */
99 
100  // Useful variables
101  const auto state = determineState();
102  const auto elem_arg = makeElemArg(_current_elem_info->elem());
103  const Real rho = _rho(elem_arg, state);
104 
105  if (_wall_bounded.find(_current_elem_info->elem()) != _wall_bounded.end())
106  {
107 
108  // Useful variables
109  const Real mu = _mu(elem_arg, state);
110  const Real TKE = _var.getElemValue(*_current_elem_info, state);
111 
112  // Assembly variables
113  Real production = 0.0;
114  Real destruction = 0.0;
115  Real tot_weight = 0.0;
116  std::vector<Real> y_plus_vec, velocity_grad_norm_vec;
117 
118  // Get near-wall element velocity vector
119  RealVectorValue velocity(_u_var(elem_arg, state));
120  if (_v_var)
121  velocity(1) = (*_v_var)(elem_arg, state);
122  if (_w_var)
123  velocity(2) = (*_w_var)(elem_arg, state);
124 
125  // Get faceInfo and distance vector for between all walls and elements
126  const auto & face_info_vec = libmesh_map_find(_face_infos, _current_elem_info->elem());
127  const auto & distance_vec = libmesh_map_find(_dist, _current_elem_info->elem());
128  mooseAssert(distance_vec.size(), "Should have found a distance vector");
129  mooseAssert(distance_vec.size() == face_info_vec.size(),
130  "Should be as many distance vectors as face info vectors");
131 
132  // Populate y+, near wall velocity gradient, and update weights
133  for (unsigned int i = 0; i < distance_vec.size(); i++)
134  {
135  const auto parallel_speed = NS::computeSpeed<Real>(
136  velocity - velocity * face_info_vec[i]->normal() * face_info_vec[i]->normal());
137  const auto distance = distance_vec[i];
138 
139  Real y_plus;
140  if (_wall_treatment == NS::WallTreatmentEnum::NEQ) // Non-equilibrium --> Non-iterative
141  y_plus = distance * std::sqrt(std::sqrt(_C_mu) * TKE) * rho / mu;
142  else // Equilibrium --> Iterative
143  y_plus = NS::findyPlus<Real>(mu, rho, std::max(parallel_speed, 1e-10), distance);
144 
145  y_plus_vec.push_back(y_plus);
146  const Real velocity_grad_norm = parallel_speed / distance;
147  velocity_grad_norm_vec.push_back(velocity_grad_norm);
148  tot_weight += 1.0;
149  }
150 
151  // Compute near-wall production and destruction
152  for (unsigned int i = 0; i < y_plus_vec.size(); i++)
153  {
154  const auto y_plus = y_plus_vec[i];
155 
156  const auto fi = face_info_vec[i];
157  const bool defined_on_elem_side = _var.hasFaceSide(*fi, true);
158  const Elem * const loc_elem = defined_on_elem_side ? &fi->elem() : fi->neighborPtr();
159  const Moose::FaceArg facearg = {
160  fi, Moose::FV::LimiterType::CentralDifference, false, false, loc_elem, nullptr};
161 
162  const Real wall_mut = _mu_t(facearg, state);
163  const Real wall_mu = _mu(facearg, state);
164  const auto tau_w = (wall_mut + wall_mu) * velocity_grad_norm_vec[i];
165  const auto destruction_visc = 2.0 * wall_mu / Utility::pow<2>(distance_vec[i]) / tot_weight;
166  const auto destruction_log = std::pow(_C_mu, 0.75) * rho * std::pow(TKE, 0.5) /
167  (NS::von_karman_constant * distance_vec[i]) / tot_weight;
168 
169  if (y_plus < 11.25)
170  destruction += destruction_visc;
171  else
172  {
173  destruction += destruction_log;
174  production += tau_w * std::pow(_C_mu, 0.25) / std::sqrt(TKE) /
175  (NS::von_karman_constant * distance_vec[i]) / tot_weight;
176  }
177  }
178 
179  // Assign terms to matrix to solve implicitly (they get multiplied by TKE)
180  return (destruction - production) * _current_elem_volume;
181  }
182  else
183  {
184  // Compute destruction
185  const auto destruction =
186  rho * _epsilon(elem_arg, state) / _var.getElemValue(*_current_elem_info, state);
187 
188  // Assign to matrix to solve implicitly
189  return destruction * _current_elem_volume;
190  }
191 }
192 
193 Real
195 {
196 
197  // Useful variables
198  const auto state = determineState();
199  const auto elem_arg = makeElemArg(_current_elem_info->elem());
200 
201  if (_wall_bounded.find(_current_elem_info->elem()) != _wall_bounded.end()) // Wall bounded
202  return 0.0; // Do nothing
203  else // Not wall bounded
204  {
205  // Compute TKE production
206  const auto symmetric_strain_tensor_sq_norm =
208 
209  auto production = _mu_t(elem_arg, state) * symmetric_strain_tensor_sq_norm;
210 
211  // k-Production limiter (needed for flows with stagnation zones)
212  const Real production_limit = _C_pl * _rho(elem_arg, state) * _epsilon(elem_arg, state);
213 
214  // Apply production limiter
215  production = std::min(production, production_limit);
216 
217  // Assign production to RHS
218  return production * _current_elem_volume;
219  }
220 }
const ElemInfo * _current_elem_info
const Moose::Functor< Real > * _w_var
z-velocity
std::map< const Elem *, std::vector< Real > > _dist
static constexpr Real von_karman_constant
Definition: NS.h:191
static const std::string mu_t
Definition: NS.h:125
SubProblem & _subproblem
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)
const Real _C_pl
Production Limiter Constant.
Moose::StateArg determineState() const
const Elem * elem() const
const Moose::Functor< Real > & _u_var
x-velocity
static InputParameters validParams()
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
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...
virtual const std::set< SubdomainID > & blockIDs() const
const NS::WallTreatmentEnum _wall_treatment
Method used for wall treatment.
const Moose::Functor< Real > & _rho
Density.
Real distance(const Point &p)
void addRequiredParam(const std::string &name, const std::string &doc_string)
const std::vector< BoundaryName > & _wall_boundary_names
Wall boundaries.
registerMooseObject("NavierStokesApp", LinearFVTKESourceSink)
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)
Kernel that adds contributions to the source and the sink of the turbulent kinetic energy discretized...
const Real _C_mu
C_mu constant.
const Moose::Functor< Real > & _epsilon
epsilon - dissipation rate of TKE
std::unordered_set< const Elem * > _wall_bounded
const unsigned int _dim
The dimension of the domain.
static const std::string mu
Definition: NS.h:123
const Moose::Functor< Real > & _mu
Dynamic viscosity.
void paramError(const std::string &param, Args... args) const
static InputParameters validParams()
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...
virtual Real computeMatrixContribution() override
LinearFVTKESourceSink(const InputParameters &params)
Class constructor.
virtual Real computeRightHandSideContribution() override
virtual void initialSetup() override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Real getElemValue(const ElemInfo &elem_info, const StateArg &state) const
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
template Real computeSpeed< Real >(const libMesh::VectorValue< Real > &velocity)
virtual void initialSetup()
std::map< const Elem *, std::vector< const FaceInfo * > > _face_infos
MooseUnits pow(const MooseUnits &, int)
const Moose::Functor< Real > & _mu_t
Turbulent dynamic viscosity.
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)
virtual bool hasFaceSide(const FaceInfo &fi, const bool fi_elem_side) const override
const Moose::Functor< Real > * _v_var
y-velocity