interpolationLookUpTable.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 
29 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
30 
31 template<class Type>
33 (
34  const List<scalar>& indices,
35  const bool lastDim
36 ) const
37 {
38  label totalIndex = 0;
39 
40  forAll(dim_, i)
41  {
42  label dim = 1;
43  for (int j = i + 1; j < dim_.size(); j++)
44  {
45  dim *= dim_[j] + 1;
46  }
47 
48  totalIndex +=
49  dim
50  *min
51  (
52  max(label((indices[i] - min_[i])/delta_[i]), 0),
53  dim_[i]
54  );
55  }
56 
57  if (lastDim)
58  {
59  label iLastdim = dim_.size() - 1;
60  totalIndex += Foam::min
61  (
62  max
63  (
64  label((indices[iLastdim] - min_[iLastdim])/delta_[iLastdim]),
65  0
66  ),
67  dim_[iLastdim]
68  );
69  }
70 
71  return totalIndex;
72 }
73 
74 
75 template<class Type>
77 (
78  const scalar indice
79 ) const
80 {
81  label i = 0;
82  label totalIndex =
83  Foam::min
84  (
85  Foam::max
86  (
87  label((indice - min_[i])/delta_[i]),
88  0
89  ),
90  dim_[i]
91  );
92 
93  return totalIndex;
94 }
95 
96 
97 template<class Type>
99 (
100  const scalar lookUpValue,
101  const label interfield
102 ) const
103 {
104  return lookUpValue >= min_[interfield] && lookUpValue <= max_[interfield];
105 }
106 
107 
108 template<class Type>
110 (
111  const label lo,
112  const label hi,
113  const scalar lookUpValue,
114  const label ofield,
115  const label interfield
116 ) const
117 {
118  if
119  (
120  List<scalarField>::operator[](interfield).operator[](hi)
121  != List<scalarField>::operator[](interfield).operator[](lo)
122  )
123  {
124  scalar output
125  (
126  List<scalarField>::operator[](ofield).operator[](lo)
127  + (
128  List<scalarField>::operator[](ofield).operator[](hi)
129  - List<scalarField>::operator[](ofield).operator[](lo)
130  )
131  *(
132  lookUpValue
133  - List<scalarField>::operator[](interfield).operator[](lo)
134  )
135  /(
136  List<scalarField>::operator[](interfield).operator[](hi)
137  - List<scalarField>::operator[](interfield).operator[](lo)
138  )
139  );
140  return output;
141  }
142  else
143  {
144  return List<scalarField>::operator[](ofield).operator[](lo);
145  }
146 }
147 
148 
149 template<class Type>
151 {
152  min_.setSize(entries_.size());
153  dim_.setSize(entries_.size());
154  delta_.setSize(entries_.size());
155  max_.setSize(entries_.size());
156  entryIndices_.setSize(entries_.size());
157  outputIndices_.setSize(output_.size());
158  label index = 0;
159  label tableDim = 1;
160 
161  forAll(entries_,i)
162  {
163  dim_[i] = entries_[i].template get<label>("N");
164  max_[i] = entries_[i].template get<scalar>("max");
165  min_[i] = entries_[i].template get<scalar>("min");
166  delta_[i] = (max_[i] - min_[i])/dim_[i];
167  tableDim *= dim_[i] + 1;
168  fieldIndices_.insert(entries_[i].template get<word>("name"), index);
169  entryIndices_[i] = index;
170  index++;
171  }
172 
173  forAll(output_,i)
174  {
175  fieldIndices_.insert(output_[i].template get<word>("name"), index);
176  outputIndices_[i] = index;
177  index++;
178  }
179 
180  List<scalarField>& internal = *this;
181 
182  internal.setSize(entries_.size() + output_.size());
183 
184  interpOutput_.setSize(entries_.size() + output_.size());
185 
186  forAll(internal, i)
187  {
188  internal[i].setSize(tableDim);
189  }
190 }
191 
192 
193 template<class Type>
195 (
196  const word& instance,
197  const objectRegistry& obr
198 )
199 {
200  IOdictionary control
201  (
202  IOobject
203  (
204  fileName_,
205  instance,
206  obr,
207  IOobject::MUST_READ_IF_MODIFIED,
208  IOobject::NO_WRITE
209  )
210  );
211 
212  control.readEntry("fields", entries_);
213  control.readEntry("output", output_);
214  control.readEntry("values", *this);
215 
216  dimensionTable();
217 
218  check();
219 
220  if (this->size() == 0)
221  {
223  << "table is empty" << nl << exit(FatalError);
224  }
225 }
226 
227 
228 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
229 
230 template<class Type>
232 :
233  List<scalarField>(),
234  fileName_("fileNameIsUndefined")
235 {}
236 
237 
238 template<class Type>
240 (
241  const fileName& fn,
242  const word& instance,
243  const objectRegistry& obr
244 )
245 :
247  fileName_(fn),
248  dim_(0),
249  min_(0),
250  delta_(0.0),
251  max_(0.0),
252  entries_(0),
253  output_(0),
254  entryIndices_(0),
255  outputIndices_(0),
256  interpOutput_(0)
257 {
258  readTable(instance, obr);
259 }
260 
261 
262 template<class Type>
264 (
265  const interpolationLookUpTable& interpTable
266 )
267 :
268  List<scalarField>(interpTable),
269  fileName_(interpTable.fileName_),
270  entryIndices_(interpTable.entryIndices_),
271  outputIndices_(interpTable.outputIndices_),
272  dim_(interpTable.dim_),
273  min_(interpTable.min_),
274  delta_(interpTable.delta_),
275  max_(interpTable.max_),
276  entries_(0),
277  output_(0),
278  interpOutput_(interpTable.interpOutput_)
279 {}
280 
281 
282 template<class Type>
284 (
285  const dictionary& dict
286 )
287 :
289  fileName_(dict.template get<fileName>("file").expand()),
290  dim_(0),
291  min_(0.0),
292  delta_(0.0),
293  max_(0.0),
294  entries_(dict.lookup("fields")),
295  output_(dict.lookup("output")),
296  entryIndices_(0),
297  outputIndices_(0),
298  fieldIndices_(0),
299  interpOutput_(0)
300 {
301  dimensionTable();
302 }
303 
304 
305 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
306 
307 template<class Type>
309 {
310  // check order in the first dimension.
311  scalar prevValue = List<scalarField>::operator[](0).operator[](0);
312  label dim = 1;
313  for (int j = 1; j < dim_.size(); j++)
314  {
315  dim *= dim_[j] + 1;
316  }
317 
318  for (label i = 1; i < dim_[0]; i++)
319  {
320  label index = i*dim;
321  const scalar currValue =
322  List<scalarField>::operator[](0).operator[](index);
323 
324  // avoid duplicate values (divide-by-zero error)
325  if (currValue <= prevValue)
326  {
328  << "out-of-order value: " << currValue
329  << " at index " << index << nl << exit(FatalError);
330  }
331  prevValue = currValue;
332  }
333 }
334 
335 
336 template<class Type>
338 (
339  Ostream& os,
340  const fileName& fn,
341  const word& instance,
342  const objectRegistry& obr
343 ) const
344 {
345  IOdictionary control
346  (
347  IOobject
348  (
349  fn,
350  instance,
351  obr,
352  IOobject::NO_READ,
353  IOobject::NO_WRITE
354  )
355  );
356 
357  control.writeHeader(os);
358 
359  os.writeEntry("fields", entries_);
360  os.writeEntry("output", output_);
361 
362  if (this->size() == 0)
363  {
365  << "table is empty" << nl << exit(FatalError);
366  }
367  os.writeEntry("values", *this);
368 }
369 
370 
371 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
372 
373 template<class Type>
376 {
377  const label n = this->size();
378 
379  if (n <= 1)
380  {
382  << "table has (" << n << ") columns" << nl << exit(FatalError);
383  }
384  else if (i < 0)
385  {
387  << "index (" << i << ") underflow" << nl << exit(FatalError);
388  }
389  else if (i >= n)
390  {
392  << "index (" << i << ") overflow" << nl << exit(FatalError);
393  }
394 
396 }
397 
398 
399 template<class Type>
400 const Foam::scalarField&
402 {
403  const label n = this->size();
404 
405  if (n <= 1)
406  {
408  << "table has (" << n << ") columns" << nl << exit(FatalError);
409  }
410  else if (i < 0)
411  {
413  << "index (" << i << ") underflow" << nl << exit(FatalError);
414  }
415  else if (i >= n)
416  {
418  << "index (" << i << ") overflow" << nl
419  << exit(FatalError);
420  }
421 
423 }
424 
425 
426 template<class Type>
428 {
429  return fieldIndices_.found(fieldName);
430 }
431 
432 
433 template<class Type>
434 const Foam::scalarList&
436 {
437  const label lo = index(retvals);
438  findHi(lo, retvals);
439  return interpOutput_;
440 }
441 
442 
443 template<class Type>
445 (
446  const label lo,
447  const scalar retvals
448 )
449 {
450  forAll(outputIndices_,j)
451  {
452  scalar tmp = 0;
453  label ofield = outputIndices_[j];
454  scalar baseValue = List<scalarField>::operator[](ofield).operator[](lo);
455 
456  forAll(entryIndices_,i)
457  {
458  if (checkRange(retvals, entryIndices_[i]))
459  {
460  label dim = 1;
461 
462  label hi = Foam::min(lo + dim, (*this)[0].size() - 1);
463 
464  tmp += interpolate(lo, hi, retvals, ofield, entryIndices_[i])
465  - baseValue;
466  }
467  interpOutput_[entryIndices_[i]] = retvals;
468  }
469 
470  tmp += baseValue;
471  interpOutput_[outputIndices_[j]] = tmp;
472  }
473 }
474 
475 
476 // ************************************************************************* //
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::output
static Ostream & output(Ostream &os, const IntRange< T > &range)
Definition: IntRanges.C:66
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::interpolationLookUpTable::found
bool found(const word &fieldName) const
Return true if the field exists in the table.
Definition: interpolationLookUpTable.C:427
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::interpolationLookUpTable::operator[]
const scalarField & operator[](const label) const
Return an element of constant List<scalar, Type>
Definition: interpolationLookUpTable.C:401
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::check
static void check(const int retVal, const char *what)
Definition: ptscotchDecomp.C:80
Foam::Field< scalar >
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
Foam::interpolationLookUpTable::interpolationLookUpTable
interpolationLookUpTable()
Construct null.
Definition: interpolationLookUpTable.C:231
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::interpolationLookUpTable::lookUp
const List< scalar > & lookUp(const scalar)
Return the output list given a single input scalar.
Definition: interpolationLookUpTable.C:435
Foam::IOobject::writeHeader
bool writeHeader(Ostream &os) const
Write header with current type()
Definition: IOobjectWriteHeader.C:277
Foam::interpolationLookUpTable
A list of lists. Interpolates based on the first dimension. The values must be positive and monotonic...
Definition: interpolationLookUpTable.H:63
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::interpolationLookUpTable::write
void write(Ostream &, const fileName &, const word &instance, const objectRegistry &) const
Write lookup table to filename.
Definition: interpolationLookUpTable.C:338
Foam::fac::interpolate
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.