nonUniformTableThermophysicalFunction.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
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
36 defineTypeNameAndDebug(nonUniformTable, 0);
37
39 (
41 nonUniformTable,
43 );
44}
45
46
47// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48
49Foam::nonUniformTable::nonUniformTable
50(
51 const word& name,
52 const dictionary& dict
53)
54:
55 name_(name),
56 values_(),
57 Trange_(),
58 deltaT_(GREAT),
59 jumpTable_()
60{
61 dict.readEntry(name_, values_);
62
63 if (values_.size() < 2)
64 {
66 << "Table" << nl
67 << " " << name_ << nl
68 << " has fewer than 2 entries." << nl
70 }
71
72 Trange_.min() = values_.first().first();
73 Trange_.max() = values_.last().first();
74
75 for (label i = 1; i < values_.size(); ++i)
76 {
77 #ifdef FULLDEBUG
78 // Check list is monotonically increasing...
79 if (values_[i].first() <= values_[i-1].first())
80 {
82 << "Table" << nl
83 << " " << name_ << nl
84 << " out-of-order value: " << values_[i].first()
85 << " at index " << i << nl
86 << exit(FatalError);
87 }
88 #endif
89
90 deltaT_ = min(deltaT_, values_[i].first() - values_[i-1].first());
91 }
92
93 deltaT_ *= 0.9;
94
95 jumpTable_.resize(Trange_.mag()/deltaT_ + 1);
96
97 label i = 0;
98 forAll(jumpTable_, j)
99 {
100 const scalar T = Trange_.min() + j*deltaT_;
101
102 if (T > values_[i+1].first())
103 {
104 ++i;
105 }
106
107 jumpTable_[j] = i;
108 }
109}
110
111
112Foam::nonUniformTable::nonUniformTable
113(
114 const dictionary& dict
115)
116:
117 nonUniformTable("values", dict)
118{}
119
120
121// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
122
123Foam::scalar Foam::nonUniformTable::f
124(
125 scalar p,
126 scalar T
127) const
128{
129 const label i = index(p, T);
130 const scalar Ti = values_[i].first();
131 const scalar lambda = (T - Ti)/(values_[i + 1].first() - Ti);
132
133 return
134 values_[i].second()
135 + lambda*(values_[i + 1].second() - values_[i].second());
136}
137
138
139Foam::scalar Foam::nonUniformTable::dfdT
140(
141 scalar p,
142 scalar T
143) const
144{
145 const label i = index(p, T);
146
147 return
148 (values_[i + 1].second() - values_[i].second())
149 /(values_[i + 1].first() - values_[i].first());
150}
151
152
153void Foam::nonUniformTable::writeData(Ostream& os) const
154{
155 os.writeEntry("values", values_);
156}
157
158
159// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
void min(const dimensioned< Type > &dt)
Use the minimum of the field and specified value.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
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
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
IOerror FatalIOError
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333