incompressible/icoFoam, solver + tutorials

This commit is contained in:
Vuko Vukcevic 2016-04-05 11:23:41 +02:00
parent 5d2fb1e6f2
commit db568490b7
5 changed files with 26 additions and 31 deletions

View file

@ -26,6 +26,10 @@ Application
Description Description
Transient solver for incompressible, laminar flow of Newtonian fluids. Transient solver for incompressible, laminar flow of Newtonian fluids.
Consistent formulation without time-step and relaxation dependence by Jasak
Author
Hrvoje Jasak, Wikki Ltd. All rights reserved
\*---------------------------------------------------------------------------*/ \*---------------------------------------------------------------------------*/
@ -53,32 +57,36 @@ int main(int argc, char *argv[])
# include "readPISOControls.H" # include "readPISOControls.H"
# include "CourantNo.H" # include "CourantNo.H"
fvVectorMatrix UEqn // Convection-diffusion matrix
fvVectorMatrix HUEqn
( (
fvm::ddt(U) fvm::div(phi, U)
+ fvm::div(phi, U)
- fvm::laplacian(nu, U) - fvm::laplacian(nu, U)
); );
solve(UEqn == -fvc::grad(p)); // Time derivative matrix
fvVectorMatrix ddtUEqn(fvm::ddt(U));
solve(ddtUEqn + HUEqn == -fvc::grad(p));
// Prepare clean Ap without time derivative contribution
// HJ, 26/Oct/2015
volScalarField aU = HUEqn.A();
// --- PISO loop // --- PISO loop
for (int corr = 0; corr < nCorr; corr++) for (int corr = 0; corr < nCorr; corr++)
{ {
volScalarField rUA = 1.0/UEqn.A(); U = HUEqn.H()/aU;
phi = (fvc::interpolate(U) & mesh.Sf());
U = rUA*UEqn.H();
phi = (fvc::interpolate(U) & mesh.Sf())
+ fvc::ddtPhiCorr(rUA, U, phi);
adjustPhi(phi, U, p); adjustPhi(phi, U, p);
for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++) for (int nonOrth = 0; nonOrth <= nNonOrthCorr; nonOrth++)
{ {
fvScalarMatrix pEqn fvScalarMatrix pEqn
( (
fvm::laplacian(rUA, p) == fvc::div(phi) fvm::laplacian(1/aU, p) == fvc::div(phi)
); );
pEqn.setReference(pRefCell, pRefValue); pEqn.setReference(pRefCell, pRefValue);
@ -92,7 +100,12 @@ int main(int argc, char *argv[])
# include "continuityErrs.H" # include "continuityErrs.H"
U -= rUA*fvc::grad(p); // Note: cannot call H(U) here because the velocity is not complete
// HJ, 22/Jan/2016
U = 1.0/(aU + ddtUEqn.A())*
(
U*aU - fvc::grad(p) + ddtUEqn.H()
);
U.correctBoundaryConditions(); U.correctBoundaryConditions();
} }

View file

@ -28,7 +28,7 @@ gradSchemes
divSchemes divSchemes
{ {
default none; default none;
div(phi,U) Gauss limitedLinearV 1; div(phi,U) Gauss linearUpwind leastSquares;
} }
laplacianSchemes laplacianSchemes

View file

@ -54,12 +54,6 @@ boundaryField
type mixingPlane; type mixingPlane;
} }
top
{
type fixedValue;
value uniform (0 0 0);
}
frontAndBack frontAndBack
{ {
type empty; type empty;

View file

@ -54,12 +54,6 @@ boundaryField
type mixingPlane; type mixingPlane;
} }
top
{
type fixedValue;
value uniform (0 0 0);
}
frontAndBack frontAndBack
{ {
type empty; type empty;

View file

@ -54,12 +54,6 @@ boundaryField
type mixingPlane; type mixingPlane;
} }
top
{
type fixedValue;
value uniform (0 0 0);
}
frontAndBack frontAndBack
{ {
type empty; type empty;