OpenFOAM logo
The Open Source CFD Toolbox
  Source Guide OpenCFD Solutions Contact OpenFOAM

bEqn.H

Go to the documentation of this file.
00001 if (ign.ignited())
00002 {
00003     // progress variable
00004     // ~~~~~~~~~~~~~~~~~
00005     volScalarField c = scalar(1) - b;
00006 
00007     // Unburnt gas density
00008     // ~~~~~~~~~~~~~~~~~~~
00009     volScalarField rhou = thermo.rhou();
00010 
00011     // Calculate flame normal etc.
00012     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
00013 
00014     volVectorField n = fvc::grad(b);
00015 
00016     volScalarField mgb = mag(n);
00017 
00018     dimensionedScalar dMgb = 1.0e-3*
00019         (b*c*mgb)().weightedAverage(mesh.V())
00020        /((b*c)().weightedAverage(mesh.V()) + SMALL)
00021       + dimensionedScalar("ddMgb", mgb.dimensions(), SMALL);
00022 
00023     mgb += dMgb;
00024 
00025     surfaceVectorField SfHat = mesh.Sf()/mesh.magSf();
00026     surfaceVectorField nfVec = fvc::interpolate(n);
00027     nfVec += SfHat*(fvc::snGrad(b) - (SfHat & nfVec));
00028     nfVec /= (mag(nfVec) + dMgb);
00029     surfaceScalarField nf = (mesh.Sf() & nfVec);
00030     n /= mgb;
00031 
00032 
00033 #   include "StCorr.H"
00034 
00035     // Calculate turbulent flame speed flux
00036     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
00037     surfaceScalarField phiSt = fvc::interpolate(rhou*StCorr*Su*Xi)*nf;
00038 
00039     scalar StCoNum = max
00040     (
00041         mesh.surfaceInterpolation::deltaCoeffs()
00042        *mag(phiSt)/(fvc::interpolate(rho)*mesh.magSf())
00043     ).value()*runTime.deltaT().value();
00044 
00045     Info<< "Max St-Courant Number = " << StCoNum << endl;
00046 
00047     // Create b equation
00048     // ~~~~~~~~~~~~~~~~~
00049     fvScalarMatrix bEqn
00050     (
00051         fvm::ddt(rho, b)
00052       + mvConvection->fvmDiv(phi, b)
00053       + fvm::div(phiSt, b, "div(phiSt,b)")
00054       - fvm::Sp(fvc::div(phiSt), b)
00055       - fvm::laplacian(turbulence->alphaEff(), b)
00056     );
00057 
00058 
00059     // Add ignition cell contribution to b-equation
00060     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
00061 #   include "ignite.H"
00062 
00063 
00064     // Solve for b
00065     // ~~~~~~~~~~~
00066     bEqn.solve();
00067 
00068     Info<< "min(b) = " << min(b).value() << endl;
00069 
00070 
00071     // Calculate coefficients for Gulder's flame speed correlation
00072     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
00073 
00074     volScalarField up = uPrimeCoef*sqrt((2.0/3.0)*turbulence->k());
00075   //volScalarField up = sqrt(mag(diag(n * n) & diag(turbulence->r())));
00076 
00077     volScalarField epsilon = pow(uPrimeCoef, 3)*turbulence->epsilon();
00078 
00079     volScalarField tauEta = sqrt(thermo.muu()/(rhou*epsilon));
00080 
00081     volScalarField Reta = up/
00082     (
00083         sqrt(epsilon*tauEta) + dimensionedScalar("1e-8", up.dimensions(), 1e-8)
00084     );
00085 
00086   //volScalarField l = 0.337*k*sqrt(k)/epsilon;
00087   //Reta *= max((l - dimensionedScalar("dl", dimLength, 1.5e-3))/l, 0.0);
00088 
00089     // Calculate Xi flux
00090     // ~~~~~~~~~~~~~~~~~
00091     surfaceScalarField phiXi =
00092         phiSt
00093       - fvc::interpolate(fvc::laplacian(turbulence->alphaEff(), b)/mgb)*nf
00094       + fvc::interpolate(rho)*fvc::interpolate(Su*(1.0/Xi - Xi))*nf;
00095 
00096 
00097     // Calculate mean and turbulent strain rates
00098     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
00099 
00100     volVectorField Ut = U + Su*Xi*n;
00101     volScalarField sigmat = (n & n)*fvc::div(Ut) - (n & fvc::grad(Ut) & n);
00102 
00103     volScalarField sigmas =
00104         ((n & n)*fvc::div(U) - (n & fvc::grad(U) & n))/Xi
00105       + (
00106             (n & n)*fvc::div(Su*n)
00107           - (n & fvc::grad(Su*n) & n)
00108         )*(Xi + scalar(1))/(2*Xi);
00109 
00110 
00111     // Calculate the unstrained laminar flame speed
00112     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
00113     volScalarField Su0 = unstrainedLaminarFlameSpeed()();
00114 
00115 
00116     // Calculate the laminar flame speed in equilibrium with the applied strain
00117     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
00118     volScalarField SuInf = Su0*max(scalar(1) - sigmas/sigmaExt, scalar(0.01));
00119 
00120     if (SuModel == "unstrained")
00121     {
00122         Su == Su0;
00123     }
00124     else if (SuModel == "equilibrium")
00125     {
00126         Su == SuInf;
00127     }
00128     else if (SuModel == "transport")
00129     {
00130         // Solve for the strained laminar flame speed
00131         // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
00132 
00133         volScalarField Rc =
00134             (sigmas*SuInf*(Su0 - SuInf) + sqr(SuMin)*sigmaExt)
00135            /(sqr(Su0 - SuInf) + sqr(SuMin));
00136 
00137         solve
00138         (
00139             fvm::ddt(rho, Su)
00140           + fvm::div(phi + phiXi, Su, "div(phiXi,Su)")
00141           - fvm::Sp(fvc::div(phiXi), Su)
00142           ==
00143           - fvm::SuSp(-rho*Rc*Su0/Su, Su)
00144           - fvm::SuSp(rho*(sigmas + Rc), Su)
00145         );
00146 
00147         // Limit the maximum Su
00148         // ~~~~~~~~~~~~~~~~~~~~
00149         Su.min(SuMax);
00150         Su.max(SuMin);
00151     }
00152     else
00153     {
00154         FatalError
00155             << args.executable() << " : Unknown Su model " << SuModel
00156             << abort(FatalError);
00157     }
00158 
00159 
00160     // Calculate Xi according to the selected flame wrinkling model
00161     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
00162 
00163     if (XiModel == "fixed")
00164     {
00165         // Do nothing, Xi is fixed!
00166     }
00167     else if (XiModel == "algebraic")
00168     {
00169         // Simple algebraic model for Xi based on Gulders correlation
00170         // with a linear correction function to give a plausible profile for Xi
00171         // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
00172         Xi == scalar(1) +
00173             (scalar(1) + (2*XiShapeCoef)*(scalar(0.5) - b))
00174            *XiCoef*sqrt(up/(Su + SuMin))*Reta;
00175     }
00176     else if (XiModel == "transport")
00177     {
00178         // Calculate Xi transport coefficients based on Gulders correlation
00179         // and DNS data for the rate of generation
00180         // with a linear correction function to give a plausible profile for Xi
00181         // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
00182 
00183         volScalarField XiEqStar =
00184             scalar(1.001) + XiCoef*sqrt(up/(Su + SuMin))*Reta;
00185 
00186         volScalarField XiEq =
00187             scalar(1.001)
00188           + (scalar(1) + (2*XiShapeCoef)*(scalar(0.5) - b))
00189            *(XiEqStar - scalar(1.001));
00190 
00191         volScalarField Gstar = 0.28/tauEta;
00192         volScalarField R = Gstar*XiEqStar/(XiEqStar - scalar(1));
00193         volScalarField G = R*(XiEq - scalar(1.001))/XiEq;
00194 
00195         //R *= (Gstar + 2*mag(dev(symm(fvc::grad(U)))))/Gstar;
00196 
00197         // Solve for the flame wrinkling
00198         // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
00199         solve
00200         (
00201             fvm::ddt(rho, Xi)
00202           + fvm::div(phi + phiXi, Xi, "div(phiXi,Xi)")
00203           - fvm::Sp(fvc::div(phiXi), Xi)
00204          ==
00205             rho*R
00206           - fvm::Sp(rho*(R - G), Xi)
00207           - fvm::Sp
00208             (
00209                 rho*max
00210                 (
00211                     sigmat - sigmas,
00212                     dimensionedScalar("0", sigmat.dimensions(), 0)
00213                 ),
00214                 Xi
00215             )
00216         );
00217 
00218 
00219         // Correct boundedness of Xi
00220         // ~~~~~~~~~~~~~~~~~~~~~~~~~
00221         Xi.max(1.0);
00222         Info<< "max(Xi) = " << max(Xi).value() << endl;
00223         Info<< "max(XiEq) = " << max(XiEq).value() << endl;
00224     }
00225     else
00226     {
00227         FatalError
00228             << args.executable() << " : Unknown Xi model " << XiModel
00229             << abort(FatalError);
00230     }
00231 
00232     Info<< "Combustion progress = "
00233         << 100*(scalar(1) - b)().weightedAverage(mesh.V()).value() << "%"
00234         << endl;
00235 
00236     St = Xi*Su;
00237 }
Copyright © 2000-2009 OpenCFD Ltd