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.schemes().div("div(phi,ft_b_ha_hau)")
9  )
10 );
11 
12 volScalarField Db("Db", rho*turbulence->nuEff());
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  (
71  betav*fvm::ddt(rho, b)
72  + mvConvection->fvmDiv(phi, b)
73  + fvm::div(phiSt, b)
74  - fvm::Sp(fvc::div(phiSt), b)
75  - fvm::laplacian(Db, b)
76  ==
77  betav*fvModels.source(rho, b)
78  );
79 
80 
81  // Add ignition cell contribution to b-equation
82  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
83  #include "ignite.H"
84 
85  // Solve for b
86  // ~~~~~~~~~~~
87  bEqn.relax();
88 
90 
91  bEqn.solve();
92 
94 
95  Info<< "min(b) = " << min(b).value() << endl;
96 
97  if (thermo.containsSpecie("ft"))
98  {
99  volScalarField& ft = thermo.Y("ft");
100 
101  Info<< "Combustion progress = "
102  << 100*(1.0 - b)().weightedAverage(mesh.V()*ft).value() << "%"
103  << endl;
104  }
105  else
106  {
107  Info<< "Combustion progress = "
108  << 100*(1.0 - b)().weightedAverage(mesh.V()).value() << "%"
109  << endl;
110  }
111 
112  // Correct the flame-wrinkling
113  flameWrinkling->correct(mvConvection);
114 
115  St = Xi*Su;
116 }
dimensionedScalar StCorr("StCorr", dimless, 1.0)
Calculates and outputs the mean and maximum Courant Numbers.
label n
volScalarField Db("Db", rho *turbulence->nuEff())
tmp< fv::convectionScheme< scalar > > mvConvection(fv::convectionScheme< scalar >::New(mesh, fields, phi, mesh.schemes().div("div(phi,ft_b_ha_hau)")))
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
tmp< fvMatrix< Type > > source(const VolField< Type > &field) const
Return source for an equation.
Foam::fvConstraints & fvConstraints(Foam::fvConstraints::New(mesh))
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
autoPtr< XiModel > flameWrinkling
Create the flame-wrinkling model.
Definition: createFields.H:202
Info<< "Creating field dpdt\n"<< endl;volScalarField dpdt(IOobject("dpdt", runTime.name(), mesh), mesh, dimensionedScalar(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))
volScalarField & b
Definition: createFields.H:25
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:228
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
const dimensionedScalar c
Speed of light in a vacuum.
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< Type > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
tmp< VolField< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvcLaplacian.C:45
tmp< VolField< typename outerProduct< vector, Type >::type > > reconstruct(const SurfaceField< Type > &ssf)
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
Definition: fvcSup.C:44
tmp< VolField< Type > > Sp(const volScalarField &sp, const VolField< Type > &vf)
Definition: fvcSup.C:67
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< SurfaceField< Type > > snGrad(const VolField< Type > &vf, const word &name)
Definition: fvcSnGrad.C:45
VolField< vector > volVectorField
Definition: volFieldsFwd.H:65
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
SurfaceField< scalar > surfaceScalarField
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
dimensioned< scalar > mag(const dimensioned< Type > &)
SurfaceField< vector > surfaceVectorField
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
autoPtr< incompressible::momentumTransportModel > turbulence(incompressible::momentumTransportModel::New(U, phi, viscosity))
fluidMulticomponentThermo & thermo
Definition: createFields.H:31