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 rUA = 1.0/UEqn.A();
00003     surfaceScalarField rUAf = fvc::interpolate(rUA);
00004 
00005     U = rUA*UEqn.H();
00006 
00007     surfaceScalarField phiU
00008     (
00009         "phiU",
00010         (fvc::interpolate(U) & mesh.Sf())
00011       + fvc::ddtPhiCorr(rUA, rho, U, phi)
00012     );
00013 
00014     adjustPhi(phiU, U, p_rgh);
00015 
00016     phi = phiU +
00017         (
00018             fvc::interpolate(interface.sigmaK())*fvc::snGrad(alpha1)
00019           - ghf*fvc::snGrad(rho)
00020         )*rUAf*mesh.magSf();
00021 
00022     Pair<tmp<volScalarField> > vDotP = twoPhaseProperties->vDotP();
00023     const volScalarField& vDotcP = vDotP[0]();
00024     const volScalarField& vDotvP = vDotP[1]();
00025 
00026     for(int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
00027     {
00028         fvScalarMatrix p_rghEqn
00029         (
00030             fvc::div(phi) - fvm::laplacian(rUAf, p_rgh)
00031           - (vDotvP - vDotcP)*(pSat - rho*gh) + fvm::Sp(vDotvP - vDotcP, p_rgh)
00032         );
00033 
00034         p_rghEqn.setReference(pRefCell, pRefValue);
00035 
00036         p_rghEqn.solve
00037         (
00038             mesh.solver
00039             (
00040                 p_rgh.select
00041                 (
00042                     finalIter
00043                  && corr == nCorr-1
00044                  && nonOrth == nNonOrthCorr
00045                 )
00046             )
00047         );
00048 
00049         if (nonOrth == nNonOrthCorr)
00050         {
00051             phi += p_rghEqn.flux();
00052         }
00053     }
00054 
00055     U += rUA*fvc::reconstruct((phi - phiU)/rUAf);
00056     U.correctBoundaryConditions();
00057 
00058     #include <finiteVolume/continuityErrs.H>
00059 
00060     p == p_rgh + rho*gh;
00061 
00062     if (p_rgh.needReference())
00063     {
00064         p += dimensionedScalar
00065         (
00066             "p",
00067             p.dimensions(),
00068             pRefValue - getRefCellValue(p, pRefCell)
00069         );
00070         p_rgh = p - rho*gh;
00071     }
00072 }
00073 
00074 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines