YEqn.H
Go to the documentation of this file.
1tmp<fv::convectionScheme<scalar>> mvConvection
2(
3 fv::convectionScheme<scalar>::New
4 (
5 mesh,
6 fields,
7 phi,
8 mesh.divScheme("div(phi,Yi_h)")
9 )
10);
11
12{
13 reaction->correct();
14 Qdot = reaction->Qdot();
15 volScalarField Yt(0.0*Y[0]);
16
18 {
19 if (i != inertIndex && composition.active(i))
20 {
21 volScalarField& Yi = Y[i];
22
23 fvScalarMatrix YiEqn
24 (
25 fvm::ddt(rho, Yi)
26 + mvConvection->fvmDiv(phi, Yi)
27 - fvm::laplacian(turbulence->muEff(), Yi)
28 ==
29 reaction->R(Yi)
30 + fvOptions(rho, Yi)
31 );
32
33 YiEqn.relax();
34
35 fvOptions.constrain(YiEqn);
36
37 YiEqn.solve(mesh.solver("Yi"));
38
39 fvOptions.correct(Yi);
40
41 Yi.max(0.0);
42 Yt += Yi;
43 }
44 }
45
46 Y[inertIndex] = scalar(1) - Yt;
47 Y[inertIndex].max(0.0);
48}
fv::options & fvOptions
surfaceScalarField & phi
basicSpecieMixture & composition
PtrList< volScalarField > & Y
Qdot
Definition: YEqn.H:14
volScalarField Yt(0.0 *Y[0])
tmp< fv::convectionScheme< scalar > > mvConvection(fv::convectionScheme< scalar >::New(mesh, fields, phi, mesh.divScheme("div(phi,Yi_h)")))
dynamicFvMesh & mesh
label inertIndex
compressible::turbulenceModel & turbulence
CombustionModel< rhoReactionThermo > & reaction
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333