https://mooseframework.inl.gov
NavierStokesMethods.h
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 
10 #pragma once
11 
12 #include <vector>
13 #include "Moose.h"
14 #include "MooseUtils.h"
15 #include "ADReal.h"
16 #include "metaphysicl/raw_type.h"
17 #include "FEProblemBase.h"
18 #include "SubProblem.h"
19 
20 namespace NS
21 {
28 int delta(unsigned int i, unsigned int j);
29 
35 int computeSign(const Real & a);
36 
44 unsigned int getIndex(const Real & p, const std::vector<Real> & bounds);
45 
67  const Real & Re, const Real & rho, const Real & mu, const Real & drho, const Real & dmu);
68 
86  const Real & cp,
87  const Real & k,
88  const Real & dmu,
89  const Real & dcp,
90  const Real & dk);
91 
102 ADReal findUStar(const ADReal & mu, const ADReal & rho, const ADReal & u, Real dist);
103 
113 ADReal findyPlus(const ADReal & mu, const ADReal & rho, const ADReal & u, Real dist);
114 
115 using MooseUtils::isZero;
116 
121 
126 void getWallBoundedElements(const std::vector<BoundaryName> & wall_boundary_name,
127  const FEProblemBase & fe_problem,
128  const SubProblem & subproblem,
129  const std::set<SubdomainID> & block_ids,
130  std::map<const Elem *, bool> & wall_bounded_map);
131 
136 void getWallDistance(const std::vector<BoundaryName> & wall_boundary_name,
137  const FEProblemBase & fe_problem,
138  const SubProblem & subproblem,
139  const std::set<SubdomainID> & block_ids,
140  std::map<const Elem *, std::vector<Real>> & dist_map);
141 
146 void getElementFaceArgs(const std::vector<BoundaryName> & wall_boundary_name,
147  const FEProblemBase & fe_problem,
148  const SubProblem & subproblem,
149  const std::set<SubdomainID> & block_ids,
150  std::map<const Elem *, std::vector<const FaceInfo *>> & face_info_map);
151 
155 template <typename T, typename VectorType, typename PointType>
156 T
157 divergence(const TensorValue<T> & gradient,
158  const VectorType & value,
159  const PointType & point,
160  const Moose::CoordinateSystemType & coord_sys,
161  const unsigned int rz_radial_coord)
162 {
163  mooseAssert((coord_sys == Moose::COORD_XYZ) || (coord_sys == Moose::COORD_RZ),
164  "This function only supports calculations of divergence in Cartesian and "
165  "axisymmetric coordinate systems");
166  auto div = gradient.tr();
167  if (coord_sys == Moose::COORD_RZ)
168  // u_r / r
169  div += value(rz_radial_coord) / point(rz_radial_coord);
170  return div;
171 }
172 
181 template <typename T1, typename T2, typename T3>
182 auto
183 wallHeatTransferCoefficient(const T1 & Nu, const T2 & k, const T3 & D_h)
184 {
185  return Nu * k / D_h;
186 }
187 }
int computeSign(const Real &a)
Sign function, returns $+1$ if $a$ is positive and $-1$ if $a$ is negative.
T divergence(const TensorValue< T > &gradient, const VectorType &value, const PointType &point, const Moose::CoordinateSystemType &coord_sys, const unsigned int rz_radial_coord)
Compute the divergence of a vector given its matrix of derivatives.
Real prandtlPropertyDerivative(const Real &mu, const Real &cp, const Real &k, const Real &dmu, const Real &dcp, const Real &dk)
Computes the derivative of the Prandtl number, $Pr{ C_p}{k}$, with respect to an arbitrary variale $$...
bool isZero(const T &value, const Real tolerance=TOLERANCE *TOLERANCE *TOLERANCE)
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...
int delta(unsigned int i, unsigned int j)
Delta function, which returns zero if $i j$ and unity if $i=j$.
DualNumber< Real, DNDerivativeType, true > ADReal
static const std::string cp
Definition: NS.h:121
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)
static const std::string mu
Definition: NS.h:123
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...
Real reynoldsPropertyDerivative(const Real &Re, const Real &rho, const Real &mu, const Real &drho, const Real &dmu)
Computes the derivative of the Reynolds number, $Re { Vd}{}$, with respect to an arbitrary variable $...
auto wallHeatTransferCoefficient(const T1 &Nu, const T2 &k, const T3 &D_h)
Compute wall heat transfer coefficient.
unsigned int getIndex(const Real &p, const std::vector< Real > &bounds)
Determines the index $i$ in a sorted array such that the input point is within the $i$-th and $i+1$-t...
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
CoordinateSystemType
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
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
ADReal findUStar(const ADReal &mu, const ADReal &rho, const ADReal &u, Real dist)
Finds the friction velocity using standard velocity wall functions formulation.
ADReal computeSpeed(const ADRealVectorValue &velocity)
Compute the speed (velocity norm) given the supplied velocity.
static const std::string k
Definition: NS.h:130