incompressible/channelFoam, solver + tutorials
This commit is contained in:
parent
1670cebef1
commit
edb598687a
2 changed files with 28 additions and 17 deletions
|
@ -26,6 +26,10 @@ Application
|
|||
|
||||
Description
|
||||
Incompressible LES solver for flow in a channel.
|
||||
Consistent formulation without time-step and relaxation dependence by Jasak
|
||||
|
||||
Author
|
||||
Hrvoje Jasak, Wikki Ltd. All rights reserved
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
|
@ -60,38 +64,41 @@ int main(int argc, char *argv[])
|
|||
|
||||
sgsModel->correct();
|
||||
|
||||
fvVectorMatrix UEqn
|
||||
// Convection-diffusion matrix
|
||||
fvVectorMatrix HUEqn
|
||||
(
|
||||
fvm::ddt(U)
|
||||
+ fvm::div(phi, U)
|
||||
fvm::div(phi, U)
|
||||
+ sgsModel->divDevBeff(U)
|
||||
==
|
||||
flowDirection*gradP
|
||||
);
|
||||
|
||||
// Time derivative matrix
|
||||
fvVectorMatrix ddtUEqn(fvm::ddt(U));
|
||||
|
||||
if (momentumPredictor)
|
||||
{
|
||||
solve(UEqn == -fvc::grad(p));
|
||||
solve(ddtUEqn + HUEqn == -fvc::grad(p));
|
||||
}
|
||||
|
||||
// Prepare clean Ap without time derivative contribution
|
||||
// HJ, 26/Oct/2015
|
||||
volScalarField aU = HUEqn.A();
|
||||
|
||||
// --- PISO loop
|
||||
|
||||
volScalarField rUA = 1.0/UEqn.A();
|
||||
|
||||
for (int corr = 0; corr < nCorr; corr++)
|
||||
{
|
||||
U = rUA*UEqn.H();
|
||||
phi = (fvc::interpolate(U) & mesh.Sf())
|
||||
+ fvc::ddtPhiCorr(rUA, U, phi);
|
||||
U = HUEqn.H()/aU;
|
||||
phi = (fvc::interpolate(U) & mesh.Sf());
|
||||
|
||||
adjustPhi(phi, U, p);
|
||||
|
||||
for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
|
||||
for (int nonOrth = 0; nonOrth <= nNonOrthCorr; nonOrth++)
|
||||
{
|
||||
fvScalarMatrix pEqn
|
||||
(
|
||||
fvm::laplacian(rUA, p) == fvc::div(phi)
|
||||
fvm::laplacian(1/aU, p) == fvc::div(phi)
|
||||
);
|
||||
|
||||
pEqn.setReference(pRefCell, pRefValue);
|
||||
|
@ -111,13 +118,17 @@ 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();
|
||||
}
|
||||
|
||||
|
||||
// Correct driving force for a constant mass flow rate
|
||||
|
||||
// Extract the velocity in the flow direction
|
||||
|
@ -127,9 +138,9 @@ int main(int argc, char *argv[])
|
|||
// Calculate the pressure gradient increment needed to
|
||||
// adjust the average flow-rate to the correct value
|
||||
dimensionedScalar gragPplus =
|
||||
(magUbar - magUbarStar)/rUA.weightedAverage(mesh.V());
|
||||
(magUbar - magUbarStar)*aU.weightedAverage(mesh.V());
|
||||
|
||||
U += flowDirection*rUA*gragPplus;
|
||||
U += gragPplus/aU*flowDirection;
|
||||
|
||||
gradP += gragPplus;
|
||||
|
||||
|
|
|
@ -29,7 +29,7 @@ gradSchemes
|
|||
divSchemes
|
||||
{
|
||||
default none;
|
||||
div(phi,U) Gauss linear;
|
||||
div(phi,U) Gauss linearUpwind Gauss linear;
|
||||
div(phi,k) Gauss limitedLinear 1;
|
||||
div(phi,B) Gauss limitedLinear 1;
|
||||
div(B) Gauss linear;
|
||||
|
|
Reference in a new issue