Loading...
Searching...
No Matches
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-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
36template<class Type>
37void Foam::uniformInterpolationTable<Type>::checkTable() const
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
52template<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
81template<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 IOobject::NO_REGISTER
98 // No register: could be used by multiple patches if used in BCs
99 ),
100 List<scalar>(2, Zero),
101 x0_(dict.get<scalar>("x0")),
102 dx_(dict.get<scalar>("dx")),
103 log10_(dict.getOrDefault<Switch>("log10", false)),
104 bound_(dict.getOrDefault<Switch>("bound", false))
105{
106 if (initialiseOnly)
107 {
108 const scalar xMax = dict.get<scalar>("xMax");
109 const label nIntervals = static_cast<label>(xMax - x0_)/dx_ + 1;
110 this->setSize(nIntervals);
111 }
112 else
113 {
114 dict.readEntry("data", *this);
116
117 checkTable();
118}
119
121template<class Type>
123(
125)
126:
127 IOobject(uit),
128 List<scalar>(uit),
129 x0_(uit.x0_),
130 dx_(uit.dx_),
131 log10_(uit.log10_),
132 bound_(uit.bound_)
133{
134 checkTable();
135}
136
137
138// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
139
140template<class Type>
142{}
143
144
145// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
146
147template<class Type>
149{
150 if (bound_)
151 {
152 x = clamp(x, x0_, (xMax() - SMALL*dx_));
153 }
154 else
155 {
156 if (x < x0_)
157 {
159 << "Supplied value is less than minimum table value:" << nl
160 << "xMin=" << x0_ << ", xMax=" << xMax() << ", x=" << x << nl
161 << exit(FatalError);
162 }
163
164 if (x > xMax())
165 {
167 << "Supplied value is greater than maximum table value:" << nl
168 << "xMin=" << x0_ << ", xMax=" << xMax() << ", x=" << x << nl
169 << exit(FatalError);
170 }
171 }
172
173 const label i = static_cast<label>((x - x0_)/dx_);
174
175 const scalar xLo = x0_ + i*dx_;
176
177 Type fx = (x - xLo)/dx_*(operator[](i+1) - operator[](i)) + operator[](i);
178
179 if (debug)
180 {
181 Info<< "Table: " << name() << ", x=" << x
182 << ", x_lo=" << xLo << ", x_hi=" << xLo + dx_
183 << ", f(x_lo)=" << operator[](i) << ", f(x_hi)=" << operator[](i+1)
184 << ", f(x)=" << fx << endl;
186
187 return fx;
188}
189
190
191template<class Type>
193(
194 scalar x
195) const
196{
197 if (log10_)
198 {
199 if (x > 0)
200 {
201 x = ::log10(x);
202 }
203 else if (bound_ && (x <= 0))
204 {
205 x = x0_;
206 }
207 else
208 {
210 << "Table " << name() << nl
211 << "Supplied value must be greater than 0 when in log10 mode"
212 << nl << "x=" << x << nl << exit(FatalError);
213 }
215
216 return interpolate(x);
217}
218
219
220template<class Type>
222{
223 IOdictionary dict(*this);
224
225 dict.add("data", static_cast<const List<scalar>&>(*this));
226 dict.add("x0", x0_);
227 dict.add("dx", dx_);
228 if (log10_)
229 {
230 dict.add("log10", log10_);
231 }
232 if (bound_)
233 {
234 dict.add("bound", bound_);
235 }
236
237 dict.regIOobject::writeObject
238 (
239 IOstreamOption(IOstreamOption::ASCII, dict.time().writeCompression()),
240 true
241 );
242}
243
244
245// ************************************************************************* //
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition IOobject.C:450
A simple container for options an IOstream can normally have.
@ ASCII
"ascii" (normal default)
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void setSize(label n)
Definition List.H:536
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition Switch.H:81
bool get(const label i) const
Definition UList.H:868
scalar & operator[](const label i)
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Registry of regIOobjects.
constant condensation/saturation model.
Type interpolateLog10(scalar x) const
Interpolate - takes log10 flag into account.
scalar xMax() const
Return the maximum x value.
uniformInterpolationTable(const IOobject &, const bool readFields)
Construct from IOobject and readFields flag.
Type interpolate(scalar x) const
Interpolate.
A class for handling words, derived from Foam::string.
Definition word.H:66
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const auto & io
auto & name
Different types of constants.
Namespace for handling debugging switches.
Definition debug.C:45
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
messageStream Info
Information stream (stdout output on master, null elsewhere).
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
dimensionedScalar log10(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition curveTools.C:75
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
const label nIntervals(pdfDictionary.get< label >("nIntervals"))
const scalar xMax