nonUniformTableThermophysicalFunction.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 
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace thermophysicalFunctions
34 {
35  defineTypeNameAndDebug(nonUniformTable, 0);
36 
38  (
39  thermophysicalFunction,
40  nonUniformTable,
41  dictionary
42  );
43 }
44 }
45 
46 template<>
48 (
49  "Tuple2<scalar,scalar>"
50 );
51 
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
56 (
57  const word& name,
58  const dictionary& dict
59 )
60 :
61  name_(name),
62  Tlow_(great),
63  Thigh_(-great),
64  values_(dict.lookup(name)),
65  deltaT_(great)
66 {
67  if (values_.size() < 2)
68  {
70  << "Table " << nl
71  << " " << name_ << nl
72  << " has less than 2 entries."
73  << exit(FatalIOError);
74  }
75  else
76  {
77  Tlow_ = values_.first().first();
78  Thigh_ = values_.last().first();
79 
80  for(label i = 1; i<values_.size(); i++)
81  {
82  deltaT_ = min(deltaT_, values_[i].first() - values_[i - 1].first());
83  }
84 
85  deltaT_ *= 0.9;
86 
87  jumpTable_.setSize((Thigh_ - Tlow_)/deltaT_ + 1);
88 
89  label i = 0;
90  forAll(jumpTable_, j)
91  {
92  const scalar T = Tlow_ + j*deltaT_;
93 
94  if (T > values_[i + 1].first())
95  {
96  i++;
97  }
98 
99  jumpTable_[j] = i;
100  }
101  }
102 }
103 
104 
106 (
107  const dictionary& dict
108 )
109 :
110  nonUniformTable("values", dict)
111 {}
112 
113 
114 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
115 
117 (
118  scalar p,
119  scalar T
120 ) const
121 {
122  const label i = index(p, T);
123  const scalar Ti = values_[i].first();
124  const scalar lambda = (T - Ti)/(values_[i + 1].first() - Ti);
125 
126  return
127  values_[i].second()
128  + lambda*(values_[i + 1].second() - values_[i].second());
129 }
130 
131 
133 (
134  scalar p,
135  scalar T
136 ) const
137 {
138  const label i = index(p, T);
139 
140  return
141  (values_[i + 1].second() - values_[i].second())
142  /(values_[i + 1].first() - values_[i].first());
143 }
144 
145 
147 {
148  writeEntry(os, "values", values_);
149 }
150 
151 
152 // ************************************************************************* //
#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
defineTypeNameAndDebug(APIdiffCoef, 0)
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:158
Macros for easy insertion into run-time selection tables.
A class for handling words, derived from string.
Definition: word.H:59
scalar dfdT(scalar p, scalar T) const
Evaluate the derivative of the function and return the result.
addToRunTimeSelectionTable(thermophysicalFunction, APIdiffCoef, dictionary)
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
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
scalar f(scalar p, scalar T) const
Evaluate the function and return the result.
nonUniformTable(const word &name, const dictionary &dict)
Construct from entry name and dictionary.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
Non-uniform tabulated property function that linearly interpolates between the values.
static const char *const typeName
Definition: Tuple2.H:73
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812
IOerror FatalIOError