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 {
35  defineTypeNameAndDebug(Antoine, 0);
36  addToRunTimeSelectionTable(saturationModel, Antoine, dictionary);
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 (
45  const dictionary& dict,
46  const phaseInterface& interface
47 )
48 :
49  saturationModel(dict, interface),
50  A_("A", dimless, dict),
51  B_("B", dimTemperature, dict),
52  C_("C", dimTemperature, dict)
53 {}
54 
55 
56 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
57 
59 {}
60 
61 
62 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
63 
66 (
67  const volScalarField& T
68 ) const
69 {
70  return
72  *exp(A_ + B_/(C_ + T));
73 }
74 
75 
78 (
79  const volScalarField& T
80 ) const
81 {
82  return - pSat(T)*B_/sqr(C_ + T);
83 }
84 
85 
88 (
89  const volScalarField& T
90 ) const
91 {
92  return A_ + B_/(C_ + T);
93 }
94 
95 
98 (
99  const volScalarField& p
100 ) const
101 {
102  return
103  B_/(log(p*dimensionedScalar(dimless/dimPressure, 1)) - A_)
104  - C_;
105 }
106 
107 
108 // ************************************************************************* //
dictionary dict
dimensionedScalar log(const dimensionedScalar &ds)
Antoine(const dictionary &dict, const phaseInterface &interface)
Construct from a dictionary and an interface.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimPressure
const dimensionSet dimless
Macros for easy insertion into run-time selection tables.
virtual tmp< volScalarField > pSat(const volScalarField &T) const
Saturation pressure.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
virtual tmp< volScalarField > lnPSat(const volScalarField &T) const
Natural log of the saturation pressure.
dimensionedScalar exp(const dimensionedScalar &ds)
virtual tmp< volScalarField > pSatPrime(const volScalarField &T) const
Saturation pressure derivative w.r.t. temperature.
virtual ~Antoine()
Destructor.
virtual tmp< volScalarField > Tsat(const volScalarField &p) const
Saturation temperature.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionSet dimTemperature
Namespace for OpenFOAM.