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     rho = thermo.rho();
00003 
00004     volScalarField rUA = 1.0/UEqn.A();
00005     U = rUA*UEqn.H();
00006 
00007     if (transonic)
00008     {
00009         surfaceScalarField phiv =
00010             (fvc::interpolate(U) & mesh.Sf())
00011           + fvc::ddtPhiCorr(rUA, rho, U, phi);
00012 
00013         phi = fvc::interpolate(rho)*phiv;
00014 
00015         surfaceScalarField phid
00016         (
00017             "phid",
00018             fvc::interpolate(thermo.psi())*phiv
00019         );
00020 
00021         fvScalarMatrix pDDtEqn
00022         (
00023             fvc::ddt(rho) + fvc::div(phi)
00024           + correction(fvm::ddt(psi, p) + fvm::div(phid, p))
00025         );
00026 
00027         // Thermodynamic density needs to be updated by psi*d(p) after the
00028         // pressure solution - done in 2 parts. Part 1:
00029         thermo.rho() -= psi*p;
00030 
00031         for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
00032         {
00033             fvScalarMatrix pEqn
00034             (
00035                 pDDtEqn - fvm::laplacian(rho*rUA, p)
00036             );
00037 
00038             if (ocorr == nOuterCorr && corr == nCorr && nonOrth == nNonOrthCorr)
00039             {
00040                 pEqn.solve(mesh.solver(p.name() + "Final"));
00041             }
00042             else
00043             {
00044                 pEqn.solve();
00045             }
00046 
00047             if (nonOrth == nNonOrthCorr)
00048             {
00049                 phi += pEqn.flux();
00050             }
00051         }
00052 
00053         // Second part of thermodynamic density update
00054         thermo.rho() += psi*p;
00055     }
00056     else
00057     {
00058         phi =
00059             fvc::interpolate(rho)
00060            *(
00061                 (fvc::interpolate(U) & mesh.Sf())
00062               + fvc::ddtPhiCorr(rUA, rho, U, phi)
00063             );
00064 
00065         fvScalarMatrix pDDtEqn
00066         (
00067             fvc::ddt(rho) + psi*correction(fvm::ddt(p))
00068           + fvc::div(phi)
00069         );
00070 
00071         // Thermodynamic density needs to be updated by psi*d(p) after the
00072         // pressure solution - done in 2 parts. Part 1:
00073         thermo.rho() -= psi*p;
00074 
00075         for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
00076         {
00077             fvScalarMatrix pEqn
00078             (
00079                 pDDtEqn - fvm::laplacian(rho*rUA, p)
00080             );
00081 
00082             if (ocorr == nOuterCorr && corr == nCorr && nonOrth == nNonOrthCorr)
00083             {
00084                 pEqn.solve(mesh.solver(p.name() + "Final"));
00085             }
00086             else
00087             {
00088                 pEqn.solve();
00089             }
00090 
00091             if (nonOrth == nNonOrthCorr)
00092             {
00093                 phi += pEqn.flux();
00094             }
00095         }
00096 
00097         // Second part of thermodynamic density update
00098         thermo.rho() += psi*p;
00099     }
00100 
00101     #include <finiteVolume/rhoEqn.H>
00102     #include <finiteVolume/compressibleContinuityErrs.H>
00103 
00104     U -= rUA*fvc::grad(p);
00105     U.correctBoundaryConditions();
00106 
00107     DpDt = fvc::DDt(surfaceScalarField("phiU", phi/fvc::interpolate(rho)), p);
00108 }
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines