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-2024 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>
86 (
87  const scalar x
88 ) const
89 {
90  const bool under = x < values_.first().first();
91  const bool over = x > values_.last().first();
92 
93  auto errorMessage = [&]()
94  {
95  return "value (" + name(x) + ") " + (under ? "under" : "over") + "flow";
96  };
97 
98  if (under || over)
99  {
100  switch (boundsHandling_)
101  {
102  case tableBase::boundsHandling::error:
103  {
105  << errorMessage() << nl << exit(FatalError);
106  break;
107  }
108  case tableBase::boundsHandling::warn:
109  {
111  << errorMessage() << nl << endl;
112  break;
113  }
114  case tableBase::boundsHandling::clamp:
115  {
116  break;
117  }
119  {
120  const scalar t0 = values_.first().first();
121  const scalar t1 = values_.last().first();
122  const scalar dt = t1 - t0;
123  const label n = floor((x - t0)/dt);
124  return x - n*dt;
125  }
126  }
127  }
128 
129  return x;
130 }
131 
132 
133 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
134 
135 template<class Type>
137 (
138  const word& name,
140  const word& interpolationScheme,
141  const autoPtr<TableReader<Type>>& reader,
142  const List<Tuple2<scalar, Type>>& table
143 )
144 :
145  FieldFunction1<Type, Table<Type>>(name),
146  boundsHandling_(boundsHandling),
147  interpolationScheme_(interpolationScheme),
148  reader_(reader, false),
149  values_(table)
150 {}
151 
152 
153 template<class Type>
155 (
156  const word& name,
157  const unitConversions& units,
158  const dictionary& dict
159 )
160 :
161  FieldFunction1<Type, Table<Type>>(name),
162  boundsHandling_
163  (
164  dict.found("outOfBounds")
165  ? tableBase::boundsHandlingNames_.read(dict.lookup("outOfBounds"))
166  : tableBase::boundsHandling::clamp
167  ),
168  interpolationScheme_
169  (
170  dict.lookupOrDefault<word>
171  (
172  "interpolationScheme",
174  )
175  ),
176  reader_(TableReader<Type>::New(name, units, dict)),
177  values_(reader_->read(units, dict))
178 {
179  check();
180 }
181 
182 
183 template<class Type>
185 (
186  const word& name,
187  const unitConversion& xUnits,
188  const unitConversion& valueUnits,
189  const dictionary& dict
190 )
191 :
192  Table(name, {xUnits, valueUnits}, dict)
193 {}
194 
195 
196 template<class Type>
198 (
199  const word& name,
201  Istream& is
202 )
203 :
204  FieldFunction1<Type, Table<Type>>(name),
205  boundsHandling_(tableBase::boundsHandling::clamp),
206  interpolationScheme_(linearInterpolationWeights::typeName),
207  reader_(new TableReaders::Embedded<Type>()),
208  values_(TableReaders::Embedded<Type>().read(units, is))
209 {
210  check();
211 }
212 
213 
214 template<class Type>
216 :
217  FieldFunction1<Type, Table<Type>>(tbl),
218  boundsHandling_(tbl.boundsHandling_),
219  interpolationScheme_(tbl.interpolationScheme_),
220  reader_(tbl.reader_, false),
221  values_(tbl.values_),
222  tableSamplesPtr_(tbl.tableSamplesPtr_),
223  interpolatorPtr_(tbl.interpolatorPtr_)
224 {}
225 
226 
227 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
228 
229 template<class Type>
231 {}
232 
233 
234 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
235 
236 template<class Type>
238 (
239  const scalar x
240 ) const
241 {
242  const scalar bx = bound(x);
243 
244  Type y = Zero;
245 
246  interpolator().valueWeights(bx, indices_, weights_);
247  forAll(indices_, i)
248  {
249  y += weights_[i]*values_[indices_[i]].second();
250  }
251 
252  return y;
253 }
254 
255 
256 template<class Type>
258 (
259  const scalar x1,
260  const scalar x2
261 ) const
262 {
263  const scalar bx1 = bound(x1);
264  const scalar bx2 = bound(x2);
265 
266  Type sumY = Zero;
267 
268  interpolator().integrationWeights(bx1, bx2, indices_, weights_);
269  forAll(indices_, i)
270  {
271  sumY += weights_[i]*values_[indices_[i]].second();
272  }
273 
274  if (boundsHandling_ == tableBase::boundsHandling::repeat)
275  {
276  const scalar t0 = values_.first().first();
277  const scalar t1 = values_.last().first();
278  const label n = floor(((x2 - x1) - (bx2 - bx1))/(t1 - t0) + 0.5);
279 
280  if (n != 0)
281  {
282  Type sumY01 = Zero;
283 
284  interpolator().integrationWeights(t0, t1, indices_, weights_);
285 
286  forAll(indices_, i)
287  {
288  sumY01 += weights_[i]*values_[indices_[i]].second();
289  }
290 
291  sumY += n*sumY01;
292  }
293  }
294 
295  return sumY;
296 }
297 
298 
299 template<class Type>
302 {
303  tmp<scalarField> tfld(new scalarField(values_.size(), 0.0));
304  scalarField& fld = tfld.ref();
305 
306  forAll(values_, i)
307  {
308  fld[i] = values_[i].first();
309  }
310 
311  return tfld;
312 }
313 
314 
315 template<class Type>
318 {
319  tmp<Field<Type>> tfld(new Field<Type>(values_.size(), Zero));
320  Field<Type>& fld = tfld.ref();
321 
322  forAll(values_, i)
323  {
324  fld[i] = values_[i].second();
325  }
326 
327  return tfld;
328 }
329 
330 
331 template<class Type>
333 (
334  Ostream& os,
335  const unitConversions& units
336 ) const
337 {
339  (
340  os,
341  "outOfBounds",
343  tableBase::boundsHandlingNames_[boundsHandling_]
344  );
345 
347  (
348  os,
349  "interpolationScheme",
350  linearInterpolationWeights::typeName,
351  interpolationScheme_
352  );
353 
354  reader_->write(os, units, values_);
355 }
356 
357 
358 // ************************************************************************* //
scalar y
bool found
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Templated interpolated tabulated data Function1.
Definition: Table.H:122
virtual tmp< scalarField > x() const
Return the reference values.
Definition: Table.C:301
virtual Type integral(const scalar x1, const scalar x2) const
Integrate between two scalars.
Definition: Table.C:258
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:137
virtual void write(Ostream &os, const unitConversions &units) const
Write data to dictionary stream.
Definition: Table.C:333
virtual tmp< Field< Type > > y() const
Return the dependent values.
Definition: Table.C:317
virtual Type value(const scalar x) const
Return Table value as a function of scalar x.
Definition: Table.C:238
virtual ~Table()
Destructor.
Definition: Table.C:230
static const NamedEnum< boundsHandling, 4 > boundsHandlingNames_
Enumeration names for handling out-of-bound values.
Definition: tableBase.H:63
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
Macros for creating standard TableReader-s.
Definition: TableReader.H:51
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:66
T & first()
Return the first element of the list.
Definition: UListI.H:114
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 keyword definitions, which are a keyword followed by any number of values (e....
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:181
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(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(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
#define WarningInFunction
Report a warning using Foam::Warning.
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
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:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const HashTable< unitConversion > & units()
Get the table of unit conversions.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
bool bound(volScalarField &, const dimensionedScalar &min)
Bound the given scalar field where it is below the specified min value.
Definition: bound.C:31
error FatalError
List< Type > repeat(const UList< Type > &a, const UList< Type > &b)
static const char nl
Definition: Ostream.H:266
dictionary dict