integratedNonUniformTableThermophysicalFunction.C
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 -------------------------------------------------------------------------------
11 License
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 
27 \*---------------------------------------------------------------------------*/
28 
30 #include "specie.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace thermophysicalFunctions
38 {
40 
42  (
46  );
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
55 (
56  const word& name,
57  const dictionary& dict
58 )
59 :
61  intf_(values().size()),
62  intfByT_(values().size())
63 {
64  intf_[0] = 0;
65  intfByT_[0] = 0;
66 
67  for (label i = 1; i < intf_.size(); ++i)
68  {
69  intf_[i] = intf_[i-1] + intfdT(0, values()[i].first());
70  intfByT_[i] = intfByT_[i-1] + intfByTdT(0, values()[i].first());
71  }
72 
73  const scalar intfStd = intfdT(Pstd, Tstd);
74  const scalar intfByTStd = intfByTdT(Pstd, Tstd);
75 
76  forAll(intf_, i)
77  {
78  intf_[i] -= intfStd;
79  intfByT_[i] -= intfByTStd;
80  }
81 }
82 
83 
86 (
87  const dictionary& dict
88 )
89 :
91 {}
92 
93 
94 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95 
97 (
98  scalar p,
99  scalar T
100 ) const
101 {
102  const label i = index(p, T);
103  const scalar Ti = values()[i].first();
104  const scalar fi = values()[i].second();
105  const scalar dT = T - Ti;
106  const scalar lambda = dT/(values()[i + 1].first() - Ti);
107 
108  return
109  intf_[i]
110  + (fi + 0.5*lambda*(values()[i + 1].second() - fi))*dT;
111 }
112 
113 
115 (
116  scalar p,
117  scalar T
118 ) const
119 {
120  const label i = index(p, T);
121  const scalar Ti = values()[i].first();
122  const scalar fi = values()[i].second();
123  const scalar gradf =
124  (values()[i + 1].second() - fi)/(values()[i + 1].first() - Ti);
125 
126  return
127  intfByT_[i] + ((fi - gradf*Ti)*log(T/Ti) + gradf*(T - Ti));
128 }
129 
130 
132 (
133  Ostream& os
134 ) const
135 {
137 }
138 
139 
140 // ************************************************************************* //
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::thermophysicalFunction
Abstract base class for thermo-physical functions.
Definition: thermophysicalFunction.H:53
Foam::constant::standard::Pstd
const dimensionedScalar Pstd
Standard pressure.
Definition: thermodynamicConstants.C:48
Foam::HashTableOps::values
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:149
writeData
const bool writeData(pdfDictionary.get< bool >("writeData"))
Foam::constant::standard::Tstd
const dimensionedScalar Tstd
Standard temperature.
Definition: thermodynamicConstants.C:49
specie.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::thermophysicalFunctions::addToRunTimeSelectionTable
addToRunTimeSelectionTable(thermophysicalFunction, integratedNonUniformTable, dictionary)
Foam::thermophysicalFunctions::integratedNonUniformTable::integratedNonUniformTable
integratedNonUniformTable(const word &name, const dictionary &dict)
Construct from entry name and dictionary.
Definition: integratedNonUniformTableThermophysicalFunction.C:55
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::thermophysicalFunctions::integratedNonUniformTable::writeData
void writeData(Ostream &os) const
Write the function coefficients.
Definition: integratedNonUniformTableThermophysicalFunction.C:132
integratedNonUniformTableThermophysicalFunction.H
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::thermophysicalFunctions::integratedNonUniformTable
Non-uniform tabulated property function that linearly interpolates between the values.
Definition: integratedNonUniformTableThermophysicalFunction.H:75
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
lambda
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:262
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
T
const volScalarField & T
Definition: createFieldRefs.H:2
Foam::thermophysicalFunctions::integratedNonUniformTable::intfByTdT
scalar intfByTdT(scalar p, scalar T) const
Integrate the function divided by T and return the result.
Definition: integratedNonUniformTableThermophysicalFunction.C:115
Foam::thermophysicalFunctions::defineTypeNameAndDebug
defineTypeNameAndDebug(integratedNonUniformTable, 0)
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::nonUniformTable
Definition: nonUniformTableThermophysicalFunction.H:76
Foam::thermophysicalFunctions::integratedNonUniformTable::intfdT
scalar intfdT(scalar p, scalar T) const
Integrate the function and return the result.
Definition: integratedNonUniformTableThermophysicalFunction.C:97