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 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92
|
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2016 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 "phaseChangeTwoPhaseMixture.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
defineTypeNameAndDebug(phaseChangeTwoPhaseMixture, 0);
defineRunTimeSelectionTable(phaseChangeTwoPhaseMixture, components);
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::phaseChangeTwoPhaseMixture::phaseChangeTwoPhaseMixture
(
const word& type,
const volVectorField& U,
const surfaceScalarField& phi
)
:
incompressibleTwoPhaseMixture(U, phi),
phaseChangeTwoPhaseMixtureCoeffs_(subDict(type + "Coeffs")),
pSat_("pSat", dimPressure, lookup("pSat"))
{}
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
Foam::Pair<Foam::tmp<Foam::volScalarField>>
Foam::phaseChangeTwoPhaseMixture::vDotAlphal() const
{
volScalarField alphalCoeff(1.0/rho1() - alpha1_*(1.0/rho1() - 1.0/rho2()));
Pair<tmp<volScalarField>> mDotAlphal = this->mDotAlphal();
return Pair<tmp<volScalarField>>
(
alphalCoeff*mDotAlphal[0],
alphalCoeff*mDotAlphal[1]
);
}
Foam::Pair<Foam::tmp<Foam::volScalarField>>
Foam::phaseChangeTwoPhaseMixture::vDotP() const
{
dimensionedScalar pCoeff(1.0/rho1() - 1.0/rho2());
Pair<tmp<volScalarField>> mDotP = this->mDotP();
return Pair<tmp<volScalarField>>(pCoeff*mDotP[0], pCoeff*mDotP[1]);
}
bool Foam::phaseChangeTwoPhaseMixture::read()
{
if (incompressibleTwoPhaseMixture::read())
{
phaseChangeTwoPhaseMixtureCoeffs_ = subDict(type() + "Coeffs");
lookup("pSat") >> pSat_;
return true;
}
else
{
return false;
}
}
// ************************************************************************* //
|