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-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 "ArdenBuck.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace saturationModels
34 {
35  defineTypeNameAndDebug(ArdenBuck, 0);
36  addToRunTimeSelectionTable(saturationModel, ArdenBuck, dictionary);
37 }
38 }
39 
40 
41 static const Foam::dimensionedScalar zeroC("", Foam::dimTemperature, 273.15);
42 static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21);
43 static const Foam::dimensionedScalar B("", Foam::dimless, 18.678);
44 static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5);
45 static const Foam::dimensionedScalar D("", Foam::dimTemperature, 257.14);
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
51 Foam::saturationModels::ArdenBuck::xByTC
52 (
53  const volScalarField& TC
54 ) const
55 {
56  return (B - TC/C)/(D + TC);
57 }
58 
59 
60 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
61 
63 (
64  const dictionary& dict,
65  const phaseInterface& interface
66 )
67 :
68  saturationModel(dict, interface)
69 {}
70 
71 
72 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
73 
75 {}
76 
77 
78 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79 
82 (
83  const volScalarField& T
84 ) const
85 {
86  volScalarField TC(T - zeroC);
87 
88  return A*exp(TC*xByTC(TC));
89 }
90 
91 
94 (
95  const volScalarField& T
96 ) const
97 {
98  volScalarField TC(T - zeroC);
99 
100  volScalarField x(xByTC(TC));
101 
102  return A*exp(TC*x)*(D*x - TC/C)/(D + TC);
103 }
104 
105 
108 (
109  const volScalarField& T
110 ) const
111 {
112  volScalarField TC(T - zeroC);
113 
114  return log(A.value()) + TC*xByTC(TC);
115 }
116 
117 
120 (
121  const volScalarField& p
122 ) const
123 {
125 
126  return volScalarField::null();
127 }
128 
129 
130 // ************************************************************************* //
dictionary dict
dimensionedScalar log(const dimensionedScalar &ds)
const dimensionSet dimPressure
virtual tmp< volScalarField > Tsat(const volScalarField &p) const
Saturation temperature.
const dimensionSet dimless
Generic dimensioned Type class.
Macros for easy insertion into run-time selection tables.
virtual tmp< volScalarField > pSatPrime(const volScalarField &T) const
Saturation pressure derivative w.r.t. temperature.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
dimensionedScalar exp(const dimensionedScalar &ds)
virtual ~ArdenBuck()
Destructor.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
static const GeometricField< scalar, fvPatchField, volMesh > & null()
Return a null geometric field.
virtual tmp< volScalarField > lnPSat(const volScalarField &T) const
Natural log of the saturation pressure.
ArdenBuck(const dictionary &dict, const phaseInterface &interface)
Construct from a dictionary and an interface.
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
virtual tmp< volScalarField > pSat(const volScalarField &T) const
Saturation pressure.
const dimensionSet dimTemperature
Namespace for OpenFOAM.