twoPhaseEulerFoam: transform to solve for p_rgh

to avoid excessive pressure/buoyancy balance errors on non-orthogonal meshes
Resolves bug-report http://openfoam.org/mantisbt/view.php?id=1379
This commit is contained in:
Henry
2015-03-17 22:40:09 +00:00
parent afcd290d75
commit 50b2edbedc
35 changed files with 485 additions and 116 deletions

View File

@ -1,3 +1,6 @@
#include "readGravitationalAcceleration.H"
#include "readhRef.H"
Info<< "Creating twoPhaseSystem\n" << endl;
twoPhaseSystem fluid(mesh, g);
@ -27,6 +30,13 @@
fluid.lookup("pMin")
);
Info<< "Calculating field g.h\n" << endl;
dimensionedScalar ghRef(g & (cmptMag(g.value())/mag(g.value()))*hRef);
volScalarField gh("gh", (g & mesh.C()) - ghRef);
surfaceScalarField ghf("ghf", (g & mesh.Cf()) - ghRef);
rhoThermo& thermo1 = phase1.thermo();
rhoThermo& thermo2 = phase2.thermo();
@ -38,6 +48,20 @@
volScalarField& rho2 = thermo2.rho();
const volScalarField& psi2 = thermo2.psi();
Info<< "Reading field p_rgh\n" << endl;
volScalarField p_rgh
(
IOobject
(
"p_rgh",
runTime.timeName(),
mesh,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
mesh
);
volVectorField U
(
IOobject
@ -99,7 +123,14 @@
label pRefCell = 0;
scalar pRefValue = 0.0;
setRefCell(p, pimple.dict(), pRefCell, pRefValue);
setRefCell
(
p,
p_rgh,
pimple.dict(),
pRefCell,
pRefValue
);
Info<< "Creating field dpdt\n" << endl;
volScalarField dpdt

View File

@ -28,49 +28,60 @@
mrfZones.makeAbsolute(phi2);
// Phase-1 pressure flux (e.g. due to particle-particle pressure)
surfaceScalarField phiP1
surfaceScalarField phiF1
(
"phiP1",
"phiF1",
fvc::interpolate(rAU1*phase1.turbulence().pPrime())
*fvc::snGrad(alpha1)*mesh.magSf()
);
phiP1.boundaryField() == 0;
phiF1.boundaryField() == 0;
// Phase-2 pressure flux (e.g. due to particle-particle pressure)
surfaceScalarField phiP2
surfaceScalarField phiF2
(
"phiP2",
"phiF2",
fvc::interpolate(rAU2*phase2.turbulence().pPrime())
*fvc::snGrad(alpha2)*mesh.magSf()
);
phiP2.boundaryField() == 0;
phiF2.boundaryField() == 0;
volScalarField rho("rho", fluid.rho());
surfaceScalarField ghSnGradRho(ghf*fvc::snGrad(rho)*mesh.magSf());
// Add the phase-1 buoyancy force
phiF1 +=
rAlphaAU1f
*(
ghSnGradRho
- alpha2f*fvc::interpolate(rho1 - rho2)*(g & mesh.Sf())
)/fvc::interpolate(rho1);
// Add the phase-2 buoyancy force
phiF2 +=
rAlphaAU2f
*(
ghSnGradRho
- alpha1f*fvc::interpolate(rho2 - rho1)*(g & mesh.Sf())
)/fvc::interpolate(rho2);
// Phase-1 predicted flux
surfaceScalarField phiHbyA1
(
IOobject::groupName("phiHbyA", phase1.name()),
(fvc::interpolate(HbyA1) & mesh.Sf())
+ rAlphaAU1f*fvc::ddtCorr(U1, phi1)
+ fvc::interpolate(rAU1*dragCoeff)*phi2
- phiF1
);
// Phase-2 predicted flux
surfaceScalarField phiHbyA2
(
IOobject::groupName("phiHbyA", phase2.name()),
(fvc::interpolate(HbyA2) & mesh.Sf())
+ rAlphaAU2f*fvc::ddtCorr(U2, phi2)
);
phiHbyA1 +=
(
fvc::interpolate(rAU1*dragCoeff)*phi2
- phiP1
+ rAlphaAU1f*(g & mesh.Sf())
);
phiHbyA2 +=
(
fvc::interpolate(rAU2*dragCoeff)*phi1
- phiP2
+ rAlphaAU2f*(g & mesh.Sf())
+ fvc::interpolate(rAU2*dragCoeff)*phi1
- phiF2
);
mrfZones.makeRelative(phiHbyA1);
@ -98,7 +109,7 @@
// Update the fixedFluxPressure BCs to ensure flux consistency
setSnGrad<fixedFluxPressureFvPatchScalarField>
(
p.boundaryField(),
p_rgh.boundaryField(),
(
phiHbyA.boundaryField()
- mrfZones.relative
@ -134,8 +145,8 @@
)/rho1
+ (alpha1/rho1)*correction
(
psi1*fvm::ddt(p)
+ fvm::div(phid1, p) - fvm::Sp(fvc::div(phid1), p)
psi1*fvm::ddt(p_rgh)
+ fvm::div(phid1, p_rgh) - fvm::Sp(fvc::div(phid1), p_rgh)
);
deleteDemandDrivenData(pEqnComp1().faceFluxCorrectionPtr());
pEqnComp1().relax();
@ -147,8 +158,8 @@
)/rho2
+ (alpha2/rho2)*correction
(
psi2*fvm::ddt(p)
+ fvm::div(phid2, p) - fvm::Sp(fvc::div(phid2), p)
psi2*fvm::ddt(p_rgh)
+ fvm::div(phid2, p_rgh) - fvm::Sp(fvc::div(phid2), p_rgh)
);
deleteDemandDrivenData(pEqnComp2().faceFluxCorrectionPtr());
pEqnComp2().relax();
@ -160,31 +171,31 @@
contErr1
- fvc::Sp(fvc::ddt(alpha1) + fvc::div(alphaPhi1), rho1)
)/rho1
+ (alpha1*psi1/rho1)*correction(fvm::ddt(p));
+ (alpha1*psi1/rho1)*correction(fvm::ddt(p_rgh));
pEqnComp2 =
(
contErr2
- fvc::Sp(fvc::ddt(alpha2) + fvc::div(alphaPhi2), rho2)
)/rho2
+ (alpha2*psi2/rho2)*correction(fvm::ddt(p));
+ (alpha2*psi2/rho2)*correction(fvm::ddt(p_rgh));
}
// Cache p prior to solve for density update
volScalarField p_0(p);
volScalarField p_rgh_0(p_rgh);
while (pimple.correctNonOrthogonal())
{
fvScalarMatrix pEqnIncomp
(
fvc::div(phiHbyA)
- fvm::laplacian(rAUf, p)
- fvm::laplacian(rAUf, p_rgh)
);
solve
(
pEqnComp1() + pEqnComp2() + pEqnIncomp,
mesh.solver(p.select(pimple.finalInnerIter()))
mesh.solver(p_rgh.select(pimple.finalInnerIter()))
);
if (pimple.finalNonOrthogonalIter())
@ -209,22 +220,22 @@
fluid.dgdt() =
(
alpha1*(pEqnComp2 & p)
- alpha2*(pEqnComp1 & p)
alpha1*(pEqnComp2 & p_rgh)
- alpha2*(pEqnComp1 & p_rgh)
);
p.relax();
p_rgh.relax();
p = max(p_rgh + rho*gh, pMin);
p_rgh = p - rho*gh;
mSfGradp = pEqnIncomp.flux()/rAUf;
U1 = HbyA1
+ fvc::reconstruct
(
rAlphaAU1f
*(
(g & mesh.Sf())
+ mSfGradp/fvc::interpolate(rho1)
)
- phiP1
rAlphaAU1f*mSfGradp/fvc::interpolate(rho1)
- phiF1
);
U1.correctBoundaryConditions();
fvOptions.correct(U1);
@ -232,12 +243,8 @@
U2 = HbyA2
+ fvc::reconstruct
(
rAlphaAU2f
*(
(g & mesh.Sf())
+ mSfGradp/fvc::interpolate(rho2)
)
- phiP2
rAlphaAU2f*mSfGradp/fvc::interpolate(rho2)
- phiF2
);
U2.correctBoundaryConditions();
fvOptions.correct(U2);
@ -246,11 +253,9 @@
}
}
p = max(p, pMin);
// Update densities from change in p
rho1 += psi1*(p - p_0);
rho2 += psi2*(p - p_0);
rho1 += psi1*(p_rgh - p_rgh_0);
rho2 += psi2*(p_rgh - p_rgh_0);
K1 = 0.5*magSqr(U1);
K2 = 0.5*magSqr(U2);

View File

@ -49,7 +49,6 @@ int main(int argc, char *argv[])
pimpleControl pimple(mesh);
#include "readGravitationalAcceleration.H"
#include "createFields.H"
#include "createMRFZones.H"
#include "createFvOptions.H"