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-2020 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  }
118  case tableBase::boundsHandling::repeat:
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 List<Tuple2<scalar, Type>>& table
142 )
143 :
145  boundsHandling_(boundsHandling),
146  interpolationScheme_(interpolationScheme),
147  values_(table)
148 {}
149 
150 
151 template<class Type>
153 (
154  const word& name,
155  const dictionary& dict
156 )
157 :
159  boundsHandling_
160  (
161  dict.found("outOfBounds")
162  ? tableBase::boundsHandlingNames_.read(dict.lookup("outOfBounds"))
163  : tableBase::boundsHandling::clamp
164  ),
165  interpolationScheme_
166  (
167  dict.lookupOrDefault<word>
168  (
169  "interpolationScheme",
170  linearInterpolationWeights::typeName
171  )
172  ),
173  values_(),
174  reader_(TableReader<Type>::New(name, dict, this->values_))
175 {
176  check();
177 }
178 
179 
180 template<class Type>
182 :
183  FieldFunction1<Type, Table<Type>>(tbl),
184  boundsHandling_(tbl.boundsHandling_),
185  interpolationScheme_(tbl.interpolationScheme_),
186  values_(tbl.values_),
187  tableSamplesPtr_(tbl.tableSamplesPtr_),
188  interpolatorPtr_(tbl.interpolatorPtr_),
189  reader_(tbl.reader_, false)
190 {}
191 
192 
193 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
194 
195 template<class Type>
197 {}
198 
199 
200 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
201 
202 template<class Type>
204 (
205  const scalar x
206 ) const
207 {
208  const scalar bx = bound(x);
209 
210  Type y = Zero;
211 
212  interpolator().valueWeights(bx, indices_, weights_);
213  forAll(indices_, i)
214  {
215  y += weights_[i]*values_[indices_[i]].second();
216  }
217 
218  return y;
219 }
220 
221 
222 template<class Type>
224 (
225  const scalar x1,
226  const scalar x2
227 ) const
228 {
229  const scalar bx1 = bound(x1), bx2 = bound(x2);
230 
231  Type sumY = Zero;
232 
233  interpolator().integrationWeights(bx1, bx2, indices_, weights_);
234  forAll(indices_, i)
235  {
236  sumY += weights_[i]*values_[indices_[i]].second();
237  }
238 
239  if (boundsHandling_ == tableBase::boundsHandling::repeat)
240  {
241  const scalar t0 = values_.first().first();
242  const scalar t1 = values_.last().first();
243  const scalar dt = t1 - t0;
244  const label n = floor((x2 - t0)/dt) - floor((x1 - t0)/dt);
245 
246  if (n != 0)
247  {
248  Type sumY01 = Zero;
249 
250  interpolator().integrationWeights(t0, t1, indices_, weights_);
251 
252  forAll(indices_, i)
253  {
254  sumY01 += weights_[i]*values_[indices_[i]].second();
255  }
256  sumY += n*sumY01;
257  }
258  }
259 
260  return sumY;
261 }
262 
263 
264 template<class Type>
267 {
268  tmp<scalarField> tfld(new scalarField(values_.size(), 0.0));
269  scalarField& fld = tfld.ref();
270 
271  forAll(values_, i)
272  {
273  fld[i] = values_[i].first();
274  }
275 
276  return tfld;
277 }
278 
279 
280 template<class Type>
283 {
284  tmp<Field<Type>> tfld(new Field<Type>(values_.size(), Zero));
285  Field<Type>& fld = tfld.ref();
286 
287  forAll(values_, i)
288  {
289  fld[i] = values_[i].second();
290  }
291 
292  return tfld;
293 }
294 
295 
296 template<class Type>
298 {
300  (
301  os,
302  "outOfBounds",
303  tableBase::boundsHandlingNames_[tableBase::boundsHandling::clamp],
304  tableBase::boundsHandlingNames_[boundsHandling_]
305  );
306 
308  (
309  os,
310  "interpolationScheme",
311  linearInterpolationWeights::typeName,
312  interpolationScheme_
313  );
314 
315  reader_->write(os, values_);
316 }
317 
318 
319 // ************************************************************************* //
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:643
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
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:323
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:204
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
Table(const word &name, const tableBase::boundsHandling boundsHandling, const word &interpolationScheme, const List< Tuple2< scalar, Type >> &table)
Construct from components.
Definition: Table.C:137
Abstract base class for interpolating in 1D.
Templated interpolated tabulated data Function1.
Definition: Table.H:119
virtual void write(Ostream &os) const
Write data to dictionary stream.
Definition: Table.C:297
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))
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
A class for handling words, derived from string.
Definition: word.H:59
virtual ~Table()
Destructor.
Definition: Table.C:196
virtual Type integral(const scalar x1, const scalar x2) const
Integrate between two scalars.
Definition: Table.C:224
static const zero Zero
Definition: zero.H:97
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:282
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
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< scalarField > x() const
Return the reference values.
Definition: Table.C:266
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844