Newer
Older
surfaceScalarField alpha1f(fvc::interpolate(alpha1));
surfaceScalarField alpha2f(scalar(1) - alpha1f);
Henry
committed
rAU1 = 1.0/U1Eqn.A();
rAU2 = 1.0/U2Eqn.A();
Henry
committed
surfaceScalarField rAlphaAU1f(fvc::interpolate(alpha1*rAU1));
surfaceScalarField rAlphaAU2f(fvc::interpolate(alpha2*rAU2));
// Update the phi BCs from Us before p BCs are updated
phi.boundaryField() =
alpha1f.boundaryField()
*(mesh.Sf().boundaryField() & U1.boundaryField())
+ alpha2f.boundaryField()
*(mesh.Sf().boundaryField() & U2.boundaryField());
Henry
committed
Henry
committed
volVectorField HbyA1
(
IOobject::groupName("HbyA", phase1.name()),
U1
);
Henry
committed
Henry
committed
volVectorField HbyA2
(
IOobject::groupName("HbyA", phase2.name()),
U2
);
mrfZones.absoluteFlux(phi1.oldTime());
mrfZones.absoluteFlux(phi1);
mrfZones.absoluteFlux(phi2.oldTime());
mrfZones.absoluteFlux(phi2);
Henry
committed
// Phase-1 pressure flux (e.g. due to particle-particle pressure)
surfaceScalarField phiP1
(
"phiP1",
fvc::interpolate((1.0/rho1)*rAU1*turbulence1().pPrime())
*fvc::snGrad(alpha1)*mesh.magSf()
);
Henry
committed
// Phase-2 pressure flux (e.g. due to particle-particle pressure)
surfaceScalarField phiP2
(
"phiP2",
fvc::interpolate((1.0/rho2)*rAU2*turbulence2().pPrime())
*fvc::snGrad(alpha2)*mesh.magSf()
);
Henry
committed
IOobject::groupName("phiHbyA", phase1.name()),
Henry
committed
+ fvc::ddtPhiCorr(rAU1, alpha1, U1, phi1)
Henry
committed
IOobject::groupName("phiHbyA", phase2.name()),
Henry
committed
+ fvc::ddtPhiCorr(rAU2, alpha2, U2, phi2)
);
phi = alpha1f*phiHbyA1 + alpha2f*phiHbyA2;
mrfZones.relativeFlux(phi);
phiHbyA1 +=
(
Henry
committed
fvc::interpolate((1.0/rho1)*rAU1*dragCoeff)*phi2
- phiP1
+ rAlphaAU1f*(g & mesh.Sf())
);
mrfZones.relativeFlux(phiHbyA1);
phiHbyA2 +=
(
Henry
committed
fvc::interpolate((1.0/rho2)*rAU2*dragCoeff)*phi1
- phiP2
+ rAlphaAU2f*(g & mesh.Sf())
mrfZones.relativeFlux(phiHbyA2);
mrfZones.relativeFlux(phi1.oldTime());
mrfZones.relativeFlux(phi1);
mrfZones.relativeFlux(phi2.oldTime());
mrfZones.relativeFlux(phi2);
surfaceScalarField phiHbyA("phiHbyA", alpha1f*phiHbyA1 + alpha2f*phiHbyA2);
Henry
committed
HbyA1 += (1.0/rho1)*rAU1*dragCoeff*U2;
HbyA2 += (1.0/rho2)*rAU2*dragCoeff*U1;
Henry
committed
mag
(
alpha1f*rAlphaAU1f/fvc::interpolate(rho1)
+ alpha2f*rAlphaAU2f/fvc::interpolate(rho2)
)
Henry
committed
tmp<fvScalarMatrix> pEqnComp1;
tmp<fvScalarMatrix> pEqnComp2;
if (pimple.transonic())
Henry
committed
surfaceScalarField phid1
(
IOobject::groupName("phid", phase1.name()),
fvc::interpolate(psi1)*phi1
);
surfaceScalarField phid2
Henry
committed
IOobject::groupName("phid", phase2.name()),
fvc::interpolate(psi2)*phi2
Henry
committed
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
pEqnComp1 =
fvc::ddt(rho1)
+ fvc::div(phi1, rho1) - fvc::Sp(fvc::div(phi1), rho1)
+ correction
(
psi1*fvm::ddt(p)
+ fvm::div(phid1, p) - fvm::Sp(fvc::div(phid1), p)
);
deleteDemandDrivenData(pEqnComp1().faceFluxCorrectionPtr());
pEqnComp1().relax();
pEqnComp2 =
fvc::ddt(rho2)
+ fvc::div(phi2, rho2) - fvc::Sp(fvc::div(phi2), rho2)
+ correction
(
psi2*fvm::ddt(p)
+ fvm::div(phid2, p) - fvm::Sp(fvc::div(phid2), p)
);
deleteDemandDrivenData(pEqnComp2().faceFluxCorrectionPtr());
pEqnComp2().relax();
}
else
{
pEqnComp1 =
fvc::ddt(rho1) + psi1*correction(fvm::ddt(p))
+ fvc::div(phi1, rho1) - fvc::Sp(fvc::div(phi1), rho1);
pEqnComp2 =
fvc::ddt(rho2) + psi2*correction(fvm::ddt(p))
+ fvc::div(phi2, rho2) - fvc::Sp(fvc::div(phi2), rho2);
}
Henry
committed
// Cache p prior to solve for density update
volScalarField p_0(p);
while (pimple.correctNonOrthogonal())
{
fvScalarMatrix pEqnIncomp
(
fvc::div(phiHbyA)
- fvm::laplacian(Dp, p)
);
solve
(
(
(alpha1/rho1)*pEqnComp1()
+ (alpha2/rho2)*pEqnComp2()
)
+ pEqnIncomp,
mesh.solver(p.select(pimple.finalInnerIter()))
);
if (pimple.finalNonOrthogonalIter())
Henry
committed
surfaceScalarField mSfGradp(pEqnIncomp.flux()/Dp);
phi1.boundaryField() ==
(mesh.Sf().boundaryField() & U1.boundaryField());
Henry
committed
phi1 = phiHbyA1 + rAlphaAU1f*mSfGradp/fvc::interpolate(rho1);
phi2.boundaryField() ==
(mesh.Sf().boundaryField() & U2.boundaryField());
Henry
committed
phi2 = phiHbyA2 + rAlphaAU2f*mSfGradp/fvc::interpolate(rho2);
Henry
committed
dgdt =
(
pos(alpha2)*(pEqnComp2 & p)/rho2
- pos(alpha1)*(pEqnComp1 & p)/rho1
);
Henry
committed
mSfGradp = pEqnIncomp.flux()/Dp;
Henry
committed
+ fvc::reconstruct
(
rAlphaAU1f
*(
(g & mesh.Sf())
+ mSfGradp/fvc::interpolate(rho1)
)
- phiP1
);
Henry
committed
+ fvc::reconstruct
(
rAlphaAU2f
*(
(g & mesh.Sf())
+ mSfGradp/fvc::interpolate(rho2)
)
- phiP2
);
Henry
committed
U = fluid.U();
Henry
committed
p = max(p, pMin);
// Update densities from change in p
rho1 += psi1*(p - p_0);
rho2 += psi2*(p - p_0);
K1 = 0.5*magSqr(U1);
K2 = 0.5*magSqr(U2);
if (thermo1.dpdt() || thermo2.dpdt())
{
dpdt = fvc::ddt(p);
}
}