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 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,
34  const dictionary& dict
35 )
36 :
38  low_(great),
39  high_(-great),
40  values_(),
41  delta_(great),
42  reader_(TableReader<Type>::New(name, dict, this->values_))
43 {
44  if (values_.size() < 2)
45  {
47  << "Table " << nl
48  << " " << name << nl
49  << " has less than 2 entries."
50  << exit(FatalIOError);
51  }
52  else
53  {
54  low_ = values_.first().first();
55  high_ = values_.last().first();
56 
57  for(label i = 1; i<values_.size(); i++)
58  {
59  delta_ = min(delta_, values_[i].first() - values_[i - 1].first());
60  }
61 
62  delta_ *= 0.9;
63 
64  jumpTable_.setSize((high_ - low_)/delta_ + 1);
65 
66  label i = 0;
67  forAll(jumpTable_, j)
68  {
69  const scalar x = low_ + j*delta_;
70 
71  if (x > values_[i + 1].first())
72  {
73  i++;
74  }
75 
76  jumpTable_[j] = i;
77  }
78  }
79 }
80 
81 
82 template<class Type>
84 (
85  const NonUniformTable<Type>& nut
86 )
87 :
89  low_(nut.low_),
90  high_(nut.high_),
91  values_(nut.values_),
92  delta_(nut.delta_),
93  jumpTable_(nut.jumpTable_),
94  reader_(nut.reader_, false)
95 {}
96 
97 
98 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
99 
100 template<class Type>
102 (
103  scalar x
104 ) const
105 {
106  const label i = index(x);
107  const scalar xi = values_[i].first();
108  const scalar lambda = (x - xi)/(values_[i + 1].first() - xi);
109 
110  return
111  values_[i].second()
112  + lambda*(values_[i + 1].second() - values_[i].second());
113 }
114 
115 
116 template<class Type>
118 (
119  const scalar x1,
120  const scalar x2
121 ) const
122 {
124  return Zero;
125 }
126 
127 
128 template<class Type>
130 (
131  scalar T
132 ) const
133 {
134  const label i = index(T);
135 
136  return
137  (values_[i + 1].second() - values_[i].second())
138  /(values_[i + 1].first() - values_[i].first());
139 }
140 
141 
142 template<class Type>
144 {
145  reader_->write(os, values_);
146 }
147 
148 
149 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
150 
151 template<class Type>
152 void Foam::Function1s::NonUniformTable<Type>::operator=
153 (
155 )
156 {
157  low_ = nut.low_;
158  high_ = nut.high_;
159  values_ = nut.values_;
160  delta_ = nut.delta_;
161  jumpTable_ = nut.jumpTable_;
162  reader_ = nut.reader_->clone();
163 }
164 
165 
166 // ************************************************************************* //
Base class to read table data for tables.
Definition: TableReader.H:49
#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
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
virtual Type integral(const scalar x1, const scalar x2) const
Integrate between two scalar values.
virtual Type value(scalar x) const
Evaluate the function and return the result.
Type dfdT(scalar T) const
Evaluate the derivative of the function and return the result.
A class for handling words, derived from string.
Definition: word.H:59
static const zero Zero
Definition: zero.H:97
void write(Ostream &os) const
Write the function coefficients.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
static const char nl
Definition: Ostream.H:260
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Non-uniform tabulated property function that linearly interpolates between the values.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:335
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:370
NonUniformTable(const word &name, const dictionary &dict)
Construct from name and dictionary.
scalar nut
IOerror FatalIOError