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-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 "saturationModels.H"
27 #include "ArdenBuck.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace saturationModels
35 {
38 
39  static const coefficient zeroC("zeroC", dimTemperature, 273.15);
40  static const coefficient A("A", dimPressure, 611.21);
41  static const coefficient B("B", dimless, 18.678);
42  static const coefficient C("C", dimTemperature, 234.5);
43  static const coefficient D("D", dimTemperature, 257.14);
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
50 template<class FieldType>
52 Foam::saturationModels::ArdenBuck::xByTC(const FieldType& TC) const
53 {
54  return (B - TC/C)/(D + TC);
55 }
56 
57 
58 template<class FieldType>
60 Foam::saturationModels::ArdenBuck::pSat(const FieldType& T) const
61 {
62  const FieldType TC(T - zeroC);
63 
64  return A*exp(TC*xByTC(TC));
65 }
66 
67 
68 template<class FieldType>
70 Foam::saturationModels::ArdenBuck::pSatPrime(const FieldType& T) const
71 {
72  const FieldType TC(T - zeroC);
73  const FieldType x(xByTC(TC));
74 
75  return A*exp(TC*x)*(D*x - TC/C)/(D + TC);
76 }
77 
78 
79 template<class FieldType>
81 Foam::saturationModels::ArdenBuck::lnPSat(const FieldType& T) const
82 {
83  const FieldType TC(T - zeroC);
84 
85  return log(A.value()) + TC*xByTC(TC);
86 }
87 
88 
89 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
90 
92 :
94 {}
95 
96 
97 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
98 
100 {}
101 
102 
103 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
104 
105 IMPLEMENT_PSAT(saturationModels::ArdenBuck, scalarField);
106 
107 
108 IMPLEMENT_PSAT(saturationModels::ArdenBuck, volScalarField::Internal);
109 
110 
111 IMPLEMENT_PSAT(saturationModels::ArdenBuck, volScalarField);
112 
113 
114 // ************************************************************************* //
IMPLEMENT_PSAT(saturationModels::ArdenBuck, scalarField)
Macros for easy insertion into run-time selection tables.
Graphite solid properties.
Definition: C.H:51
DimensionedField< Type, GeoMesh, PrimitiveField > Internal
Type of the internal field from which this GeometricField is derived.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
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:91
virtual ~ArdenBuck()
Destructor.
Definition: ArdenBuck.C:99
Model to describe the dependence of saturation pressure on temperature.
A class for managing temporary objects.
Definition: tmp.H:55
static const coefficient zeroC("zeroC", dimTemperature, 273.15)
addToRunTimeSelectionTable(saturationPressureModel, Antoine, dictionary)
static const coefficient C("C", dimTemperature, 234.5)
static const coefficient D("D", dimTemperature, 257.14)
static const coefficient B("B", dimless, 18.678)
defineTypeNameAndDebug(Antoine, 0)
static const coefficient A("A", dimPressure, 611.21)
Namespace for OpenFOAM.
dimensionedScalar exp(const dimensionedScalar &ds)
const dimensionSet dimPressure
const dimensionSet dimless
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimTemperature
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
dictionary dict
Structure to store a dimensioned coefficient of the saturation model.