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-2018 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  const dictionary& dict,
63  const objectRegistry& db
64 )
65 :
66  saturationModel(db)
67 {}
68 
69 
70 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
71 
73 {}
74 
75 
76 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
77 
80 (
81  const volScalarField& T
82 ) const
83 {
84  volScalarField TC(T - zeroC);
85 
86  return A*exp(TC*xByTC(TC));
87 }
88 
89 
92 (
93  const volScalarField& T
94 ) const
95 {
96  volScalarField TC(T - zeroC);
97 
98  volScalarField x(xByTC(TC));
99 
100  return A*exp(TC*x)*(D*x - TC/C)/(D + TC);
101 }
102 
103 
106 (
107  const volScalarField& T
108 ) const
109 {
110  volScalarField TC(T - zeroC);
111 
112  return log(A.value()) + TC*xByTC(TC);
113 }
114 
115 
118 (
119  const volScalarField& p
120 ) const
121 {
123 
124  return volScalarField::null();
125 }
126 
127 
128 // ************************************************************************* //
dictionary dict
dimensionedScalar log(const dimensionedScalar &ds)
virtual tmp< volScalarField > Tsat(const volScalarField &p) const
Saturation temperature.
Generic dimensioned Type class.
Macros for easy insertion into run-time selection tables.
virtual tmp< volScalarField > pSatPrime(const volScalarField &T) const
Saturation pressure derivetive w.r.t. temperature.
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
virtual ~ArdenBuck()
Destructor.
const dimensionSet dimPressure
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
volScalarField & C
ArdenBuck(const dictionary &dict, const objectRegistry &db)
Construct from a dictionary.
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 > lnPSat(const volScalarField &T) const
Natural log of the saturation pressure.
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:366
virtual tmp< volScalarField > pSat(const volScalarField &T) const
Saturation pressure.
Namespace for OpenFOAM.