Update sonicLiquidFoam to PIMPLE
This commit is contained in:
parent
ad6bce6a29
commit
91e9838649
4 changed files with 62 additions and 57 deletions
|
@ -1,6 +1,3 @@
|
|||
{
|
||||
# include "rhoEqn.H"
|
||||
}
|
||||
{
|
||||
scalar sumLocalContErr =
|
||||
(sum(mag(rho - rho0 - psi*(p - p0)))/sum(rho)).value();
|
||||
|
|
|
@ -52,58 +52,71 @@ int main(int argc, char *argv[])
|
|||
{
|
||||
Info<< "Time = " << runTime.timeName() << nl << endl;
|
||||
|
||||
# include "readPISOControls.H"
|
||||
# include "readPIMPLEControls.H"
|
||||
# include "compressibleCourantNo.H"
|
||||
|
||||
# include "rhoEqn.H"
|
||||
|
||||
fvVectorMatrix UEqn
|
||||
(
|
||||
fvm::ddt(rho, U)
|
||||
+ fvm::div(phi, U)
|
||||
- fvm::laplacian(mu, U)
|
||||
);
|
||||
|
||||
solve(UEqn == -fvc::grad(p));
|
||||
|
||||
|
||||
// --- PISO loop
|
||||
|
||||
for (int corr = 0; corr < nCorr; corr++)
|
||||
// --- PIMPLE loop
|
||||
label oCorr = 0;
|
||||
do
|
||||
{
|
||||
volScalarField rUA = 1.0/UEqn.A();
|
||||
U = rUA*UEqn.H();
|
||||
|
||||
surfaceScalarField phid
|
||||
fvVectorMatrix UEqn
|
||||
(
|
||||
"phid",
|
||||
psi*
|
||||
fvm::ddt(rho, U)
|
||||
+ fvm::div(phi, U)
|
||||
- fvm::laplacian(mu, U)
|
||||
);
|
||||
|
||||
solve(UEqn == -fvc::grad(p));
|
||||
|
||||
// --- PISO loop
|
||||
for (int corr = 0; corr < nCorr; corr++)
|
||||
{
|
||||
volScalarField rAU("rAU", 1.0/UEqn.A());
|
||||
surfaceScalarField rhorAUf
|
||||
(
|
||||
(fvc::interpolate(U) & mesh.Sf())
|
||||
+ fvc::ddtPhiCorr(rUA, rho, U, phi)
|
||||
)
|
||||
);
|
||||
"rhorAUf",
|
||||
fvc::interpolate(rho*rAU)
|
||||
);
|
||||
|
||||
phi = (rhoO/psi)*phid;
|
||||
U = rAU*UEqn.H();
|
||||
|
||||
fvScalarMatrix pEqn
|
||||
(
|
||||
fvm::ddt(psi, p)
|
||||
+ fvc::div(phi)
|
||||
+ fvm::div(phid, p)
|
||||
- fvm::laplacian(rho*rUA, p)
|
||||
);
|
||||
surfaceScalarField phid
|
||||
(
|
||||
"phid",
|
||||
psi*
|
||||
(
|
||||
(fvc::interpolate(U) & mesh.Sf())
|
||||
+ fvc::ddtPhiCorr(rAU, rho, U, phi)
|
||||
)
|
||||
);
|
||||
|
||||
pEqn.solve();
|
||||
phi = (rhoO/psi)*phid;
|
||||
|
||||
phi += pEqn.flux();
|
||||
fvScalarMatrix pEqn
|
||||
(
|
||||
fvm::ddt(psi, p)
|
||||
+ fvc::div(phi)
|
||||
+ fvm::div(phid, p)
|
||||
- fvm::laplacian(rhorAUf, p)
|
||||
);
|
||||
|
||||
# include "compressibleContinuityErrs.H"
|
||||
pEqn.solve();
|
||||
|
||||
U -= rUA*fvc::grad(p);
|
||||
U.correctBoundaryConditions();
|
||||
}
|
||||
phi += pEqn.flux();
|
||||
|
||||
# include "rhoEqn.H"
|
||||
# include "compressibleContinuityErrs.H"
|
||||
|
||||
// Correct velocity
|
||||
U -= rAU*fvc::grad(p);
|
||||
U.correctBoundaryConditions();
|
||||
}
|
||||
} while (++oCorr < nOuterCorr);
|
||||
|
||||
|
||||
// Correct density
|
||||
rho = rhoO + psi*p;
|
||||
|
||||
runTime.write();
|
||||
|
@ -115,7 +128,7 @@ int main(int argc, char *argv[])
|
|||
|
||||
Info<< "End\n" << endl;
|
||||
|
||||
return 0;
|
||||
return(0);
|
||||
}
|
||||
|
||||
|
||||
|
|
|
@ -37,7 +37,7 @@ laplacianSchemes
|
|||
{
|
||||
default none;
|
||||
laplacian(mu,U) Gauss linear corrected;
|
||||
laplacian((rho*(1|A(U))),p) Gauss linear corrected;
|
||||
laplacian(rhorAUf,p) Gauss linear corrected;
|
||||
}
|
||||
|
||||
interpolationSchemes
|
||||
|
|
|
@ -19,32 +19,27 @@ solvers
|
|||
{
|
||||
p
|
||||
{
|
||||
solver PBiCG;
|
||||
solver BiCGStab;
|
||||
preconditioner DILU;
|
||||
tolerance 1e-06;
|
||||
minIter 1;
|
||||
tolerance 1e-08;
|
||||
relTol 0;
|
||||
}
|
||||
|
||||
U
|
||||
{
|
||||
solver PBiCG;
|
||||
solver BiCGStab;
|
||||
preconditioner DILU;
|
||||
tolerance 1e-05;
|
||||
tolerance 1e-08;
|
||||
relTol 0;
|
||||
}
|
||||
|
||||
rho
|
||||
{
|
||||
solver PCG;
|
||||
preconditioner DIC;
|
||||
tolerance 1e-05;
|
||||
relTol 0;
|
||||
}
|
||||
{}
|
||||
}
|
||||
|
||||
PISO
|
||||
PIMPLE
|
||||
{
|
||||
nCorrectors 2;
|
||||
nOuterCorrectors 2;
|
||||
nCorrectors 2;
|
||||
nNonOrthogonalCorrectors 0;
|
||||
}
|
||||
|
||||
|
|
Reference in a new issue