TableBase.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-2019 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"
29 
30 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
31 
32 template<class Type, class Function1Type>
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 template<class Type, class Function1Type>
58 {
59  if (!table_.size())
60  {
62  << "Table for entry " << this->name_ << " is invalid (empty)"
63  << nl << exit(FatalError);
64  }
65 
66  label n = table_.size();
67  scalar prevValue = table_[0].first();
68 
69  for (label i = 1; i < n; ++i)
70  {
71  const scalar currValue = table_[i].first();
72 
73  // avoid duplicate values (divide-by-zero error)
74  if (currValue <= prevValue)
75  {
77  << "out-of-order value: " << currValue << " at index " << i
78  << exit(FatalError);
79  }
80  prevValue = currValue;
81  }
82 }
83 
84 
85 template<class Type, class Function1Type>
87 (
88  const scalar x
89 ) const
90 {
91  const bool under = x < table_.first().first();
92  const bool over = x > table_.last().first();
93 
94  auto errorMessage = [&]()
95  {
96  return "value (" + name(x) + ") " + (under ? "under" : "over") + "flow";
97  };
98 
99  if (under || over)
100  {
101  switch (boundsHandling_)
102  {
103  case tableBase::boundsHandling::error:
104  {
106  << errorMessage() << nl << exit(FatalError);
107  break;
108  }
109  case tableBase::boundsHandling::warn:
110  {
112  << errorMessage() << nl << endl;
113  break;
114  }
115  case tableBase::boundsHandling::clamp:
116  {
117  break;
118  }
119  case tableBase::boundsHandling::repeat:
120  {
121  const scalar t0 = table_.first().first();
122  const scalar t1 = table_.last().first();
123  const scalar dt = t1 - t0;
124  const label n = floor((x - t0)/dt);
125  return x - n*dt;
126  }
127  }
128  }
129 
130  return x;
131 }
132 
133 
134 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
135 
136 template<class Type, class Function1Type>
138 (
139  const word& name,
140  const dictionary& dict
141 )
142 :
143  tableBase(),
145  name_(name),
146  boundsHandling_
147  (
148  dict.found("outOfBounds")
149  ? tableBase::boundsHandlingNames_.read(dict.lookup("outOfBounds"))
150  : tableBase::boundsHandling::clamp
151  ),
152  interpolationScheme_
153  (
154  dict.lookupOrDefault<word>
155  (
156  "interpolationScheme",
157  linearInterpolationWeights::typeName
158  )
159  ),
160  table_()
161 {}
162 
163 
164 template<class Type, class Function1Type>
166 (
167  const word& name,
168  const tableBase::boundsHandling boundsHandling,
169  const word& interpolationScheme,
170  const List<Tuple2<scalar, Type>>& table
171 )
172 :
173  tableBase(),
175  name_(name),
176  boundsHandling_(boundsHandling),
177  interpolationScheme_(interpolationScheme),
178  table_(table)
179 {}
180 
181 
182 template<class Type, class Function1Type>
184 (
186 )
187 :
188  tableBase(),
190  name_(tbl.name_),
191  boundsHandling_(tbl.boundsHandling_),
192  interpolationScheme_(tbl.interpolationScheme_),
193  table_(tbl.table_),
194  tableSamplesPtr_(tbl.tableSamplesPtr_),
195  interpolatorPtr_(tbl.interpolatorPtr_)
196 {}
197 
198 
199 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
200 
201 template<class Type, class Function1Type>
203 {}
204 
205 
206 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
207 
208 template<class Type, class Function1Type>
210 (
211  const scalar x
212 ) const
213 {
214  const scalar bx = bound(x);
215 
216  Type y = Zero;
217 
218  interpolator().valueWeights(bx, indices_, weights_);
219  forAll(indices_, i)
220  {
221  y += weights_[i]*table_[indices_[i]].second();
222  }
223 
224  return y;
225 }
226 
227 
228 template<class Type, class Function1Type>
230 (
231  const scalar x1,
232  const scalar x2
233 ) const
234 {
235  const scalar bx1 = bound(x1), bx2 = bound(x2);
236 
237  Type sumY = Zero;
238 
239  interpolator().integrationWeights(bx1, bx2, indices_, weights_);
240  forAll(indices_, i)
241  {
242  sumY += weights_[i]*table_[indices_[i]].second();
243  }
244 
245  if (boundsHandling_ == tableBase::boundsHandling::repeat)
246  {
247  const scalar t0 = table_.first().first();
248  const scalar t1 = table_.last().first();
249  const scalar dt = t1 - t0;
250  const label n = floor((x2 - t0)/dt) - floor((x1 - t0)/dt);
251 
252  if (n != 0)
253  {
254  Type sumY01 = Zero;
255 
256  interpolator().integrationWeights(t0, t1, indices_, weights_);
257 
258  forAll(indices_, i)
259  {
260  sumY01 += weights_[i]*table_[indices_[i]].second();
261  }
262  sumY += n*sumY01;
263  }
264  }
265 
266  return sumY;
267 }
268 
269 
270 template<class Type, class Function1Type>
273 {
274  tmp<scalarField> tfld(new scalarField(table_.size(), 0.0));
275  scalarField& fld = tfld.ref();
276 
277  forAll(table_, i)
278  {
279  fld[i] = table_[i].first();
280  }
281 
282  return tfld;
283 }
284 
285 
286 template<class Type, class Function1Type>
289 {
290  tmp<Field<Type>> tfld(new Field<Type>(table_.size(), Zero));
291  Field<Type>& fld = tfld.ref();
292 
293  forAll(table_, i)
294  {
295  fld[i] = table_[i].second();
296  }
297 
298  return tfld;
299 }
300 
301 
302 template<class Type, class Function1Type>
304 (
305  Ostream& os
306 ) const
307 {
309  (
310  os,
311  "outOfBounds",
312  tableBase::boundsHandlingNames_[tableBase::boundsHandling::clamp],
313  tableBase::boundsHandlingNames_[boundsHandling_]
314  );
315 
317  (
318  os,
319  "interpolationScheme",
320  linearInterpolationWeights::typeName,
321  interpolationScheme_
322  );
323 }
324 
325 
326 template<class Type, class Function1Type>
328 (
329  Ostream& os
330 ) const
331 {
333  os << token::END_STATEMENT << nl;
334 
335  os << indent << word(this->name() + "Coeffs") << nl;
336  os << indent << token::BEGIN_BLOCK << incrIndent << nl;
337  writeEntries(os);
338  os << decrIndent << indent << token::END_BLOCK << endl;
339 }
340 
341 
342 // ************************************************************************* //
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:667
#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
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
List< Tuple2< scalar, Type > > table_
Table data.
Definition: TableBase.H:75
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
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:174
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
virtual Type value(const scalar x) const
Return Table value.
Definition: TableBase.C:210
T & first()
Return the first element of the list.
Definition: UListI.H:114
Abstract base class for interpolating in 1D.
const word interpolationScheme_
Interpolation type.
Definition: TableBase.H:72
scalar y
const tableBase::boundsHandling boundsHandling_
Enumeration for handling out-of-bound values.
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))
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
scalar bound(const scalar x) const
Bound the argument to the table. Errors or warns, or shifts the.
Definition: TableBase.C:87
virtual ~TableBase()
Destructor.
Definition: TableBase.C:202
static const zero Zero
Definition: zero.H:97
virtual Type integrate(const scalar x1, const scalar x2) const
Integrate between two (scalar) values.
Definition: TableBase.C:230
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
const bool writeData(readBool(pdfDictionary.lookup("writeData")))
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const word name_
Table name.
Definition: TableBase.H:66
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
virtual tmp< Field< Type > > y() const
Return the dependent values.
Definition: TableBase.C:288
TableBase(const word &name, const dictionary &dict)
Construct from dictionary. Table is not populated.
Definition: TableBase.C:138
virtual tmp< scalarField > x() const
Return the reference values.
Definition: TableBase.C:272
A class for managing temporary objects.
Definition: PtrList.H:53
const interpolationWeights & interpolator() const
Return (demand driven) interpolator.
Definition: TableBase.C:34
autoPtr< interpolationWeights > interpolatorPtr_
Interpolator method.
Definition: TableBase.H:81
virtual void writeEntries(Ostream &os) const
Write entries only in dictionary format.
Definition: TableBase.C:304
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:228
autoPtr< scalarField > tableSamplesPtr_
Extracted values.
Definition: TableBase.H:78
Base class for table with bounds handling, interpolation and integration.
Definition: TableBase.H:56
virtual void writeData(Ostream &os) const
Write all table data in dictionary format.
Definition: TableBase.C:328
void check() const
Check the table for size and consistency.
Definition: TableBase.C:57
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812