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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 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 
29 #include "Time.H"
30 #include "IOstream.H"
31 #include "IOdictionary.H"
32 
33 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
34 
35 template<class Type>
37 {
38  if (size() < 2)
39  {
41  << "Table " << name() << ": must have at least 2 values." << nl
42  << "Table size = " << size() << nl
43  << " min, interval width = " << x0_ << ", " << dx_ << nl
44  << exit(FatalError);
45  }
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
51 template<class Type>
53 (
54  const IOobject& io,
55  bool readFields
56 )
57 :
58  IOobject(io),
59  List<scalar>(2, Zero),
60  x0_(0.0),
61  dx_(1.0),
62  log10_(false),
63  bound_(false)
64 {
65  if (readFields)
66  {
67  IOdictionary dict(io);
68 
69  dict.readEntry("data", *this);
70  dict.readEntry("x0", x0_);
71  dict.readEntry("dx", dx_);
72  dict.readIfPresent("log10", log10_);
73  dict.readIfPresent("bound", bound_);
74  }
75 
76  checkTable();
77 }
78 
79 
80 template<class Type>
82 (
83  const word& tableName,
84  const objectRegistry& db,
85  const dictionary& dict,
86  const bool initialiseOnly
87 )
88 :
89  IOobject
90  (
91  tableName,
92  db.time().constant(),
93  db,
94  IOobject::NO_READ,
95  IOobject::NO_WRITE,
96  false // if used in BCs, could be used by multiple patches
97  ),
98  List<scalar>(2, Zero),
99  x0_(dict.get<scalar>("x0")),
100  dx_(dict.get<scalar>("dx")),
101  log10_(dict.lookupOrDefault<Switch>("log10", false)),
102  bound_(dict.lookupOrDefault<Switch>("bound", false))
103 {
104  if (initialiseOnly)
105  {
106  const scalar xMax = dict.get<scalar>("xMax");
107  const label nIntervals = static_cast<label>(xMax - x0_)/dx_ + 1;
108  this->setSize(nIntervals);
109  }
110  else
111  {
112  dict.readEntry("data", *this);
113  }
114 
115  checkTable();
116 }
117 
118 
119 template<class Type>
121 (
122  const uniformInterpolationTable& uit
123 )
124 :
125  IOobject(uit),
126  List<scalar>(uit),
127  x0_(uit.x0_),
128  dx_(uit.dx_),
129  log10_(uit.log10_),
130  bound_(uit.bound_)
131 {
132  checkTable();
133 }
134 
135 
136 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
137 
138 template<class Type>
140 {}
141 
142 
143 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
144 
145 template<class Type>
147 {
148  if (bound_)
149  {
150  x = max(min(xMax() - SMALL*dx_, x), x0_);
151  }
152  else
153  {
154  if (x < x0_)
155  {
157  << "Supplied value is less than minimum table value:" << nl
158  << "xMin=" << x0_ << ", xMax=" << xMax() << ", x=" << x << nl
159  << exit(FatalError);
160  }
161 
162  if (x > xMax())
163  {
165  << "Supplied value is greater than maximum table value:" << nl
166  << "xMin=" << x0_ << ", xMax=" << xMax() << ", x=" << x << nl
167  << exit(FatalError);
168  }
169  }
170 
171  const label i = static_cast<label>((x - x0_)/dx_);
172 
173  const scalar xLo = x0_ + i*dx_;
174 
175  Type fx = (x - xLo)/dx_*(operator[](i+1) - operator[](i)) + operator[](i);
176 
177  if (debug)
178  {
179  Info<< "Table: " << name() << ", x=" << x
180  << ", x_lo=" << xLo << ", x_hi=" << xLo + dx_
181  << ", f(x_lo)=" << operator[](i) << ", f(x_hi)=" << operator[](i+1)
182  << ", f(x)=" << fx << endl;
183  }
184 
185  return fx;
186 }
187 
188 
189 template<class Type>
191 (
192  scalar x
193 ) const
194 {
195  if (log10_)
196  {
197  if (x > 0)
198  {
199  x = ::log10(x);
200  }
201  else if (bound_ && (x <= 0))
202  {
203  x = x0_;
204  }
205  else
206  {
208  << "Table " << name() << nl
209  << "Supplied value must be greater than 0 when in log10 mode"
210  << nl << "x=" << x << nl << exit(FatalError);
211  }
212  }
213 
214  return interpolate(x);
215 }
216 
217 
218 template<class Type>
220 {
221  IOdictionary dict(*this);
222 
223  dict.add("data", static_cast<const List<scalar>&>(*this));
224  dict.add("x0", x0_);
225  dict.add("dx", dx_);
226  if (log10_)
227  {
228  dict.add("log10", log10_);
229  }
230  if (bound_)
231  {
232  dict.add("bound", bound_);
233  }
234 
235  dict.regIOobject::writeObject
236  (
237  IOstream::ASCII,
238  IOstream::currentVersion,
239  dict.time().writeCompression(),
240  true
241  );
242 }
243 
244 
245 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
setSize
points setSize(newPointi)
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:70
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::uniformInterpolationTable
Table with uniform interval in independant variable, with linear interpolation.
Definition: uniformInterpolationTable.H:72
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::objectRegistry::time
const Time & time() const
Return time.
Definition: objectRegistry.H:186
uniformInterpolationTable.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::uniformInterpolationTable::interpolateLog10
Type interpolateLog10(scalar x) const
Interpolate - takes log10 flag into account.
Definition: uniformInterpolationTable.C:191
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
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:62
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::log10
dimensionedScalar log10(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:263
IOstream.H
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
nIntervals
const label nIntervals(pdfDictionary.get< label >("nIntervals"))
Foam::uniformInterpolationTable::~uniformInterpolationTable
~uniformInterpolationTable()
Destructor.
Definition: uniformInterpolationTable.C:139
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::readFields
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
Definition: ReadFieldsTemplates.C:312
IOdictionary.H
Time.H
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::nl
constexpr char nl
Definition: Ostream.H:372
xMax
const scalar xMax
Definition: createFields.H:35
Foam::List< scalar >
Foam::uniformInterpolationTable::write
void write() const
Write.
Definition: uniformInterpolationTable.C:219
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::uniformInterpolationTable::uniformInterpolationTable
uniformInterpolationTable(const IOobject &, const bool readFields)
Construct from IOobject and readFields flag.
Definition: uniformInterpolationTable.C:53
Foam::uniformInterpolationTable::interpolate
Type interpolate(scalar x) const
Interpolate.
Definition: uniformInterpolationTable.C:146
Foam::TimePaths::constant
const word & constant() const
Return constant name.
Definition: TimePathsI.H:88
Foam::fac::interpolate
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.