Go to the documentation of this file.00001 if (ign.ignited())
00002 {
00003
00004
00005 volScalarField c = scalar(1) - b;
00006
00007
00008
00009 volScalarField rhou = thermo.rhou();
00010
00011
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 <engine/StCorr.H>
00034
00035
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
00048
00049 fvScalarMatrix bEqn
00050 (
00051 fvm::ddt(rho, b)
00052 + mvConvection->fvmDiv(phi, b)
00053 + fvm::div(phiSt, b, "div(phiSt,bprog)")
00054 - fvm::Sp(fvc::div(phiSt), b)
00055 - fvm::laplacian(turbulence->alphaEff(), b)
00056 );
00057
00058
00059
00060
00061 #include <engine/ignite.H>
00062
00063
00064
00065
00066 bEqn.solve();
00067
00068 Info<< "min(bprog) = " << min(b).value() << endl;
00069
00070
00071
00072
00073
00074 volScalarField up = uPrimeCoef*sqrt((2.0/3.0)*turbulence->k());
00075
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
00087
00088
00089
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
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
00112
00113 volScalarField Su0 = unstrainedLaminarFlameSpeed()();
00114
00115
00116
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
00131
00132
00133 volScalarField Rc =
00134 (sigmas*SuInf*(Su0 - SuInf) + sqr(SuMin)*sigmaExt)
00135 /(sqr(Su0 - SuInf) + sqr(SuMin));
00136
00137 fvScalarMatrix SuEqn
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 SuEqn.relax();
00148 SuEqn.solve();
00149
00150
00151
00152 Su.min(SuMax);
00153 Su.max(SuMin);
00154 }
00155 else
00156 {
00157 FatalError
00158 << args.executable() << " : Unknown Su model " << SuModel
00159 << abort(FatalError);
00160 }
00161
00162
00163
00164
00165
00166 if (XiModel == "fixed")
00167 {
00168
00169 }
00170 else if (XiModel == "algebraic")
00171 {
00172
00173
00174
00175 Xi == scalar(1) +
00176 (scalar(1) + (2*XiShapeCoef)*(scalar(0.5) - b))
00177 *XiCoef*sqrt(up/(Su + SuMin))*Reta;
00178 }
00179 else if (XiModel == "transport")
00180 {
00181
00182
00183
00184
00185
00186 volScalarField XiEqStar =
00187 scalar(1.001) + XiCoef*sqrt(up/(Su + SuMin))*Reta;
00188
00189 volScalarField XiEq =
00190 scalar(1.001)
00191 + (scalar(1) + (2*XiShapeCoef)*(scalar(0.5) - b))
00192 *(XiEqStar - scalar(1.001));
00193
00194 volScalarField Gstar = 0.28/tauEta;
00195 volScalarField R = Gstar*XiEqStar/(XiEqStar - scalar(1));
00196 volScalarField G = R*(XiEq - scalar(1.001))/XiEq;
00197
00198
00199
00200
00201
00202 fvScalarMatrix XiEqn
00203 (
00204 fvm::ddt(rho, Xi)
00205 + fvm::div(phi + phiXi, Xi, "div(phiXi,Xi)")
00206 - fvm::Sp(fvc::div(phiXi), Xi)
00207 ==
00208 rho*R
00209 - fvm::Sp(rho*(R - G), Xi)
00210 - fvm::Sp
00211 (
00212 rho*max
00213 (
00214 sigmat - sigmas,
00215 dimensionedScalar("0", sigmat.dimensions(), 0)
00216 ),
00217 Xi
00218 )
00219 );
00220
00221 XiEqn.relax();
00222 XiEqn.solve();
00223
00224
00225
00226 Xi.max(1.0);
00227 Info<< "max(Xi) = " << max(Xi).value() << endl;
00228 Info<< "max(XiEq) = " << max(XiEq).value() << endl;
00229 }
00230 else
00231 {
00232 FatalError
00233 << args.executable() << " : Unknown Xi model " << XiModel
00234 << abort(FatalError);
00235 }
00236
00237 Info<< "Combustion progress = "
00238 << 100*(scalar(1) - b)().weightedAverage(mesh.V()).value() << "%"
00239 << endl;
00240
00241 St = Xi*Su;
00242 }
00243
00244