ArdenBuck.C 3.38 KB
Newer Older
1 2 3 4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
OpenFOAM bot's avatar
OpenFOAM bot committed
5
    \\  /    A nd           | www.openfoam.com
6
     \\/     M anipulation  |
OpenFOAM bot's avatar
OpenFOAM bot committed
7
-------------------------------------------------------------------------------
OpenFOAM bot's avatar
OpenFOAM bot committed
8
    Copyright (C) 2015-2018 OpenFOAM Foundation
9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34
-------------------------------------------------------------------------------
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
{
35
namespace saturationModels
36 37
{
    defineTypeNameAndDebug(ArdenBuck, 0);
38
    addToRunTimeSelectionTable(saturationModel, ArdenBuck, dictionary);
39 40 41 42 43 44 45 46 47 48 49 50
}
}

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>
51
Foam::saturationModels::ArdenBuck::xByTC
52 53 54 55 56 57 58 59 60 61
(
    const volScalarField& TC
) const
{
    return (B - TC/C)/(D + TC);
}


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

62 63 64 65 66
Foam::saturationModels::ArdenBuck::ArdenBuck
(
    const dictionary& dict,
    const objectRegistry& db
)
67
:
68
    saturationModel(db)
69 70 71 72 73
{}


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

74
Foam::saturationModels::ArdenBuck::~ArdenBuck()
75 76 77 78 79 80
{}


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

Foam::tmp<Foam::volScalarField>
81
Foam::saturationModels::ArdenBuck::pSat
82 83 84 85 86 87 88 89 90 91 92
(
    const volScalarField& T
) const
{
    volScalarField TC(T - zeroC);

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


Foam::tmp<Foam::volScalarField>
93
Foam::saturationModels::ArdenBuck::pSatPrime
94 95 96 97 98 99 100 101 102 103 104 105 106
(
    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>
107
Foam::saturationModels::ArdenBuck::lnPSat
108 109 110 111 112 113 114 115 116 117
(
    const volScalarField& T
) const
{
    volScalarField TC(T - zeroC);

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


118 119 120 121 122 123
Foam::tmp<Foam::volScalarField>
Foam::saturationModels::ArdenBuck::Tsat
(
    const volScalarField& p
) const
{
124
    NotImplemented;
125 126 127 128 129

    return volScalarField::null();
}


130
// ************************************************************************* //