Table.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2025 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 "Table.H"
28 
29 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
30 
31 template<class Type>
34 {
35  if (interpolatorPtr_.empty())
36  {
37  // Re-work table into linear list
38  tableSamplesPtr_.reset(new scalarField(values_.size()));
39  scalarField& tableSamples = tableSamplesPtr_();
40  forAll(values_, i)
41  {
42  tableSamples[i] = values_[i].first();
43  }
44  interpolatorPtr_ = interpolationWeights::New
45  (
46  interpolationScheme_,
47  tableSamples
48  );
49  }
50 
51  return interpolatorPtr_();
52 }
53 
54 
55 template<class Type>
57 {
58  if (!values_.size())
59  {
61  << "Table for entry " << this->name() << " is invalid (empty)"
62  << nl << exit(FatalError);
63  }
64 
65  label n = values_.size();
66  scalar prevValue = values_[0].first();
67 
68  for (label i = 1; i < n; ++i)
69  {
70  const scalar currValue = values_[i].first();
71 
72  // avoid duplicate values (divide-by-zero error)
73  if (currValue <= prevValue)
74  {
76  << "out-of-order value: " << currValue << " at index " << i
77  << exit(FatalError);
78  }
79  prevValue = currValue;
80  }
81 }
82 
83 
84 template<class Type>
85 void Foam::Function1s::Table<Type>::checkX(const scalar x) const
86 {
87  const bool under = x < values_.first().first();
88  const bool over = x > values_.last().first();
89 
90  auto errorMessage = [&]()
91  {
92  return "value (" + name(x) + ") " + (under ? "under" : "over") + "flow";
93  };
94 
95  if (under || over)
96  {
97  switch (boundsHandling_)
98  {
99  case tableBase::boundsHandling::error:
100  {
102  << errorMessage() << nl << exit(FatalError);
103  break;
104  }
105  case tableBase::boundsHandling::warn:
106  {
108  << errorMessage() << nl << endl;
109  break;
110  }
111  default:
112  {
113  break;
114  }
115  }
116  }
117 }
118 
119 
120 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
121 
122 template<class Type>
124 (
125  const word& name,
127  const word& interpolationScheme,
128  const autoPtr<TableReader<Type>>& reader,
129  const List<Tuple2<scalar, Type>>& table
130 )
131 :
132  FieldFunction1<Type, Table<Type>>(name),
133  boundsHandling_(boundsHandling),
134  interpolationScheme_(interpolationScheme),
135  reader_(reader, false),
136  values_(table)
137 {}
138 
139 
140 template<class Type>
142 (
143  const word& name,
144  const unitConversions& units,
145  const dictionary& dict
146 )
147 :
148  FieldFunction1<Type, Table<Type>>(name),
149  boundsHandling_
150  (
151  dict.found("outOfBounds")
152  ? tableBase::boundsHandlingNames.read(dict.lookup("outOfBounds"))
153  : tableBase::boundsHandling::clamp
154  ),
155  interpolationScheme_
156  (
157  dict.lookupOrDefault<word>
158  (
159  "interpolationScheme",
161  )
162  ),
163  reader_(TableReader<Type>::New(name, units, dict)),
164  values_(reader_->read(units, dict))
165 {
166  check();
167 }
168 
169 
170 template<class Type>
172 (
173  const word& name,
174  const unitConversion& xUnits,
175  const unitConversion& valueUnits,
176  const dictionary& dict
177 )
178 :
179  Table(name, {xUnits, valueUnits}, dict)
180 {}
181 
182 
183 template<class Type>
185 (
186  const word& name,
187  const unitConversions& units,
188  Istream& is
189 )
190 :
191  FieldFunction1<Type, Table<Type>>(name),
192  boundsHandling_(tableBase::boundsHandling::clamp),
193  interpolationScheme_(linearInterpolationWeights::typeName),
194  reader_(new TableReaders::Embedded<Type>()),
195  values_(TableReaders::Embedded<Type>().read(units, is))
196 {
197  check();
198 }
199 
200 
201 template<class Type>
203 :
204  FieldFunction1<Type, Table<Type>>(tbl),
205  boundsHandling_(tbl.boundsHandling_),
206  interpolationScheme_(tbl.interpolationScheme_),
207  reader_(tbl.reader_, false),
208  values_(tbl.values_),
209  tableSamplesPtr_(tbl.tableSamplesPtr_),
210  interpolatorPtr_(tbl.interpolatorPtr_)
211 {}
212 
213 
214 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
215 
216 template<class Type>
218 {}
219 
220 
221 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
222 
223 template<class Type>
225 (
226  const scalar xArg
227 ) const
228 {
229  scalar x(xArg);
230 
231  checkX(x);
232 
233  const scalar x0 = values_.first().first();
234  const scalar x1 = values_.last().first();
235 
236  Type y = Zero;
237 
238  switch (boundsHandling_)
239  {
240  case tableBase::boundsHandling::clamp:
241  {
242  // Integration weights clamp by default, so do nothing
243  break;
244  }
245  case tableBase::boundsHandling::zero:
246  {
247  // Return zero if outside the range
248  if (x <= x0 || x >= x1)
249  {
250  return Zero;
251  }
252  break;
253  }
255  {
256  // Shift the value into the range of the table
257  const scalar n = floor((x - x0)/(x1 - x0));
258  x -= n*(x1 - x0);
259  break;
260  }
261  default:
262  {
263  break;
264  }
265  }
266 
267  // Evaluate
268  interpolator().valueWeights(x, indices_, weights_);
269  forAll(indices_, i)
270  {
271  y += weights_[i]*values_[indices_[i]].second();
272  }
273 
274  return y;
275 }
276 
277 
278 template<class Type>
280 (
281  const scalar xAarg,
282  const scalar xBarg
283 ) const
284 {
285  scalar xA(xAarg), xB(xBarg);
286 
287  checkX(xA);
288  checkX(xB);
289 
290  const scalar x0 = values_.first().first();
291  const scalar x1 = values_.last().first();
292 
293  Type sumY = Zero;
294 
295  switch (boundsHandling_)
296  {
297  case tableBase::boundsHandling::clamp:
298  {
299  // Integration weights clamp by default, so do nothing
300  break;
301  }
302  case tableBase::boundsHandling::zero:
303  {
304  // Actually clamp to remove any integral components from the ends
305  xA = min(max(xA, x0), x1);
306  xB = min(max(xB, x0), x1);
307  break;
308  }
310  {
311  // Shift the values into the range of the table and add any
312  // repetitions necessary to the integral
313  const scalar nA = floor((xA - x0)/(x1 - x0));
314  const scalar nB = floor((xB - x0)/(x1 - x0));
315  if (nA != nB)
316  {
317  interpolator().integrationWeights(x0, x1, indices_, weights_);
318  forAll(indices_, i)
319  {
320  sumY += (nB - nA)*weights_[i]*values_[indices_[i]].second();
321  }
322  }
323  xA -= nA*(x1 - x0);
324  xB -= nB*(x1 - x0);
325  break;
326  }
327  default:
328  {
329  break;
330  }
331  }
332 
333  // Integrate between the bounds
334  interpolator().integrationWeights(xA, xB, indices_, weights_);
335  forAll(indices_, i)
336  {
337  sumY += weights_[i]*values_[indices_[i]].second();
338  }
339 
340  return sumY;
341 }
342 
343 
344 template<class Type>
347 {
348  tmp<scalarField> tfld(new scalarField(values_.size(), 0.0));
349  scalarField& fld = tfld.ref();
350 
351  forAll(values_, i)
352  {
353  fld[i] = values_[i].first();
354  }
355 
356  return tfld;
357 }
358 
359 
360 template<class Type>
363 {
364  tmp<Field<Type>> tfld(new Field<Type>(values_.size(), Zero));
365  Field<Type>& fld = tfld.ref();
366 
367  forAll(values_, i)
368  {
369  fld[i] = values_[i].second();
370  }
371 
372  return tfld;
373 }
374 
375 
376 template<class Type>
378 (
379  Ostream& os,
380  const unitConversions& units
381 ) const
382 {
384  (
385  os,
386  "outOfBounds",
387  tableBase::boundsHandlingNames[tableBase::boundsHandling::clamp],
388  tableBase::boundsHandlingNames[boundsHandling_]
389  );
390 
392  (
393  os,
394  "interpolationScheme",
395  linearInterpolationWeights::typeName,
396  interpolationScheme_
397  );
398 
399  reader_->write(os, units, values_);
400 }
401 
402 
403 // ************************************************************************* //
scalar y
bool found
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Templated interpolated tabulated data Function1.
Definition: Table.H:122
virtual Type integral(const scalar xA, const scalar xB) const
Integrate between two scalars.
Definition: Table.C:280
virtual tmp< scalarField > x() const
Return the reference values.
Definition: Table.C:346
Table(const word &name, const tableBase::boundsHandling boundsHandling, const word &interpolationScheme, const autoPtr< TableReader< Type >> &reader, const List< Tuple2< scalar, Type >> &table)
Construct from components.
Definition: Table.C:124
virtual void write(Ostream &os, const unitConversions &units) const
Write data to dictionary stream.
Definition: Table.C:378
virtual tmp< Field< Type > > y() const
Return the dependent values.
Definition: Table.C:362
virtual Type value(const scalar x) const
Return Table value as a function of scalar x.
Definition: Table.C:225
virtual ~Table()
Destructor.
Definition: Table.C:217
boundsHandling
Enumeration for handling out-of-bound values.
Definition: tableBase.H:55
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
Base class to read table data for tables.
Definition: TableReader.H:51
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:66
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Abstract base class for interpolating in 1D.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
Unit conversion structure. Contains the associated dimensions and the multiplier with which to conver...
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
volScalarField scalarField(fieldObject, mesh)
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter) { positions[particlei++]=iter().position(mesh);} for(i=0;i< pTraits< point >::nComponents;i++) { forAll(positions, particlei) { gmvFile<< component(positions[particlei], i)<< ' ';} gmvFile<< nl;}}forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
#define WarningInFunction
Report a warning using Foam::Warning.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
bool read(const char *, int32_t &)
Definition: int32IO.C:85
void writeEntryIfDifferent(Ostream &os, const word &entryName, const EntryType &value1, const EntryType &value2)
Helper function to write the keyword and entry only if the.
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
labelList second(const UList< labelPair > &p)
Definition: patchToPatch.C:49
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const HashTable< unitConversion > & units()
Get the table of unit conversions.
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
List< Type > repeat(const UList< Type > &a, const UList< Type > &b)
tmp< DimensionedField< TypeR, GeoMesh, Field > > New(const tmp< DimensionedField< TypeR, GeoMesh, Field >> &tdf1, const word &name, const dimensionSet &dimensions)
static const char nl
Definition: Ostream.H:267
dictionary dict