general.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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "general.H"
28 #include "Tuple2.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34  namespace tabulatedWallFunctions
35  {
38  (
39  tabulatedWallFunction,
40  general,
41  dictionary
42  );
43  }
44 
45  template<>
46  const char* Foam::NamedEnum
47  <
49  1
50  >::names[] =
51  {
52  "linear"
53  };
54 
55 }
56 
57 const
60 
61 
62 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
63 
65 {
66  scalarList Rey(uPlus_.size(), 0.0);
67 
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  }
77 
78  // Populate the U+ vs Re table
80  {
81  scalar Re = i*invertedTable_.dx() + invertedTable_.x0();
83  }
84 }
85 
86 
88 (
89  const scalar xi,
90  const scalarList& x,
91  const scalarList& fx
92 ) const
93 {
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;
114 
115  return (xi - x[i1])/(x[i2] - x[i1])*(fx[i2] - fx[i1]) + fx[i1];
116  }
117 
118  break;
119  }
120  default:
121  {
123  << "Unknown interpolation method" << nl
124  << abort(FatalError);
125  }
126  }
127 
128  return 0.0;
129 }
130 
131 
132 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
133 
135 (
136  const dictionary& dict,
137  const polyMesh& mesh
138 )
139 :
140  tabulatedWallFunction(dict, mesh, typeName),
141  interpType_(interpolationTypeNames_[coeffDict_.lookup("interpType")]),
142  yPlus_(),
143  uPlus_(),
144  log10YPlus_(coeffDict_.lookup("log10YPlus")),
145  log10UPlus_(coeffDict_.lookup("log10UPlus"))
146 {
147  List<Tuple2<scalar, scalar> > inputTable = coeffDict_.lookup("inputTable");
148  if (inputTable.size() < 2)
149  {
151  << "Input table must have at least 2 values" << nl
152  << exit(FatalError);
153  }
154 
155  yPlus_.setSize(inputTable.size());
156  uPlus_.setSize(inputTable.size());
157 
158  forAll(inputTable, i)
159  {
160  if (log10YPlus_)
161  {
162  yPlus_[i] = pow(10, inputTable[i].first());
163  }
164  else
165  {
166  yPlus_[i] = inputTable[i].first();
167  }
168 
169  if (log10UPlus_)
170  {
171  uPlus_[i] = pow(10, inputTable[i].second());
172  }
173  else
174  {
175  uPlus_[i] = inputTable[i].second();
176  }
177  }
178 
179  invertTable();
180 
181  if (debug)
182  {
183  writeData(Info);
184  }
185 }
186 
187 
188 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
189 
191 {}
192 
193 
194 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
195 
197 (
198  const scalar uPlus
199 ) const
200 {
201  return interpolate(uPlus, uPlus_, yPlus_);
202 }
203 
204 
206 (
207  const scalar uPlus
208 ) const
209 {
210  return uPlus*yPlus(uPlus);
211 }
212 
213 
214 // * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * * //
215 
217 {
218  if (invertedTable_.log10())
219  {
220  os << "log10(Re), y+, u+:" << endl;
221  forAll(invertedTable_, i)
222  {
223  scalar uPlus = invertedTable_[i];
224  scalar Re = ::log10(this->Re(uPlus));
225  scalar yPlus = this->yPlus(uPlus);
226  os << Re << ", " << yPlus << ", " << uPlus << endl;
227  }
228  }
229  else
230  {
231  os << "Re, y+, u+:" << endl;
232  forAll(invertedTable_, i)
233  {
234  scalar uPlus = invertedTable_[i];
235  scalar Re = this->Re(uPlus);
236  scalar yPlus = this->yPlus(uPlus);
237  os << Re << ", " << yPlus << ", " << uPlus << endl;
238  }
239  }
240 }
241 
242 
243 // ************************************************************************* //
Foam::scalarList
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
Foam::uniformInterpolationTable::log10
const Switch & log10() const
Return the log10(x) flag.
Definition: uniformInterpolationTableI.H:41
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Tuple2.H
Foam::tabulatedWallFunctions::general::yPlus_
List< scalar > yPlus_
Input y+ values.
Definition: general.H:104
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
Foam::Re
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
Foam::uniformInterpolationTable::x0
scalar x0() const
Return the lower limit.
Definition: uniformInterpolationTableI.H:27
uPlus
scalar uPlus
Definition: evaluateNearWall.H:18
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::tabulatedWallFunctions::general::~general
virtual ~general()
Destructor.
Foam::tabulatedWallFunctions::general::general
general(const dictionary &dict, const polyMesh &mesh)
Foam::uniformInterpolationTable::dx
scalar dx() const
Return the fixed interval.
Definition: uniformInterpolationTableI.H:34
List::size
int size() const
Definition: ListI.H:83
Foam::interpolate
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::tabulatedWallFunctions::general::interpolationType
interpolationType
Enumeration listing available interpolation types.
Definition: general.H:88
general
general distribution model
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::tabulatedWallFunctions::general::invertTable
virtual void invertTable()
Invert the table.
Foam::log10
dimensionedScalar log10(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:254
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
Foam::tabulatedWallFunctions::general::writeData
virtual void writeData(Ostream &os) const
Write to Ostream.
dict
dictionary dict
Definition: searchingEngine.H:14
writeData
const bool writeData(readBool(pdfDictionary.lookup("writeData")))
Foam::FatalError
error FatalError
Foam::tabulatedWallFunctions::general::uPlus_
List< scalar > uPlus_
Input U+ values.
Definition: general.H:107
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::tabulatedWallFunctions::general::yPlus
virtual scalar yPlus(const scalar uPlus) const
Return y+ as a function of u+.
Foam::tabulatedWallFunctions::general::interpolate
virtual scalar interpolate(const scalar xi, const scalarList &x, const scalarList &fx) const
Interpolate.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::tabulatedWallFunctions::general::interpolationTypeNames_
static const NamedEnum< interpolationType, 1 > interpolationTypeNames_
Definition: general.H:93
x
x
Definition: LISASMDCalcMethod2.H:52
Rey
scalar Rey
Definition: evaluateNearWall.H:28
yPlus
scalar yPlus
Definition: evaluateNearWall.H:16
List
Definition: Test.C:19
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::tabulatedWallFunctions::general::Re
virtual scalar Re(const scalar uPlus) const
Return Reynolds number as a function of u+.
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::tabulatedWallFunctions::tabulatedWallFunction::invertedTable_
uniformInterpolationTable< scalar > invertedTable_
Inverted table.
Definition: tabulatedWallFunction.H:73
Foam::NamedEnum
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:52
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress