FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

pEqn.H

Go to the documentation of this file.
00001 {
00002     volScalarField rAU = 1.0/UEqn.A();
00003     surfaceScalarField rAUf = fvc::interpolate(rAU);
00004 
00005     U = rAU*UEqn.H();
00006     surfaceScalarField phiU
00007     (
00008         "phiU",
00009         (fvc::interpolate(U) & mesh.Sf())
00010       + fvc::ddtPhiCorr(rAU, rho, U, phi)
00011     );
00012 
00013     adjustPhi(phiU, U, p_rgh);
00014 
00015     phi = phiU - ghf*fvc::snGrad(rho)*rAUf*mesh.magSf();
00016 
00017     for(int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
00018     {
00019         fvScalarMatrix p_rghEqn
00020         (
00021             fvm::laplacian(rAUf, p_rgh) == fvc::div(phi)
00022         );
00023 
00024         p_rghEqn.setReference(pRefCell, getRefCellValue(p_rgh, pRefCell));
00025 
00026         p_rghEqn.solve
00027         (
00028             mesh.solver
00029             (
00030                 p_rgh.select
00031                 (
00032                     finalIter
00033                  && corr == nCorr-1
00034                  && nonOrth == nNonOrthCorr
00035                 )
00036             )
00037         );
00038 
00039         if (nonOrth == nNonOrthCorr)
00040         {
00041             phi -= p_rghEqn.flux();
00042         }
00043     }
00044 
00045     U += rAU*fvc::reconstruct((phi - phiU)/rAUf);
00046     U.correctBoundaryConditions();
00047 
00048     #include <finiteVolume/continuityErrs.H>
00049 
00050     p == p_rgh + rho*gh;
00051 
00052     if (p_rgh.needReference())
00053     {
00054         p += dimensionedScalar
00055         (
00056             "p",
00057             p.dimensions(),
00058             pRefValue - getRefCellValue(p, pRefCell)
00059         );
00060         p_rgh = p - rho*gh;
00061     }
00062 }
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines