uniformInterpolationTable.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 
27 #include "Time.H"
28 
29 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
30 
31 template<class Type>
33 {
34  if (size() < 2)
35  {
37  << "Table " << name() << ": must have at least 2 values." << nl
38  << "Table size = " << size() << nl
39  << " min, interval width = " << x0_ << ", " << dx_ << nl
40  << exit(FatalError);
41  }
42 }
43 
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
47 template<class Type>
49 (
50  const IOobject& io,
51  bool readFields
52 )
53 :
54  IOobject(io),
55  List<scalar>(2, 0.0),
56  x0_(0.0),
57  dx_(1.0),
58  log10_(false),
59  bound_(false)
60 {
61  if (readFields)
62  {
63  IOdictionary dict(io);
64 
65  dict.lookup("data") >> *this;
66  dict.lookup("x0") >> x0_;
67  dict.lookup("dx") >> dx_;
68  dict.readIfPresent("log10", log10_);
69  dict.readIfPresent("bound", bound_);
70  }
71 
72  checkTable();
73 }
74 
75 
76 template<class Type>
78 (
79  const word& tableName,
80  const objectRegistry& db,
81  const dictionary& dict,
82  const bool initialiseOnly
83 )
84 :
85  IOobject
86  (
87  tableName,
88  db.time().constant(),
89  db,
90  IOobject::NO_READ,
91  IOobject::NO_WRITE,
92  false // if used in BCs, could be used by multiple patches
93  ),
94  List<scalar>(2, 0.0),
95  x0_(readScalar(dict.lookup("x0"))),
96  dx_(readScalar(dict.lookup("dx"))),
97  log10_(dict.lookupOrDefault<Switch>("log10", false)),
98  bound_(dict.lookupOrDefault<Switch>("bound", false))
99 {
100  if (initialiseOnly)
101  {
102  const scalar xMax = readScalar(dict.lookup("xMax"));
103  const label nIntervals = static_cast<label>(xMax - x0_)/dx_ + 1;
104  this->setSize(nIntervals);
105  }
106  else
107  {
108  dict.lookup("data") >> *this;
109  }
110 
111  checkTable();
112 }
113 
114 
115 template<class Type>
117 (
118  const uniformInterpolationTable& uit
119 )
120 :
121  IOobject(uit),
122  List<scalar>(uit),
123  x0_(uit.x0_),
124  dx_(uit.dx_),
125  log10_(uit.log10_),
126  bound_(uit.bound_)
127 {
128  checkTable();
129 }
130 
131 
132 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
133 
134 template<class Type>
136 {}
137 
138 
139 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
140 
141 template<class Type>
143 {
144  if (bound_)
145  {
146  x = max(min(xMax() - SMALL*dx_, x), x0_);
147  }
148  else
149  {
150  if (x < x0_)
151  {
153  << "Supplied value is less than minimum table value:" << nl
154  << "xMin=" << x0_ << ", xMax=" << xMax() << ", x=" << x << nl
155  << exit(FatalError);
156  }
157 
158  if (x > xMax())
159  {
161  << "Supplied value is greater than maximum table value:" << nl
162  << "xMin=" << x0_ << ", xMax=" << xMax() << ", x=" << x << nl
163  << exit(FatalError);
164  }
165  }
166 
167  const label i = static_cast<label>((x - x0_)/dx_);
168 
169  const scalar xLo = x0_ + i*dx_;
170 
171  Type fx = (x - xLo)/dx_*(operator[](i+1) - operator[](i)) + operator[](i);
172 
173  if (debug)
174  {
175  Info<< "Table: " << name() << ", x=" << x
176  << ", x_lo=" << xLo << ", x_hi=" << xLo + dx_
177  << ", f(x_lo)=" << operator[](i) << ", f(x_hi)=" << operator[](i+1)
178  << ", f(x)=" << fx << endl;
179  }
180 
181  return fx;
182 }
183 
184 
185 template<class Type>
187 (
188  scalar x
189 ) const
190 {
191  if (log10_)
192  {
193  if (x > 0)
194  {
195  x = ::log10(x);
196  }
197  else if (bound_ && (x <= 0))
198  {
199  x = x0_;
200  }
201  else
202  {
204  << "Table " << name() << nl
205  << "Supplied value must be greater than 0 when in log10 mode"
206  << nl << "x=" << x << nl << exit(FatalError);
207  }
208  }
209 
210  return interpolate(x);
211 }
212 
213 
214 template<class Type>
216 {
217  IOdictionary dict(*this);
218 
219  dict.add("data", static_cast<const List<scalar>&>(*this));
220  dict.add("x0", x0_);
221  dict.add("dx", dx_);
222  if (log10_)
223  {
224  dict.add("log10", log10_);
225  }
226  if (bound_)
227  {
228  dict.add("bound", bound_);
229  }
230 
231  dict.regIOobject::write();
232 }
233 
234 
235 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
const Time & time() const
Return time.
const label nIntervals(readLabel(pdfDictionary.lookup("nIntervals")))
dictionary dict
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none.
Definition: Switch.H:60
Table with uniform interval in independant variable, with linear interpolation.
tmp< surfaceScalarField > interpolate(const RhoType &rho)
uniformInterpolationTable(const IOobject &, const bool readFields)
Construct from IOobject and readFields flag.
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:737
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const HashSet< word > &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the specified type.
Definition: ReadFields.C:244
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
points setSize(newPointi)
A class for handling words, derived from string.
Definition: word.H:59
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
const word & constant() const
Return constant name.
Definition: TimePaths.H:124
static const char nl
Definition: Ostream.H:262
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const scalar xMax
Definition: createFields.H:35
Type interpolateLog10(scalar x) const
Interpolate - takes log10 flag into account.
messageStream Info
Registry of regIOobjects.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Type interpolate(scalar x) const
Interpolate.
dimensionedScalar log10(const dimensionedScalar &ds)
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451