openfoam/applications/solvers/multiphase/multiphaseEulerFoam/multiphaseSystem/phaseModel/phaseModel.C

248 lines
5.7 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / 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 "phaseModel.H"
#include "diameterModel.H"
#include "fixedValueFvPatchFields.H"
#include "slipFvPatchFields.H"
#include "partialSlipFvPatchFields.H"
#include "surfaceInterpolate.H"
#include "fvcFlux.H"
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::phaseModel::phaseModel
(
const word& phaseName,
const dictionary& phaseDict,
const fvMesh& mesh
)
:
volScalarField
(
IOobject
(
IOobject::groupName("alpha", phaseName),
mesh.time().timeName(),
mesh,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
mesh
),
name_(phaseName),
phaseDict_(phaseDict),
nu_
(
"nu",
dimensionSet(0, 2, -1, 0, 0),
phaseDict_
),
kappa_
(
"kappa",
dimensionSet(1, 1, -3, -1, 0),
phaseDict_
),
Cp_
(
"Cp",
dimensionSet(0, 2, -2, -1, 0),
phaseDict_
),
rho_
(
"rho",
dimDensity,
phaseDict_
),
U_
(
IOobject
(
IOobject::groupName("U", phaseName),
mesh.time().timeName(),
mesh,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
mesh
),
DDtU_
(
IOobject
(
IOobject::groupName("DDtU", phaseName),
mesh.time().timeName(),
mesh
),
mesh,
dimensionedVector("0", dimVelocity/dimTime, Zero)
),
alphaPhi_
(
IOobject
(
IOobject::groupName("alphaPhi", phaseName),
mesh.time().timeName(),
mesh
),
mesh,
dimensionedScalar("0", dimensionSet(0, 3, -1, 0, 0), 0)
)
{
alphaPhi_.setOriented();
const word phiName = IOobject::groupName("phi", name_);
IOobject phiHeader
(
phiName,
mesh.time().timeName(),
mesh,
IOobject::NO_READ
);
if (phiHeader.typeHeaderOk<surfaceScalarField>(true))
{
Info<< "Reading face flux field " << phiName << endl;
phiPtr_.reset
(
new surfaceScalarField
(
IOobject
(
phiName,
mesh.time().timeName(),
mesh,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
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,
mesh.time().timeName(),
mesh,
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
fvc::flux(U_),
phiTypes
)
);
}
dPtr_ = diameterModel::New
(
phaseDict_,
*this
);
}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam::phaseModel::~phaseModel()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
Foam::autoPtr<Foam::phaseModel> Foam::phaseModel::clone() const
{
NotImplemented;
return autoPtr<phaseModel>(nullptr);
}
void Foam::phaseModel::correct()
{
//nuModel_->correct();
}
bool Foam::phaseModel::read(const dictionary& phaseDict)
{
phaseDict_ = phaseDict;
//if (nuModel_->read(phaseDict_))
{
phaseDict_.lookup("nu") >> nu_.value();
phaseDict_.lookup("kappa") >> kappa_.value();
phaseDict_.lookup("Cp") >> Cp_.value();
phaseDict_.lookup("rho") >> rho_.value();
return true;
}
// else
// {
// return false;
// }
return true;
}
Foam::tmp<Foam::volScalarField> Foam::phaseModel::d() const
{
return dPtr_().d();
}
// ************************************************************************* //