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-2021 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,
139  const tableBase::boundsHandling boundsHandling,
140  const word& interpolationScheme,
141  const autoPtr<TableReader<Type>>& reader,
142  const List<Tuple2<scalar, Type>>& table
143 )
144 :
146  boundsHandling_(boundsHandling),
147  interpolationScheme_(interpolationScheme),
148  values_(table),
149  reader_(reader, false)
150 {}
151 
152 
153 template<class Type>
155 (
156  const word& name,
157  const dictionary& dict
158 )
159 :
161  boundsHandling_
162  (
163  dict.found("outOfBounds")
164  ? tableBase::boundsHandlingNames_.read(dict.lookup("outOfBounds"))
165  : tableBase::boundsHandling::clamp
166  ),
167  interpolationScheme_
168  (
169  dict.lookupOrDefault<word>
170  (
171  "interpolationScheme",
172  linearInterpolationWeights::typeName
173  )
174  ),
175  values_(),
176  reader_(TableReader<Type>::New(name, dict, this->values_))
177 {
178  check();
179 }
180 
181 
182 template<class Type>
184 :
185  FieldFunction1<Type, Table<Type>>(tbl),
186  boundsHandling_(tbl.boundsHandling_),
187  interpolationScheme_(tbl.interpolationScheme_),
188  values_(tbl.values_),
189  tableSamplesPtr_(tbl.tableSamplesPtr_),
190  interpolatorPtr_(tbl.interpolatorPtr_),
191  reader_(tbl.reader_, false)
192 {}
193 
194 
195 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
196 
197 template<class Type>
199 {}
200 
201 
202 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
203 
204 template<class Type>
206 (
207  const scalar x
208 ) const
209 {
210  const scalar bx = bound(x);
211 
212  Type y = Zero;
213 
214  interpolator().valueWeights(bx, indices_, weights_);
215  forAll(indices_, i)
216  {
217  y += weights_[i]*values_[indices_[i]].second();
218  }
219 
220  return y;
221 }
222 
223 
224 template<class Type>
226 (
227  const scalar x1,
228  const scalar x2
229 ) const
230 {
231  const scalar bx1 = bound(x1), bx2 = bound(x2);
232 
233  Type sumY = Zero;
234 
235  interpolator().integrationWeights(bx1, bx2, indices_, weights_);
236  forAll(indices_, i)
237  {
238  sumY += weights_[i]*values_[indices_[i]].second();
239  }
240 
241  if (boundsHandling_ == tableBase::boundsHandling::repeat)
242  {
243  const scalar t0 = values_.first().first();
244  const scalar t1 = values_.last().first();
245  const scalar dt = t1 - t0;
246  const label n = floor((x2 - t0)/dt) - floor((x1 - t0)/dt);
247 
248  if (n != 0)
249  {
250  Type sumY01 = Zero;
251 
252  interpolator().integrationWeights(t0, t1, indices_, weights_);
253 
254  forAll(indices_, i)
255  {
256  sumY01 += weights_[i]*values_[indices_[i]].second();
257  }
258  sumY += n*sumY01;
259  }
260  }
261 
262  return sumY;
263 }
264 
265 
266 template<class Type>
269 {
270  tmp<scalarField> tfld(new scalarField(values_.size(), 0.0));
271  scalarField& fld = tfld.ref();
272 
273  forAll(values_, i)
274  {
275  fld[i] = values_[i].first();
276  }
277 
278  return tfld;
279 }
280 
281 
282 template<class Type>
285 {
286  tmp<Field<Type>> tfld(new Field<Type>(values_.size(), Zero));
287  Field<Type>& fld = tfld.ref();
288 
289  forAll(values_, i)
290  {
291  fld[i] = values_[i].second();
292  }
293 
294  return tfld;
295 }
296 
297 
298 template<class Type>
300 {
302  (
303  os,
304  "outOfBounds",
305  tableBase::boundsHandlingNames_[tableBase::boundsHandling::clamp],
306  tableBase::boundsHandlingNames_[boundsHandling_]
307  );
308 
310  (
311  os,
312  "interpolationScheme",
313  linearInterpolationWeights::typeName,
314  interpolationScheme_
315  );
316 
317  reader_->write(os, values_);
318 }
319 
320 
321 // ************************************************************************* //
Base class to read table data for tables.
Definition: TableReader.H:49
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
List< Type > repeat(const UList< Type > &a, const UList< Type > &b)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
A 2-tuple for storing two objects of different types.
Definition: HashTable.H:65
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
virtual Type value(const scalar x) const
Return Table value as a function of scalar x.
Definition: Table.C:206
void writeEntryIfDifferent(Ostream &os, const word &entryName, const EntryType &value1, const EntryType &value2)
Helper function to write the keyword and entry only if the.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
T & first()
Return the first element of the list.
Definition: UListI.H:114
Abstract base class for interpolating in 1D.
Templated interpolated tabulated data Function1.
Definition: Table.H:118
virtual void write(Ostream &os) const
Write data to dictionary stream.
Definition: Table.C:299
scalar y
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))
A class for handling words, derived from string.
Definition: word.H:59
virtual ~Table()
Destructor.
Definition: Table.C:198
virtual Type integral(const scalar x1, const scalar x2) const
Integrate between two scalars.
Definition: Table.C:226
static const zero Zero
Definition: zero.H:97
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
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
volScalarField scalarField(fieldObject, mesh)
boundsHandling
Enumeration for handling out-of-bound values.
Definition: tableBase.H:54
static const char nl
Definition: Ostream.H:260
virtual tmp< Field< Type > > y() const
Return the dependent values.
Definition: Table.C:284
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define WarningInFunction
Report a warning using Foam::Warning.
label n
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< scalarField > x() const
Return the reference values.
Definition: Table.C:268
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864