ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / applications / utilities / preProcessing / wallFunctionTable / tabulatedWallFunction / general / general.C
blob3ef4cdfd24f69578bedf74f12af14b85c5be9cf1
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM is free software: you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by
13     the Free Software Foundation, either version 3 of the License, or
14     (at your option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "general.H"
27 #include "addToRunTimeSelectionTable.H"
28 #include "Tuple2.H"
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 namespace Foam
34     namespace tabulatedWallFunctions
35     {
36         defineTypeNameAndDebug(general, 0);
37         addToRunTimeSelectionTable
38         (
39             tabulatedWallFunction,
40             general,
41             dictionary
42         );
43     }
45     template<>
46     const char* Foam::NamedEnum
47     <
48         Foam::tabulatedWallFunctions::general::interpolationType,
49         1
50     >::names[] =
51     {
52         "linear"
53     };
57 const
58 Foam::NamedEnum<Foam::tabulatedWallFunctions::general::interpolationType, 1>
59     Foam::tabulatedWallFunctions::general::interpolationTypeNames_;
62 // * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * * //
64 void Foam::tabulatedWallFunctions::general::invertTable()
66     scalarList Rey(uPlus_.size(), 0.0);
68     // Calculate Reynolds number
69     forAll(uPlus_, i)
70     {
71         Rey[i] = yPlus_[i]*uPlus_[i];
72         if (invertedTable_.log10())
73         {
74             Rey[i] = ::log10(max(ROOTVSMALL, Rey[i]));
75         }
76     }
78     // Populate the U+ vs Re table
79     forAll(invertedTable_, i)
80     {
81         scalar Re = i*invertedTable_.dx() + invertedTable_.x0();
82         invertedTable_[i] = max(0, interpolate(Re, Rey, uPlus_));
83     }
87 Foam::scalar Foam::tabulatedWallFunctions::general::interpolate
89     const scalar xi,
90     const scalarList& x,
91     const scalarList& fx
92 ) const
94     switch (interpType_)
95     {
96         case itLinear:
97         {
98             if (xi <= x[0])
99             {
100                 return fx[0];
101             }
102             else if (xi >= x.last())
103             {
104                 return fx.last();
105             }
106             else
107             {
108                 label i2 = 0;
109                 while (x[i2] < xi)
110                 {
111                     i2++;
112                 }
113                 label i1 = i2 - 1;
115                 return (xi - x[i1])/(x[i2] - x[i1])*(fx[i2] - fx[i1]) + fx[i1];
116             }
118             break;
119         }
120         default:
121         {
122             FatalErrorIn
123             (
124                 "tabulatedWallFunctions::general::interpolate"
125                 "("
126                     "const scalar, "
127                     "const scalarList&, "
128                     "const scalarList&"
129                 ")"
130             )   << "Unknown interpolation method" << nl
131                 << abort(FatalError);
132         }
133     }
135     return 0.0;
139 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
141 Foam::tabulatedWallFunctions::general::general
143     const dictionary& dict,
144     const polyMesh& mesh
147     tabulatedWallFunction(dict, mesh, typeName),
148     interpType_(interpolationTypeNames_[coeffDict_.lookup("interpType")]),
149     yPlus_(),
150     uPlus_(),
151     log10YPlus_(coeffDict_.lookup("log10YPlus")),
152     log10UPlus_(coeffDict_.lookup("log10UPlus"))
154     List<Tuple2<scalar, scalar> > inputTable = coeffDict_.lookup("inputTable");
155     if (inputTable.size() < 2)
156     {
157         FatalErrorIn
158         (
159             "tabulatedWallFunctions::general::general"
160             "("
161                 "const dictionary&, "
162                 "const polyMesh&"
163             ")"
164         )   << "Input table must have at least 2 values" << nl
165             << exit(FatalError);
166     }
168     yPlus_.setSize(inputTable.size());
169     uPlus_.setSize(inputTable.size());
171     forAll(inputTable, i)
172     {
173         if (log10YPlus_)
174         {
175             yPlus_[i] = pow(10, inputTable[i].first());
176         }
177         else
178         {
179             yPlus_[i] = inputTable[i].first();
180         }
182         if (log10UPlus_)
183         {
184             uPlus_[i] = pow(10, inputTable[i].second());
185         }
186         else
187         {
188             uPlus_[i] = inputTable[i].second();
189         }
190     }
192     invertTable();
194     if (debug)
195     {
196         writeData(Info);
197     }
201 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
203 Foam::tabulatedWallFunctions::general::~general()
207 // * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * * //
209 Foam::scalar Foam::tabulatedWallFunctions::general::yPlus
211     const scalar uPlus
212 ) const
214     return interpolate(uPlus, uPlus_, yPlus_);
218 Foam::scalar Foam::tabulatedWallFunctions::general::Re
220     const scalar uPlus
221 ) const
223     return uPlus*yPlus(uPlus);
227 // * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * * //
229 void Foam::tabulatedWallFunctions::general::writeData(Ostream& os) const
231     if (invertedTable_.log10())
232     {
233         os  << "log10(Re), y+, u+:" << endl;
234         forAll(invertedTable_, i)
235         {
236             scalar uPlus = invertedTable_[i];
237             scalar Re = ::log10(this->Re(uPlus));
238             scalar yPlus = this->yPlus(uPlus);
239             os  << Re << ", " << yPlus << ", " << uPlus << endl;
240         }
241     }
242     else
243     {
244         os  << "Re, y+, u+:" << endl;
245         forAll(invertedTable_, i)
246         {
247             scalar uPlus = invertedTable_[i];
248             scalar Re = this->Re(uPlus);
249             scalar yPlus = this->yPlus(uPlus);
250             os  << Re << ", " << yPlus << ", " << uPlus << endl;
251         }
252     }
256 // ************************************************************************* //