SpaldingsLaw.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2015 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "SpaldingsLaw.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  namespace tabulatedWallFunctions
36  {
37  defineTypeNameAndDebug(SpaldingsLaw, 0);
39  (
40  tabulatedWallFunction,
41  SpaldingsLaw,
42  dictionary
43  );
44  }
45 }
46 
48 
49 const Foam::scalar
51 
52 
53 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
54 
56 {
57  // Initialise u+ and R
58  scalar Re = 0.0;
59  scalar uPlus = 1;
60 
61  // Populate the table
63  {
64  if (invertedTable_.log10())
65  {
66  Re = pow(10, (i*invertedTable_.dx() + invertedTable_.x0()));
67  }
68  else
69  {
71  }
72 
73  // Use latest available u+ estimate
74  if (i > 0)
75  {
76  uPlus = invertedTable_[i-1];
77  }
78 
79  // Newton iterations to determine u+
80  label iter = 0;
81  scalar error = GREAT;
82  do
83  {
84  scalar kUPlus = min(kappa_*uPlus, 50);
85  scalar A =
86  E_*sqr(uPlus)
87  + uPlus
88  *(exp(kUPlus) - pow3(kUPlus)/6 - 0.5*sqr(kUPlus) - kUPlus - 1);
89 
90  scalar f = - Re + A/E_;
91 
92  scalar df =
93  (
94  2*E_*uPlus
95  + exp(kUPlus)*(kUPlus + 1)
96  - 2.0/3.0*pow3(kUPlus)
97  - 1.5*sqr(kUPlus)
98  - 2*kUPlus
99  - 1
100  )/E_;
101 
102  scalar uPlusNew = uPlus - f/(df + ROOTVSMALL);
103  error = mag((uPlus - uPlusNew)/uPlusNew);
104  uPlus = uPlusNew;
105  } while (error > tolerance_ && ++iter < maxIters_);
106 
107  if (iter == maxIters_)
108  {
110  << "Newton iterations not converged:" << nl
111  << " iters = " << iter << ", error = " << error << endl;
112  }
113 
114  // Set new values - constrain u+ >= 0
115  invertedTable_[i] = max(0, uPlus);
116  }
117 }
118 
119 
120 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
121 
123 (
124  const dictionary& dict,
125  const polyMesh& mesh
126 )
127 :
128  tabulatedWallFunction(dict, mesh, typeName),
129  kappa_(coeffDict_.get<scalar>("kappa")),
130  E_(coeffDict_.get<scalar>("E"))
131 {
132  invertFunction();
133 
134  if (debug)
135  {
136  writeData(Info);
137  }
138 }
139 
140 
141 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
142 
144 (
145  const scalar uPlus
146 ) const
147 {
148  scalar kUPlus = min(kappa_*uPlus, 50);
149 
150  return
151  uPlus
152  + 1/E_*(exp(kUPlus) - pow3(kUPlus)/6 - 0.5*sqr(kUPlus) - kUPlus - 1);
153 }
154 
155 
157 (
158  const scalar uPlus
159 ) const
160 {
161  return uPlus*yPlus(uPlus);
162 }
163 
164 
166 {
167  if (invertedTable_.log10())
168  {
169  os << "log10(Re), y+, u+:" << endl;
170  forAll(invertedTable_, i)
171  {
172  scalar uPlus = invertedTable_[i];
173  scalar Re = ::log10(this->Re(uPlus));
174  scalar yPlus = this->yPlus(uPlus);
175  os << Re << ", " << yPlus << ", " << uPlus << endl;
176  }
177  }
178  else
179  {
180  os << "Re, y+, u+:" << endl;
181  forAll(invertedTable_, i)
182  {
183  scalar uPlus = invertedTable_[i];
184  scalar Re = this->Re(uPlus);
185  scalar yPlus = this->yPlus(uPlus);
186  os << Re << ", " << yPlus << ", " << uPlus << endl;
187  }
188  }
189 }
190 
191 
192 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::tabulatedWallFunctions::SpaldingsLaw::tolerance_
static const scalar tolerance_
Definition: SpaldingsLaw.H:87
Foam::uniformInterpolationTable::log10
const Switch & log10() const
Definition: uniformInterpolationTableI.H:36
writeData
const bool writeData(pdfDictionary.get< bool >("writeData"))
Foam::uniformInterpolationTable::x0
scalar x0() const
Definition: uniformInterpolationTableI.H:22
A
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
uPlus
scalar uPlus
Definition: evaluateNearWall.H:18
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:254
Foam::tabulatedWallFunctions::SpaldingsLaw::maxIters_
static const label maxIters_
Definition: SpaldingsLaw.H:84
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Definition: hashSets.C:26
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
SpaldingsLaw.H
Foam::uniformInterpolationTable::dx
scalar dx() const
Definition: uniformInterpolationTableI.H:29
Foam::tabulatedWallFunctions::SpaldingsLaw::kappa_
scalar kappa_
Definition: SpaldingsLaw.H:75
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:82
Foam::Info
messageStream Info
Foam::log10
dimensionedScalar log10(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:256
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:68
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Definition: hashSets.C:40
Foam::tabulatedWallFunctions::SpaldingsLaw::SpaldingsLaw
SpaldingsLaw(const dictionary &dict, const polyMesh &mesh)
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::tabulatedWallFunctions::SpaldingsLaw::E_
scalar E_
Definition: SpaldingsLaw.H:78
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Definition: atmBoundaryLayer.C:26
Foam::tabulatedWallFunctions::SpaldingsLaw::writeData
virtual void writeData(Ostream &os) const
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:44
Foam::nl
constexpr char nl
Definition: Ostream.H:424
Foam::tabulatedWallFunctions::SpaldingsLaw::invertFunction
virtual void invertFunction()
f
labelList f(nPoints)
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Definition: createFields.H:11
Foam::Re
scalarField Re(const UList< complex > &cf)
Definition: complexField.C:152
Foam::PtrListOps::get
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
yPlus
scalar yPlus
Definition: evaluateNearWall.H:16
Foam::tabulatedWallFunctions::SpaldingsLaw::Re
virtual scalar Re(const scalar uPlus) const
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::tabulatedWallFunctions::tabulatedWallFunction::invertedTable_
uniformInterpolationTable< scalar > invertedTable_
Definition: tabulatedWallFunction.H:71
WarningInFunction
#define WarningInFunction
Definition: messageStream.H:353
Foam::tabulatedWallFunctions::SpaldingsLaw::yPlus
virtual scalar yPlus(const scalar uPlus) const