All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
tableThermophysicalFunction.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) 2019-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 
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace thermophysicalFunctions
34 {
35  defineTypeNameAndDebug(table, 0);
36  addToRunTimeSelectionTable(thermophysicalFunction, table, dictionary);
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 :
45  dictName_(dict.name()),
46  Tlow_(dict.lookup<scalar>("Tlow")),
47  Thigh_(dict.lookup<scalar>("Thigh")),
48  values_(dict.lookup("values"))
49 {
50  if (values_.size() < 2)
51  {
53  << "Table " << nl
54  << " " << dictName_ << nl
55  << " has less than 2 entries."
56  << exit(FatalError);
57  }
58  else
59  {
60  deltaT_ = (Thigh_ - Tlow_)/(values_.size() - 1);
61  }
62 }
63 
64 
65 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
66 
67 Foam::scalar Foam::thermophysicalFunctions::table::f(scalar p, scalar T) const
68 {
69  const scalar nd = (T - Tlow_)/deltaT_;
70  const label i = nd;
71 
72  if (nd < 0 || i > values_.size() - 2)
73  {
75  << "Temperature " << T << " out of range "
76  << Tlow_ << " to " << Thigh_ << nl
77  << " of table " << dictName_
78  << exit(FatalError);
79  }
80 
81  const scalar Ti = Tlow_ + i*deltaT_;
82  const scalar lambda = (T - Ti)/deltaT_;
83 
84  return values_[i] + lambda*(values_[i + 1] - values_[i]);
85 }
86 
87 
89 {
90  writeEntry(os, "Tlow", Tlow_);
91  writeEntry(os, "Thigh", Thigh_);
92  writeEntry(os, "values", values_);
93 }
94 
95 
96 // ************************************************************************* //
dictionary dict
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)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
table(const dictionary &dict)
Construct from dictionary.
void write(Ostream &os) const
Write the function coefficients.
Macros for easy insertion into run-time selection tables.
stressControl lookup("compactNormalStress") >> compactNormalStress
addToRunTimeSelectionTable(thermophysicalFunction, APIdiffCoef, dictionary)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
static const char nl
Definition: Ostream.H:260
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
volScalarField & p
scalar f(scalar p, scalar T) const
Evaluate the function and return the result.
Namespace for OpenFOAM.