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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "TableBase.H"
27 #include "Time.H"
28 #include "interpolationWeights.H"
29 
30 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
31 
32 template<class Type>
34 {
35  if (interpolatorPtr_.empty())
36  {
37  // Re-work table into linear list
38  tableSamplesPtr_.reset(new scalarField(table_.size()));
39  scalarField& tableSamples = tableSamplesPtr_();
40  forAll(table_, i)
41  {
42  tableSamples[i] = table_[i].first();
43  }
44  interpolatorPtr_ = interpolationWeights::New
45  (
46  interpolationScheme_,
47  tableSamples
48  );
49  }
50 
51  return interpolatorPtr_();
52 }
53 
54 
55 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
56 
57 template<class Type>
59 :
60  DataEntry<Type>(name),
61  name_(name),
62  boundsHandling_
63  (
64  wordToBoundsHandling
65  (
66  dict.lookupOrDefault<word>("outOfBounds", "clamp")
67  )
68  ),
69  interpolationScheme_
70  (
71  dict.lookupOrDefault<word>("interpolationScheme", "linear")
72  ),
73  table_(),
74  dimensions_(dimless)
75 {}
76 
77 
78 template<class Type>
80 :
81  DataEntry<Type>(tbl),
82  name_(tbl.name_),
83  boundsHandling_(tbl.boundsHandling_),
84  interpolationScheme_(tbl.interpolationScheme_),
85  table_(tbl.table_),
86  dimensions_(tbl.dimensions_),
87  tableSamplesPtr_(tbl.tableSamplesPtr_),
88  interpolatorPtr_(tbl.interpolatorPtr_)
89 {}
90 
91 
92 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
93 
94 template<class Type>
96 {}
97 
98 
99 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
100 
101 template<class Type>
103 (
104  const boundsHandling& bound
105 ) const
106 {
107  word enumName("warn");
108 
109  switch (bound)
110  {
111  case ERROR:
112  {
113  enumName = "error";
114  break;
115  }
116  case WARN:
117  {
118  enumName = "warn";
119  break;
120  }
121  case CLAMP:
122  {
123  enumName = "clamp";
124  break;
125  }
126  case REPEAT:
127  {
128  enumName = "repeat";
129  break;
130  }
131  }
132 
133  return enumName;
134 }
135 
136 
137 template<class Type>
140 (
141  const word& bound
142 ) const
143 {
144  if (bound == "error")
145  {
146  return ERROR;
147  }
148  else if (bound == "warn")
149  {
150  return WARN;
151  }
152  else if (bound == "clamp")
153  {
154  return CLAMP;
155  }
156  else if (bound == "repeat")
157  {
158  return REPEAT;
159  }
160  else
161  {
163  << "bad outOfBounds specifier " << bound << " using 'warn'"
164  << endl;
165 
166  return WARN;
167  }
168 }
169 
170 
171 template<class Type>
174 (
175  const boundsHandling& bound
176 )
177 {
178  boundsHandling prev = boundsHandling_;
179  boundsHandling_ = bound;
180 
181  return prev;
182 }
183 
184 
185 template<class Type>
187 {
188  if (!table_.size())
189  {
191  << "Table for entry " << this->name_ << " is invalid (empty)"
192  << nl << exit(FatalError);
193  }
194 
195  label n = table_.size();
196  scalar prevValue = table_[0].first();
197 
198  for (label i = 1; i < n; ++i)
199  {
200  const scalar currValue = table_[i].first();
201 
202  // avoid duplicate values (divide-by-zero error)
203  if (currValue <= prevValue)
204  {
206  << "out-of-order value: " << currValue << " at index " << i
207  << exit(FatalError);
208  }
209  prevValue = currValue;
210  }
211 }
212 
213 
214 template<class Type>
216 (
217  const scalar x,
218  scalar& xDash
219 ) const
220 {
221  if (x < table_[0].first())
222  {
223  switch (boundsHandling_)
224  {
225  case ERROR:
226  {
228  << "value (" << x << ") underflow"
229  << exit(FatalError);
230  break;
231  }
232  case WARN:
233  {
235  << "value (" << x << ") underflow" << nl
236  << endl;
237 
238  // fall-through to 'CLAMP'
239  }
240  case CLAMP:
241  {
242  xDash = table_[0].first();
243  return true;
244  break;
245  }
246  case REPEAT:
247  {
248  // adjust x to >= minX
249  scalar span = table_.last().first() - table_[0].first();
250  xDash = fmod(x - table_[0].first(), span) + table_[0].first();
251  break;
252  }
253  }
254  }
255  else
256  {
257  xDash = x;
258  }
259 
260  return false;
261 }
262 
263 
264 template<class Type>
266 (
267  const scalar x,
268  scalar& xDash
269 ) const
270 {
271  if (x > table_.last().first())
272  {
273  switch (boundsHandling_)
274  {
275  case ERROR:
276  {
278  << "value (" << x << ") overflow"
279  << exit(FatalError);
280  break;
281  }
282  case WARN:
283  {
285  << "value (" << x << ") overflow" << nl
286  << endl;
287 
288  // fall-through to 'CLAMP'
289  }
290  case CLAMP:
291  {
292  xDash = table_.last().first();
293  return true;
294  break;
295  }
296  case REPEAT:
297  {
298  // adjust x to >= minX
299  scalar span = table_.last().first() - table_[0].first();
300  xDash = fmod(x - table_[0].first(), span) + table_[0].first();
301  break;
302  }
303  }
304  }
305  else
306  {
307  xDash = x;
308  }
309 
310  return false;
311 }
312 
313 
314 template<class Type>
316 {
317  forAll(table_, i)
318  {
319  scalar value = table_[i].first();
320  table_[i].first() = t.userTimeToTime(value);
321  }
322 
323  tableSamplesPtr_.clear();
324  interpolatorPtr_.clear();
325 }
326 
327 
328 template<class Type>
329 Type Foam::TableBase<Type>::value(const scalar x) const
330 {
331  scalar xDash = x;
332 
333  if (checkMinBounds(x, xDash))
334  {
335  return table_[0].second();
336  }
337 
338  if (checkMaxBounds(xDash, xDash))
339  {
340  return table_.last().second();
341  }
342 
343  // Use interpolator
344  interpolator().valueWeights(xDash, currentIndices_, currentWeights_);
345 
346  Type t = currentWeights_[0]*table_[currentIndices_[0]].second();
347  for (label i = 1; i < currentIndices_.size(); i++)
348  {
349  t += currentWeights_[i]*table_[currentIndices_[i]].second();
350  }
351 
352  return t;
353 }
354 
355 
356 template<class Type>
357 Type Foam::TableBase<Type>::integrate(const scalar x1, const scalar x2) const
358 {
359  // Use interpolator
360  interpolator().integrationWeights(x1, x2, currentIndices_, currentWeights_);
361 
362  Type sum = currentWeights_[0]*table_[currentIndices_[0]].second();
363  for (label i = 1; i < currentIndices_.size(); i++)
364  {
365  sum += currentWeights_[i]*table_[currentIndices_[i]].second();
366  }
367 
368  return sum;
369 }
370 
371 
372 template<class Type>
374 dimValue(const scalar x) const
375 {
376  return dimensioned<Type>("dimensionedValue", dimensions_, this->value(x));
377 }
378 
379 
380 template<class Type>
382 (
383  const scalar x1, const scalar x2
384 ) const
385 {
386  return dimensioned<Type>
387  (
388  "dimensionedValue",
389  dimensions_,
390  this->integrate(x2, x1)
391  );
392 }
393 
394 
395 template<class Type>
397 {
398  tmp<scalarField> tfld(new scalarField(table_.size(), 0.0));
399  scalarField& fld = tfld();
400 
401  forAll(table_, i)
402  {
403  fld[i] = table_[i].first();
404  }
405 
406  return tfld;
407 }
408 
409 
410 template<class Type>
412 {
413  tmp<Field<Type> > tfld(new Field<Type>(table_.size(), pTraits<Type>::zero));
414  Field<Type>& fld = tfld();
415 
416  forAll(table_, i)
417  {
418  fld[i] = table_[i].second();
419  }
420 
421  return tfld;
422 }
423 
424 
425 // * * * * * * * * * * * * * * IOStream operators * * * * * * * * * * * * * //
426 
427 #include "TableBaseIO.C"
428 
429 // ************************************************************************* //
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:48
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Foam::TableBase::x
virtual tmp< scalarField > x() const
Return the reference values.
Definition: TableBase.C:396
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::TableBase::y
virtual tmp< Field< Type > > y() const
Return the dependent values.
Definition: TableBase.C:411
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::TableBase::dimValue
virtual dimensioned< Type > dimValue(const scalar x) const
Return dimensioned constant value.
Definition: TableBase.C:374
Foam::bound
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
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:55
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::TableBase::boundsHandling
boundsHandling
Enumeration for handling out-of-bound values.
Definition: TableBase.H:72
Foam::TableBase::checkMinBounds
bool checkMinBounds(const scalar x, scalar &xDash) const
Check minimum table bounds.
Definition: TableBase.C:216
Foam::TableBase::integrate
virtual Type integrate(const scalar x1, const scalar x2) const
Integrate between two (scalar) values.
Definition: TableBase.C:357
n
label n
Definition: TABSMDCalcMethod2.H:31
TableBase.H
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::nl
static const char nl
Definition: Ostream.H:260
interpolationWeights.H
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::TableBase::TableBase
TableBase(const word &name, const dictionary &dict)
Construct from dictionary - note table is not populated.
Definition: TableBase.C:58
Foam::TableBase::~TableBase
virtual ~TableBase()
Destructor.
Definition: TableBase.C:95
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::interpolationWeights
Abstract base class for interpolating in 1D.
Definition: interpolationWeights.H:55
fld
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::dimensioned< Type >
Foam::TableBase::checkMaxBounds
bool checkMaxBounds(const scalar x, scalar &xDash) const
Check maximum table bounds.
Definition: TableBase.C:266
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::TableBase::outOfBounds
boundsHandling outOfBounds(const boundsHandling &bound)
Set the out-of-bounds handling from enum, return previous setting.
Definition: TableBase.C:174
TableBaseIO.C
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:50
Foam::TableBase::boundsHandlingToWord
word boundsHandlingToWord(const boundsHandling &bound) const
Return the out-of-bounds handling as a word.
Definition: TableBase.C:103
scalarField
volScalarField scalarField(fieldObject, mesh)
Foam::TableBase::value
virtual Type value(const scalar x) const
Return Table value.
Definition: TableBase.C:329
Foam::TableBase::check
void check() const
Check the table for size and consistency.
Definition: TableBase.C:186
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:333
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::TableBase::wordToBoundsHandling
boundsHandling wordToBoundsHandling(const word &bound) const
Return the out-of-bounds handling as an enumeration.
Definition: TableBase.C:140
Foam::TableBase::convertTimeBase
virtual void convertTimeBase(const Time &t)
Convert time.
Definition: TableBase.C:315
Foam::TableBase::dimIntegrate
virtual dimensioned< Type > dimIntegrate(const scalar x1, const scalar x2) const
Integrate between two values and return dimensioned type.
Definition: TableBase.C:382
Foam::TableBase::interpolator
const interpolationWeights & interpolator() const
Return (demand driven) interpolator.
Definition: TableBase.C:33
Foam::DataEntry
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: DataEntry.H:52
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::TableBase
Base class for table with bounds handling, interpolation and integration.
Definition: TableBase.H:47