All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
nonUniformTableThermophysicalFunction.H
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 Class
25  Foam::thermophysicalFunctions::nonUniformTable
26 
27 Description
28  Non-uniform tabulated property function that linearly interpolates between
29  the values.
30 
31  To speed-up the search of the non-uniform table a uniform jump-table is
32  created on construction which is used for fast indirect addressing into
33  the table.
34 
35 Usage
36  \nonUniformTable
37  Property | Description
38  values | List of (temperature property) value pairs
39  \endnonUniformTable
40 
41  Example for the density of water between 280 and 350K
42  \verbatim
43  rho
44  {
45  type nonUniformTable;
46 
47  values
48  (
49  (280 999.87)
50  (300 995.1)
51  (350 973.7)
52  );
53  }
54  \endverbatim
55 
56 \*---------------------------------------------------------------------------*/
57 
58 #ifndef nonUniformTableThermophysicalFunction_H
59 #define nonUniformTableThermophysicalFunction_H
60 
61 #include "thermophysicalFunction.H"
62 #include "Tuple2.H"
63 
64 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65 
66 namespace Foam
67 {
68 namespace thermophysicalFunctions
69 {
70 
71 /*---------------------------------------------------------------------------*\
72  Class nonUniformTable Declaration
73 \*---------------------------------------------------------------------------*/
74 
75 class nonUniformTable
76 :
78 {
79  // Private member data
80 
81  //- Name of dictionary from which this function is instantiated
82  word name_;
83 
84  //- Lowest temperature in the nonUniformTable
85  scalar Tlow_;
86 
87  //- Highest temperature in the nonUniformTable
88  scalar Thigh_;
89 
90  //- Table values
92 
93  //- Temperature increment derived from Tlow_, Thigh_ and values_.size()
94  scalar deltaT_;
95 
96  List<label> jumpTable_;
97 
98 
99 protected:
100 
101  //- Return the lower index of the interval in the table
102  // corresponding to the given temperature
103  inline label index(scalar p, scalar T) const;
104 
105 
106 public:
107 
108  //- Runtime type information
109  TypeName("nonUniformTable");
110 
111 
112  // Constructors
113 
114  //- Construct from entry name and dictionary
115  nonUniformTable(const word& name, const dictionary& dict);
116 
117  //- Construct from dictionary
118  nonUniformTable(const dictionary& dict);
119 
120 
121  // Member Functions
122 
123  //- Return the non-uniform table of values
124  const List<Tuple2<scalar, scalar>>& values() const
125  {
126  return values_;
127  }
128 
129  //- Evaluate the function and return the result
130  scalar f(scalar p, scalar T) const;
131 
132  //- Evaluate the derivative of the function and return the result
133  scalar dfdT(scalar p, scalar T) const;
134 
135  //- Write the function coefficients
136  void write(Ostream& os) const;
137 };
138 
139 
140 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
141 
142 } // End namespace thermophysicalFunctions
143 } // End namespace Foam
144 
145 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
146 
148 
149 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
150 
151 #endif
152 
153 // ************************************************************************* //
dictionary dict
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
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
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
TypeName("nonUniformTable")
Runtime type information.
Abstract base class for thermo-physical functions.
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.
const List< Tuple2< scalar, scalar > > & values() const
Return the non-uniform table of values.
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
label index(scalar p, scalar T) const
Return the lower index of the interval in the table.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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.
Non-uniform tabulated property function that linearly interpolates between the values.
volScalarField & p
Namespace for OpenFOAM.