nonUniformTableThermophysicalFunction.H
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) 2020 OpenFOAM Foundation
9 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::thermophysicalFunctions::nonUniformTable
29
30Description
31 Non-uniform tabulated property function that linearly interpolates between
32 the values.
33
34 To speed-up the search of the non-uniform table a uniform jump-table is
35 created on construction which is used for fast indirect addressing into
36 the table.
37
38Usage
39 \table
40 Property | Description
41 values | List of (temperature property) value pairs
42 \endTable
43
44 Example for the density of water between 280 and 350K
45 \verbatim
46 rho
47 {
48 type nonUniformTable;
49
50 values
51 (
52 (280 999.87)
53 (300 995.1)
54 (350 973.7)
55 );
56 }
57 \endverbatim
58
59\*---------------------------------------------------------------------------*/
60
61#ifndef nonUniformTableThermophysicalFunction_H
62#define nonUniformTableThermophysicalFunction_H
63
65#include "MinMax.H"
66#include "Tuple2.H"
67
68// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
69
70namespace Foam
71{
72
73/*---------------------------------------------------------------------------*\
74 Class nonUniformTable Declaration
75\*---------------------------------------------------------------------------*/
76
77class nonUniformTable
78:
79 public thermophysicalFunction
80{
81 // Private Data
82
83 //- Name of dictionary from which this function is instantiated
84 word name_;
85
86 //- Table values
87 List<Tuple2<scalar, scalar>> values_;
88
89 //- Lowest/highest temperature in the table
90 MinMax<scalar> Trange_;
91
92 //- Temperature increment, based on Trange_ and values_.size()
93 scalar deltaT_;
94
95 //- Lookup indices into values_ table based on deltaT_
96 List<label> jumpTable_;
97
98
99protected:
100
101 //- Return the lower index of the interval in the table
102 //- corresponding to the given temperature.
103 // Checks the temperature range
104 inline label index(scalar p, scalar T) const;
105
106
107public:
108
109 //- Runtime type information
110 TypeName("nonUniformTable");
111
112
113 // Constructors
114
115 //- Construct from entry name and dictionary
116 nonUniformTable(const word& name, const dictionary& dict);
117
118 //- Construct from dictionary, using "values" for the lookup
119 explicit nonUniformTable(const dictionary& dict);
120
121
122 // Member Functions
123
124 //- Return the non-uniform table of values
125 const List<Tuple2<scalar, scalar>>& values() const
126 {
127 return values_;
128 }
129
130 //- Evaluate the function and return the result
131 scalar f(scalar p, scalar T) const;
132
133 //- Evaluate the derivative of the function and return the result
134 scalar dfdT(scalar p, scalar T) const;
135
136 //- Write the function coefficients
137 void writeData(Ostream& os) const;
138};
139
140
141// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
142
143} // End namespace Foam
144
145// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
146
148
149// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
150
151#endif
152
153// ************************************************************************* //
volScalarField & p
const volScalarField & T
OBJstream os(runTime.globalPath()/outputName)
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:149
Namespace for OpenFOAM.
labelList f(nPoints)
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73
const bool writeData(pdfDictionary.get< bool >("writeData"))