/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2013-2017 OpenFOAM Foundation
Copyright (C) 2020 OpenCFD Ltd.
-------------------------------------------------------------------------------
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 "multiphaseMixtureThermo.H"
#include "alphaContactAngleFvPatchScalarField.H"
#include "Time.H"
#include "subCycle.H"
#include "MULES.H"
#include "fvcDiv.H"
#include "fvcGrad.H"
#include "fvcSnGrad.H"
#include "fvcFlux.H"
#include "fvcMeshPhi.H"
#include "surfaceInterpolate.H"
#include "unitConversion.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
defineTypeNameAndDebug(multiphaseMixtureThermo, 0);
}
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
void Foam::multiphaseMixtureThermo::calcAlphas()
{
scalar level = 0.0;
alphas_ == 0.0;
for (const phaseModel& phase : phases_)
{
alphas_ += level * phase;
level += 1.0;
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::multiphaseMixtureThermo::multiphaseMixtureThermo
(
const volVectorField& U,
const surfaceScalarField& phi
)
:
psiThermo(U.mesh(), word::null),
phases_(lookup("phases"), phaseModel::iNew(p_, T_)),
mesh_(U.mesh()),
U_(U),
phi_(phi),
rhoPhi_
(
IOobject
(
"rhoPhi",
mesh_.time().timeName(),
mesh_,
IOobject::NO_READ,
IOobject::NO_WRITE
),
mesh_,
dimensionedScalar(dimMass/dimTime, Zero)
),
alphas_
(
IOobject
(
"alphas",
mesh_.time().timeName(),
mesh_,
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
mesh_,
dimensionedScalar(dimless, Zero)
),
sigmas_(lookup("sigmas")),
dimSigma_(1, 0, -2, 0, 0),
deltaN_
(
"deltaN",
1e-8/cbrt(average(mesh_.V()))
)
{
rhoPhi_.setOriented();
calcAlphas();
alphas_.write();
correct();
}
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
void Foam::multiphaseMixtureThermo::correct()
{
for (phaseModel& phase : phases_)
{
phase.correct();
}
auto phasei = phases_.cbegin();
psi_ = phasei()*phasei().thermo().psi();
mu_ = phasei()*phasei().thermo().mu();
alpha_ = phasei()*phasei().thermo().alpha();
for (++phasei; phasei != phases_.cend(); ++phasei)
{
psi_ += phasei()*phasei().thermo().psi();
mu_ += phasei()*phasei().thermo().mu();
alpha_ += phasei()*phasei().thermo().alpha();
}
}
void Foam::multiphaseMixtureThermo::correctRho(const volScalarField& dp)
{
for (phaseModel& phase : phases_)
{
phase.thermo().rho() += phase.thermo().psi()*dp;
}
}
Foam::word Foam::multiphaseMixtureThermo::thermoName() const
{
auto phasei = phases_.cbegin();
word name = phasei().thermo().thermoName();
for (++phasei; phasei != phases_.cend(); ++phasei)
{
name += ',' + phasei().thermo().thermoName();
}
return name;
}
bool Foam::multiphaseMixtureThermo::incompressible() const
{
for (const phaseModel& phase : phases_)
{
if (!phase.thermo().incompressible())
{
return false;
}
}
return true;
}
bool Foam::multiphaseMixtureThermo::isochoric() const
{
for (const phaseModel& phase : phases_)
{
if (!phase.thermo().isochoric())
{
return false;
}
}
return true;
}
Foam::tmp Foam::multiphaseMixtureThermo::he
(
const volScalarField& p,
const volScalarField& T
) const
{
auto phasei = phases_.cbegin();
tmp the(phasei()*phasei().thermo().he(p, T));
for (++phasei; phasei != phases_.cend(); ++phasei)
{
the.ref() += phasei()*phasei().thermo().he(p, T);
}
return the;
}
Foam::tmp Foam::multiphaseMixtureThermo::he
(
const scalarField& p,
const scalarField& T,
const labelList& cells
) const
{
auto phasei = phases_.cbegin();
tmp the
(
scalarField(phasei(), cells)*phasei().thermo().he(p, T, cells)
);
for (++phasei; phasei != phases_.cend(); ++phasei)
{
the.ref() +=
scalarField(phasei(), cells)*phasei().thermo().he(p, T, cells);
}
return the;
}
Foam::tmp Foam::multiphaseMixtureThermo::he
(
const scalarField& p,
const scalarField& T,
const label patchi
) const
{
auto phasei = phases_.cbegin();
tmp the
(
phasei().boundaryField()[patchi]*phasei().thermo().he(p, T, patchi)
);
for (++phasei; phasei != phases_.cend(); ++phasei)
{
the.ref() +=
phasei().boundaryField()[patchi]*phasei().thermo().he(p, T, patchi);
}
return the;
}
Foam::tmp Foam::multiphaseMixtureThermo::hc() const
{
auto phasei = phases_.cbegin();
tmp thc(phasei()*phasei().thermo().hc());
for (++phasei; phasei != phases_.cend(); ++phasei)
{
thc.ref() += phasei()*phasei().thermo().hc();
}
return thc;
}
Foam::tmp Foam::multiphaseMixtureThermo::THE
(
const scalarField& h,
const scalarField& p,
const scalarField& T0,
const labelList& cells
) const
{
NotImplemented;
return T0;
}
Foam::tmp Foam::multiphaseMixtureThermo::THE
(
const scalarField& h,
const scalarField& p,
const scalarField& T0,
const label patchi
) const
{
NotImplemented;
return T0;
}
Foam::tmp Foam::multiphaseMixtureThermo::rho() const
{
auto phasei = phases_.cbegin();
tmp trho(phasei()*phasei().thermo().rho());
for (++phasei; phasei != phases_.cend(); ++phasei)
{
trho.ref() += phasei()*phasei().thermo().rho();
}
return trho;
}
Foam::tmp Foam::multiphaseMixtureThermo::rho
(
const label patchi
) const
{
auto phasei = phases_.cbegin();
tmp trho
(
phasei().boundaryField()[patchi]*phasei().thermo().rho(patchi)
);
for (++phasei; phasei != phases_.cend(); ++phasei)
{
trho.ref() +=
phasei().boundaryField()[patchi]*phasei().thermo().rho(patchi);
}
return trho;
}
Foam::tmp Foam::multiphaseMixtureThermo::Cp() const
{
auto phasei = phases_.cbegin();
tmp tCp(phasei()*phasei().thermo().Cp());
for (++phasei; phasei != phases_.cend(); ++phasei)
{
tCp.ref() += phasei()*phasei().thermo().Cp();
}
return tCp;
}
Foam::tmp Foam::multiphaseMixtureThermo::Cp
(
const scalarField& p,
const scalarField& T,
const label patchi
) const
{
auto phasei = phases_.cbegin();
tmp tCp
(
phasei().boundaryField()[patchi]*phasei().thermo().Cp(p, T, patchi)
);
for (++phasei; phasei != phases_.cend(); ++phasei)
{
tCp.ref() +=
phasei().boundaryField()[patchi]*phasei().thermo().Cp(p, T, patchi);
}
return tCp;
}
Foam::tmp Foam::multiphaseMixtureThermo::Cv() const
{
auto phasei = phases_.cbegin();
tmp tCv(phasei()*phasei().thermo().Cv());
for (++phasei; phasei != phases_.cend(); ++phasei)
{
tCv.ref() += phasei()*phasei().thermo().Cv();
}
return tCv;
}
Foam::tmp Foam::multiphaseMixtureThermo::Cv
(
const scalarField& p,
const scalarField& T,
const label patchi
) const
{
auto phasei = phases_.cbegin();
tmp tCv
(
phasei().boundaryField()[patchi]*phasei().thermo().Cv(p, T, patchi)
);
for (++phasei; phasei != phases_.cend(); ++phasei)
{
tCv.ref() +=
phasei().boundaryField()[patchi]*phasei().thermo().Cv(p, T, patchi);
}
return tCv;
}
Foam::tmp Foam::multiphaseMixtureThermo::gamma() const
{
auto phasei = phases_.cbegin();
tmp tgamma(phasei()*phasei().thermo().gamma());
for (++phasei; phasei != phases_.cend(); ++phasei)
{
tgamma.ref() += phasei()*phasei().thermo().gamma();
}
return tgamma;
}
Foam::tmp Foam::multiphaseMixtureThermo::gamma
(
const scalarField& p,
const scalarField& T,
const label patchi
) const
{
auto phasei = phases_.cbegin();
tmp tgamma
(
phasei().boundaryField()[patchi]*phasei().thermo().gamma(p, T, patchi)
);
for (++phasei; phasei != phases_.cend(); ++phasei)
{
tgamma.ref() +=
phasei().boundaryField()[patchi]
*phasei().thermo().gamma(p, T, patchi);
}
return tgamma;
}
Foam::tmp Foam::multiphaseMixtureThermo::Cpv() const
{
auto phasei = phases_.cbegin();
tmp tCpv(phasei()*phasei().thermo().Cpv());
for (++phasei; phasei != phases_.cend(); ++phasei)
{
tCpv.ref() += phasei()*phasei().thermo().Cpv();
}
return tCpv;
}
Foam::tmp Foam::multiphaseMixtureThermo::Cpv
(
const scalarField& p,
const scalarField& T,
const label patchi
) const
{
auto phasei = phases_.cbegin();
tmp tCpv
(
phasei().boundaryField()[patchi]*phasei().thermo().Cpv(p, T, patchi)
);
for (++phasei; phasei != phases_.cend(); ++phasei)
{
tCpv.ref() +=
phasei().boundaryField()[patchi]
*phasei().thermo().Cpv(p, T, patchi);
}
return tCpv;
}
Foam::tmp Foam::multiphaseMixtureThermo::CpByCpv() const
{
auto phasei = phases_.cbegin();
tmp tCpByCpv(phasei()*phasei().thermo().CpByCpv());
for (++phasei; phasei != phases_.cend(); ++phasei)
{
tCpByCpv.ref() += phasei()*phasei().thermo().CpByCpv();
}
return tCpByCpv;
}
Foam::tmp Foam::multiphaseMixtureThermo::CpByCpv
(
const scalarField& p,
const scalarField& T,
const label patchi
) const
{
auto phasei = phases_.cbegin();
tmp tCpByCpv
(
phasei().boundaryField()[patchi]*phasei().thermo().CpByCpv(p, T, patchi)
);
for (++phasei; phasei != phases_.cend(); ++phasei)
{
tCpByCpv.ref() +=
phasei().boundaryField()[patchi]
*phasei().thermo().CpByCpv(p, T, patchi);
}
return tCpByCpv;
}
Foam::tmp Foam::multiphaseMixtureThermo::W() const
{
auto phasei = phases_.cbegin();
tmp tW(phasei()*phasei().thermo().W());
for (++phasei; phasei != phases_.cend(); ++phasei)
{
tW.ref() += phasei()*phasei().thermo().W();
}
return tW;
}
Foam::tmp Foam::multiphaseMixtureThermo::nu() const
{
return mu()/rho();
}
Foam::tmp Foam::multiphaseMixtureThermo::nu
(
const label patchi
) const
{
return mu(patchi)/rho(patchi);
}
Foam::tmp Foam::multiphaseMixtureThermo::kappa() const
{
auto phasei = phases_.cbegin();
tmp tkappa(phasei()*phasei().thermo().kappa());
for (++phasei; phasei != phases_.cend(); ++phasei)
{
tkappa.ref() += phasei()*phasei().thermo().kappa();
}
return tkappa;
}
Foam::tmp Foam::multiphaseMixtureThermo::kappa
(
const label patchi
) const
{
auto phasei = phases_.cbegin();
tmp tkappa
(
phasei().boundaryField()[patchi]*phasei().thermo().kappa(patchi)
);
for (++phasei; phasei != phases_.cend(); ++phasei)
{
tkappa.ref() +=
phasei().boundaryField()[patchi]*phasei().thermo().kappa(patchi);
}
return tkappa;
}
Foam::tmp Foam::multiphaseMixtureThermo::alphahe() const
{
PtrDictionary::const_iterator phasei = phases_.begin();
tmp talphaEff(phasei()*phasei().thermo().alphahe());
for (++phasei; phasei != phases_.end(); ++phasei)
{
talphaEff.ref() += phasei()*phasei().thermo().alphahe();
}
return talphaEff;
}
Foam::tmp Foam::multiphaseMixtureThermo::alphahe
(
const label patchi
) const
{
PtrDictionary::const_iterator phasei = phases_.begin();
tmp talphaEff
(
phasei().boundaryField()[patchi]
*phasei().thermo().alphahe(patchi)
);
for (++phasei; phasei != phases_.end(); ++phasei)
{
talphaEff.ref() +=
phasei().boundaryField()[patchi]
*phasei().thermo().alphahe(patchi);
}
return talphaEff;
}
Foam::tmp Foam::multiphaseMixtureThermo::kappaEff
(
const volScalarField& alphat
) const
{
auto phasei = phases_.cbegin();
tmp tkappaEff(phasei()*phasei().thermo().kappaEff(alphat));
for (++phasei; phasei != phases_.cend(); ++phasei)
{
tkappaEff.ref() += phasei()*phasei().thermo().kappaEff(alphat);
}
return tkappaEff;
}
Foam::tmp Foam::multiphaseMixtureThermo::kappaEff
(
const scalarField& alphat,
const label patchi
) const
{
auto phasei = phases_.cbegin();
tmp tkappaEff
(
phasei().boundaryField()[patchi]
*phasei().thermo().kappaEff(alphat, patchi)
);
for (++phasei; phasei != phases_.cend(); ++phasei)
{
tkappaEff.ref() +=
phasei().boundaryField()[patchi]
*phasei().thermo().kappaEff(alphat, patchi);
}
return tkappaEff;
}
Foam::tmp Foam::multiphaseMixtureThermo::alphaEff
(
const volScalarField& alphat
) const
{
auto phasei = phases_.cbegin();
tmp talphaEff(phasei()*phasei().thermo().alphaEff(alphat));
for (++phasei; phasei != phases_.cend(); ++phasei)
{
talphaEff.ref() += phasei()*phasei().thermo().alphaEff(alphat);
}
return talphaEff;
}
Foam::tmp Foam::multiphaseMixtureThermo::alphaEff
(
const scalarField& alphat,
const label patchi
) const
{
auto phasei = phases_.cbegin();
tmp talphaEff
(
phasei().boundaryField()[patchi]
*phasei().thermo().alphaEff(alphat, patchi)
);
for (++phasei; phasei != phases_.cend(); ++phasei)
{
talphaEff.ref() +=
phasei().boundaryField()[patchi]
*phasei().thermo().alphaEff(alphat, patchi);
}
return talphaEff;
}
Foam::tmp Foam::multiphaseMixtureThermo::rCv() const
{
auto phasei = phases_.cbegin();
tmp trCv(phasei()/phasei().thermo().Cv());
for (++phasei; phasei != phases_.cend(); ++phasei)
{
trCv.ref() += phasei()/phasei().thermo().Cv();
}
return trCv;
}
Foam::tmp
Foam::multiphaseMixtureThermo::surfaceTensionForce() const
{
tmp tstf
(
new surfaceScalarField
(
IOobject
(
"surfaceTensionForce",
mesh_.time().timeName(),
mesh_
),
mesh_,
dimensionedScalar(dimensionSet(1, -2, -2, 0, 0), Zero)
)
);
surfaceScalarField& stf = tstf.ref();
stf.setOriented();
forAllConstIters(phases_, phase1)
{
const phaseModel& alpha1 = *phase1;
auto phase2 = phase1;
for (++phase2; phase2 != phases_.cend(); ++phase2)
{
const phaseModel& alpha2 = *phase2;
auto sigma = sigmas_.cfind(interfacePair(alpha1, alpha2));
if (!sigma.found())
{
FatalErrorInFunction
<< "Cannot find interface " << interfacePair(alpha1, alpha2)
<< " in list of sigma values"
<< exit(FatalError);
}
stf += dimensionedScalar("sigma", dimSigma_, *sigma)
*fvc::interpolate(K(alpha1, alpha2))*
(
fvc::interpolate(alpha2)*fvc::snGrad(alpha1)
- fvc::interpolate(alpha1)*fvc::snGrad(alpha2)
);
}
}
return tstf;
}
void Foam::multiphaseMixtureThermo::solve()
{
const Time& runTime = mesh_.time();
const dictionary& alphaControls = mesh_.solverDict("alpha");
label nAlphaSubCycles(alphaControls.get