25#ifndef DUMUX_TURBULENCE_PROPERTIES_HH
26#define DUMUX_TURBULENCE_PROPERTIES_HH
30#include<dune/common/fvector.hh>
37template<
class Scalar,
unsigned dim,
bool verbose = false>
46 const Dune::FieldVector<Scalar, dim> position,
47 const Scalar kinematicViscosity,
50 bool print=verbose)
const
55 const Scalar re_x =
reynoldsNumber(velocity, position[0], kinematicViscosity,
false);
56 const Scalar c_f = pow((2.0 * log10(re_x) - 0.65), -2.3);
57 const Scalar wallShearStress = 0.5 * c_f *
density * velocity * velocity;
58 const Scalar frictionVelocity = sqrt(wallShearStress /
density);
59 const Scalar yPlus = position[yCoordDim] * frictionVelocity / kinematicViscosity;
62 std::cout <<
"turbulence properties at (";
63 for (
unsigned int dimIdx = 0; dimIdx < dim; ++dimIdx)
64 std::cout << position[dimIdx] <<
",";
65 std::cout <<
")" << std::endl;
66 std::cout <<
"estimated Re_x : " << re_x <<
" [-]" << std::endl;
67 std::cout <<
"estimated c_f : " << c_f <<
" [-]" << std::endl;
68 std::cout <<
"estimated tau_w : " << wallShearStress <<
" [kg/(m*s^2)]" << std::endl;
69 std::cout <<
"estimated UStar : " << frictionVelocity <<
" [m/s]" << std::endl;
70 std::cout <<
"estimated yPlus : " << yPlus <<
" [-]" << std::endl;
71 std::cout << std::endl;
81 const Scalar kinematicViscosity,
82 bool print=verbose)
const
89 std::cout <<
"estimated Re_d : " << re_d <<
" [-]" << std::endl;
90 std::cout <<
"estimated l_ent : " <<
entranceLength <<
" [m]"<< std::endl;
91 std::cout << std::endl;
100 const Scalar charLengthScale,
101 const Scalar kinematicViscosity,
102 bool print=verbose)
const
105 return abs(velocity * charLengthScale / kinematicViscosity);
113 bool print=verbose)
const
129 bool print=verbose)
const
145 const Scalar kinematicViscosity,
146 bool print=verbose)
const
149 const Scalar k = 1.5 * velocity * velocity
153 std::cout <<
"estimated k : " << k <<
" [m^2/s^2]" << std::endl;
164 const Scalar kinematicViscosity,
165 bool print=verbose)
const
169 const Scalar factor = 0.1643;
173 std::cout <<
"estimated eps. : " << epsilon <<
" [m^2/s^3]" << std::endl;
185 const Scalar kinematicViscosity,
186 bool print=verbose)
const
190 const Scalar factor = 0.54772;
192 const Scalar omega = pow(k, 0.5) / (L*factor);
195 std::cout <<
"estimated omega : " << omega <<
" [1/s]" << std::endl;
207 const Scalar kinematicViscosity,
208 bool print=verbose)
const
218 std::cout <<
"estimated nu~ : " <<
viscosityTilde <<
" [m^2/s]" << std::endl;
Geometry::ctype diameter(const Geometry &geo)
Computes the longest distance between points of a geometry.
Definition: diameter.hh:36
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
This class contains different functions for estimating turbulence properties.
Definition: turbulenceproperties.hh:39
Scalar viscosityTilde(const Scalar velocity, const Scalar diameter, const Scalar kinematicViscosity, bool print=verbose) const
Estimates the viscosity tilde based on a formula given in in the ANSYS Fluent user guide .
Definition: turbulenceproperties.hh:205
Scalar yPlusEstimation(const Scalar velocity, const Dune::FieldVector< Scalar, dim > position, const Scalar kinematicViscosity, const Scalar density, int yCoordDim=dim-1, bool print=verbose) const
Estimates dimensionless wall distance based on a formula given in http://www.cfd-online....
Definition: turbulenceproperties.hh:45
Scalar entranceLength(const Scalar velocity, const Scalar diameter, const Scalar kinematicViscosity, bool print=verbose) const
Estimates the entrance length for this pipe.
Definition: turbulenceproperties.hh:79
Scalar turbulenceLengthScale(const Scalar charLengthScale, bool print=verbose) const
Estimates the turbulence length scale based on a formula given in the ANSYS Fluent user guide .
Definition: turbulenceproperties.hh:128
Scalar turbulentKineticEnergy(const Scalar velocity, const Scalar diameter, const Scalar kinematicViscosity, bool print=verbose) const
Estimates the turbulent kinetic energy based on a formula given in the ANSYS Fluent user guide .
Definition: turbulenceproperties.hh:143
Scalar dissipationRate(const Scalar velocity, const Scalar diameter, const Scalar kinematicViscosity, bool print=verbose) const
Estimates the dissipation rate based on a formula given in the ANSYS Fluent user guide .
Definition: turbulenceproperties.hh:183
Scalar turbulenceIntensity(const Scalar reynoldsNumber, bool print=verbose) const
Estimates the turbulence intensity based on a formula given in the ANSYS Fluent user guide .
Definition: turbulenceproperties.hh:112
Scalar dissipation(const Scalar velocity, const Scalar diameter, const Scalar kinematicViscosity, bool print=verbose) const
Estimates the dissipation based on a formula given in the ANSYS Fluent user guide .
Definition: turbulenceproperties.hh:162
Scalar reynoldsNumber(const Scalar velocity, const Scalar charLengthScale, const Scalar kinematicViscosity, bool print=verbose) const
Calculates the Reynolds number.
Definition: turbulenceproperties.hh:99