NonUniformTable1.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) 2020-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 "NonUniformTable1.H"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class Type>
32 (
33  const word& name,
35  const dictionary& dict
36 )
37 :
38  FieldFunction1<Type, NonUniformTable<Type>>(name),
39  low_(great),
40  high_(-great),
41  values_(),
42  delta_(great),
43  reader_(TableReader<Type>::New(name, units, dict))
44 {
45  assertNoConvertUnits(typeName, units, dict);
46 
47  values_ = reader_->read(units, dict);
48 
49  if (values_.size() < 2)
50  {
52  << "Table " << nl
53  << " " << name << nl
54  << " has less than 2 entries."
55  << exit(FatalIOError);
56  }
57  else
58  {
59  low_ = values_.first().first();
60  high_ = values_.last().first();
61 
62  for(label i = 1; i<values_.size(); i++)
63  {
64  delta_ = min(delta_, values_[i].first() - values_[i - 1].first());
65  }
66 
67  delta_ *= 0.9;
68 
69  jumpTable_.setSize((high_ - low_)/delta_ + 1);
70 
71  label i = 0;
72  forAll(jumpTable_, j)
73  {
74  const scalar x = low_ + j*delta_;
75 
76  if (x > values_[i + 1].first())
77  {
78  i++;
79  }
80 
81  jumpTable_[j] = i;
82  }
83  }
84 }
85 
86 
87 template<class Type>
89 (
91 )
92 :
93  FieldFunction1<Type, NonUniformTable<Type>>(nut),
94  low_(nut.low_),
95  high_(nut.high_),
96  values_(nut.values_),
97  delta_(nut.delta_),
98  jumpTable_(nut.jumpTable_),
99  reader_(nut.reader_, false)
100 {}
101 
102 
103 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
104 
105 template<class Type>
107 (
108  scalar x
109 ) const
110 {
111  const label i = index(x);
112  const scalar xi = values_[i].first();
113  const scalar lambda = (x - xi)/(values_[i + 1].first() - xi);
114 
115  return
116  values_[i].second()
117  + lambda*(values_[i + 1].second() - values_[i].second());
118 }
119 
120 
121 template<class Type>
123 (
124  const scalar x1,
125  const scalar x2
126 ) const
127 {
129  return Zero;
130 }
131 
132 
133 template<class Type>
135 (
136  scalar T
137 ) const
138 {
139  const label i = index(T);
140 
141  return
142  (values_[i + 1].second() - values_[i].second())
143  /(values_[i + 1].first() - values_[i].first());
144 }
145 
146 
147 template<class Type>
149 (
150  Ostream& os,
151  const unitConversions& units
152 ) const
153 {
154  reader_->write(os, units, values_);
155 }
156 
157 
158 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
159 
160 template<class Type>
162 (
164 )
165 {
166  low_ = nut.low_;
167  high_ = nut.high_;
168  values_ = nut.values_;
169  delta_ = nut.delta_;
170  jumpTable_ = nut.jumpTable_;
171  reader_ = nut.reader_->clone();
172 }
173 
174 
175 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const word & name() const
Return the name of the entry.
Definition: Function1.C:78
Non-uniform tabulated property function that linearly interpolates between the values.
NonUniformTable(const word &name, const Function1s::unitConversions &units, const dictionary &dict)
Construct from name and dictionary.
Type dfdT(scalar T) const
Evaluate the derivative of the function and return the result.
virtual Type integral(const scalar x1, const scalar x2) const
Integrate between two scalar values.
void write(Ostream &os, const unitConversions &units) const
Write the function coefficients.
virtual Type value(scalar x) const
Evaluate the function and return the result.
void setSize(const label)
Reset size of List.
Definition: List.C:281
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 list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
A class for handling words, derived from string.
Definition: word.H:62
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
const scalar nut
dimensionedScalar lambda(viscosity->lookup("lambda"))
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
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
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
labelList second(const UList< labelPair > &p)
Definition: patchToPatch.C:49
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const HashTable< unitConversion > & units()
Get the table of unit conversions.
IOerror FatalIOError
void assertNoConvertUnits(const word &typeName, const Function1s::unitConversions &units, const dictionary &dict)
Generate an error in an context where unit conversions are not supported.
static const char nl
Definition: Ostream.H:266
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict