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-------------------------------------------------------------------------------
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
27\*---------------------------------------------------------------------------*/
28
30#include "specie.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
36{
37namespace thermophysicalFunctions
38{
39 defineTypeNameAndDebug(integratedNonUniformTable, 0);
40
42 (
44 integratedNonUniformTable,
46 );
47}
48}
49
50
51// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52
53Foam::thermophysicalFunctions::integratedNonUniformTable::
54integratedNonUniformTable
55(
56 const word& name,
57 const dictionary& dict
58)
59:
60 nonUniformTable(name, dict),
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
84Foam::thermophysicalFunctions::integratedNonUniformTable::
85integratedNonUniformTable
86(
87 const dictionary& dict
88)
89:
90 integratedNonUniformTable("values", dict)
91{}
92
93
94// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95
96Foam::scalar Foam::thermophysicalFunctions::integratedNonUniformTable::intfdT
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
114Foam::scalar Foam::thermophysicalFunctions::integratedNonUniformTable::intfByTdT
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
131void Foam::thermophysicalFunctions::integratedNonUniformTable::writeData
132(
133 Ostream& os
134) const
135{
136 nonUniformTable::writeData(os);
137}
138
139
140// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Abstract base class for thermo-physical functions.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
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.
dimensionedScalar log(const dimensionedScalar &ds)
dictionary dict
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333