Antoine.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2015 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  saturationModel(),
46  A_("A", dimless, dict),
47  B_("B", dimTemperature, dict),
48  C_("C", dimTemperature, dict)
49 {}
50 
51 
52 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
53 
55 {}
56 
57 
58 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
59 
62 (
63  const volScalarField& T
64 ) const
65 {
66  return
68  *exp(A_ + B_/(C_ + T));
69 }
70 
71 
74 (
75  const volScalarField& T
76 ) const
77 {
78  return - pSat(T)*B_/sqr(C_ + T);
79 }
80 
81 
84 (
85  const volScalarField& T
86 ) const
87 {
88  return A_ + B_/(C_ + T);
89 }
90 
91 
94 (
95  const volScalarField& p
96 ) const
97 {
98  return
99  B_/(log(p*dimensionedScalar("one", dimless/dimPressure, 1)) - A_)
100  - C_;
101 }
102 
103 
104 // ************************************************************************* //
dictionary dict
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual tmp< volScalarField > pSatPrime(const volScalarField &T) const
Saturation pressure derivetive w.r.t. temperature.
virtual tmp< volScalarField > lnPSat(const volScalarField &T) const
Natural log of the saturation pressure.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dimensionedScalar exp(const dimensionedScalar &ds)
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
Antoine(const dictionary &dict)
Construct from a dictionary.
const dimensionSet dimPressure
virtual ~Antoine()
Destructor.
virtual tmp< volScalarField > Tsat(const volScalarField &p) const
Saturation temperature.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
virtual tmp< volScalarField > pSat(const volScalarField &T) const
Saturation pressure.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: PtrList.H:54
Namespace for OpenFOAM.