FSD.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) 2011-2024 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 "FSD.H"
28 #include "LESModel.H"
29 #include "fvcGrad.H"
30 #include "fvcDiv.H"
32 
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace combustionModels
39 {
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
49 (
50  const word& modelType,
53  const word& combustionProperties
54 )
55 :
57  (
58  modelType,
59  thermo,
60  turb,
61  combustionProperties
62  ),
63  reactionRateFlameArea_
64  (
66  (
67  this->coeffs(),
68  this->mesh(),
69  *this
70  )
71  ),
72  ft_
73  (
74  IOobject
75  (
76  this->thermo().phasePropertyName("ft"),
77  this->mesh().time().name(),
78  this->mesh(),
79  IOobject::NO_READ,
80  IOobject::AUTO_WRITE
81  ),
82  this->mesh(),
84  ),
85  YFuelFuelStream_(dimensionedScalar(dimless, 1.0)),
86  YO2OxiStream_(dimensionedScalar(dimless, 0.23)),
87  Cv_(this->coeffs().template lookup<scalar>("Cv")),
88  C_(5),
89  ftMin_(0),
90  ftMax_(1),
91  ftDim_(300),
92  ftVarMin_(this->coeffs().template lookup<scalar>("ftVarMin"))
93 {}
94 
95 
96 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
97 
99 {}
100 
101 
102 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
103 
104 void Foam::combustionModels::FSD::calculateSourceNorm()
105 {
106  this->fresCorrect();
107 
108  const label fuelI = this->fuelIndex();
109 
110  const volScalarField& YFuel = this->thermo().Y()[fuelI];
111 
112  const volScalarField& YO2 = this->thermo().Y("O2");
113 
114  const dimensionedScalar s = this->s();
115 
116  ft_ =
117  (s*YFuel - (YO2 - YO2OxiStream_))/(s*YFuelFuelStream_ + YO2OxiStream_);
118 
119 
120  volVectorField nft(fvc::grad(ft_));
121 
122  volScalarField mgft(mag(nft));
123 
124  const surfaceVectorField SfHat(this->mesh().Sf()/this->mesh().magSf());
125 
126  const volScalarField cAux(scalar(1) - ft_);
127 
128  const dimensionedScalar dMgft = 1e-3*
129  (ft_*cAux*mgft)().weightedAverage(this->mesh().V())
130  /((ft_*cAux)().weightedAverage(this->mesh().V()) + small)
131  + dimensionedScalar(mgft.dimensions(), small);
132 
133  mgft += dMgft;
134 
135  nft /= mgft;
136 
137  const volVectorField& U = YO2.db().lookupObject<volVectorField>("U");
138 
139  const volScalarField sigma
140  (
141  (nft & nft)*fvc::div(U) - (nft & fvc::grad(U) & nft)
142  );
143 
144  reactionRateFlameArea_->correct(sigma);
145 
146  const volScalarField& omegaFuel = reactionRateFlameArea_->omega();
147 
148  const scalar ftStoich =
149  YO2OxiStream_.value()
150  /(
151  s.value()*YFuelFuelStream_.value() + YO2OxiStream_.value()
152  );
153 
155  (
157  (
158  this->thermo().phasePropertyName("Pc"),
159  U.mesh(),
161  )
162  );
163 
164  volScalarField& pc = tPc.ref();
165 
166  tmp<volScalarField> tomegaFuel
167  (
169  (
170  this->thermo().phasePropertyName("omegaFuelBar"),
171  U.mesh(),
172  dimensionedScalar(omegaFuel.dimensions(), 0)
173  )
174  );
175 
176  volScalarField& omegaFuelBar = tomegaFuel.ref();
177 
178  // Calculation of the mixture fraction variance (ftVar)
179  const compressible::LESModel& lesModel =
181  (
182  momentumTransportModel::typeName
183  );
184 
185  const volScalarField& delta = lesModel.delta();
186  const volScalarField ftVar(Cv_*sqr(delta)*sqr(mgft));
187 
188  // Thickened flame (average flame thickness for counterflow configuration
189  // is 1.5 mm)
190 
191  const volScalarField deltaF
192  (
193  lesModel.delta()/dimensionedScalar(dimLength, 1.5e-3)
194  );
195 
196  // Linear correlation between delta and flame thickness
197  const volScalarField omegaF(max(deltaF*(4.0/3.0) + (2.0/3.0), scalar(1)));
198 
199  const scalar deltaFt = 1/ftDim_;
200 
201  forAll(ft_, celli)
202  {
203  if (ft_[celli] > ftMin_ && ft_[celli] < ftMax_)
204  {
205  const scalar ftCell = ft_[celli];
206 
207  if (ftVar[celli] > ftVarMin_) // sub-grid beta pdf of ft_
208  {
209  const scalar ftVarc = ftVar[celli];
210  const scalar a =
211  max(ftCell*(ftCell*(1 - ftCell)/ftVarc - 1), 0);
212  const scalar b = max(a/ftCell - a, 0);
213 
214  for (int i=1; i<ftDim_; i++)
215  {
216  const scalar ft = i*deltaFt;
217  pc[celli] += pow(ft, a - 1)*pow(1 - ft, b - 1)*deltaFt;
218  }
219 
220  for (int i=1; i<ftDim_; i++)
221  {
222  const scalar ft = i*deltaFt;
223  omegaFuelBar[celli] +=
224  omegaFuel[celli]/omegaF[celli]
225  *exp
226  (
227  -sqr(ft - ftStoich)
228  /(2*sqr(0.01*omegaF[celli]))
229  )
230  *pow(ft, a - 1)
231  *pow(1 - ft, b - 1)
232  *deltaFt;
233  }
234  omegaFuelBar[celli] /= max(pc[celli], 1e-4);
235  }
236  else
237  {
238  omegaFuelBar[celli] =
239  omegaFuel[celli]/omegaF[celli]
240  *exp(-sqr(ftCell - ftStoich)/(2*sqr(0.01*omegaF[celli])));
241  }
242  }
243  else
244  {
245  omegaFuelBar[celli] = 0;
246  }
247  }
248 
249 
250  // Combustion progress variable, c
251 
252  List<label> productsIndex(2, label(-1));
253  {
254  label i = 0;
255  forAll(this->specieProd(), specieI)
256  {
257  if (this->specieProd()[specieI] < 0)
258  {
259  productsIndex[i] = specieI;
260  i++;
261  }
262  }
263  }
264 
265 
266  // Flamelet probability of the progress c based on IFC (reuse pc)
267  scalar YprodTotal = 0;
268  forAll(productsIndex, j)
269  {
270  YprodTotal += this->Yprod0()[productsIndex[j]];
271  }
272 
273  forAll(ft_, celli)
274  {
275  if (ft_[celli] < ftStoich)
276  {
277  pc[celli] = ft_[celli]*(YprodTotal/ftStoich);
278  }
279  else
280  {
281  pc[celli] = (1 - ft_[celli])*(YprodTotal/(1 - ftStoich));
282  }
283  }
284 
285  tmp<volScalarField> tproducts
286  (
288  (
289  this->thermo().phasePropertyName("products"),
290  U.mesh(),
292  )
293  );
294 
295  volScalarField& products = tproducts.ref();
296 
297  forAll(productsIndex, j)
298  {
299  const label specieI = productsIndex[j];
300  const volScalarField& Yp = this->thermo().Y()[specieI];
301  products += Yp;
302  }
303 
305  (
306  max(scalar(1) - products/max(pc, scalar(1e-5)), scalar(0))
307  );
308 
309  pc = min(C_*c, scalar(1));
310 
311  const volScalarField fres(this->fres(fuelI));
312 
313  this->wFuel_ == mgft*pc*omegaFuelBar;
314 }
315 
316 
318 {
319  this->wFuel_ ==
321 
322  calculateSourceNorm();
323 }
324 
325 
327 {
329  {
330  this->coeffs().lookup("Cv") >> Cv_ ;
331  this->coeffs().lookup("ftVarMin") >> ftVarMin_;
332  reactionRateFlameArea_->read(this->coeffs());
333  return true;
334  }
335  else
336  {
337  return false;
338  }
339 }
340 
341 
342 // ************************************************************************* //
scalar delta
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
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,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:309
Base class for combustion models.
Flame Surface Density (FDS) combustion model.
Definition: FSD.H:81
FSD(const word &modelType, const fluidMulticomponentThermo &thermo, const compressibleMomentumTransportModel &turb, const word &combustionProperties)
Construct from components.
Definition: FSD.C:49
virtual void correct()
Correct combustion rate.
Definition: FSD.C:317
virtual ~FSD()
Destructor.
Definition: FSD.C:98
virtual bool read()
Update properties.
Definition: FSD.C:326
Base class for single-step combustion models.
virtual bool read()
Update properties from given dictionary.
Base class for single-phase compressible turbulence models.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Base-class for multi-component fluid thermodynamic properties.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
Abstract class for reaction rate per flame area unit.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
Calculate the divergence of the given field.
Calculate the gradient of the given field.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
volScalarField & b
Definition: createFields.H:27
U
Definition: pEqn.H:72
defineTypeNameAndDebug(diffusion, 0)
addToRunTimeSelectionTable(combustionModel, diffusion, dictionary)
LESModel< momentumTransportModel > LESModel
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
const dimensionedScalar c
Speed of light in a vacuum.
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
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:106
dimensionedScalar exp(const dimensionedScalar &ds)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const dimensionSet dimless
const dimensionSet dimLength
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void pow3(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimMass
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fluidMulticomponentThermo & thermo
Definition: createFields.H:31