openfoam/applications/solvers/multiphase/twoPhaseEulerFoam/twoPhaseSystem/phaseModel/phaseModel.C

174 lines
4.6 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2013 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 "phaseModel.H"
#include "twoPhaseSystem.H"
#include "diameterModel.H"
#include "dragModel.H"
#include "heatTransferModel.H"
#include "fixedValueFvPatchFields.H"
#include "slipFvPatchFields.H"
#include "partialSlipFvPatchFields.H"
#include "surfaceInterpolate.H"
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::phaseModel::phaseModel
(
const twoPhaseSystem& fluid,
const dictionary& phaseProperties,
const word& phaseName
)
:
volScalarField
(
IOobject
(
IOobject::groupName("alpha", phaseName),
fluid.mesh().time().timeName(),
fluid.mesh(),
IOobject::READ_IF_PRESENT,
IOobject::AUTO_WRITE
),
fluid.mesh(),
dimensionedScalar("alpha", dimless, 0)
),
fluid_(fluid),
name_(phaseName),
phaseDict_
(
phaseProperties.subDict(name_)
),
thermo_(rhoThermo::New(fluid.mesh(), name_)),
U_
(
IOobject
(
IOobject::groupName("U", name_),
fluid.mesh().time().timeName(),
fluid.mesh(),
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
fluid.mesh()
)
{
thermo_->validate("phaseModel " + name_, "h", "e");
const word phiName = IOobject::groupName("phi", name_);
IOobject phiHeader
(
phiName,
fluid_.mesh().time().timeName(),
fluid_.mesh(),
IOobject::NO_READ
);
if (phiHeader.headerOk())
{
Info<< "Reading face flux field " << phiName << endl;
phiPtr_.reset
(
new surfaceScalarField
(
IOobject
(
phiName,
fluid_.mesh().time().timeName(),
fluid_.mesh(),
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
fluid_.mesh()
)
);
}
else
{
Info<< "Calculating face flux field " << phiName << endl;
wordList phiTypes
(
U_.boundaryField().size(),
calculatedFvPatchScalarField::typeName
);
forAll(U_.boundaryField(), i)
{
if
(
isA<fixedValueFvPatchVectorField>(U_.boundaryField()[i])
|| isA<slipFvPatchVectorField>(U_.boundaryField()[i])
|| isA<partialSlipFvPatchVectorField>(U_.boundaryField()[i])
)
{
phiTypes[i] = fixedValueFvPatchScalarField::typeName;
}
}
phiPtr_.reset
(
new surfaceScalarField
(
IOobject
(
phiName,
fluid_.mesh().time().timeName(),
fluid_.mesh(),
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
fvc::interpolate(U_) & fluid_.mesh().Sf(),
phiTypes
)
);
}
dPtr_ = diameterModel::New
(
phaseDict_,
*this
);
}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam::phaseModel::~phaseModel()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
Foam::tmp<Foam::volScalarField> Foam::phaseModel::d() const
{
return dPtr_().d();
}
// ************************************************************************* //