Antoine.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 "Antoine.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace saturationModels
34 {
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 template<class FieldType>
46 Foam::saturationModels::Antoine::pSat(const FieldType& T) const
47 {
48  return dimensionedScalar(dimPressure, 1)*exp(A_ + B_/(C_ + T));
49 }
50 
51 
52 template<class FieldType>
54 Foam::saturationModels::Antoine::pSatPrime(const FieldType& T) const
55 {
56  return - pSat(T)*B_/sqr(C_ + T);
57 }
58 
59 
60 template<class FieldType>
62 Foam::saturationModels::Antoine::lnPSat(const FieldType& T) const
63 {
64  return A_ + B_/(C_ + T);
65 }
66 
67 
68 template<class FieldType>
70 Foam::saturationModels::Antoine::Tsat(const FieldType& p) const
71 {
72  return B_/(log(p*dimensionedScalar(dimless/dimPressure, 1)) - A_) - C_;
73 }
74 
75 
76 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
77 
79 :
82  A_("A", dimless, dict),
83  B_("B", dimTemperature, dict),
84  C_("C", dimTemperature, dict)
85 {}
86 
87 
88 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
89 
91 {}
92 
93 
94 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95 
97 
98 
100 
101 
103 
104 
106 
107 
108 // ************************************************************************* //
IMPLEMENT_TSAT(Antoine, volScalarField::Internal)
IMPLEMENT_PSAT(Antoine, volScalarField::Internal)
Macros for easy insertion into run-time selection tables.
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
Antoine equation for the vapour pressure.
Definition: Antoine.H:62
Antoine(const dictionary &dict)
Construct from a dictionary.
Definition: Antoine.C:78
virtual ~Antoine()
Destructor.
Definition: Antoine.C:90
Model to describe the dependence of saturation pressure on temperature, and vice versa.
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
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)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict
volScalarField & p