AntoineExtended.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2015-2018 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "AntoineExtended.H"
30
31// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32
33namespace Foam
34{
35namespace saturationModels
36{
39 (
43 );
44}
45}
46
47
48// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49
51(
52 const dictionary& dict,
53 const objectRegistry& db
54)
55:
56 Antoine(dict, db),
57 D_("D", dimless, dict),
58 F_("F", dimless, dict),
59 E_("E", dimless/pow(dimTemperature, F_), dict)
60{}
61
62
63// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
64
66{}
67
68
69// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
70
73(
74 const volScalarField& T
75) const
76{
77 return
79 *exp(A_ + B_/(C_ + T) + E_*pow(T, F_))
80 *pow(T, D_);
81}
82
83
86(
87 const volScalarField& T
88) const
89{
90 return pSat(T)*((D_ + E_*F_*pow(T, F_))/T - B_/sqr(C_ + T));
91}
92
93
96(
97 const volScalarField& T
98) const
99{
100 return
101 A_
102 + B_/(C_ + T)
104 + E_*pow(T, F_);
105}
106
107
110(
111 const volScalarField& p
112) const
113{
115
116 return volScalarField::null();
117}
118
119
120// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
static const GeometricField< scalar, fvPatchField, volMesh > & null()
Return a null geometric field.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Registry of regIOobjects.
Extended Antoine equation for the vapour pressure.
virtual tmp< volScalarField > lnPSat(const volScalarField &T) const
Natural log of the saturation pressure.
virtual tmp< volScalarField > Tsat(const volScalarField &p) const
Saturation temperature.
virtual tmp< volScalarField > pSatPrime(const volScalarField &T) const
Saturation pressure derivetive w.r.t. temperature.
virtual tmp< volScalarField > pSat(const volScalarField &T) const
Saturation pressure.
Antoine equation for the vapour pressure.
Definition: Antoine.H:62
A class for managing temporary objects.
Definition: tmp.H:65
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
const volScalarField & T
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:517
Namespace for OpenFOAM.
const dimensionSet dimPressure
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar exp(const dimensionedScalar &ds)
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:54
dictionary dict