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-2016 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>
35 {
36  if (interpolatorPtr_.empty())
37  {
38  // Re-work table into linear list
39  tableSamplesPtr_.reset(new scalarField(table_.size()));
40  scalarField& tableSamples = tableSamplesPtr_();
41  forAll(table_, i)
42  {
43  tableSamples[i] = table_[i].first();
44  }
45  interpolatorPtr_ = interpolationWeights::New
46  (
47  interpolationScheme_,
48  tableSamples
49  );
50  }
51 
52  return interpolatorPtr_();
53 }
54 
55 
56 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
57 
58 template<class Type>
60 (
61  const word& name,
62  const dictionary& dict
63 )
64 :
66  name_(name),
67  boundsHandling_
68  (
69  wordToBoundsHandling
70  (
71  dict.lookupOrDefault<word>("outOfBounds", "clamp")
72  )
73  ),
74  interpolationScheme_
75  (
76  dict.lookupOrDefault<word>("interpolationScheme", "linear")
77  ),
78  table_()
79 {}
80 
81 
82 template<class Type>
84 :
85  Function1<Type>(tbl),
86  name_(tbl.name_),
87  boundsHandling_(tbl.boundsHandling_),
88  interpolationScheme_(tbl.interpolationScheme_),
89  table_(tbl.table_),
90  tableSamplesPtr_(tbl.tableSamplesPtr_),
91  interpolatorPtr_(tbl.interpolatorPtr_)
92 {}
93 
94 
95 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
96 
97 template<class Type>
99 {}
100 
101 
102 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103 
104 template<class Type>
106 (
107  const boundsHandling& bound
108 ) const
109 {
110  word enumName("warn");
111 
112  switch (bound)
113  {
114  case ERROR:
115  {
116  enumName = "error";
117  break;
118  }
119  case WARN:
120  {
121  enumName = "warn";
122  break;
123  }
124  case CLAMP:
125  {
126  enumName = "clamp";
127  break;
128  }
129  case REPEAT:
130  {
131  enumName = "repeat";
132  break;
133  }
134  }
135 
136  return enumName;
137 }
138 
139 
140 template<class Type>
143 (
144  const word& bound
145 ) const
146 {
147  if (bound == "error")
148  {
149  return ERROR;
150  }
151  else if (bound == "warn")
152  {
153  return WARN;
154  }
155  else if (bound == "clamp")
156  {
157  return CLAMP;
158  }
159  else if (bound == "repeat")
160  {
161  return REPEAT;
162  }
163  else
164  {
166  << "bad outOfBounds specifier " << bound << " using 'warn'"
167  << endl;
168 
169  return WARN;
170  }
171 }
172 
173 
174 template<class Type>
177 (
178  const boundsHandling& bound
179 )
180 {
183 
184  return prev;
185 }
186 
187 
188 template<class Type>
190 {
191  if (!table_.size())
192  {
194  << "Table for entry " << this->name_ << " is invalid (empty)"
195  << nl << exit(FatalError);
196  }
197 
198  label n = table_.size();
199  scalar prevValue = table_[0].first();
200 
201  for (label i = 1; i < n; ++i)
202  {
203  const scalar currValue = table_[i].first();
204 
205  // avoid duplicate values (divide-by-zero error)
206  if (currValue <= prevValue)
207  {
209  << "out-of-order value: " << currValue << " at index " << i
210  << exit(FatalError);
211  }
212  prevValue = currValue;
213  }
214 }
215 
216 
217 template<class Type>
219 (
220  const scalar x,
221  scalar& xDash
222 ) const
223 {
224  if (x < table_[0].first())
225  {
226  switch (boundsHandling_)
227  {
228  case ERROR:
229  {
231  << "value (" << x << ") underflow"
232  << exit(FatalError);
233  break;
234  }
235  case WARN:
236  {
238  << "value (" << x << ") underflow" << nl
239  << endl;
240 
241  // fall-through to 'CLAMP'
242  }
243  case CLAMP:
244  {
245  xDash = table_[0].first();
246  return true;
247  break;
248  }
249  case REPEAT:
250  {
251  // adjust x to >= minX
252  scalar span = table_.last().first() - table_[0].first();
253  xDash = fmod(x - table_[0].first(), span) + table_[0].first();
254  break;
255  }
256  }
257  }
258  else
259  {
260  xDash = x;
261  }
262 
263  return false;
264 }
265 
266 
267 template<class Type>
269 (
270  const scalar x,
271  scalar& xDash
272 ) const
273 {
274  if (x > table_.last().first())
275  {
276  switch (boundsHandling_)
277  {
278  case ERROR:
279  {
281  << "value (" << x << ") overflow"
282  << exit(FatalError);
283  break;
284  }
285  case WARN:
286  {
288  << "value (" << x << ") overflow" << nl
289  << endl;
290 
291  // fall-through to 'CLAMP'
292  }
293  case CLAMP:
294  {
295  xDash = table_.last().first();
296  return true;
297  break;
298  }
299  case REPEAT:
300  {
301  // adjust x to >= minX
302  scalar span = table_.last().first() - table_[0].first();
303  xDash = fmod(x - table_[0].first(), span) + table_[0].first();
304  break;
305  }
306  }
307  }
308  else
309  {
310  xDash = x;
311  }
312 
313  return false;
314 }
315 
316 
317 template<class Type>
319 {
320  forAll(table_, i)
321  {
322  scalar value = table_[i].first();
323  table_[i].first() = t.userTimeToTime(value);
324  }
325 
326  tableSamplesPtr_.clear();
327  interpolatorPtr_.clear();
328 }
329 
330 
331 template<class Type>
333 {
334  scalar xDash = x;
335 
336  if (checkMinBounds(x, xDash))
337  {
338  return table_[0].second();
339  }
340 
341  if (checkMaxBounds(xDash, xDash))
342  {
343  return table_.last().second();
344  }
345 
346  // Use interpolator
348 
349  Type t = currentWeights_[0]*table_[currentIndices_[0]].second();
350  for (label i = 1; i < currentIndices_.size(); i++)
351  {
352  t += currentWeights_[i]*table_[currentIndices_[i]].second();
353  }
354 
355  return t;
356 }
357 
358 
359 template<class Type>
361 (
362  const scalar x1,
363  const scalar x2
364 ) const
365 {
366  // Use interpolator
368 
369  Type sum = currentWeights_[0]*table_[currentIndices_[0]].second();
370  for (label i = 1; i < currentIndices_.size(); i++)
371  {
372  sum += currentWeights_[i]*table_[currentIndices_[i]].second();
373  }
374 
375  return sum;
376 }
377 
378 
379 template<class Type>
381 {
382  tmp<scalarField> tfld(new scalarField(table_.size(), 0.0));
383  scalarField& fld = tfld.ref();
384 
385  forAll(table_, i)
386  {
387  fld[i] = table_[i].first();
388  }
389 
390  return tfld;
391 }
392 
393 
394 template<class Type>
396 {
397  tmp<Field<Type>> tfld(new Field<Type>(table_.size(), Zero));
398  Field<Type>& fld = tfld.ref();
399 
400  forAll(table_, i)
401  {
402  fld[i] = table_[i].second();
403  }
404 
405  return tfld;
406 }
407 
408 
409 template<class Type>
411 {
412  if (boundsHandling_ != CLAMP)
413  {
415  << token::END_STATEMENT << nl;
416  }
417  if (interpolationScheme_ != "linear")
418  {
419  os.writeKeyword("interpolationScheme") << interpolationScheme_
420  << token::END_STATEMENT << nl;
421  }
422 }
423 
424 
425 template<class Type>
427 {
429  os << nl << indent << table_ << token::END_STATEMENT << nl;
430  writeEntries(os);
431 }
432 
433 
434 // ************************************************************************* //
Issue warning and clamp value (default)
Definition: TableBase.H:67
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:53
virtual Type value(const scalar x) const
Return Table value.
Definition: TableBase.C:332
bool checkMinBounds(const scalar x, scalar &xDash) const
Check minimum table bounds.
Definition: TableBase.C:219
boundsHandling wordToBoundsHandling(const word &bound) const
Return the out-of-bounds handling as an enumeration.
Definition: TableBase.C:143
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
boundsHandling outOfBounds(const boundsHandling &bound)
Set the out-of-bounds handling from enum, return previous setting.
Definition: TableBase.C:177
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:223
TableBase(const word &name, const dictionary &dict)
Construct from dictionary - note table is not populated.
Definition: TableBase.C:60
virtual void writeEntries(Ostream &os) const
Write keywords only in dictionary format. Used for non-inline.
Definition: TableBase.C:410
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
List< Tuple2< scalar, Type > > table_
Table data.
Definition: TableBase.H:87
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
virtual scalar userTimeToTime(const scalar theta) const
Convert the user-time (e.g. CA deg) to real-time (s).
Definition: TimeState.C:52
Exit with a FatalError.
Definition: TableBase.H:66
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
autoPtr< interpolationWeights > interpolatorPtr_
Interpolator method.
Definition: TableBase.H:93
virtual void writeData(Ostream &os) const
Write in dictionary format.
Definition: Function1.C:121
T & first()
Return the first element of the list.
Definition: UListI.H:114
Abstract base class for interpolating in 1D.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
word boundsHandlingToWord(const boundsHandling &bound) const
Return the out-of-bounds handling as a word.
Definition: TableBase.C:106
virtual ~TableBase()
Destructor.
Definition: TableBase.C:98
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Treat as a repeating list.
Definition: TableBase.H:69
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))
Base class for table with bounds handling, interpolation and integration.
Definition: TableBase.H:55
bool checkMaxBounds(const scalar x, scalar &xDash) const
Check maximum table bounds.
Definition: TableBase.C:269
const word name_
Table name.
Definition: TableBase.H:78
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const interpolationWeights & interpolator() const
Return (demand driven) interpolator.
Definition: TableBase.C:34
static const zero Zero
Definition: zero.H:91
virtual tmp< scalarField > x() const
Return the reference values.
Definition: TableBase.C:380
virtual bool integrationWeights(const scalar t1, const scalar t2, labelList &indices, scalarField &weights) const =0
Calculate weights and indices to calculate integrand of t1..t2.
virtual Type integrate(const scalar x1, const scalar x2) const
Integrate between two (scalar) values.
Definition: TableBase.C:361
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
volScalarField scalarField(fieldObject, mesh)
static const char nl
Definition: Ostream.H:262
void check() const
Check the table for size and consistency.
Definition: TableBase.C:189
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
virtual tmp< Field< Type > > y() const
Return the dependent values.
Definition: TableBase.C:395
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const word interpolationScheme_
Interpolation type.
Definition: TableBase.H:84
virtual void convertTimeBase(const Time &t)
Convert time.
Definition: TableBase.C:318
autoPtr< scalarField > tableSamplesPtr_
Extracted values.
Definition: TableBase.H:90
const boundsHandling boundsHandling_
Enumeration for handling out-of-bound values.
Definition: TableBase.H:81
labelList currentIndices_
Cached indices and weights.
Definition: TableBase.H:96
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
#define WarningInFunction
Report a warning using Foam::Warning.
virtual void writeData(Ostream &os) const
Write all table data in dictionary format.
Definition: TableBase.C:426
label n
A class for managing temporary objects.
Definition: PtrList.H:54
T & last()
Return the last element of the list.
Definition: UListI.H:128
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
virtual bool valueWeights(const scalar t, labelList &indices, scalarField &weights) const =0
Calculate weights and indices to calculate t from samples.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Clamp value to the start/end value.
Definition: TableBase.H:68
boundsHandling
Enumeration for handling out-of-bound values.
Definition: TableBase.H:64