ArdenBuck.C 3.24 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
    \\  /    A nd           | Copyright (C) 2015 OpenFOAM Foundation
     \\/     M anipulation  |
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

    OpenFOAM is free software: you can redistribute it and/or modify it
    under the terms of the GNU General Public License as published by
    the Free Software Foundation, either version 3 of the License, or
    (at your option) any later version.

    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
    for more details.

    You should have received a copy of the GNU General Public License
    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.

\*---------------------------------------------------------------------------*/

#include "ArdenBuck.H"
#include "addToRunTimeSelectionTable.H"

// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //

namespace Foam
{
33
namespace saturationModels
34
35
{
    defineTypeNameAndDebug(ArdenBuck, 0);
36
    addToRunTimeSelectionTable(saturationModel, ArdenBuck, dictionary);
37
38
39
40
41
42
43
44
45
46
47
48
}
}

static const Foam::dimensionedScalar zeroC("", Foam::dimTemperature, 273.15);
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21);
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678);
static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5);
static const Foam::dimensionedScalar D("", Foam::dimTemperature, 257.14);

// * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //

Foam::tmp<Foam::volScalarField>
49
Foam::saturationModels::ArdenBuck::xByTC
50
51
52
53
54
55
56
57
58
59
(
    const volScalarField& TC
) const
{
    return (B - TC/C)/(D + TC);
}


// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //

60
Foam::saturationModels::ArdenBuck::ArdenBuck(const dictionary& dict)
61
:
62
    saturationModel()
63
64
65
66
67
{}


// * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //

68
Foam::saturationModels::ArdenBuck::~ArdenBuck()
69
70
71
72
73
74
{}


// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //

Foam::tmp<Foam::volScalarField>
75
Foam::saturationModels::ArdenBuck::pSat
76
77
78
79
80
81
82
83
84
85
86
(
    const volScalarField& T
) const
{
    volScalarField TC(T - zeroC);

    return A*exp(TC*xByTC(TC));
}


Foam::tmp<Foam::volScalarField>
87
Foam::saturationModels::ArdenBuck::pSatPrime
88
89
90
91
92
93
94
95
96
97
98
99
100
(
    const volScalarField& T
) const
{
    volScalarField TC(T - zeroC);

    volScalarField x(xByTC(TC));

    return A*exp(TC*x)*(D*x - TC/C)/(D + TC);
}


Foam::tmp<Foam::volScalarField>
101
Foam::saturationModels::ArdenBuck::lnPSat
102
103
104
105
106
107
108
109
110
111
(
    const volScalarField& T
) const
{
    volScalarField TC(T - zeroC);

    return log(A.value()) + TC*xByTC(TC);
}


112
113
114
115
116
117
Foam::tmp<Foam::volScalarField>
Foam::saturationModels::ArdenBuck::Tsat
(
    const volScalarField& p
) const
{
118
    NotImplemented;
119
120
121
122
123

    return volScalarField::null();
}


124
// ************************************************************************* //