YEqn.H
Go to the documentation of this file.
1 tmp<fv::convectionScheme<scalar>> mvConvection(nullptr);
2 
3 if (Y.size())
4 {
5  mvConvection = tmp<fv::convectionScheme<scalar>>
6  (
8  (
9  mesh,
10  fields,
11  phi,
12  mesh.divScheme("div(phi,Yi_h)")
13  )
14  );
15 }
16 
17 {
18  reaction.correct();
19 
20  tmp<volScalarField> Yt(nullptr);
21 
22  if (Y.size())
23  {
24  Yt = tmp<volScalarField>
25  (
26  new volScalarField
27  (
28  IOobject("Yt", runTime.timeName(), mesh),
29  mesh,
31  )
32  );
33  }
34 
35  forAll(Y, i)
36  {
37  if (i != inertIndex && composition.active(i))
38  {
39  volScalarField& Yi = Y[i];
40 
41  fvScalarMatrix YiEqn
42  (
43  fvm::ddt(rho, Yi)
44  + mvConvection->fvmDiv(phi, Yi)
45  - fvm::laplacian(turbulence.muEff(), Yi)
46  ==
47  reaction.R(Yi)
48  + fvOptions(rho, Yi)
49  );
50 
51  YiEqn.relax();
52 
53  fvOptions.constrain(YiEqn);
54 
55  YiEqn.solve("Yi");
56 
57  fvOptions.correct(Yi);
58 
59  Yi.max(0.0);
60  Yt.ref() += Yi;
61  }
62  }
63 
64  if (Y.size())
65  {
66  Y[inertIndex] = scalar(1) - Yt;
67  Y[inertIndex].max(0.0);
68  }
69 }
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
fv::options & fvOptions
surfaceScalarField & phi
basicSpecieMixture & composition
engineTime & runTime
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
label inertIndex
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dynamicFvMesh & mesh
Y[inertIndex]
Definition: YEqn.H:45
CombustionModel< rhoReactionThermo > & reaction
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
forAll(Y, i)
Definition: YEqn.H:16
tmp< fv::convectionScheme< scalar > > mvConvection(fv::convectionScheme< scalar >::New(mesh, fields, phi, mesh.divScheme("div(phi,Yi_h)")))
volScalarField Yt(0.0 *Y[0])