FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

bEqn.H

Go to the documentation of this file.
00001 tmp<fv::convectionScheme<scalar> > mvConvection
00002 (
00003     fv::convectionScheme<scalar>::New
00004     (
00005         mesh,
00006         fields,
00007         phi,
00008         mesh.divScheme("div(phi,ft_b_h_hu)")
00009     )
00010 );
00011 
00012 volScalarField Db("Db", turbulence->muEff());
00013 
00014 if (ign.ignited())
00015 {
00016     // Calculate the unstrained laminar flame speed
00017     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
00018     Su = unstrainedLaminarFlameSpeed()();
00019 
00020     const volScalarField& Xi = flameWrinkling->Xi();
00021 
00022     // progress variable
00023     // ~~~~~~~~~~~~~~~~~
00024     volScalarField c("c", 1.0 - b);
00025 
00026     // Unburnt gas density
00027     // ~~~~~~~~~~~~~~~~~~~
00028     volScalarField rhou = thermo.rhou();
00029 
00030     // Calculate flame normal etc.
00031     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
00032 
00033     //volVectorField n = fvc::grad(b);
00034     volVectorField n = fvc::reconstruct(fvc::snGrad(b,"snGrad(bprog)")*mesh.magSf());
00035 
00036     volScalarField mgb("mgb", mag(n));
00037 
00038     dimensionedScalar dMgb("dMgb", mgb.dimensions(), SMALL);
00039 
00040     {
00041         volScalarField bc = b*c;
00042 
00043         dMgb += 1.0e-3*
00044             (bc*mgb)().weightedAverage(mesh.V())
00045            /(bc.weightedAverage(mesh.V()) + SMALL);
00046     }
00047 
00048     mgb += dMgb;
00049 
00050     surfaceVectorField Sfhat = mesh.Sf()/mesh.magSf();
00051     surfaceVectorField nfVec = fvc::interpolate(n);
00052     nfVec += Sfhat*(fvc::snGrad(b,"snGrad(bprog)") - (Sfhat & nfVec));
00053     nfVec /= (mag(nfVec) + dMgb);
00054     surfaceScalarField nf("nf", mesh.Sf() & nfVec);
00055     n /= mgb;
00056 
00057 #   include <engine/StCorr.H>
00058 
00059     // Calculate turbulent flame speed flux
00060     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
00061     surfaceScalarField phiSt("phiSt", fvc::interpolate(rhou*StCorr*St)*nf);
00062 
00063 #   include "StCourantNo.H"
00064 
00065     Db = flameWrinkling->Db();
00066 
00067     // Create b equation
00068     // ~~~~~~~~~~~~~~~~~
00069     fvScalarMatrix bEqn
00070     (
00071         betav*fvm::ddt(rho, b)
00072       + mvConvection->fvmDiv(phi, b)
00073       + fvm::div(phiSt, b, "div(phiSt,bprog)")
00074       - fvm::Sp(fvc::div(phiSt), b)
00075       - fvm::laplacian(Db, b, "laplacian(DB,bprog)")
00076     );
00077 
00078 
00079     // Add ignition cell contribution to b-equation
00080     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
00081 #   include <engine/ignite.H>
00082 
00083     // Solve for b
00084     // ~~~~~~~~~~~
00085     bEqn.solve();
00086 
00087     Info<< "min(bprog) = " << min(b).value() << endl;
00088 
00089     if (composition.contains("ft"))
00090     {
00091         volScalarField& ft = composition.Y("ft");
00092 
00093         Info<< "Combustion progress = "
00094             << 100*(1.0 - b)().weightedAverage(mesh.V()*ft).value() << "%"
00095             << endl;
00096     }
00097     else
00098     {
00099         Info<< "Combustion progress = "
00100             << 100*(1.0 - b)().weightedAverage(mesh.V()).value() << "%"
00101             << endl;
00102     }
00103 
00104     // Correct the flame-wrinkling
00105     flameWrinkling->correct(mvConvection);
00106 
00107     St = Xi*Su;
00108 }
00109 
00110 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines