AntoineExtended.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 "AntoineExtended.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace saturationModels
34 {
35  defineTypeNameAndDebug(AntoineExtended, 0);
37  (
38  saturationModel,
39  AntoineExtended,
40  dictionary
41  );
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
49 (
50  const dictionary& dict,
51  const phaseInterface& interface
52 )
53 :
54  Antoine(dict, interface),
55  D_("D", dimless, dict),
56  F_("F", dimless, dict),
57  E_("E", dimless/pow(dimTemperature, F_), dict)
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
62 
64 {}
65 
66 
67 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
68 
71 (
72  const volScalarField& T
73 ) const
74 {
75  return
77  *exp(A_ + B_/(C_ + T) + E_*pow(T, F_))
78  *pow(T, D_);
79 }
80 
81 
84 (
85  const volScalarField& T
86 ) const
87 {
88  return pSat(T)*((D_ + E_*F_*pow(T, F_))/T - B_/sqr(C_ + T));
89 }
90 
91 
94 (
95  const volScalarField& T
96 ) const
97 {
98  return
99  A_
100  + B_/(C_ + T)
102  + E_*pow(T, F_);
103 }
104 
105 
108 (
109  const volScalarField& p
110 ) const
111 {
113 
114  return volScalarField::null();
115 }
116 
117 
118 // ************************************************************************* //
dictionary dict
dimensionedScalar log(const dimensionedScalar &ds)
virtual ~AntoineExtended()
Destructor.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimPressure
dimensionedScalar C_
Constant C.
Definition: Antoine.H:72
const dimensionSet dimless
virtual tmp< volScalarField > pSat(const volScalarField &T) const
Saturation pressure.
virtual tmp< volScalarField > pSatPrime(const volScalarField &T) const
Saturation pressure derivative w.r.t. temperature.
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
AntoineExtended(const dictionary &dict, const phaseInterface &interface)
Construct from a dictionary and an interface.
dimensionedScalar exp(const dimensionedScalar &ds)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
virtual tmp< volScalarField > lnPSat(const volScalarField &T) const
Natural log of the saturation pressure.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
static const GeometricField< scalar, fvPatchField, volMesh > & null()
Return a null geometric field.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual tmp< volScalarField > Tsat(const volScalarField &p) const
Saturation temperature.
dimensionedScalar B_
Constant B.
Definition: Antoine.H:69
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:353
dimensionedScalar A_
Constant A.
Definition: Antoine.H:66
const dimensionSet dimTemperature
Namespace for OpenFOAM.