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     tmp<fvScalarMatrix> p_rghEqnComp;
00006 
00007     if (transonic)
00008     {
00009         p_rghEqnComp =
00010         (
00011             fvm::ddt(p_rgh)
00012           + fvm::div(phi, p_rgh)
00013           - fvm::Sp(fvc::div(phi), p_rgh)
00014         );
00015     }
00016     else
00017     {
00018         p_rghEqnComp =
00019         (
00020             fvm::ddt(p_rgh)
00021           + fvc::div(phi, p_rgh)
00022           - fvc::Sp(fvc::div(phi), p_rgh)
00023         );
00024     }
00025 
00026 
00027     U = rUA*UEqn.H();
00028 
00029     surfaceScalarField phiU
00030     (
00031         "phiU",
00032         (fvc::interpolate(U) & mesh.Sf())
00033       + fvc::ddtPhiCorr(rUA, rho, U, phi)
00034     );
00035 
00036     phi = phiU +
00037         (
00038             fvc::interpolate(interface.sigmaK())*fvc::snGrad(alpha1)
00039           - ghf*fvc::snGrad(rho)
00040         )*rUAf*mesh.magSf();
00041 
00042     for(int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
00043     {
00044         fvScalarMatrix p_rghEqnIncomp
00045         (
00046             fvc::div(phi)
00047           - fvm::laplacian(rUAf, p_rgh)
00048         );
00049 
00050         solve
00051         (
00052             (
00053                 max(alpha1, scalar(0))*(psi1/rho1)
00054               + max(alpha2, scalar(0))*(psi2/rho2)
00055             )
00056            *p_rghEqnComp()
00057           + p_rghEqnIncomp,
00058             mesh.solver
00059             (
00060                 p_rgh.select
00061                 (
00062                     oCorr == nOuterCorr-1
00063                  && corr == nCorr-1
00064                  && nonOrth == nNonOrthCorr
00065                 )
00066             )
00067         );
00068 
00069         if (nonOrth == nNonOrthCorr)
00070         {
00071             dgdt =
00072                 (pos(alpha2)*(psi2/rho2) - pos(alpha1)*(psi1/rho1))
00073                *(p_rghEqnComp & p_rgh);
00074             phi += p_rghEqnIncomp.flux();
00075         }
00076     }
00077 
00078     U += rUA*fvc::reconstruct((phi - phiU)/rUAf);
00079     U.correctBoundaryConditions();
00080 
00081     p = max
00082     (
00083         (p_rgh + gh*(alpha1*rho10 + alpha2*rho20))
00084        /(1.0 - gh*(alpha1*psi1 + alpha2*psi2)),
00085         pMin
00086     );
00087 
00088     rho1 = rho10 + psi1*p;
00089     rho2 = rho20 + psi2*p;
00090 
00091     Info<< "max(U) " << max(mag(U)).value() << endl;
00092     Info<< "min(p_rgh) " << min(p_rgh).value() << endl;
00093 }
00094 
00095 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines