3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_TURBULENCE_PROPERTIES_HH
26#define DUMUX_TURBULENCE_PROPERTIES_HH
27
28#include <iostream>
29
30#include<dune/common/fvector.hh>
31
32namespace Dumux {
33
37template<class Scalar, unsigned dim, bool verbose = false>
39{
40public:
45 Scalar yPlusEstimation(const Scalar velocity,
46 const Dune::FieldVector<Scalar, dim> position,
47 const Scalar kinematicViscosity,
48 const Scalar density,
49 int yCoordDim=dim-1,
50 bool print=verbose) const
51 {
52 using std::pow;
53 using std::log10;
54 using std::sqrt;
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); // for re_x < 10^9
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;
60 if (print)
61 {
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;
72 }
73 return yPlus;
74 }
75
79 Scalar entranceLength(const Scalar velocity,
80 const Scalar diameter,
81 const Scalar kinematicViscosity,
82 bool print=verbose) const
83 {
84 using std::pow;
85 const Scalar re_d = reynoldsNumber(velocity, diameter, kinematicViscosity, false);
86 const Scalar entranceLength = 4.4 * pow(re_d, 1.0/6.0) * diameter;
87 if (print)
88 {
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;
92 }
93 return entranceLength;
94 }
95
99 Scalar reynoldsNumber(const Scalar velocity,
100 const Scalar charLengthScale/*e.g. diameter*/,
101 const Scalar kinematicViscosity,
102 bool print=verbose) const
103 {
104 using std::abs;
105 return abs(velocity * charLengthScale / kinematicViscosity);
106 }
107
113 bool print=verbose) const
114 {
115 using std::pow;
116 const Scalar turbulenceIntensity = 0.16 * pow(reynoldsNumber, -0.125);
117 if (print)
118 {
119 std::cout << "estimated I : " << turbulenceIntensity << " [-]" << std::endl;
120 }
121 return turbulenceIntensity;
122 }
123
128 Scalar turbulenceLengthScale(const Scalar charLengthScale/*e.g. diameter*/,
129 bool print=verbose) const
130 {
131 const Scalar turbulenceLengthScale = 0.07 * charLengthScale;
132 if (print)
133 {
134 std::cout << "estimated l_turb: " << turbulenceLengthScale << " [m]" << std::endl;
135 }
137 }
138
143 Scalar turbulentKineticEnergy(const Scalar velocity,
144 const Scalar diameter,
145 const Scalar kinematicViscosity,
146 bool print=verbose) const
147 {
148 const Scalar re_d = reynoldsNumber(velocity, diameter, kinematicViscosity, false);
149 const Scalar k = 1.5 * velocity * velocity
150 * turbulenceIntensity(re_d, false) * turbulenceIntensity(re_d, false);
151 if (print)
152 {
153 std::cout << "estimated k : " << k << " [m^2/s^2]" << std::endl;
154 }
155 return k;
156 }
157
162 Scalar dissipation(const Scalar velocity,
163 const Scalar diameter,
164 const Scalar kinematicViscosity,
165 bool print=verbose) const
166 {
167 using std::pow;
168 const Scalar k = turbulentKineticEnergy(velocity, diameter, kinematicViscosity, false);
169 const Scalar factor = 0.1643; // = cMu^(3/4) = 0.09^(3/4)
170 const Scalar epsilon = factor * pow(k, 1.5) / turbulenceLengthScale(diameter, false);
171 if (print)
172 {
173 std::cout << "estimated eps. : " << epsilon << " [m^2/s^3]" << std::endl;
174 }
175 return epsilon;
176 }
177
183 Scalar dissipationRate(const Scalar velocity,
184 const Scalar diameter,
185 const Scalar kinematicViscosity,
186 bool print=verbose) const
187 {
188 using std::pow;
189 const Scalar k = turbulentKineticEnergy(velocity, diameter, kinematicViscosity, false);
190 const Scalar factor = 0.54772; // = cMu^(1/4) = 0.09^(1/4)
191 const Scalar L = turbulenceLengthScale(diameter , false);
192 const Scalar omega = pow(k, 0.5) / (L*factor);
193 if (print)
194 {
195 std::cout << "estimated omega : " << omega << " [1/s]" << std::endl;
196
197 }
198 return omega;
199 }
200
205 Scalar viscosityTilde(const Scalar velocity,
206 const Scalar diameter,
207 const Scalar kinematicViscosity,
208 bool print=verbose) const
209 {
210 using std::abs;
211 using std::sqrt;
212 const Scalar re_d = reynoldsNumber(velocity, diameter, kinematicViscosity, false);
213 const Scalar viscosityTilde = sqrt(1.5) * abs(velocity)
214 * turbulenceIntensity(re_d)
216 if (print)
217 {
218 std::cout << "estimated nu~ : " << viscosityTilde << " [m^2/s]" << std::endl;
219 }
220 return viscosityTilde;
221 }
222};
223} // end namespace Dumux
224
225#endif // DUMUX_TURBULENCE_PROPERTIES_HH
Geometry::ctype diameter(const Geometry &geo)
Computes the longest distance between points of a geometry.
Definition: diameter.hh:36
Definition: adapt.hh:29
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