polynomialTemperature.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) 2015-2025 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 
26 #include "polynomialTemperature.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace saturationModels
34 {
37  (
41  );
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
49 (
50  const dictionary& dict
51 )
52 :
54  C_(dict.lookup("C<8>"))
55 {}
56 
57 
58 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
59 
61 {}
62 
63 
64 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
65 
68 (
69  const scalarField& p
70 ) const
71 {
72  tmp<scalarField> tTsat(new scalarField(p.size(), scalar(0)));
73  scalarField& Tsat = tTsat.ref();
74 
75  forAll(Tsat, celli)
76  {
77  Tsat[celli] = C_.value(p[celli]);
78  }
79 
80  return tTsat;
81 }
82 
83 
86 (
87  const scalarField& p
88 ) const
89 {
90  tmp<scalarField> tTsatPrime(new scalarField(p.size(), scalar(0)));
91  scalarField& TsatPrime = tTsatPrime.ref();
92 
93  forAll(TsatPrime, celli)
94  {
95  TsatPrime[celli] = C_.derivative(p[celli]);
96  }
97 
98  return tTsatPrime;
99 }
100 
101 
102 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Polynomial equation for the saturation vapour temperature in terms of the vapour pressure (in Pa).
virtual tmp< scalarField > Tsat(const scalarField &p) const
Saturation temperature for scalarField.
virtual tmp< scalarField > TsatPrime(const scalarField &p) const
Saturation temperature derivative w.r.t. pressure for scalarField.
polynomialTemperature(const dictionary &dict)
Construct from a dictionary.
Model to describe the dependence of saturation temperature on pressure.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
addToRunTimeSelectionTable(saturationPressureModel, Antoine, dictionary)
defineTypeNameAndDebug(Antoine, 0)
Namespace for OpenFOAM.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dictionary dict
volScalarField & p