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  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 #include "Time.H"
31 #include "IOstream.H"
32 #include "IOdictionary.H"
33 
34 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
35 
36 template<class Type>
38 {
39  if (size() < 2)
40  {
42  << "Table " << name() << ": must have at least 2 values." << nl
43  << "Table size = " << size() << nl
44  << " min, interval width = " << x0_ << ", " << dx_ << nl
45  << exit(FatalError);
46  }
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
52 template<class Type>
54 (
55  const IOobject& io,
56  bool readFields
57 )
58 :
59  IOobject(io),
60  List<scalar>(2, Zero),
61  x0_(0.0),
62  dx_(1.0),
63  log10_(false),
64  bound_(false)
65 {
66  if (readFields)
67  {
68  IOdictionary dict(io);
69 
70  dict.readEntry("data", *this);
71  dict.readEntry("x0", x0_);
72  dict.readEntry("dx", dx_);
73  dict.readIfPresent("log10", log10_);
74  dict.readIfPresent("bound", bound_);
75  }
76 
77  checkTable();
78 }
79 
80 
81 template<class Type>
83 (
84  const word& tableName,
85  const objectRegistry& db,
86  const dictionary& dict,
87  const bool initialiseOnly
88 )
89 :
90  IOobject
91  (
92  tableName,
93  db.time().constant(),
94  db,
95  IOobject::NO_READ,
96  IOobject::NO_WRITE,
97  false // if used in BCs, could be used by multiple patches
98  ),
99  List<scalar>(2, Zero),
100  x0_(dict.get<scalar>("x0")),
101  dx_(dict.get<scalar>("dx")),
102  log10_(dict.getOrDefault<Switch>("log10", false)),
103  bound_(dict.getOrDefault<Switch>("bound", false))
104 {
105  if (initialiseOnly)
106  {
107  const scalar xMax = dict.get<scalar>("xMax");
108  const label nIntervals = static_cast<label>(xMax - x0_)/dx_ + 1;
109  this->setSize(nIntervals);
110  }
111  else
112  {
113  dict.readEntry("data", *this);
114  }
115 
116  checkTable();
117 }
118 
119 
120 template<class Type>
122 (
123  const uniformInterpolationTable& uit
124 )
125 :
126  IOobject(uit),
127  List<scalar>(uit),
128  x0_(uit.x0_),
129  dx_(uit.dx_),
130  log10_(uit.log10_),
131  bound_(uit.bound_)
132 {
133  checkTable();
134 }
135 
136 
137 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
138 
139 template<class Type>
141 {}
142 
143 
144 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
145 
146 template<class Type>
148 {
149  if (bound_)
150  {
151  x = max(min(xMax() - SMALL*dx_, x), x0_);
152  }
153  else
154  {
155  if (x < x0_)
156  {
158  << "Supplied value is less than minimum table value:" << nl
159  << "xMin=" << x0_ << ", xMax=" << xMax() << ", x=" << x << nl
160  << exit(FatalError);
161  }
162 
163  if (x > xMax())
164  {
166  << "Supplied value is greater than maximum table value:" << nl
167  << "xMin=" << x0_ << ", xMax=" << xMax() << ", x=" << x << nl
168  << exit(FatalError);
169  }
170  }
171 
172  const label i = static_cast<label>((x - x0_)/dx_);
173 
174  const scalar xLo = x0_ + i*dx_;
175 
176  Type fx = (x - xLo)/dx_*(operator[](i+1) - operator[](i)) + operator[](i);
177 
178  if (debug)
179  {
180  Info<< "Table: " << name() << ", x=" << x
181  << ", x_lo=" << xLo << ", x_hi=" << xLo + dx_
182  << ", f(x_lo)=" << operator[](i) << ", f(x_hi)=" << operator[](i+1)
183  << ", f(x)=" << fx << endl;
184  }
185 
186  return fx;
187 }
188 
189 
190 template<class Type>
192 (
193  scalar x
194 ) const
195 {
196  if (log10_)
197  {
198  if (x > 0)
199  {
200  x = ::log10(x);
201  }
202  else if (bound_ && (x <= 0))
203  {
204  x = x0_;
205  }
206  else
207  {
209  << "Table " << name() << nl
210  << "Supplied value must be greater than 0 when in log10 mode"
211  << nl << "x=" << x << nl << exit(FatalError);
212  }
213  }
214 
215  return interpolate(x);
216 }
217 
218 
219 template<class Type>
221 {
222  IOdictionary dict(*this);
223 
224  dict.add("data", static_cast<const List<scalar>&>(*this));
225  dict.add("x0", x0_);
226  dict.add("dx", dx_);
227  if (log10_)
228  {
229  dict.add("log10", log10_);
230  }
231  if (bound_)
232  {
233  dict.add("bound", bound_);
234  }
235 
236  dict.regIOobject::writeObject
237  (
238  IOstreamOption(IOstream::ASCII, dict.time().writeCompression()),
239  true
240  );
241 }
242 
243 
244 // ************************************************************************* //
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:169
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:77
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::uniformInterpolationTable
Table with uniform interval in independent variable, with linear interpolation.
Definition: uniformInterpolationTable.H:72
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
uniformInterpolationTable.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::uniformInterpolationTable::interpolateLog10
Type interpolateLog10(scalar x) const
Interpolate - takes log10 flag into account.
Definition: uniformInterpolationTable.C:192
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::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::log10
dimensionedScalar log10(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:263
IOstream.H
Foam::IOstreamOption
The IOstreamOption is a simple container for options an IOstream can normally have.
Definition: IOstreamOption.H:63
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:123
nIntervals
const label nIntervals(pdfDictionary.get< label >("nIntervals"))
Foam::uniformInterpolationTable::~uniformInterpolationTable
~uniformInterpolationTable()
Destructor.
Definition: uniformInterpolationTable.C:140
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:453
Foam::nl
constexpr char nl
Definition: Ostream.H:404
xMax
const scalar xMax
Definition: createFields.H:35
Foam::List< scalar >
Foam::uniformInterpolationTable::write
void write() const
Write.
Definition: uniformInterpolationTable.C:220
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::uniformInterpolationTable::uniformInterpolationTable
uniformInterpolationTable(const IOobject &, const bool readFields)
Construct from IOobject and readFields flag.
Definition: uniformInterpolationTable.C:54
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::uniformInterpolationTable::interpolate
Type interpolate(scalar x) const
Interpolate.
Definition: uniformInterpolationTable.C:147
Foam::TimePaths::constant
const word & constant() const
Return constant name.
Definition: TimePathsI.H:96
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.
Foam::objectRegistry::time
const Time & time() const noexcept
Return time registry.
Definition: objectRegistry.H:178