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 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  name_(name),
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>("interpolationScheme", "linear")
83  ),
84  table_(),
85  tableSamplesPtr_(nullptr),
86  interpolatorPtr_(nullptr)
87 {}
88 
89 
90 template<class Type>
92 :
93  Function1<Type>(tbl),
94  name_(tbl.name_),
95  bounding_(tbl.bounding_),
96  interpolationScheme_(tbl.interpolationScheme_),
97  table_(tbl.table_),
98  tableSamplesPtr_(tbl.tableSamplesPtr_.clone()),
99  interpolatorPtr_(tbl.interpolatorPtr_) // steal/reuse (missing clone!)
100 {}
101 
102 
103 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
104 
105 template<class Type>
107 {}
108 
109 
110 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111 
112 template<class Type>
114 {
115  if (!table_.size())
116  {
118  << "Table for entry " << this->name_ << " is invalid (empty)"
119  << nl << exit(FatalError);
120  }
121 
122  scalar prevValue(0);
123 
124  label i = 0;
125  for (const auto& item : table_)
126  {
127  const scalar& currValue = item.first();
128 
129  // Avoid duplicate values (divide-by-zero error)
130  if (i && currValue <= prevValue)
131  {
133  << "out-of-order value: "
134  << currValue << " at index " << i << nl
135  << exit(FatalError);
136  }
137  prevValue = currValue;
138  ++i;
139  }
140 }
141 
142 
143 template<class Type>
145 (
146  const scalar x,
147  scalar& xDash
148 ) const
149 {
150  const scalar minLimit = table_.first().first();
151  const scalar maxLimit = table_.last().first();
152 
153  if (x < minLimit)
154  {
155  switch (bounding_)
156  {
158  {
160  << "value (" << x << ") less than lower "
161  << "bound (" << minLimit << ")" << nl
162  << exit(FatalError);
163  break;
164  }
166  {
168  << "value (" << x << ") less than lower "
169  << "bound (" << minLimit << ")" << nl
170  << " Continuing with the first entry" << endl;
171 
172  // Behaviour as per CLAMP
173  xDash = minLimit;
174  return true;
175  break;
176  }
178  {
179  xDash = minLimit;
180  return true;
181  break;
182  }
184  {
185  // Adjust x to >= minX
186  const scalar span = maxLimit - minLimit;
187  xDash = fmod(x - minLimit, span) + minLimit;
188  break;
189  }
190  }
191  }
192  else
193  {
194  xDash = x;
195  }
196 
197  return false;
198 }
199 
200 
201 template<class Type>
203 (
204  const scalar x,
205  scalar& xDash
206 ) const
207 {
208  const scalar minLimit = table_.first().first();
209  const scalar maxLimit = table_.last().first();
210 
211  if (x > maxLimit)
212  {
213  switch (bounding_)
214  {
216  {
218  << "value (" << x << ") greater than upper "
219  << "bound (" << maxLimit << ")" << nl
220  << exit(FatalError);
221  break;
222  }
224  {
226  << "value (" << x << ") greater than upper "
227  << "bound (" << maxLimit << ")" << nl
228  << " Continuing with the last entry" << endl;
229 
230  // Behaviour as per CLAMP
231  xDash = maxLimit;
232  return true;
233  break;
234  }
236  {
237  xDash = maxLimit;
238  return true;
239  break;
240  }
242  {
243  // Adjust x to >= minX
244  const scalar span = maxLimit - minLimit;
245  xDash = fmod(x - minLimit, span) + minLimit;
246  break;
247  }
248  }
249  }
250  else
251  {
252  xDash = x;
253  }
254 
255  return false;
256 }
257 
258 
259 template<class Type>
261 {
262  for (auto& item : table_)
263  {
264  item.first() = t.userTimeToTime(item.first());
265  }
266 
267  tableSamplesPtr_.clear();
268  interpolatorPtr_.clear();
269 }
270 
271 
272 template<class Type>
274 {
275  scalar xDash = x;
276 
277  if (checkMinBounds(x, xDash))
278  {
279  return table_.first().second();
280  }
281 
282  if (checkMaxBounds(xDash, xDash))
283  {
284  return table_.last().second();
285  }
286 
287  // Use interpolator
288  interpolator().valueWeights(xDash, currentIndices_, currentWeights_);
289 
290  Type t = currentWeights_[0]*table_[currentIndices_[0]].second();
291  for (label i = 1; i < currentIndices_.size(); i++)
292  {
293  t += currentWeights_[i]*table_[currentIndices_[i]].second();
294  }
295 
296  return t;
297 }
298 
299 
300 template<class Type>
302 (
303  const scalar x1,
304  const scalar x2
305 ) const
306 {
307  // Use interpolator
308  interpolator().integrationWeights(x1, x2, currentIndices_, currentWeights_);
309 
310  Type sum = currentWeights_[0]*table_[currentIndices_[0]].second();
311  for (label i = 1; i < currentIndices_.size(); i++)
312  {
313  sum += currentWeights_[i]*table_[currentIndices_[i]].second();
314  }
315 
316  return sum;
317 }
318 
319 
320 template<class Type>
322 {
323  auto tfld = tmp<scalarField>::New(table_.size(), Zero);
324  auto& fld = tfld.ref();
325 
326  forAll(table_, i)
327  {
328  fld[i] = table_[i].first();
329  }
330 
331  return tfld;
332 }
333 
334 
335 template<class Type>
337 {
338  auto tfld = tmp<Field<Type>>::New(table_.size(), Zero);
339  auto& fld = tfld.ref();
340 
341  forAll(table_, i)
342  {
343  fld[i] = table_[i].second();
344  }
345 
346  return tfld;
347 }
348 
349 
350 template<class Type>
352 {
353  if (bounds::repeatableBounding::CLAMP != bounding_)
354  {
355  os.writeEntry
356  (
357  "outOfBounds",
359  );
360  }
361 
363  (
364  "interpolationScheme",
365  "linear",
366  interpolationScheme_
367  );
368 }
369 
370 
371 template<class Type>
373 {
375  os << nl << indent << table_ << token::END_STATEMENT << nl;
376  writeEntries(os);
377 }
378 
379 
380 // ************************************************************************* //
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:203
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:59
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:302
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:273
Foam::Function1Types::TableBase::convertTimeBase
virtual void convertTimeBase(const Time &t)
Convert time.
Definition: TableBase.C:260
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:321
Foam::Function1Types::TableBase::check
void check() const
Check the table for size and consistency.
Definition: TableBase.C:113
Foam::Function1::writeData
virtual void writeData(Ostream &os) const
Write in dictionary format.
Definition: Function1.C:165
samples
scalarField samples(nIntervals, Zero)
Foam::token::END_STATEMENT
End entry [isseparator].
Definition: token.H:122
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:106
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:372
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:372
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:351
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:59
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:336
Foam::Function1Types::TableBase::checkMinBounds
bool checkMinBounds(const scalar x, scalar &xDash) const
Check minimum table bounds.
Definition: TableBase.C:145
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:298