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 rAU = 1.0/UEqn.A();
00005     U = rAU*UEqn.H();
00006 
00007     if (pZones.size() > 0)
00008     {
00009         // ddtPhiCorr not well defined for cases with porosity
00010         phi = fvc::interpolate(rho)*(fvc::interpolate(U) & mesh.Sf());
00011     }
00012     else
00013     {
00014         phi =
00015             fvc::interpolate(rho)
00016            *(
00017                 (fvc::interpolate(U) & mesh.Sf())
00018               + fvc::ddtPhiCorr(rAU, rho, U, phi)
00019             );
00020     }
00021 
00022     {
00023         fvScalarMatrix pDDtEqn
00024         (
00025             fvc::ddt(rho) + psi*correction(fvm::ddt(p))
00026           + fvc::div(phi)
00027         );
00028 
00029         // Thermodynamic density needs to be updated by psi*d(p) after the
00030         // pressure solution - done in 2 parts. Part 1:
00031         thermo.rho() -= psi*p;
00032 
00033         for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
00034         {
00035             fvScalarMatrix pEqn
00036             (
00037                 pDDtEqn - fvm::laplacian(rho*rAU, p)
00038              ==
00039                 parcels.Srho()
00040               + massSource.SuTot()
00041             );
00042 
00043             if (corr == nCorr-1 && nonOrth == nNonOrthCorr)
00044             {
00045                 pEqn.solve(mesh.solver("pFinal"));
00046             }
00047             else
00048             {
00049                 pEqn.solve();
00050             }
00051 
00052             if (nonOrth == nNonOrthCorr)
00053             {
00054                 phi += pEqn.flux();
00055             }
00056         }
00057 
00058         // Second part of thermodynamic density update
00059         thermo.rho() += psi*p;
00060     }
00061 
00062     #include "rhoEqn.H"
00063     #include <finiteVolume/compressibleContinuityErrs.H>
00064 
00065     U -= rAU*fvc::grad(p);
00066     U.correctBoundaryConditions();
00067 }
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines