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-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 
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 )
67 :
69  bounding_
70  (
72  (
73  "outOfBounds",
74  dict,
75  bounds::repeatableBounding::CLAMP,
76  true // Failsafe behaviour
77  )
78  ),
79  interpolationScheme_
80  (
81  dict.getOrDefault<word>("interpolationScheme", "linear")
82  ),
83  table_(),
84  tableSamplesPtr_(nullptr),
85  interpolatorPtr_(nullptr)
86 {}
87 
88 
89 template<class Type>
91 :
92  Function1<Type>(tbl),
93  bounding_(tbl.bounding_),
94  interpolationScheme_(tbl.interpolationScheme_),
95  table_(tbl.table_),
96  tableSamplesPtr_(nullptr),
97  interpolatorPtr_(nullptr)
98 {}
99 
100 
101 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
102 
103 template<class Type>
105 {}
106 
107 
108 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
109 
110 template<class Type>
112 {
113  if (!table_.size())
114  {
116  << "Table for entry " << this->name() << " is invalid (empty)"
117  << nl << exit(FatalError);
118  }
119 
120  scalar prevValue(0);
121 
122  label i = 0;
123  for (const auto& item : table_)
124  {
125  const scalar& currValue = item.first();
126 
127  // Avoid duplicate values (divide-by-zero error)
128  if (i && currValue <= prevValue)
129  {
131  << "out-of-order value: "
132  << currValue << " at index " << i << nl
133  << exit(FatalError);
134  }
135  prevValue = currValue;
136  ++i;
137  }
138 }
139 
140 
141 template<class Type>
143 (
144  const scalar x,
145  scalar& xDash
146 ) const
147 {
148  const scalar minLimit = table_.first().first();
149  const scalar maxLimit = table_.last().first();
150 
151  if (x < minLimit)
152  {
153  switch (bounding_)
154  {
156  {
158  << "value (" << x << ") less than lower "
159  << "bound (" << minLimit << ")" << nl
160  << exit(FatalError);
161  break;
162  }
164  {
166  << "value (" << x << ") less than lower "
167  << "bound (" << minLimit << ")" << nl
168  << " Continuing with the first entry" << endl;
169 
170  // Behaviour as per CLAMP
171  xDash = minLimit;
172  return true;
173  break;
174  }
176  {
177  xDash = minLimit;
178  return true;
179  break;
180  }
182  {
183  // Adjust x to >= minX
184  const scalar span = maxLimit - minLimit;
185  xDash = fmod(x - minLimit, span) + minLimit;
186  break;
187  }
188  }
189  }
190  else
191  {
192  xDash = x;
193  }
194 
195  return false;
196 }
197 
198 
199 template<class Type>
201 (
202  const scalar x,
203  scalar& xDash
204 ) const
205 {
206  const scalar minLimit = table_.first().first();
207  const scalar maxLimit = table_.last().first();
208 
209  if (x > maxLimit)
210  {
211  switch (bounding_)
212  {
214  {
216  << "value (" << x << ") greater than upper "
217  << "bound (" << maxLimit << ")" << nl
218  << exit(FatalError);
219  break;
220  }
222  {
224  << "value (" << x << ") greater than upper "
225  << "bound (" << maxLimit << ")" << nl
226  << " Continuing with the last entry" << endl;
227 
228  // Behaviour as per CLAMP
229  xDash = maxLimit;
230  return true;
231  break;
232  }
234  {
235  xDash = maxLimit;
236  return true;
237  break;
238  }
240  {
241  // Adjust x to >= minX
242  const scalar span = maxLimit - minLimit;
243  xDash = fmod(x - minLimit, span) + minLimit;
244  break;
245  }
246  }
247  }
248  else
249  {
250  xDash = x;
251  }
252 
253  return false;
254 }
255 
256 
257 template<class Type>
259 {
260  for (auto& item : table_)
261  {
262  item.first() = t.userTimeToTime(item.first());
263  }
264 
265  tableSamplesPtr_.clear();
266  interpolatorPtr_.clear();
267 }
268 
269 
270 template<class Type>
272 {
273  scalar xDash = x;
274 
275  if (checkMinBounds(x, xDash))
276  {
277  return table_.first().second();
278  }
279 
280  if (checkMaxBounds(xDash, xDash))
281  {
282  return table_.last().second();
283  }
284 
285  // Use interpolator
286  interpolator().valueWeights(xDash, currentIndices_, currentWeights_);
287 
288  Type t = currentWeights_[0]*table_[currentIndices_[0]].second();
289  for (label i = 1; i < currentIndices_.size(); i++)
290  {
291  t += currentWeights_[i]*table_[currentIndices_[i]].second();
292  }
293 
294  return t;
295 }
296 
297 
298 template<class Type>
300 (
301  const scalar x1,
302  const scalar x2
303 ) const
304 {
305  // Use interpolator
306  interpolator().integrationWeights(x1, x2, currentIndices_, currentWeights_);
307 
308  Type sum = currentWeights_[0]*table_[currentIndices_[0]].second();
309  for (label i = 1; i < currentIndices_.size(); i++)
310  {
311  sum += currentWeights_[i]*table_[currentIndices_[i]].second();
312  }
313 
314  return sum;
315 }
316 
317 
318 template<class Type>
320 {
321  auto tfld = tmp<scalarField>::New(table_.size(), Zero);
322  auto& fld = tfld.ref();
323 
324  forAll(table_, i)
325  {
326  fld[i] = table_[i].first();
327  }
328 
329  return tfld;
330 }
331 
332 
333 template<class Type>
335 {
336  auto tfld = tmp<Field<Type>>::New(table_.size(), Zero);
337  auto& fld = tfld.ref();
338 
339  forAll(table_, i)
340  {
341  fld[i] = table_[i].second();
342  }
343 
344  return tfld;
345 }
346 
347 
348 template<class Type>
350 {
351  if (bounds::repeatableBounding::CLAMP != bounding_)
352  {
353  os.writeEntry
354  (
355  "outOfBounds",
357  );
358  }
359 
361  (
362  "interpolationScheme",
363  "linear",
364  interpolationScheme_
365  );
366 }
367 
368 
369 template<class Type>
371 {
373  os << nl << indent << table_ << token::END_STATEMENT << nl;
374  writeEntries(os);
375 }
376 
377 
378 // ************************************************************************* //
Foam::Ostream::writeEntryIfDifferent
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:244
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:201
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
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:350
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:300
Foam::Function1
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:86
Foam::Function1Types::TableBase::value
virtual Type value(const scalar x) const
Return Table value.
Definition: TableBase.C:271
Foam::Function1Types::TableBase::convertTimeBase
virtual void convertTimeBase(const Time &t)
Convert time.
Definition: TableBase.C:258
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
TableBase.H
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::Function1Types::TableBase::x
virtual tmp< scalarField > x() const
Return the reference values.
Definition: TableBase.C:319
Foam::Function1Types::TableBase::check
void check() const
Check the table for size and consistency.
Definition: TableBase.C:111
Foam::Function1::writeData
virtual void writeData(Ostream &os) const
Write in dictionary format.
Definition: Function1.C:166
samples
scalarField samples(nIntervals, Zero)
Foam::token::END_STATEMENT
End entry [isseparator].
Definition: token.H:121
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:121
Foam::interpolationWeights
Abstract base class for interpolating in 1D.
Definition: interpolationWeights.H:58
Foam::Function1Types::TableBase::~TableBase
virtual ~TableBase()
Destructor.
Definition: TableBase.C:104
Foam::indent
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:320
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:381
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::Function1Types::TableBase::writeData
virtual void writeData(Ostream &os) const
Write all table data in dictionary format.
Definition: TableBase.C:370
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:349
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:232
Foam::Function1Types::TableBase::TableBase
TableBase(const word &name, const dictionary &dict)
Construct from dictionary - note table is not populated.
Definition: TableBase.C:63
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:334
Foam::Function1Types::TableBase::checkMinBounds
bool checkMinBounds(const scalar x, scalar &xDash) const
Check minimum table bounds.
Definition: TableBase.C:143
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:303