/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2015-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 .
\*---------------------------------------------------------------------------*/
#include "MovingPhaseModel.H"
#include "phaseSystem.H"
#include "phaseCompressibleTurbulenceModel.H"
#include "fixedValueFvPatchFields.H"
#include "slipFvPatchFields.H"
#include "partialSlipFvPatchFields.H"
#include "fvmDdt.H"
#include "fvmDiv.H"
#include "fvmSup.H"
#include "fvcDdt.H"
#include "fvcDiv.H"
#include "fvcFlux.H"
#include "surfaceInterpolate.H"
#include "fvMatrix.H"
// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
template
Foam::tmp
Foam::MovingPhaseModel::phi(const volVectorField& U) const
{
word phiName(IOobject::groupName("phi", this->name()));
IOobject phiHeader
(
phiName,
U.mesh().time().timeName(),
U.mesh(),
IOobject::NO_READ
);
if (phiHeader.headerOk())
{
Info<< "Reading face flux field " << phiName << endl;
return tmp
(
new surfaceScalarField
(
IOobject
(
phiName,
U.mesh().time().timeName(),
U.mesh(),
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
U.mesh()
)
);
}
else
{
Info<< "Calculating face flux field " << phiName << endl;
wordList phiTypes
(
U.boundaryField().size(),
calculatedFvPatchScalarField::typeName
);
forAll(U.boundaryField(), i)
{
if
(
isA(U.boundaryField()[i])
|| isA(U.boundaryField()[i])
|| isA(U.boundaryField()[i])
)
{
phiTypes[i] = fixedValueFvPatchScalarField::typeName;
}
}
return tmp
(
new surfaceScalarField
(
IOobject
(
phiName,
U.mesh().time().timeName(),
U.mesh(),
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
fvc::flux(U),
phiTypes
)
);
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
template
Foam::MovingPhaseModel::MovingPhaseModel
(
const phaseSystem& fluid,
const word& phaseName,
const label index
)
:
BasePhaseModel(fluid, phaseName, index),
U_
(
IOobject
(
IOobject::groupName("U", this->name()),
fluid.mesh().time().timeName(),
fluid.mesh(),
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
fluid.mesh()
),
phi_(phi(U_)),
alphaPhi_
(
IOobject
(
IOobject::groupName("alphaPhi", this->name()),
fluid.mesh().time().timeName(),
fluid.mesh()
),
fluid.mesh(),
dimensionedScalar("0", dimensionSet(0, 3, -1, 0, 0), 0)
),
alphaRhoPhi_
(
IOobject
(
IOobject::groupName("alphaRhoPhi", this->name()),
fluid.mesh().time().timeName(),
fluid.mesh()
),
fluid.mesh(),
dimensionedScalar("0", dimensionSet(1, 0, -1, 0, 0), 0)
),
DUDt_
(
IOobject
(
IOobject::groupName("DUDt", this->name()),
fluid.mesh().time().timeName(),
fluid.mesh()
),
fluid.mesh(),
dimensionedVector("0", dimAcceleration, Zero)
),
divU_(NULL),
turbulence_
(
phaseCompressibleTurbulenceModel::New
(
*this,
this->thermo().rho(),
U_,
alphaRhoPhi_,
phi_,
*this
)
),
continuityError_
(
IOobject
(
IOobject::groupName("continuityError", this->name()),
fluid.mesh().time().timeName(),
fluid.mesh()
),
fluid.mesh(),
dimensionedScalar("0", dimDensity/dimTime, 0)
)
{
phi_.writeOpt() = IOobject::AUTO_WRITE;
correctKinematics();
}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
template
Foam::MovingPhaseModel::~MovingPhaseModel()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template
void Foam::MovingPhaseModel::correct()
{
BasePhaseModel::correct();
this->fluid().MRF().correctBoundaryVelocity(U_);
volScalarField& rho = this->thermo().rho();
continuityError_ =
fvc::ddt(*this, rho) + fvc::div(alphaRhoPhi_)
- (this->fluid().fvOptions()(*this, rho)&rho);
}
template
void Foam::MovingPhaseModel::correctKinematics()
{
BasePhaseModel::correctKinematics();
DUDt_ = fvc::ddt(U_) + fvc::div(phi_, U_) - fvc::div(phi_)*U_;
}
template
void Foam::MovingPhaseModel::correctTurbulence()
{
BasePhaseModel::correctTurbulence();
turbulence_->correct();
}
template
void Foam::MovingPhaseModel::correctEnergyTransport()
{
BasePhaseModel::correctEnergyTransport();
turbulence_->correctEnergyTransport();
}
template
Foam::tmp
Foam::MovingPhaseModel::UEqn()
{
return
(
fvm::ddt(*this, this->thermo().rho(), U_)
+ fvm::div(alphaRhoPhi_, U_)
- fvm::Sp(continuityError_, U_)
+ this->fluid().MRF().DDt(*this*this->thermo().rho(), U_)
+ turbulence_->divDevRhoReff(U_)
);
}
template
const Foam::surfaceScalarField&
Foam::MovingPhaseModel::DbyA() const
{
if (DbyA_.valid())
{
return DbyA_;
}
else
{
return surfaceScalarField::null();
}
}
template
void Foam::MovingPhaseModel::DbyA
(
const tmp& DbyA
)
{
DbyA_ = DbyA;
}
template
Foam::tmp
Foam::MovingPhaseModel::U() const
{
return U_;
}
template
Foam::volVectorField&
Foam::MovingPhaseModel::U()
{
return U_;
}
template
Foam::tmp
Foam::MovingPhaseModel::DUDt() const
{
return DUDt_;
}
template
const Foam::tmp&
Foam::MovingPhaseModel::divU() const
{
return divU_;
}
template
void
Foam::MovingPhaseModel::divU
(
const tmp& divU
)
{
divU_ = divU;
}
template
Foam::tmp
Foam::MovingPhaseModel::continuityError() const
{
return continuityError_;
}
template
Foam::tmp
Foam::MovingPhaseModel::phi() const
{
return phi_;
}
template
Foam::surfaceScalarField&
Foam::MovingPhaseModel::phi()
{
return phi_;
}
template
Foam::tmp
Foam::MovingPhaseModel::alphaPhi() const
{
return alphaPhi_;
}
template
Foam::surfaceScalarField&
Foam::MovingPhaseModel::alphaPhi()
{
return alphaPhi_;
}
template
Foam::tmp
Foam::MovingPhaseModel::alphaRhoPhi() const
{
return alphaRhoPhi_;
}
template
Foam::surfaceScalarField&
Foam::MovingPhaseModel::alphaRhoPhi()
{
return alphaRhoPhi_;
}
template
const Foam::phaseCompressibleTurbulenceModel&
Foam::MovingPhaseModel::turbulence() const
{
return turbulence_;
}
// ************************************************************************* //