933 lines
22 KiB
C
933 lines
22 KiB
C
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
\\ / O peration |
|
|
\\ / A nd | Copyright (C) 2011 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 "multiphaseSystem.H"
|
|
#include "alphaContactAngleFvPatchScalarField.H"
|
|
#include "fixedValueFvsPatchFields.H"
|
|
#include "Time.H"
|
|
#include "subCycle.H"
|
|
#include "MULES.H"
|
|
#include "fvcSnGrad.H"
|
|
#include "fvcFlux.H"
|
|
#include "fvcAverage.H"
|
|
|
|
// * * * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * //
|
|
|
|
const Foam::scalar Foam::multiphaseSystem::convertToRad =
|
|
Foam::constant::mathematical::pi/180.0;
|
|
|
|
|
|
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
|
|
|
|
void Foam::multiphaseSystem::calcAlphas()
|
|
{
|
|
scalar level = 0.0;
|
|
alphas_ == 0.0;
|
|
|
|
forAllIter(PtrDictionary<phaseModel>, phases_, iter)
|
|
{
|
|
alphas_ += level*iter();
|
|
level += 1.0;
|
|
}
|
|
|
|
alphas_.correctBoundaryConditions();
|
|
}
|
|
|
|
|
|
void Foam::multiphaseSystem::solveAlphas()
|
|
{
|
|
surfaceScalarField phic(mag(phi_/mesh_.magSf()));
|
|
|
|
PtrList<surfaceScalarField> phiAlphaCorrs(phases_.size());
|
|
int phasei = 0;
|
|
|
|
forAllIter(PtrDictionary<phaseModel>, phases_, iter)
|
|
{
|
|
phaseModel& phase1 = iter();
|
|
volScalarField& alpha1 = phase1;
|
|
|
|
phase1.phiAlpha() =
|
|
dimensionedScalar("0", dimensionSet(0, 3, -1, 0, 0), 0);
|
|
|
|
phiAlphaCorrs.set
|
|
(
|
|
phasei,
|
|
new surfaceScalarField
|
|
(
|
|
fvc::flux
|
|
(
|
|
phi_,
|
|
phase1,
|
|
"div(phi," + alpha1.name() + ')'
|
|
)
|
|
)
|
|
);
|
|
|
|
surfaceScalarField& phiAlphaCorr = phiAlphaCorrs[phasei];
|
|
|
|
forAllIter(PtrDictionary<phaseModel>, phases_, iter2)
|
|
{
|
|
phaseModel& phase2 = iter2();
|
|
volScalarField& alpha2 = phase2;
|
|
|
|
if (&phase2 == &phase1) continue;
|
|
|
|
surfaceScalarField phir
|
|
(
|
|
(phase1.phi() - phase2.phi())
|
|
+ min(cAlpha(phase1, phase2)*phic, max(phic))
|
|
*nHatf(phase1, phase2)
|
|
);
|
|
|
|
word phirScheme
|
|
(
|
|
"div(phir," + alpha2.name() + ',' + alpha1.name() + ')'
|
|
);
|
|
|
|
phiAlphaCorr += fvc::flux
|
|
(
|
|
-fvc::flux(-phir, phase2, phirScheme),
|
|
phase1,
|
|
phirScheme
|
|
);
|
|
}
|
|
|
|
MULES::limit
|
|
(
|
|
geometricOneField(),
|
|
phase1,
|
|
phi_,
|
|
phiAlphaCorr,
|
|
zeroField(),
|
|
zeroField(),
|
|
1,
|
|
0,
|
|
3,
|
|
true
|
|
);
|
|
|
|
phasei++;
|
|
}
|
|
|
|
MULES::limitSum(phiAlphaCorrs);
|
|
|
|
volScalarField sumAlpha
|
|
(
|
|
IOobject
|
|
(
|
|
"sumAlpha",
|
|
mesh_.time().timeName(),
|
|
mesh_
|
|
),
|
|
mesh_,
|
|
dimensionedScalar("sumAlpha", dimless, 0)
|
|
);
|
|
|
|
phasei = 0;
|
|
|
|
forAllIter(PtrDictionary<phaseModel>, phases_, iter)
|
|
{
|
|
phaseModel& phase1 = iter();
|
|
|
|
surfaceScalarField& phiAlpha = phiAlphaCorrs[phasei];
|
|
phiAlpha += upwind<scalar>(mesh_, phi_).flux(phase1);
|
|
|
|
MULES::explicitSolve
|
|
(
|
|
geometricOneField(),
|
|
phase1,
|
|
phiAlpha,
|
|
zeroField(),
|
|
zeroField()
|
|
);
|
|
|
|
phase1.phiAlpha() += phiAlpha;
|
|
|
|
Info<< phase1.name() << " volume fraction, min, max = "
|
|
<< phase1.weightedAverage(mesh_.V()).value()
|
|
<< ' ' << min(phase1).value()
|
|
<< ' ' << max(phase1).value()
|
|
<< endl;
|
|
|
|
sumAlpha += phase1;
|
|
|
|
phasei++;
|
|
}
|
|
|
|
Info<< "Phase-sum volume fraction, min, max = "
|
|
<< sumAlpha.weightedAverage(mesh_.V()).value()
|
|
<< ' ' << min(sumAlpha).value()
|
|
<< ' ' << max(sumAlpha).value()
|
|
<< endl;
|
|
|
|
calcAlphas();
|
|
}
|
|
|
|
|
|
Foam::dimensionedScalar Foam::multiphaseSystem::sigma
|
|
(
|
|
const phaseModel& phase1,
|
|
const phaseModel& phase2
|
|
) const
|
|
{
|
|
scalarCoeffTable::const_iterator sigma
|
|
(
|
|
sigmas_.find(interfacePair(phase1, phase2))
|
|
);
|
|
|
|
if (sigma == sigmas_.end())
|
|
{
|
|
FatalErrorIn
|
|
(
|
|
"multiphaseSystem::sigma(const phaseModel& phase1,"
|
|
"const phaseModel& phase2) const"
|
|
) << "Cannot find interface " << interfacePair(phase1, phase2)
|
|
<< " in list of sigma values"
|
|
<< exit(FatalError);
|
|
}
|
|
|
|
return dimensionedScalar("sigma", dimSigma_, sigma());
|
|
}
|
|
|
|
|
|
Foam::scalar Foam::multiphaseSystem::cAlpha
|
|
(
|
|
const phaseModel& phase1,
|
|
const phaseModel& phase2
|
|
) const
|
|
{
|
|
scalarCoeffTable::const_iterator cAlpha
|
|
(
|
|
cAlphas_.find(interfacePair(phase1, phase2))
|
|
);
|
|
|
|
if (cAlpha == cAlphas_.end())
|
|
{
|
|
FatalErrorIn
|
|
(
|
|
"multiphaseSystem::cAlpha"
|
|
"(const phaseModel& phase1, const phaseModel& phase2) const"
|
|
) << "Cannot find interface " << interfacePair(phase1, phase2)
|
|
<< " in list of cAlpha values"
|
|
<< exit(FatalError);
|
|
}
|
|
|
|
return cAlpha();
|
|
}
|
|
|
|
|
|
Foam::dimensionedScalar Foam::multiphaseSystem::Cvm
|
|
(
|
|
const phaseModel& phase1,
|
|
const phaseModel& phase2
|
|
) const
|
|
{
|
|
scalarCoeffTable::const_iterator Cvm
|
|
(
|
|
Cvms_.find(interfacePair(phase1, phase2))
|
|
);
|
|
|
|
if (Cvm != Cvms_.end())
|
|
{
|
|
return Cvm()*phase2.rho();
|
|
}
|
|
|
|
Cvm = Cvms_.find(interfacePair(phase2, phase1));
|
|
|
|
if (Cvm != Cvms_.end())
|
|
{
|
|
return Cvm()*phase1.rho();
|
|
}
|
|
|
|
FatalErrorIn
|
|
(
|
|
"multiphaseSystem::sigma"
|
|
"(const phaseModel& phase1, const phaseModel& phase2) const"
|
|
) << "Cannot find interface " << interfacePair(phase1, phase2)
|
|
<< " in list of sigma values"
|
|
<< exit(FatalError);
|
|
|
|
return Cvm()*phase2.rho();
|
|
}
|
|
|
|
|
|
Foam::tmp<Foam::surfaceVectorField> Foam::multiphaseSystem::nHatfv
|
|
(
|
|
const volScalarField& alpha1,
|
|
const volScalarField& alpha2
|
|
) const
|
|
{
|
|
/*
|
|
// Cell gradient of alpha
|
|
volVectorField gradAlpha =
|
|
alpha2*fvc::grad(alpha1) - alpha1*fvc::grad(alpha2);
|
|
|
|
// Interpolated face-gradient of alpha
|
|
surfaceVectorField gradAlphaf = fvc::interpolate(gradAlpha);
|
|
*/
|
|
|
|
surfaceVectorField gradAlphaf
|
|
(
|
|
fvc::interpolate(alpha2)*fvc::interpolate(fvc::grad(alpha1))
|
|
- fvc::interpolate(alpha1)*fvc::interpolate(fvc::grad(alpha2))
|
|
);
|
|
|
|
// Face unit interface normal
|
|
return gradAlphaf/(mag(gradAlphaf) + deltaN_);
|
|
}
|
|
|
|
|
|
Foam::tmp<Foam::surfaceScalarField> Foam::multiphaseSystem::nHatf
|
|
(
|
|
const volScalarField& alpha1,
|
|
const volScalarField& alpha2
|
|
) const
|
|
{
|
|
// Face unit interface normal flux
|
|
return nHatfv(alpha1, alpha2) & mesh_.Sf();
|
|
}
|
|
|
|
|
|
// Correction for the boundary condition on the unit normal nHat on
|
|
// walls to produce the correct contact angle.
|
|
|
|
// The dynamic contact angle is calculated from the component of the
|
|
// velocity on the direction of the interface, parallel to the wall.
|
|
|
|
void Foam::multiphaseSystem::correctContactAngle
|
|
(
|
|
const phaseModel& phase1,
|
|
const phaseModel& phase2,
|
|
surfaceVectorField::GeometricBoundaryField& nHatb
|
|
) const
|
|
{
|
|
const volScalarField::GeometricBoundaryField& gbf
|
|
= phase1.boundaryField();
|
|
|
|
const fvBoundaryMesh& boundary = mesh_.boundary();
|
|
|
|
forAll(boundary, patchi)
|
|
{
|
|
if (isA<alphaContactAngleFvPatchScalarField>(gbf[patchi]))
|
|
{
|
|
const alphaContactAngleFvPatchScalarField& acap =
|
|
refCast<const alphaContactAngleFvPatchScalarField>(gbf[patchi]);
|
|
|
|
vectorField& nHatPatch = nHatb[patchi];
|
|
|
|
vectorField AfHatPatch
|
|
(
|
|
mesh_.Sf().boundaryField()[patchi]
|
|
/mesh_.magSf().boundaryField()[patchi]
|
|
);
|
|
|
|
alphaContactAngleFvPatchScalarField::thetaPropsTable::
|
|
const_iterator tp =
|
|
acap.thetaProps().find(interfacePair(phase1, phase2));
|
|
|
|
if (tp == acap.thetaProps().end())
|
|
{
|
|
FatalErrorIn
|
|
(
|
|
"multiphaseSystem::correctContactAngle"
|
|
"(const phaseModel& phase1, const phaseModel& phase2, "
|
|
"fvPatchVectorFieldField& nHatb) const"
|
|
) << "Cannot find interface " << interfacePair(phase1, phase2)
|
|
<< "\n in table of theta properties for patch "
|
|
<< acap.patch().name()
|
|
<< exit(FatalError);
|
|
}
|
|
|
|
bool matched = (tp.key().first() == phase1.name());
|
|
|
|
scalar theta0 = convertToRad*tp().theta0(matched);
|
|
scalarField theta(boundary[patchi].size(), theta0);
|
|
|
|
scalar uTheta = tp().uTheta();
|
|
|
|
// Calculate the dynamic contact angle if required
|
|
if (uTheta > SMALL)
|
|
{
|
|
scalar thetaA = convertToRad*tp().thetaA(matched);
|
|
scalar thetaR = convertToRad*tp().thetaR(matched);
|
|
|
|
// Calculated the component of the velocity parallel to the wall
|
|
vectorField Uwall
|
|
(
|
|
phase1.U().boundaryField()[patchi].patchInternalField()
|
|
- phase1.U().boundaryField()[patchi]
|
|
);
|
|
Uwall -= (AfHatPatch & Uwall)*AfHatPatch;
|
|
|
|
// Find the direction of the interface parallel to the wall
|
|
vectorField nWall
|
|
(
|
|
nHatPatch - (AfHatPatch & nHatPatch)*AfHatPatch
|
|
);
|
|
|
|
// Normalise nWall
|
|
nWall /= (mag(nWall) + SMALL);
|
|
|
|
// Calculate Uwall resolved normal to the interface parallel to
|
|
// the interface
|
|
scalarField uwall(nWall & Uwall);
|
|
|
|
theta += (thetaA - thetaR)*tanh(uwall/uTheta);
|
|
}
|
|
|
|
|
|
// Reset nHatPatch to correspond to the contact angle
|
|
|
|
scalarField a12(nHatPatch & AfHatPatch);
|
|
|
|
scalarField b1(cos(theta));
|
|
|
|
scalarField b2(nHatPatch.size());
|
|
|
|
forAll(b2, facei)
|
|
{
|
|
b2[facei] = cos(acos(a12[facei]) - theta[facei]);
|
|
}
|
|
|
|
scalarField det(1.0 - a12*a12);
|
|
|
|
scalarField a((b1 - a12*b2)/det);
|
|
scalarField b((b2 - a12*b1)/det);
|
|
|
|
nHatPatch = a*AfHatPatch + b*nHatPatch;
|
|
|
|
nHatPatch /= (mag(nHatPatch) + deltaN_.value());
|
|
}
|
|
}
|
|
}
|
|
|
|
|
|
Foam::tmp<Foam::volScalarField> Foam::multiphaseSystem::K
|
|
(
|
|
const phaseModel& phase1,
|
|
const phaseModel& phase2
|
|
) const
|
|
{
|
|
tmp<surfaceVectorField> tnHatfv = nHatfv(phase1, phase2);
|
|
|
|
correctContactAngle(phase1, phase2, tnHatfv().boundaryField());
|
|
|
|
// Simple expression for curvature
|
|
return -fvc::div(tnHatfv & mesh_.Sf());
|
|
}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
|
|
|
Foam::multiphaseSystem::multiphaseSystem
|
|
(
|
|
const volVectorField& U,
|
|
const surfaceScalarField& phi
|
|
)
|
|
:
|
|
transportModel(U, phi),
|
|
|
|
phases_(lookup("phases"), phaseModel::iNew(U.mesh())),
|
|
|
|
mesh_(U.mesh()),
|
|
phi_(phi),
|
|
|
|
alphas_
|
|
(
|
|
IOobject
|
|
(
|
|
"alphas",
|
|
mesh_.time().timeName(),
|
|
mesh_,
|
|
IOobject::NO_READ,
|
|
IOobject::AUTO_WRITE
|
|
),
|
|
mesh_,
|
|
dimensionedScalar("alphas", dimless, 0.0),
|
|
zeroGradientFvPatchScalarField::typeName
|
|
),
|
|
|
|
sigmas_(lookup("sigmas")),
|
|
dimSigma_(1, 0, -2, 0, 0),
|
|
cAlphas_(lookup("interfaceCompression")),
|
|
Cvms_(lookup("virtualMass")),
|
|
deltaN_
|
|
(
|
|
"deltaN",
|
|
1e-8/pow(average(mesh_.V()), 1.0/3.0)
|
|
)
|
|
{
|
|
calcAlphas();
|
|
alphas_.write();
|
|
|
|
forAllIter(PtrDictionary<phaseModel>, phases_, iter)
|
|
{
|
|
phaseModelTable_.add(iter());
|
|
}
|
|
|
|
interfaceDictTable dragModelsDict(lookup("drag"));
|
|
|
|
forAllConstIter(interfaceDictTable, dragModelsDict, iter)
|
|
{
|
|
dragModels_.insert
|
|
(
|
|
iter.key(),
|
|
dragModel::New
|
|
(
|
|
iter(),
|
|
*phaseModelTable_.find(iter.key().first())(),
|
|
*phaseModelTable_.find(iter.key().second())()
|
|
).ptr()
|
|
);
|
|
}
|
|
}
|
|
|
|
|
|
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
|
|
|
|
Foam::tmp<Foam::volScalarField> Foam::multiphaseSystem::rho() const
|
|
{
|
|
PtrDictionary<phaseModel>::const_iterator iter = phases_.begin();
|
|
|
|
tmp<volScalarField> trho = iter()*iter().rho();
|
|
|
|
for (++iter; iter != phases_.end(); ++iter)
|
|
{
|
|
trho() += iter()*iter().rho();
|
|
}
|
|
|
|
return trho;
|
|
}
|
|
|
|
|
|
Foam::tmp<Foam::volScalarField> Foam::multiphaseSystem::nu() const
|
|
{
|
|
PtrDictionary<phaseModel>::const_iterator iter = phases_.begin();
|
|
|
|
tmp<volScalarField> tnu = iter()*iter().nu();
|
|
|
|
for (++iter; iter != phases_.end(); ++iter)
|
|
{
|
|
tnu() += iter()*iter().nu();
|
|
}
|
|
|
|
return tnu;
|
|
}
|
|
|
|
|
|
Foam::tmp<Foam::volScalarField> Foam::multiphaseSystem::Cvm
|
|
(
|
|
const phaseModel& phase
|
|
) const
|
|
{
|
|
tmp<volScalarField> tCvm
|
|
(
|
|
new volScalarField
|
|
(
|
|
IOobject
|
|
(
|
|
"Cvm",
|
|
mesh_.time().timeName(),
|
|
mesh_
|
|
),
|
|
mesh_,
|
|
dimensionedScalar
|
|
(
|
|
"Cvm",
|
|
dimensionSet(1, -3, 0, 0, 0),
|
|
0
|
|
)
|
|
)
|
|
);
|
|
|
|
forAllConstIter(PtrDictionary<phaseModel>, phases_, iter)
|
|
{
|
|
const phaseModel& phase2 = iter();
|
|
|
|
if (&phase2 != &phase)
|
|
{
|
|
tCvm() += Cvm(phase, phase2)*phase2;
|
|
}
|
|
}
|
|
|
|
return tCvm;
|
|
}
|
|
|
|
|
|
Foam::tmp<Foam::volVectorField> Foam::multiphaseSystem::Svm
|
|
(
|
|
const phaseModel& phase
|
|
) const
|
|
{
|
|
tmp<volVectorField> tSvm
|
|
(
|
|
new volVectorField
|
|
(
|
|
IOobject
|
|
(
|
|
"Svm",
|
|
mesh_.time().timeName(),
|
|
mesh_
|
|
),
|
|
mesh_,
|
|
dimensionedVector
|
|
(
|
|
"Svm",
|
|
dimensionSet(1, -2, -2, 0, 0),
|
|
vector::zero
|
|
)
|
|
)
|
|
);
|
|
|
|
forAllConstIter(PtrDictionary<phaseModel>, phases_, iter)
|
|
{
|
|
const phaseModel& phase2 = iter();
|
|
|
|
if (&phase2 != &phase)
|
|
{
|
|
tSvm() += Cvm(phase, phase2)*phase2*phase2.DDtU();
|
|
}
|
|
}
|
|
|
|
// Remove lift at fixed-flux boundaries
|
|
forAll(phase.phi().boundaryField(), patchi)
|
|
{
|
|
if
|
|
(
|
|
isA<fixedValueFvsPatchScalarField>
|
|
(
|
|
phase.phi().boundaryField()[patchi]
|
|
)
|
|
)
|
|
{
|
|
tSvm().boundaryField()[patchi] = vector::zero;
|
|
}
|
|
}
|
|
|
|
return tSvm;
|
|
}
|
|
|
|
|
|
Foam::autoPtr<Foam::multiphaseSystem::dragCoeffFields>
|
|
Foam::multiphaseSystem::dragCoeffs() const
|
|
{
|
|
autoPtr<dragCoeffFields> dragCoeffsPtr(new dragCoeffFields);
|
|
|
|
forAllConstIter(dragModelTable, dragModels_, iter)
|
|
{
|
|
const dragModel& dm = *iter();
|
|
|
|
volScalarField* Kptr =
|
|
(
|
|
max
|
|
(
|
|
//fvc::average(dm.phase1()*dm.phase2()),
|
|
//fvc::average(dm.phase1())*fvc::average(dm.phase2()),
|
|
dm.phase1()*dm.phase2(),
|
|
dm.residualPhaseFraction()
|
|
)
|
|
*dm.K
|
|
(
|
|
max
|
|
(
|
|
mag(dm.phase1().U() - dm.phase2().U()),
|
|
dm.residualSlip()
|
|
)
|
|
)
|
|
).ptr();
|
|
|
|
// Remove drag at fixed-flux boundaries
|
|
forAll(dm.phase1().phi().boundaryField(), patchi)
|
|
{
|
|
if
|
|
(
|
|
isA<fixedValueFvsPatchScalarField>
|
|
(
|
|
dm.phase1().phi().boundaryField()[patchi]
|
|
)
|
|
)
|
|
{
|
|
Kptr->boundaryField()[patchi] = 0.0;
|
|
}
|
|
}
|
|
|
|
dragCoeffsPtr().insert(iter.key(), Kptr);
|
|
}
|
|
|
|
return dragCoeffsPtr;
|
|
}
|
|
|
|
|
|
Foam::tmp<Foam::volScalarField> Foam::multiphaseSystem::dragCoeff
|
|
(
|
|
const phaseModel& phase,
|
|
const dragCoeffFields& dragCoeffs
|
|
) const
|
|
{
|
|
tmp<volScalarField> tdragCoeff
|
|
(
|
|
new volScalarField
|
|
(
|
|
IOobject
|
|
(
|
|
"dragCoeff",
|
|
mesh_.time().timeName(),
|
|
mesh_
|
|
),
|
|
mesh_,
|
|
dimensionedScalar
|
|
(
|
|
"dragCoeff",
|
|
dimensionSet(1, -3, -1, 0, 0),
|
|
0
|
|
)
|
|
)
|
|
);
|
|
|
|
dragModelTable::const_iterator dmIter = dragModels_.begin();
|
|
dragCoeffFields::const_iterator dcIter = dragCoeffs.begin();
|
|
for
|
|
(
|
|
;
|
|
dmIter != dragModels_.end() && dcIter != dragCoeffs.end();
|
|
++dmIter, ++dcIter
|
|
)
|
|
{
|
|
if
|
|
(
|
|
&phase == &dmIter()->phase1()
|
|
|| &phase == &dmIter()->phase2()
|
|
)
|
|
{
|
|
tdragCoeff() += *dcIter();
|
|
}
|
|
}
|
|
|
|
return tdragCoeff;
|
|
}
|
|
|
|
|
|
Foam::tmp<Foam::surfaceScalarField>
|
|
Foam::multiphaseSystem::surfaceTensionForce() const
|
|
{
|
|
tmp<surfaceScalarField> tstf
|
|
(
|
|
new surfaceScalarField
|
|
(
|
|
IOobject
|
|
(
|
|
"surfaceTensionForce",
|
|
mesh_.time().timeName(),
|
|
mesh_
|
|
),
|
|
mesh_,
|
|
dimensionedScalar
|
|
(
|
|
"surfaceTensionForce",
|
|
dimensionSet(1, -2, -2, 0, 0),
|
|
0.0
|
|
)
|
|
)
|
|
);
|
|
|
|
surfaceScalarField& stf = tstf();
|
|
|
|
forAllConstIter(PtrDictionary<phaseModel>, phases_, iter1)
|
|
{
|
|
const phaseModel& phase1 = iter1();
|
|
|
|
PtrDictionary<phaseModel>::const_iterator iter2 = iter1;
|
|
++iter2;
|
|
|
|
for (; iter2 != phases_.end(); ++iter2)
|
|
{
|
|
const phaseModel& phase2 = iter2();
|
|
|
|
stf += sigma(phase1, phase2)
|
|
*fvc::interpolate(K(phase1, phase2))*
|
|
(
|
|
fvc::interpolate(phase2)*fvc::snGrad(phase1)
|
|
- fvc::interpolate(phase1)*fvc::snGrad(phase2)
|
|
);
|
|
}
|
|
}
|
|
|
|
return tstf;
|
|
}
|
|
|
|
|
|
Foam::tmp<Foam::volScalarField>
|
|
Foam::multiphaseSystem::nearInterface() const
|
|
{
|
|
tmp<volScalarField> tnearInt
|
|
(
|
|
new volScalarField
|
|
(
|
|
IOobject
|
|
(
|
|
"nearInterface",
|
|
mesh_.time().timeName(),
|
|
mesh_
|
|
),
|
|
mesh_,
|
|
dimensionedScalar("nearInterface", dimless, 0.0)
|
|
)
|
|
);
|
|
|
|
forAllConstIter(PtrDictionary<phaseModel>, phases_, iter)
|
|
{
|
|
tnearInt() = max(tnearInt(), pos(iter() - 0.01)*pos(0.99 - iter()));
|
|
}
|
|
|
|
return tnearInt;
|
|
}
|
|
|
|
|
|
void Foam::multiphaseSystem::solve()
|
|
{
|
|
forAllIter(PtrDictionary<phaseModel>, phases_, iter)
|
|
{
|
|
iter().correct();
|
|
}
|
|
|
|
const Time& runTime = mesh_.time();
|
|
|
|
const dictionary& pimpleDict = mesh_.solutionDict().subDict("PIMPLE");
|
|
|
|
label nAlphaSubCycles(readLabel(pimpleDict.lookup("nAlphaSubCycles")));
|
|
|
|
if (nAlphaSubCycles > 1)
|
|
{
|
|
dimensionedScalar totalDeltaT = runTime.deltaT();
|
|
|
|
PtrList<volScalarField> alpha0s(phases_.size());
|
|
PtrList<surfaceScalarField> phiSums(phases_.size());
|
|
|
|
int phasei = 0;
|
|
forAllIter(PtrDictionary<phaseModel>, phases_, iter)
|
|
{
|
|
phaseModel& phase = iter();
|
|
volScalarField& alpha = phase;
|
|
|
|
alpha0s.set
|
|
(
|
|
phasei,
|
|
new volScalarField(alpha.oldTime())
|
|
);
|
|
|
|
phiSums.set
|
|
(
|
|
phasei,
|
|
new surfaceScalarField
|
|
(
|
|
IOobject
|
|
(
|
|
"phiSum" + alpha.name(),
|
|
runTime.timeName(),
|
|
mesh_
|
|
),
|
|
mesh_,
|
|
dimensionedScalar("0", dimensionSet(0, 3, -1, 0, 0), 0)
|
|
)
|
|
);
|
|
|
|
phasei++;
|
|
}
|
|
|
|
for
|
|
(
|
|
subCycleTime alphaSubCycle
|
|
(
|
|
const_cast<Time&>(runTime),
|
|
nAlphaSubCycles
|
|
);
|
|
!(++alphaSubCycle).end();
|
|
)
|
|
{
|
|
solveAlphas();
|
|
|
|
int phasei = 0;
|
|
forAllIter(PtrDictionary<phaseModel>, phases_, iter)
|
|
{
|
|
phiSums[phasei] += (runTime.deltaT()/totalDeltaT)*iter().phi();
|
|
phasei++;
|
|
}
|
|
}
|
|
|
|
phasei = 0;
|
|
forAllIter(PtrDictionary<phaseModel>, phases_, iter)
|
|
{
|
|
phaseModel& phase = iter();
|
|
volScalarField& alpha = phase;
|
|
|
|
phase.phi() = phiSums[phasei];
|
|
|
|
// Correct the time index of the field
|
|
// to correspond to the global time
|
|
alpha.timeIndex() = runTime.timeIndex();
|
|
|
|
// Reset the old-time field value
|
|
alpha.oldTime() = alpha0s[phasei];
|
|
alpha.oldTime().timeIndex() = runTime.timeIndex();
|
|
|
|
phasei++;
|
|
}
|
|
}
|
|
else
|
|
{
|
|
solveAlphas();
|
|
}
|
|
}
|
|
|
|
|
|
bool Foam::multiphaseSystem::read()
|
|
{
|
|
if (regIOobject::read())
|
|
{
|
|
bool readOK = true;
|
|
|
|
PtrList<entry> phaseData(lookup("phases"));
|
|
label phasei = 0;
|
|
|
|
forAllIter(PtrDictionary<phaseModel>, phases_, iter)
|
|
{
|
|
readOK &= iter().read(phaseData[phasei++].dict());
|
|
}
|
|
|
|
lookup("sigmas") >> sigmas_;
|
|
lookup("interfaceCompression") >> cAlphas_;
|
|
lookup("virtualMass") >> Cvms_;
|
|
|
|
return readOK;
|
|
}
|
|
else
|
|
{
|
|
return false;
|
|
}
|
|
}
|
|
|
|
|
|
// ************************************************************************* //
|