thermophysicalPredictor.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2022-2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "XiFluid.H"
27 #include "bXiIgnition.H"
28 #include "fvcDdt.H"
29 #include "fvmDiv.H"
30 #include "fvcSup.H"
31 
32 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
33 
35 (
37  const volScalarField& Db
38 )
39 {
40  volScalarField& ft = thermo_.Y("ft");
41 
42  fvScalarMatrix ftEqn
43  (
44  fvm::ddt(rho, ft)
45  + mvConvection.fvmDiv(phi, ft)
46  - fvm::laplacian(Db, ft)
47  ==
48  fvModels().source(rho, ft)
49  );
50 
51  ftEqn.relax();
52  fvConstraints().constrain(ftEqn);
53  ftEqn.solve();
55 }
56 
57 
59 (
61  const volScalarField& Db,
62  const volScalarField& bSource
63 )
64 {
65  volScalarField& fu = thermo_.Y("fu");
66  const volScalarField& ft = thermo_.Y("ft");
67  const volScalarField& b(b_);
68 
69  // Progress variable
70  const volScalarField c("c", scalar(1) - b);
71 
72  // Unburnt gas density
73  const volScalarField rhou("rhou", thermo.rhou());
74 
75  const volScalarField fres(thermo_.fres());
76  const volScalarField fuFres(max(fu - fres, scalar(0)));
77 
78  const volScalarField fuDot
79  (
80  bSource/(b + 0.001)
81  - 2*Db*c
82  *mag(fvc::grad(ft))/max(ft, scalar(1e-6))
83  *mag(fvc::grad(fuFres))/max(fuFres, scalar(1e-6))
84  );
85 
86  fvScalarMatrix fuEqn
87  (
88  fvm::ddt(rho, fu)
89  + mvConvection.fvmDiv(phi, fu)
90  - fvm::laplacian(Db, fu)
91  ==
92  fvm::Sp(fuDot, fu)
93  - fuDot*fres
94  + fvModels().source(rho, fu)
95  );
96 
97  fuEqn.relax();
98  fvConstraints().constrain(fuEqn);
99  fuEqn.solve();
100  fvConstraints().constrain(fu);
101 }
102 
103 
105 (
107  const volScalarField& Db
108 )
109 {
110  volScalarField& egr = thermo_.Y("egr");
111 
112  fvScalarMatrix egrEqn
113  (
114  fvm::ddt(rho, egr)
115  + mvConvection.fvmDiv(phi, egr)
116  - fvm::laplacian(Db, egr)
117  ==
118  fvModels().source(rho, egr)
119  );
120 
121  egrEqn.relax();
122  fvConstraints().constrain(egrEqn);
123  egrEqn.solve();
124  fvConstraints().constrain(egr);
125 }
126 
127 
128 
130 (
131  const volScalarField& Xi,
132  const surfaceScalarField& nf,
133  const dimensionedScalar& dMgb
134 ) const
135 {
136  const UPtrListDictionary<fv::bXiIgnition> ignitionModels
137  (
138  fvModels().lookupType<fv::bXiIgnition>()
139  );
140 
141  bool igniting = false;
142 
143  forAll(ignitionModels, i)
144  {
145  if (ignitionModels[i].igniting())
146  {
147  igniting = true;
148  break;
149  }
150  }
151 
152  if (igniting)
153  {
154  tmp<volScalarField> tXi(volScalarField::New("XiCorrected", Xi));
155 
156  // Calculate kernel area from b field consistent with the
157  // discretisation of the b equation.
158  const volScalarField mgb
159  (
160  fvc::div(nf, b, "div(phiSt,b)") - b*fvc::div(nf) + dMgb
161  );
162 
163  forAll(ignitionModels, i)
164  {
165  ignitionModels[i].XiCorr(tXi.ref(), b, mgb);
166  }
167 
168  return tXi;
169  }
170  else
171  {
172  return Xi;
173  }
174 }
175 
176 
178 (
180  const volScalarField& Db
181 )
182 {
183  volScalarField& b(b_);
184 
185  // Progress variable
186  const volScalarField c("c", scalar(1) - b);
187 
188  // Unburnt gas density
189  const volScalarField rhou("rhou", thermo.rhou());
190 
191  // Calculate flame normal etc.
192  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
193 
194  volVectorField n("n", fvc::grad(b));
195 
196  volScalarField mgb("mgb", mag(n));
197 
198  const dimensionedScalar dMgb
199  (
200  "dMgb",
201  1.0e-3*
202  (b*c*mgb)().weightedAverage(mesh.V())
203  /((b*c)().weightedAverage(mesh.V()) + small)
204  + dimensionedScalar(mgb.dimensions(), small)
205  );
206 
207  mgb += dMgb;
208 
209  const surfaceVectorField SfHat(mesh.Sf()/mesh.magSf());
211  nfVec += SfHat*(fvc::snGrad(b) - (SfHat & nfVec));
212  nfVec /= (mag(nfVec) + dMgb);
213  surfaceScalarField nf("nf", mesh.Sf() & nfVec);
214  n /= mgb;
215 
216  // Turbulent flame speed flux
217  const surfaceScalarField phiSt
218  (
219  "phiSt",
220  fvc::interpolate(rhou*Su*XiCorr(Xi, nf, dMgb))*nf
221  );
222 
223  fvScalarMatrix bSourceEqn
224  (
225  fvModels().source(rho, b)
226  - fvm::div(phiSt, b)
227  + fvm::Sp(fvc::div(phiSt), b)
228  );
229 
230  // Create b equation
231  fvScalarMatrix bEqn
232  (
233  fvm::ddt(rho, b)
234  + mvConvection.fvmDiv(phi, b)
235  - fvm::laplacian(Db, b)
236  ==
237  bSourceEqn
238  );
239 
240  // Solve for b and constrain
241  bEqn.relax();
242  fvConstraints().constrain(bEqn);
243  bEqn.solve();
245 
246  // Correct the flame wrinkling
247  XiModel_->correct();
248 
249  // Correct the laminar flame speed
250  SuModel_->correct();
251 
252  if (thermo_.containsSpecie("fu"))
253  {
254  fuSolve(mvConvection, Db, (bSourceEqn & b));
255  }
256 }
257 
258 
260 (
262  const volScalarField& Db
263 )
264 {
265  volScalarField& heau = thermo_.heu();
266 
267  const volScalarField::Internal rhoByRhou(rho()/thermo.rhou()());
268 
269  fvScalarMatrix EauEqn
270  (
271  fvm::ddt(rho, heau) + mvConvection.fvmDiv(phi, heau)
272  + rhoByRhou
273  *(
274  (fvc::ddt(rho, K) + fvc::div(phi, K))()
275  + pressureWork
276  (
277  heau.name() == "eau"
278  ? mvConvection.fvcDiv(phi, p/rho)()
279  : -dpdt
280  )
281  )
282  - fvm::laplacian(Db, heau)
283 
284  // These terms cannot be used in partially-premixed combustion due to
285  // the resultant inconsistency between ft and heau transport.
286  // A possible solution would be to solve for ftu as well as ft.
287  //- fvm::div(muEff*fvc::grad(b)/(b + 0.001), heau)
288  //+ fvm::Sp(fvc::div(muEff*fvc::grad(b)/(b + 0.001)), heau)
289 
290  ==
291  fvModels().source(rho, heau)
292  );
293 
294  EauEqn.relax();
295  fvConstraints().constrain(EauEqn);
296  EauEqn.solve();
297  fvConstraints().constrain(heau);
298 }
299 
300 
302 (
304  const volScalarField& Db
305 )
306 {
307  volScalarField& hea = thermo_.he();
308 
310  (
311  fvm::ddt(rho, hea) + mvConvection.fvmDiv(phi, hea)
312  + fvc::ddt(rho, K) + fvc::div(phi, K)
313  + pressureWork
314  (
315  hea.name() == "ea"
316  ? mvConvection.fvcDiv(phi, p/rho)()
317  : -dpdt
318  )
319  - fvm::laplacian(Db, hea)
320  ==
321  (
322  buoyancy.valid()
323  ? fvModels().source(rho, hea) + rho*(U & buoyancy->g)
324  : fvModels().source(rho, hea)
325  )
326  );
327 
328  EaEqn.relax();
330  EaEqn.solve();
331  fvConstraints().constrain(hea);
332 }
333 
334 
336 {
337  const UPtrListDictionary<fv::bXiIgnition> ignitionModels
338  (
339  fvModels().lookupType<fv::bXiIgnition>()
340  );
341 
342  bool ignited = false;
343 
344  forAll(ignitionModels, i)
345  {
346  if (ignitionModels[i].ignited())
347  {
348  ignited = true;
349  break;
350  }
351  }
352 
354  (
356  (
357  mesh,
358  fields,
359  phi,
360  mesh.schemes().div("div(phi,ft_b_ha_hau)")
361  )
362  );
363 
364  const volScalarField Db
365  (
366  "Db",
367  ignited ? XiModel_->Db() : thermophysicalTransport.DEff(b)
368  );
369 
370  if (thermo_.containsSpecie("ft"))
371  {
372  ftSolve(mvConvection(), Db);
373  }
374 
375  if (ignited)
376  {
377  bSolve(mvConvection(), Db);
378  EauSolve(mvConvection(), Db);
379  }
380  else
381  {
382  if (thermo_.containsSpecie("fu"))
383  {
384  volScalarField& fu = thermo_.Y("fu");
385  const volScalarField& ft = thermo_.Y("ft");
386  fu = ft;
387  }
388  }
389 
390  if (thermo_.containsSpecie("egr"))
391  {
392  egrSolve(mvConvection(), Db);
393  }
394 
395  EaSolve(mvConvection(), Db);
396 
397  if (!ignited)
398  {
399  thermo_.heu() == thermo.he();
400  }
401 
402  thermo_.correct();
403 }
404 
405 
406 // ************************************************************************* //
fvScalarMatrix EaEqn(betav *fvm::ddt(rho, hea)+mvConvection->fvmDiv(phi, hea)+betav *fvc::ddt(rho, K)+fvc::div(phi, K)+(hea.name()=="ea" ? fvc::div(fvc::absolute(phi, rho, U), p/rho) :-betav *dpdt) - fvm::laplacian(Db, hea)+betav *fvModels.source(rho, hea))
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
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)")))
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
const word & name() const
Return name.
Definition: IOobject.H:307
Template dictionary class which manages the storage associated with it.
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:603
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:58
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1799
const surfaceVectorField & Sf() const
Return cell face area vectors.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
Finite volume models.
Definition: fvModels.H:65
ITstream & div(const word &name) const
Definition: fvSchemes.C:423
Abstract base class for convection schemes.
virtual PtrList< volScalarField > & Y()=0
Access the mass-fraction fields.
Foam::fvConstraints & fvConstraints() const
Return the fvConstraints that are created on demand.
Definition: solver.C:107
virtual void thermophysicalPredictor()
Construct and solve the energy equation,.
psiuMulticomponentThermo & thermo_
Definition: XiFluid.H:149
void EaSolve(const fv::convectionScheme< scalar > &mvConvection, const volScalarField &Db)
Solve the energy equation.
void EauSolve(const fv::convectionScheme< scalar > &mvConvection, const volScalarField &Db)
Solve the unburnt energy equation.
tmp< volScalarField > XiCorr(const volScalarField &Xi, const surfaceScalarField &nf, const dimensionedScalar &dMgb) const
Apply the early kernel growth correction to the flame-wrinkling Xi.
void egrSolve(const fv::convectionScheme< scalar > &mvConvection, const volScalarField &Db)
Solve the egr mass-fraction equation.
void fuSolve(const fv::convectionScheme< scalar > &mvConvection, const volScalarField &Db, const volScalarField &bSource)
Solve the fu equation for partially- and non- premixed mixtures.
void ftSolve(const fv::convectionScheme< scalar > &mvConvection, const volScalarField &Db)
Solve the ft equation for partially-premixed mixtures.
void bSolve(const fv::convectionScheme< scalar > &mvConvection, const volScalarField &Db)
Solve the Xi and regress variable equations.
Buoyancy related data for the Foam::solvers::isothermalFluid solver module when solving buoyant cases...
Definition: buoyancy.H:70
uniformDimensionedVectorField g
Gravitational acceleration.
Definition: buoyancy.H:83
const surfaceScalarField & phi
Mass-flux field.
const volScalarField & rho
Reference to the continuity density field.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
Foam::fvConstraints & fvConstraints(Foam::fvConstraints::New(mesh))
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
Calculate the first temporal derivative.
Calculate the field for explicit evaluation of implicit and explicit sources.
Calculate the matrix for the divergence of the given field and flux.
volScalarField & b
Definition: createFields.H:27
Info<< "Creating thermophysical transport model\n"<< endl;turbulenceThermophysicalTransportModels::unityLewisEddyDiffusivity< RASThermophysicalTransportModel< ThermophysicalTransportModel< compressibleMomentumTransportModel, fluidThermo > >> thermophysicalTransport(turbulence(), thermo, true)
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:234
K
Definition: pEqn.H:75
U
Definition: pEqn.H:72
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 > > Su(const VolField< Type > &su, const VolField< Type > &vf)
Definition: fvcSup.C:44
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
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
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
const doubleScalar e
Definition: doubleScalar.H:106
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
volScalarField & p
fluidMulticomponentThermo & thermo
Definition: createFields.H:31