ArdenBuck.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 "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 static const Foam::dimensionedScalar zeroC("", Foam::dimTemperature, 273.15);
41 static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21);
42 static const Foam::dimensionedScalar B("", Foam::dimless, 18.678);
43 static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5);
44 static const Foam::dimensionedScalar D("", Foam::dimTemperature, 257.14);
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
49 Foam::saturationModels::ArdenBuck::xByTC
50 (
51  const volScalarField& TC
52 ) const
53 {
54  return (B - TC/C)/(D + TC);
55 }
56 
57 
58 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
59 
61 :
62  saturationModel()
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
67 
69 {}
70 
71 
72 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
73 
76 (
77  const volScalarField& T
78 ) const
79 {
80  volScalarField TC(T - zeroC);
81 
82  return A*exp(TC*xByTC(TC));
83 }
84 
85 
88 (
89  const volScalarField& T
90 ) const
91 {
92  volScalarField TC(T - zeroC);
93 
94  volScalarField x(xByTC(TC));
95 
96  return A*exp(TC*x)*(D*x - TC/C)/(D + TC);
97 }
98 
99 
102 (
103  const volScalarField& T
104 ) const
105 {
106  volScalarField TC(T - zeroC);
107 
108  return log(A.value()) + TC*xByTC(TC);
109 }
110 
111 
114 (
115  const volScalarField& p
116 ) const
117 {
119 
120  return volScalarField::null();
121 }
122 
123 
124 // ************************************************************************* //
dictionary dict
dimensionedScalar log(const dimensionedScalar &ds)
ArdenBuck(const dictionary &dict)
Construct from a dictionary.
Generic dimensioned Type class.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dimensionedScalar exp(const dimensionedScalar &ds)
virtual tmp< volScalarField > lnPSat(const volScalarField &T) const
Natural log of the saturation pressure.
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
virtual ~ArdenBuck()
Destructor.
virtual tmp< volScalarField > pSat(const volScalarField &T) const
Saturation pressure.
const dimensionSet dimPressure
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
volScalarField & C
static const GeometricField< scalar, fvPatchField, volMesh > & null()
Return a null geometric field.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
virtual tmp< volScalarField > Tsat(const volScalarField &p) const
Saturation temperature.
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:54
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:366
Namespace for OpenFOAM.
virtual tmp< volScalarField > pSatPrime(const volScalarField &T) const
Saturation pressure derivetive w.r.t. temperature.