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("rUA", 1.0/UEqn().A());
00003     surfaceScalarField rUAf("(1|A(U))", fvc::interpolate(rUA));
00004 
00005     U = rUA*UEqn().H();
00006     UEqn.clear();
00007 
00008     phi = fvc::interpolate(U) & mesh.Sf();
00009     adjustPhi(phi, U, p_rgh);
00010 
00011     surfaceScalarField buoyancyPhi = rUAf*ghf*fvc::snGrad(rhok)*mesh.magSf();
00012     phi -= buoyancyPhi;
00013 
00014     for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
00015     {
00016         fvScalarMatrix p_rghEqn
00017         (
00018             fvm::laplacian(rUAf, p_rgh) == fvc::div(phi)
00019         );
00020 
00021         p_rghEqn.setReference(pRefCell, getRefCellValue(p_rgh, pRefCell));
00022 
00023         // retain the residual from the first iteration
00024         if (nonOrth == 0)
00025         {
00026             eqnResidual = p_rghEqn.solve().initialResidual();
00027             maxResidual = max(eqnResidual, maxResidual);
00028         }
00029         else
00030         {
00031             p_rghEqn.solve();
00032         }
00033 
00034         if (nonOrth == nNonOrthCorr)
00035         {
00036             // Calculate the conservative fluxes
00037             phi -= p_rghEqn.flux();
00038 
00039             // Explicitly relax pressure for momentum corrector
00040             p_rgh.relax();
00041 
00042             // Correct the momentum source with the pressure gradient flux
00043             // calculated from the relaxed pressure
00044             U -= rUA*fvc::reconstruct((buoyancyPhi + p_rghEqn.flux())/rUAf);
00045             U.correctBoundaryConditions();
00046         }
00047     }
00048 
00049     #include <finiteVolume/continuityErrs.H>
00050 
00051     p = p_rgh + rhok*gh;
00052 
00053     if (p_rgh.needReference())
00054     {
00055         p += dimensionedScalar
00056         (
00057             "p",
00058             p.dimensions(),
00059             pRefValue - getRefCellValue(p, pRefCell)
00060         );
00061         p_rgh = p - rhok*gh;
00062     }
00063 }
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines