version 3.8
turbulenceproperties.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_TURBULENCE_PROPERTIES_HH
14#define DUMUX_TURBULENCE_PROPERTIES_HH
15
16#include <iostream>
17
18#include<dune/common/fvector.hh>
19
20namespace Dumux {
21
25template<class Scalar, unsigned dim, bool verbose = false>
27{
28public:
33 Scalar yPlusEstimation(const Scalar velocity,
34 const Dune::FieldVector<Scalar, dim> position,
35 const Scalar kinematicViscosity,
36 const Scalar density,
37 int yCoordDim=dim-1,
38 bool print=verbose) const
39 {
40 using std::pow;
41 using std::log10;
42 using std::sqrt;
43 const Scalar re_x = reynoldsNumber(velocity, position[0], kinematicViscosity, false);
44 const Scalar c_f = pow((2.0 * log10(re_x) - 0.65), -2.3); // for re_x < 10^9
45 const Scalar wallShearStress = 0.5 * c_f * density * velocity * velocity;
46 const Scalar frictionVelocity = sqrt(wallShearStress / density);
47 const Scalar yPlus = position[yCoordDim] * frictionVelocity / kinematicViscosity;
48 if (print)
49 {
50 std::cout << "turbulence properties at (";
51 for (unsigned int dimIdx = 0; dimIdx < dim; ++dimIdx)
52 std::cout << position[dimIdx] << ",";
53 std::cout << ")" << std::endl;
54 std::cout << "estimated Re_x : " << re_x << " [-]" << std::endl;
55 std::cout << "estimated c_f : " << c_f << " [-]" << std::endl;
56 std::cout << "estimated tau_w : " << wallShearStress << " [kg/(m*s^2)]" << std::endl;
57 std::cout << "estimated UStar : " << frictionVelocity << " [m/s]" << std::endl;
58 std::cout << "estimated yPlus : " << yPlus << " [-]" << std::endl;
59 std::cout << std::endl;
60 }
61 return yPlus;
62 }
63
67 Scalar entranceLength(const Scalar velocity,
68 const Scalar diameter,
69 const Scalar kinematicViscosity,
70 bool print=verbose) const
71 {
72 using std::pow;
73 const Scalar re_d = reynoldsNumber(velocity, diameter, kinematicViscosity, false);
74 const Scalar entranceLength = 4.4 * pow(re_d, 1.0/6.0) * diameter;
75 if (print)
76 {
77 std::cout << "estimated Re_d : " << re_d << " [-]" << std::endl;
78 std::cout << "estimated l_ent : " << entranceLength << " [m]"<< std::endl;
79 std::cout << std::endl;
80 }
81 return entranceLength;
82 }
83
87 Scalar reynoldsNumber(const Scalar velocity,
88 const Scalar charLengthScale/*e.g. diameter*/,
89 const Scalar kinematicViscosity,
90 bool print=verbose) const
91 {
92 using std::abs;
93 return abs(velocity * charLengthScale / kinematicViscosity);
94 }
95
101 bool print=verbose) const
102 {
103 using std::pow;
104 const Scalar turbulenceIntensity = 0.16 * pow(reynoldsNumber, -0.125);
105 if (print)
106 {
107 std::cout << "estimated I : " << turbulenceIntensity << " [-]" << std::endl;
108 }
109 return turbulenceIntensity;
110 }
111
116 Scalar turbulenceLengthScale(const Scalar charLengthScale/*e.g. diameter*/,
117 bool print=verbose) const
118 {
119 const Scalar turbulenceLengthScale = 0.07 * charLengthScale;
120 if (print)
121 {
122 std::cout << "estimated l_turb: " << turbulenceLengthScale << " [m]" << std::endl;
123 }
125 }
126
131 Scalar turbulentKineticEnergy(const Scalar velocity,
132 const Scalar diameter,
133 const Scalar kinematicViscosity,
134 bool print=verbose) const
135 {
136 const Scalar re_d = reynoldsNumber(velocity, diameter, kinematicViscosity, false);
137 const Scalar k = 1.5 * velocity * velocity
138 * turbulenceIntensity(re_d, false) * turbulenceIntensity(re_d, false);
139 if (print)
140 {
141 std::cout << "estimated k : " << k << " [m^2/s^2]" << std::endl;
142 }
143 return k;
144 }
145
150 Scalar dissipation(const Scalar velocity,
151 const Scalar diameter,
152 const Scalar kinematicViscosity,
153 bool print=verbose) const
154 {
155 using std::pow;
156 const Scalar k = turbulentKineticEnergy(velocity, diameter, kinematicViscosity, false);
157 const Scalar factor = 0.1643; // = cMu^(3/4) = 0.09^(3/4)
158 const Scalar epsilon = factor * pow(k, 1.5) / turbulenceLengthScale(diameter, false);
159 if (print)
160 {
161 std::cout << "estimated eps. : " << epsilon << " [m^2/s^3]" << std::endl;
162 }
163 return epsilon;
164 }
165
171 Scalar dissipationRate(const Scalar velocity,
172 const Scalar diameter,
173 const Scalar kinematicViscosity,
174 bool print=verbose) const
175 {
176 using std::pow;
177 const Scalar k = turbulentKineticEnergy(velocity, diameter, kinematicViscosity, false);
178 const Scalar factor = 0.54772; // = cMu^(1/4) = 0.09^(1/4)
179 const Scalar L = turbulenceLengthScale(diameter , false);
180 const Scalar omega = pow(k, 0.5) / (L*factor);
181 if (print)
182 {
183 std::cout << "estimated omega : " << omega << " [1/s]" << std::endl;
184
185 }
186 return omega;
187 }
188
193 Scalar viscosityTilde(const Scalar velocity,
194 const Scalar diameter,
195 const Scalar kinematicViscosity,
196 bool print=verbose) const
197 {
198 using std::abs;
199 using std::sqrt;
200 const Scalar re_d = reynoldsNumber(velocity, diameter, kinematicViscosity, false);
201 const Scalar viscosityTilde = sqrt(1.5) * abs(velocity)
202 * turbulenceIntensity(re_d)
204 if (print)
205 {
206 std::cout << "estimated nu~ : " << viscosityTilde << " [m^2/s]" << std::endl;
207 }
208 return viscosityTilde;
209 }
210};
211} // end namespace Dumux
212
213#endif
This class contains different functions for estimating turbulence properties.
Definition: turbulenceproperties.hh:27
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:193
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:33
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:67
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:116
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:131
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:171
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:100
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:150
Scalar reynoldsNumber(const Scalar velocity, const Scalar charLengthScale, const Scalar kinematicViscosity, bool print=verbose) const
Calculates the Reynolds number.
Definition: turbulenceproperties.hh:87
Geometry::ctype diameter(const Geometry &geo)
Computes the longest distance between points of a geometry.
Definition: diameter.hh:26
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:53
void print(const Collector &collector)
prints json tree
Definition: metadata.hh:237
Definition: adapt.hh:17