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 reaction.correct();
18 
19 forAll(Y, i)
20 {
21  if (composition.solve(i))
22  {
23  volScalarField& Yi = Y[i];
24 
25  fvScalarMatrix YiEqn
26  (
27  fvm::ddt(rho, Yi) + mvConvection->fvmDiv(phi, Yi)
28  + thermophysicalTransport.divj(Yi)
29  ==
30  reaction.R(Yi)
31  + fvModels.source(rho, Yi)
32  );
33 
34  YiEqn.relax();
35 
36  fvConstraints.constrain(YiEqn);
37 
38  YiEqn.solve("Yi");
39 
41  }
42 }
43 
44 composition.normalise();
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
basicSpecieMixture & composition
combustionModel & reaction
Info<< "Calculating turbulent flame speed field St\"<< endl;volScalarField St(IOobject("St", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:228
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:57
dynamicFvMesh & mesh
Foam::fvConstraints & fvConstraints
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
phi
Definition: correctPhi.H:3
fluidReactionThermophysicalTransportModel & thermophysicalTransport
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
Foam::fvModels & fvModels
PtrList< volScalarField > & Y
forAll(Y, i)
Definition: YEqn.H:14
tmp< fv::convectionScheme< scalar > > mvConvection(fv::convectionScheme< scalar >::New(mesh, fields, phi, mesh.divScheme("div(phi,Yi_h)")))