bEqn.H
Go to the documentation of this file.
1 tmp<fv::convectionScheme<scalar> > mvConvection
2 (
4  (
5  mesh,
6  fields,
7  phi,
8  mesh.divScheme("div(phi,ft_b_ha_hau)")
9  )
10 );
11 
12 volScalarField Db("Db", turbulence->muEff());
13 
14 if (ign.ignited())
15 {
16  // Calculate the unstrained laminar flame speed
17  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
19 
20  const volScalarField& Xi = flameWrinkling->Xi();
21 
22  // progress variable
23  // ~~~~~~~~~~~~~~~~~
24  volScalarField c("c", 1.0 - b);
25 
26  // Unburnt gas density
27  // ~~~~~~~~~~~~~~~~~~~
28  volScalarField rhou(thermo.rhou());
29 
30  // Calculate flame normal etc.
31  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
32 
33  //volVectorField n(fvc::grad(b));
35 
36  volScalarField mgb("mgb", mag(n));
37 
38  dimensionedScalar dMgb("dMgb", mgb.dimensions(), SMALL);
39 
40  {
41  volScalarField bc(b*c);
42 
43  dMgb += 1.0e-3*
44  (bc*mgb)().weightedAverage(mesh.V())
45  /(bc.weightedAverage(mesh.V()) + SMALL);
46  }
47 
48  mgb += dMgb;
49 
50  surfaceVectorField Sfhat(mesh.Sf()/mesh.magSf());
52  nfVec += Sfhat*(fvc::snGrad(b) - (Sfhat & nfVec));
53  nfVec /= (mag(nfVec) + dMgb);
54  surfaceScalarField nf("nf", mesh.Sf() & nfVec);
55  n /= mgb;
56 
57  #include "StCorr.H"
58 
59  // Calculate turbulent flame speed flux
60  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
61  surfaceScalarField phiSt("phiSt", fvc::interpolate(rhou*StCorr*St)*nf);
62 
63  #include "StCourantNo.H"
64 
65  Db = flameWrinkling->Db();
66 
67  // Create b equation
68  // ~~~~~~~~~~~~~~~~~
69  fvScalarMatrix bEqn
70  (
72  + mvConvection->fvmDiv(phi, b)
73  + fvm::div(phiSt, b)
74  - fvm::Sp(fvc::div(phiSt), b)
75  - fvm::laplacian(Db, b)
76  );
77 
78 
79  // Add ignition cell contribution to b-equation
80  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
81  #include "ignite.H"
82 
83  // Solve for b
84  // ~~~~~~~~~~~
85  bEqn.solve();
86 
87  Info<< "min(b) = " << min(b).value() << endl;
88 
89  if (composition.contains("ft"))
90  {
91  volScalarField& ft = composition.Y("ft");
92 
93  Info<< "Combustion progress = "
94  << 100*(1.0 - b)().weightedAverage(mesh.V()*ft).value() << "%"
95  << endl;
96  }
97  else
98  {
99  Info<< "Combustion progress = "
100  << 100*(1.0 - b)().weightedAverage(mesh.V()).value() << "%"
101  << endl;
102  }
103 
104  // Correct the flame-wrinkling
105  flameWrinkling->correct(mvConvection);
106 
107  St = Xi*Su;
108 }
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:28
Info<< "Creating field dpdt\n"<< endl;volScalarField dpdt(IOobject("dpdt", runTime.timeName(), mesh), mesh, dimensionedScalar("dpdt", p.dimensions()/dimTime, 0));Info<< "Creating field kinetic energy K\n"<< endl;volScalarField K("K", 0.5 *magSqr(U));volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);p_rgh=p-rho *gh;mesh.setFluxRequired(p_rgh.name());multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:127
dimensioned< scalar > mag(const dimensioned< Type > &)
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< GeometricField< Type, fvPatchField, volMesh > > Sp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:67
const volScalarField & betav
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
messageStream Info
tmp< surfaceScalarField > interpolate(const RhoType &rho)
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
dynamicFvMesh & mesh
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:44
label n
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
basicMultiComponentMixture & composition
Definition: createFields.H:35
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Info<< "Creating field dpdt\n"<< endl;volScalarField dpdt(IOobject("dpdt", runTime.timeName(), mesh), mesh, dimensionedScalar("dpdt", p.dimensions()/dimTime, 0));Info<< "Creating field kinetic energy K\n"<< endl;volScalarField K("K", 0.5 *magSqr(U));Info<< "Creating the unstrained laminar flame speed\n"<< endl;autoPtr< laminarFlameSpeed > unstrainedLaminarFlameSpeed(laminarFlameSpeed::New(thermo))
Definition: createFields.H:81
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Calculates and outputs the mean and maximum Courant Numbers.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
volScalarField Db("Db", turbulence->muEff())
surfaceScalarField & phi
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
const dimensionedScalar c
Speed of light in a vacuum.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
psiReactionThermo & thermo
Definition: createFields.H:32
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedScalar StCorr("StCorr", dimless, 1.0)
tmp< fv::convectionScheme< scalar > > mvConvection(fv::convectionScheme< scalar >::New( mesh, fields, phi, mesh.divScheme("div(phi,ft_b_ha_hau)") ))
autoPtr< XiModel > flameWrinkling
Create the flame-wrinkling model.
Definition: createFields.H:193
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.