integratedNonUniformTableThermophysicalFunction.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 
27 #include "specie.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace thermophysicalFunctions
35 {
36  defineTypeNameAndDebug(integratedNonUniformTable, 0);
37 
39  (
40  thermophysicalFunction,
41  integratedNonUniformTable,
42  dictionary
43  );
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
52 (
53  const word& name,
54  const dictionary& dict
55 )
56 :
57  nonUniformTable(name, dict),
58  intf_(values().size()),
59  intfByT_(values().size())
60 {
61  intf_[0] = 0;
62  intfByT_[0] = 0;
63 
64  for(label i = 0; i<intf_.size() - 1; i++)
65  {
66  intf_[i + 1] = intf_[i] + intfdT(0, values()[i + 1].first());
67  intfByT_[i + 1] = intfByT_[i] + intfByTdT(0, values()[i + 1].first());
68  }
69 
70  const scalar intfStd = intfdT(Pstd, Tstd);
71  const scalar intfByTStd = intfByTdT(Pstd, Tstd);
72 
73  forAll(intf_, i)
74  {
75  intf_[i] -= intfStd;
76  intfByT_[i] -= intfByTStd;
77  }
78 }
79 
80 
83 (
84  const dictionary& dict
85 )
86 :
87  integratedNonUniformTable("values", dict)
88 {}
89 
90 
91 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
92 
94 (
95  scalar p,
96  scalar T
97 ) const
98 {
99  const label i = index(p, T);
100  const scalar Ti = values()[i].first();
101  const scalar fi = values()[i].second();
102  const scalar dT = T - Ti;
103  const scalar lambda = dT/(values()[i + 1].first() - Ti);
104 
105  return
106  intf_[i]
107  + (fi + 0.5*lambda*(values()[i + 1].second() - fi))*dT;
108 }
109 
110 
112 (
113  scalar p,
114  scalar T
115 ) const
116 {
117  const label i = index(p, T);
118  const scalar Ti = values()[i].first();
119  const scalar fi = values()[i].second();
120  const scalar gradf =
121  (values()[i + 1].second() - fi)/(values()[i + 1].first() - Ti);
122 
123  return
124  intfByT_[i] + ((fi - gradf*Ti)*log(T/Ti) + gradf*(T - Ti));
125 }
126 
127 
129 (
130  Ostream& os
131 ) const
132 {
134 }
135 
136 
137 // ************************************************************************* //
scalar intfdT(scalar p, scalar T) const
Integrate the function and return the result.
#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)
dimensionedScalar log(const dimensionedScalar &ds)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
const dimensionedScalar & Tstd
Standard temperature.
integratedNonUniformTable(const word &name, const dictionary &dict)
Construct from entry name and dictionary.
Macros for easy insertion into run-time selection tables.
scalar intfByTdT(scalar p, scalar T) const
Integrate the function divided by T and return the result.
A class for handling words, derived from string.
Definition: word.H:59
addToRunTimeSelectionTable(thermophysicalFunction, APIdiffCoef, dictionary)
Non-uniform tabulated property function that linearly interpolates between the 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
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
Non-uniform tabulated property function that linearly interpolates between the values.
const dimensionedScalar & Pstd
Standard pressure.
Namespace for OpenFOAM.