ArdenBuck.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-2022 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 "ArdenBuck.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace saturationModels
34 {
37 }
38 }
39 
40 
42 (
43  "zeroC",
45  273.15
46 );
47 
48 
49 static const Foam::dimensionedScalar A("A", Foam::dimPressure, 611.21);
50 static const Foam::dimensionedScalar B("B", Foam::dimless, 18.678);
52 static const Foam::dimensionedScalar D("D", Foam::dimTemperature, 257.14);
53 
54 
55 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
56 
57 template<class FieldType>
59 Foam::saturationModels::ArdenBuck::xByTC(const FieldType& TC) const
60 {
61  return (B - TC/C)/(D + TC);
62 }
63 
64 
65 template<class FieldType>
67 Foam::saturationModels::ArdenBuck::pSat(const FieldType& T) const
68 {
69  const FieldType TC(T - zeroC);
70 
71  return A*exp(TC*xByTC(TC));
72 }
73 
74 
75 template<class FieldType>
77 Foam::saturationModels::ArdenBuck::pSatPrime(const FieldType& T) const
78 {
79  const FieldType TC(T - zeroC);
80  const FieldType x(xByTC(TC));
81 
82  return A*exp(TC*x)*(D*x - TC/C)/(D + TC);
83 }
84 
85 
86 template<class FieldType>
88 Foam::saturationModels::ArdenBuck::lnPSat(const FieldType& T) const
89 {
90  const FieldType TC(T - zeroC);
91 
92  return log(A.value()) + TC*xByTC(TC);
93 }
94 
95 
96 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
97 
99 :
101 {}
102 
103 
104 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
105 
107 {}
108 
109 
110 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111 
113 
114 
116 
117 
118 // ************************************************************************* //
static const Foam::dimensionedScalar A("A", Foam::dimPressure, 611.21)
IMPLEMENT_PSAT(ArdenBuck, volScalarField::Internal)
static const Foam::dimensionedScalar B("B", Foam::dimless, 18.678)
static const Foam::dimensionedScalar D("D", Foam::dimTemperature, 257.14)
static const Foam::dimensionedScalar zeroC("zeroC", Foam::dimTemperature, 273.15)
static const Foam::dimensionedScalar C("C", Foam::dimTemperature, 234.5)
Macros for easy insertion into run-time selection tables.
Graphite solid properties.
Definition: C.H:51
DimensionedField< Type, GeoMesh > Internal
Type of the internal field from which this GeometricField is derived.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const Type & value() const
Return const reference to value.
ArdenBuck equation for the vapour pressure of moist air.
Definition: ArdenBuck.H:53
ArdenBuck(const dictionary &dict)
Construct from a dictionary.
Definition: ArdenBuck.C:98
virtual ~ArdenBuck()
Destructor.
Definition: ArdenBuck.C:106
Model to describe the dependence of saturation pressure on temperature, and vice versa.
addToRunTimeSelectionTable(saturationPressureModel, Antoine, dictionary)
defineTypeNameAndDebug(Antoine, 0)
Namespace for OpenFOAM.
dimensionedScalar exp(const dimensionedScalar &ds)
const dimensionSet dimPressure
const dimensionSet dimless
const dimensionSet dimTemperature
dimensionedScalar log(const dimensionedScalar &ds)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict