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 -------------------------------------------------------------------------------
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 
29 #include "TableBase.H"
30 #include "Time.H"
31 #include "interpolationWeights.H"
32 
33 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
34 
35 template<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_,
51  samples
52  );
53  }
54 
55  return *interpolatorPtr_;
56 }
57 
58 
59 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
60 
61 template<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  (
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 
95 template<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 
109 template<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 
140 template<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 
198 template<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 
256 template<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 
269 template<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 
297 template<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 
317 template<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 
332 template<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 
347 template<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 
368 template<class Type>
370 {
372 
373  os << nl << indent << table_;
374  os.endEntry();
375 
376  writeEntries(os);
377 }
378 
379 
380 // ************************************************************************* //
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::Function1Types::TableBase::checkMaxBounds
bool checkMaxBounds(const scalar x, scalar &xDash) const
Check maximum table bounds.
Definition: TableBase.C:200
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::Function1Types::TableBase::interpolator
const interpolationWeights & interpolator() const
Return (demand driven) interpolator.
Definition: TableBase.C:37
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::bounds::repeatableBounding::ERROR
Exit with a FatalError.
Foam::TimeState::userTimeToTime
virtual scalar userTimeToTime(const scalar theta) const
Convert the user-time (e.g. CA deg) to real-time (s).
Definition: TimeState.C:49
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::bounds::repeatableBounding::REPEAT
Treat as a repeating list.
Foam::Function1Types::TableBase::integrate
virtual Type integrate(const scalar x1, const scalar x2) const
Integrate between two (scalar) values.
Definition: TableBase.C:299
Foam::Function1
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: propellerInfo.H:291
Foam::Function1Types::TableBase::value
virtual Type value(const scalar x) const
Return Table value.
Definition: TableBase.C:270
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
TableBase.H
Foam::Function1Types::TableBase::x
virtual tmp< scalarField > x() const
Return the reference values.
Definition: TableBase.C:318
Foam::Function1::writeData
virtual void writeData(Ostream &os) const
Write in dictionary format.
Definition: Function1.C:174
samples
scalarField samples(nIntervals, Zero)
interpolationWeights.H
fld
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputLagrangian.H:23
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::bounds::normalBounding::CLAMP
Clamp value to the start/end value.
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
Foam::interpolationWeights
Abstract base class for interpolating in 1D.
Definition: interpolationWeights.H:58
os
OBJstream os(runTime.globalPath()/outputName)
Foam::Function1Types::TableBase::initialise
void initialise()
Check the table for size and consistency.
Definition: TableBase.C:110
Foam::indent
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:339
Foam::bounds::repeatableBounding::WARN
Issue warning and clamp value (this is a good default)
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
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
Foam::Function1Types::TableBase::writeData
virtual void writeData(Ostream &os) const
Write all table data in dictionary format.
Definition: TableBase.C:369
Foam::Function1Types::TableBase::TableBase
TableBase(const word &name, const dictionary &dict, const objectRegistry *obrPtr=nullptr)
Construct from dictionary - note table is not populated.
Definition: TableBase.C:63
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::Function1Types::TableBase::writeEntries
virtual void writeEntries(Ostream &os) const
Write keywords only in dictionary format.
Definition: TableBase.C:348
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::bounds::repeatableBoundingNames
const Foam::Enum< repeatableBounding > repeatableBoundingNames
Strings corresponding to the repeatableBounding.
Foam::Function1Types::TableBase
Base class for table with bounds handling, interpolation and integration.
Definition: TableBase.H:60
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::Function1Types::TableBase::y
virtual tmp< Field< Type > > y() const
Return the dependent values.
Definition: TableBase.C:333
Foam::Function1Types::TableBase::checkMinBounds
bool checkMinBounds(const scalar x, scalar &xDash) const
Check minimum table bounds.
Definition: TableBase.C:142
Foam::Function1Types::TableBase::userTimeToTime
virtual void userTimeToTime(const Time &t)
Convert time.
Definition: TableBase.C:257
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328