TableBase.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) 2019-2021 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
29#include "TableBase.H"
30#include "Time.H"
32
33// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
34
35template<class Type>
38{
39 if (!interpolatorPtr_)
40 {
41 // Re-work table into linear list
42 tableSamplesPtr_.reset(new scalarField(table_.size()));
43 auto& samples = *tableSamplesPtr_;
44 forAll(table_, i)
45 {
46 samples[i] = table_[i].first();
47 }
48 interpolatorPtr_ = interpolationWeights::New
49 (
50 interpolationScheme_,
52 );
53 }
54
55 return *interpolatorPtr_;
56}
57
58
59// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
60
61template<class Type>
63(
64 const word& name,
65 const dictionary& dict,
66 const objectRegistry* obrPtr
67)
68:
69 Function1<Type>(name, dict, obrPtr),
70 bounding_
71 (
72 bounds::repeatableBoundingNames.getOrDefault
73 (
74 "outOfBounds",
75 dict,
76 bounds::repeatableBounding::CLAMP,
77 true // Failsafe behaviour
78 )
79 ),
80 interpolationScheme_
81 (
82 dict.getOrDefault<word>
83 (
84 "interpolationScheme",
85 "linear",
86 keyType::LITERAL
87 )
88 ),
89 table_(),
90 tableSamplesPtr_(nullptr),
91 interpolatorPtr_(nullptr)
92{}
93
94
95template<class Type>
97:
98 Function1<Type>(tbl),
99 bounding_(tbl.bounding_),
100 interpolationScheme_(tbl.interpolationScheme_),
101 table_(tbl.table_),
102 tableSamplesPtr_(nullptr),
103 interpolatorPtr_(nullptr)
104{}
105
106
107// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
108
109template<class Type>
111{
112 if (!table_.size())
113 {
115 << "Table for entry " << this->name() << " is invalid (empty)"
116 << nl << exit(FatalError);
117 }
118
119 scalar prevValue(0);
120
121 label i = 0;
122 for (const auto& item : table_)
123 {
124 const scalar& currValue = item.first();
125
126 // Avoid duplicate values (divide-by-zero error)
127 if (i && currValue <= prevValue)
128 {
130 << "out-of-order value: "
131 << currValue << " at index " << i << nl
132 << exit(FatalError);
133 }
134 prevValue = currValue;
135 ++i;
136 }
137}
138
139
140template<class Type>
142(
143 const scalar x,
144 scalar& xDash
145) const
146{
147 const scalar minLimit = table_.first().first();
148 const scalar maxLimit = table_.last().first();
149
150 if (x < minLimit)
151 {
152 switch (bounding_)
153 {
155 {
157 << "value (" << x << ") less than lower "
158 << "bound (" << minLimit << ")" << nl
159 << exit(FatalError);
160 break;
161 }
163 {
165 << "value (" << x << ") less than lower "
166 << "bound (" << minLimit << ")" << nl
167 << " Continuing with the first entry" << endl;
168
169 // Behaviour as per CLAMP
170 xDash = minLimit;
171 return true;
172 break;
173 }
175 {
176 xDash = minLimit;
177 return true;
178 break;
179 }
181 {
182 // Adjust x to >= minX
183 const scalar span = maxLimit - minLimit;
184 xDash = fmod(x - minLimit, span) + minLimit;
185 break;
186 }
187 }
188 }
189 else
190 {
191 xDash = x;
192 }
193
194 return false;
195}
196
197
198template<class Type>
200(
201 const scalar x,
202 scalar& xDash
203) const
204{
205 const scalar minLimit = table_.first().first();
206 const scalar maxLimit = table_.last().first();
207
208 if (x > maxLimit)
209 {
210 switch (bounding_)
211 {
213 {
215 << "value (" << x << ") greater than upper "
216 << "bound (" << maxLimit << ")" << nl
217 << exit(FatalError);
218 break;
219 }
221 {
223 << "value (" << x << ") greater than upper "
224 << "bound (" << maxLimit << ")" << nl
225 << " Continuing with the last entry" << endl;
226
227 // Behaviour as per CLAMP
228 xDash = maxLimit;
229 return true;
230 break;
231 }
233 {
234 xDash = maxLimit;
235 return true;
236 break;
237 }
239 {
240 // Adjust x to >= minX
241 const scalar span = maxLimit - minLimit;
242 xDash = fmod(x - minLimit, span) + minLimit;
243 break;
244 }
245 }
246 }
247 else
248 {
249 xDash = x;
250 }
251
252 return false;
253}
254
255
256template<class Type>
258{
259 for (auto& item : table_)
260 {
261 item.first() = t.userTimeToTime(item.first());
262 }
263
264 tableSamplesPtr_.clear();
265 interpolatorPtr_.clear();
266}
267
268
269template<class Type>
271{
272 scalar xDash = x;
273
274 if (checkMinBounds(x, xDash))
275 {
276 return table_.first().second();
277 }
278
279 if (checkMaxBounds(xDash, xDash))
280 {
281 return table_.last().second();
282 }
283
284 // Use interpolator
285 interpolator().valueWeights(xDash, currentIndices_, currentWeights_);
286
287 Type t(currentWeights_[0]*table_[currentIndices_[0]].second());
288 for (label i = 1; i < currentIndices_.size(); i++)
289 {
290 t += currentWeights_[i]*table_[currentIndices_[i]].second();
291 }
292
293 return t;
294}
295
296
297template<class Type>
299(
300 const scalar x1,
301 const scalar x2
302) const
303{
304 // Use interpolator
305 interpolator().integrationWeights(x1, x2, currentIndices_, currentWeights_);
306
307 Type sum(currentWeights_[0]*table_[currentIndices_[0]].second());
308 for (label i = 1; i < currentIndices_.size(); i++)
309 {
310 sum += currentWeights_[i]*table_[currentIndices_[i]].second();
311 }
312
313 return sum;
314}
315
316
317template<class Type>
319{
320 auto tfld = tmp<scalarField>::New(table_.size(), Zero);
321 auto& fld = tfld.ref();
322
323 forAll(table_, i)
324 {
325 fld[i] = table_[i].first();
326 }
327
328 return tfld;
329}
330
331
332template<class Type>
334{
335 auto tfld = tmp<Field<Type>>::New(table_.size(), Zero);
336 auto& fld = tfld.ref();
337
338 forAll(table_, i)
339 {
340 fld[i] = table_[i].second();
341 }
342
343 return tfld;
344}
345
346
347template<class Type>
349{
350 if (bounds::repeatableBounding::CLAMP != bounding_)
351 {
352 os.writeEntry
353 (
354 "outOfBounds",
356 );
357 }
358
359 os.writeEntryIfDifferent<word>
360 (
361 "interpolationScheme",
362 "linear",
363 interpolationScheme_
364 );
365}
366
367
368template<class Type>
370{
372
373 os << nl << indent << table_;
374 os.endEntry();
375
376 writeEntries(os);
377}
378
379
380// ************************************************************************* //
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Base class for table with bounds handling, interpolation and integration.
Definition: TableBase.H:63
void initialise()
Check the table for size and consistency.
Definition: TableBase.C:110
const interpolationWeights & interpolator() const
Return (demand driven) interpolator.
Definition: TableBase.C:37
bool checkMinBounds(const scalar x, scalar &xDash) const
Check minimum table bounds.
Definition: TableBase.C:142
virtual tmp< scalarField > x() const
Return the reference values.
Definition: TableBase.C:318
bool checkMaxBounds(const scalar x, scalar &xDash) const
Check maximum table bounds.
Definition: TableBase.C:200
virtual void writeEntries(Ostream &os) const
Write keywords only in dictionary format.
Definition: TableBase.C:348
virtual tmp< Field< Type > > y() const
Return the dependent values.
Definition: TableBase.C:333
virtual Type integrate(const scalar x1, const scalar x2) const
Integrate between two (scalar) values.
Definition: TableBase.C:299
virtual void userTimeToTime(const Time &t)
Convert time.
Definition: TableBase.C:257
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:96
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
virtual scalar userTimeToTime(const scalar theta) const
Convert the user-time (e.g. CA deg) to real-time (s).
Definition: TimeState.C:49
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Abstract base class for interpolating in 1D.
A class for handling keywords in dictionaries.
Definition: keyType.H:71
Registry of regIOobjects.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
#define WarningInFunction
Report a warning using Foam::Warning.
@ WARN
Issue warning and clamp value (this is a good default)
@ REPEAT
Treat as a repeating list.
@ ERROR
Exit with a FatalError.
@ CLAMP
Clamp value to the start/end value.
const Foam::Enum< repeatableBounding > repeatableBoundingNames
Strings corresponding to the repeatableBounding.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:342
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
scalarField samples(nIntervals, Zero)