integratedNonUniformTableThermophysicalFunction.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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::thermophysicalFunctions::integratedNonUniformTable
28 
29 Description
30  Non-uniform tabulated property function that linearly interpolates between
31  the values.
32 
33  To speed-up the search of the non-uniform table a uniform jump-table is
34  created on construction which is used for fast indirect addressing into
35  the table.
36 
37 Usage
38  \table
39  Property | Description
40  values | List of (temperature property) value pairs
41  \endTable
42 
43  Example for the density of water between 280 and 350K
44  \verbatim
45  rho
46  {
47  type integratedNonUniformTable;
48 
49  values
50  (
51  (280 999.87)
52  (300 995.1)
53  (350 973.7)
54  );
55  }
56  \endverbatim
57 
58 \*---------------------------------------------------------------------------*/
59 
60 #ifndef integratedNonUniformTableThermophysicalFunction_H
61 #define integratedNonUniformTableThermophysicalFunction_H
62 
64 
65 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66 
67 namespace Foam
68 {
69 namespace thermophysicalFunctions
70 {
71 
72 /*---------------------------------------------------------------------------*\
73  Class integratedNonUniformTable Declaration
74 \*---------------------------------------------------------------------------*/
75 
76 class integratedNonUniformTable
77 :
78  public nonUniformTable
79 {
80  // Private Data
81 
82  List<scalar> intf_;
83 
84  List<scalar> intfByT_;
85 
86 
87 public:
88 
89  //- Runtime type information
90  TypeName("integratedNonUniformTable");
91 
92 
93  // Constructors
94 
95  //- Construct from entry name and dictionary
96  integratedNonUniformTable(const word& name, const dictionary& dict);
97 
98  //- Construct from dictionary, using "values" for the lookup
99  explicit integratedNonUniformTable(const dictionary& dict);
100 
101 
102  // Member Functions
103 
104  //- Integrate the function and return the result
105  scalar intfdT(scalar p, scalar T) const;
106 
107  //- Integrate the function divided by T and return the result
108  scalar intfByTdT(scalar p, scalar T) const;
109 
110  //- Write the function coefficients
111  void writeData(Ostream& os) const;
112 };
113 
114 
115 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
116 
117 } // End namespace thermophysicalFunctions
118 } // End namespace Foam
119 
120 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
121 
122 #endif
123 
124 // ************************************************************************* //
dictionary dict
static void writeData(Ostream &os, const Type &val)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:69
OBJstream os(runTime.globalPath()/outputName)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
volScalarField & p
Namespace for OpenFOAM.