FSD.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 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"
31 
32 namespace Foam
33 {
34 namespace combustionModels
35 {
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
39 template<class CombThermoType, class ThermoType>
40 FSD<CombThermoType, ThermoType>::FSD
41 (
42  const word& modelType,
43  const fvMesh& mesh,
44  const word& combustionProperties,
45  const word& phaseName
46 )
47 :
49  (
50  modelType,
51  mesh,
52  combustionProperties,
53  phaseName
54  ),
55  reactionRateFlameArea_
56  (
58  (
59  this->coeffs(),
60  this->mesh(),
61  *this
62  )
63  ),
64  ft_
65  (
66  IOobject
67  (
68  IOobject::groupName("ft", phaseName),
69  this->mesh().time().timeName(),
70  this->mesh(),
73  ),
74  this->mesh(),
75  dimensionedScalar("zero", dimless, 0.0)
76  ),
77  YFuelFuelStream_(dimensionedScalar("YFuelStream", dimless, 1.0)),
78  YO2OxiStream_(dimensionedScalar("YOxiStream", dimless, 0.23)),
79  Cv_(readScalar(this->coeffs().lookup("Cv"))),
80  C_(5.0),
81  ftMin_(0.0),
82  ftMax_(1.0),
83  ftDim_(300),
84  ftVarMin_(readScalar(this->coeffs().lookup("ftVarMin")))
85 {}
86 
87 
88 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
89 
90 template<class CombThermoType, class ThermoType>
92 {}
93 
94 
95 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
96 
97 template<class CombThermoType, class ThermoType>
99 {
100  this->singleMixturePtr_->fresCorrect();
101 
102  const label fuelI = this->singleMixturePtr_->fuelIndex();
103 
104  const volScalarField& YFuel = this->thermoPtr_->composition().Y()[fuelI];
105 
106  const volScalarField& YO2 = this->thermoPtr_->composition().Y("O2");
107 
108  const dimensionedScalar s = this->singleMixturePtr_->s();
109 
110  ft_ =
111  (s*YFuel - (YO2 - YO2OxiStream_))/(s*YFuelFuelStream_ + YO2OxiStream_);
112 
113 
114  volVectorField nft(fvc::grad(ft_));
115 
116  volScalarField mgft(mag(nft));
117 
118  surfaceVectorField SfHat(this->mesh().Sf()/this->mesh().magSf());
119 
120  volScalarField cAux(scalar(1) - ft_);
121 
122  dimensionedScalar dMgft = 1.0e-3*
123  (ft_*cAux*mgft)().weightedAverage(this->mesh().V())
124  /((ft_*cAux)().weightedAverage(this->mesh().V()) + SMALL)
125  + dimensionedScalar("ddMgft", mgft.dimensions(), SMALL);
126 
127  mgft += dMgft;
128 
129  nft /= mgft;
130 
131  const volVectorField& U = YO2.db().lookupObject<volVectorField>("U");
132 
133  const volScalarField sigma
134  (
135  (nft & nft)*fvc::div(U) - (nft & fvc::grad(U) & nft)
136  );
137 
138  reactionRateFlameArea_->correct(sigma);
139 
140  const volScalarField& omegaFuel = reactionRateFlameArea_->omega();
141 
142 
143  const scalar ftStoich =
144  YO2OxiStream_.value()
145  /(
146  s.value()*YFuelFuelStream_.value() + YO2OxiStream_.value()
147  );
148 
150  (
151  new volScalarField
152  (
153  IOobject
154  (
155  IOobject::groupName("Pc", this->phaseName_),
156  U.time().timeName(),
157  U.db(),
160  ),
161  U.mesh(),
162  dimensionedScalar("Pc", dimless, 0)
163  )
164  );
165 
166  volScalarField& pc = tPc.ref();
167 
168  tmp<volScalarField> tomegaFuel
169  (
170  new volScalarField
171  (
172  IOobject
173  (
174  IOobject::groupName("omegaFuelBar", this->phaseName_),
175  U.time().timeName(),
176  U.db(),
179  ),
180  U.mesh(),
182  (
183  "omegaFuelBar",
184  omegaFuel.dimensions(),
185  0
186  )
187  )
188  );
189 
190  volScalarField& omegaFuelBar = tomegaFuel.ref();
191 
192  // Calculation of the mixture fraction variance (ftVar)
193  const compressible::LESModel& lesModel =
195  (
197  );
198 
199  const volScalarField& delta = lesModel.delta();
200  const volScalarField ftVar(Cv_*sqr(delta)*sqr(mgft));
201 
202  // Thickened flame (average flame thickness for counterflow configuration
203  // is 1.5 mm)
204 
205  volScalarField deltaF
206  (
207  lesModel.delta()/dimensionedScalar("flame", dimLength, 1.5e-3)
208  );
209 
210  // Linear correlation between delta and flame thickness
211  volScalarField omegaF(max(deltaF*(4.0/3.0) + (2.0/3.0), scalar(1)));
212 
213  scalar deltaFt = 1.0/ftDim_;
214 
215  forAll(ft_, celli)
216  {
217  if (ft_[celli] > ftMin_ && ft_[celli] < ftMax_)
218  {
219  scalar ftCell = ft_[celli];
220 
221  if (ftVar[celli] > ftVarMin_) //sub-grid beta pdf of ft_
222  {
223  scalar ftVarc = ftVar[celli];
224  scalar a =
225  max(ftCell*(ftCell*(1.0 - ftCell)/ftVarc - 1.0), 0.0);
226  scalar b = max(a/ftCell - a, 0.0);
227 
228  for (int i=1; i<ftDim_; i++)
229  {
230  scalar ft = i*deltaFt;
231  pc[celli] += pow(ft, a-1.0)*pow(1.0 - ft, b - 1.0)*deltaFt;
232  }
233 
234  for (int i=1; i<ftDim_; i++)
235  {
236  scalar ft = i*deltaFt;
237  omegaFuelBar[celli] +=
238  omegaFuel[celli]/omegaF[celli]
239  *exp
240  (
241  -sqr(ft - ftStoich)
242  /(2.0*sqr(0.01*omegaF[celli]))
243  )
244  *pow(ft, a - 1.0)
245  *pow(1.0 - ft, b - 1.0)
246  *deltaFt;
247  }
248  omegaFuelBar[celli] /= max(pc[celli], 1e-4);
249  }
250  else
251  {
252  omegaFuelBar[celli] =
253  omegaFuel[celli]/omegaF[celli]
254  *exp(-sqr(ftCell - ftStoich)/(2.0*sqr(0.01*omegaF[celli])));
255  }
256  }
257  else
258  {
259  omegaFuelBar[celli] = 0.0;
260  }
261  }
262 
263 
264  // Combustion progress variable, c
265 
266  List<label> productsIndex(2, label(-1));
267  {
268  label i = 0;
269  forAll(this->singleMixturePtr_->specieProd(), specieI)
270  {
271  if (this->singleMixturePtr_->specieProd()[specieI] < 0)
272  {
273  productsIndex[i] = specieI;
274  i++;
275  }
276  }
277  }
278 
279 
280  // Flamelet probability of the progress c based on IFC (reuse pc)
281  scalar YprodTotal = 0;
282  forAll(productsIndex, j)
283  {
284  YprodTotal += this->singleMixturePtr_->Yprod0()[productsIndex[j]];
285  }
286 
287  forAll(ft_, celli)
288  {
289  if (ft_[celli] < ftStoich)
290  {
291  pc[celli] = ft_[celli]*(YprodTotal/ftStoich);
292  }
293  else
294  {
295  pc[celli] = (1.0 - ft_[celli])*(YprodTotal/(1.0 - ftStoich));
296  }
297  }
298 
299  tmp<volScalarField> tproducts
300  (
301  new volScalarField
302  (
303  IOobject
304  (
305  IOobject::groupName("products", this->phaseName_),
306  U.time().timeName(),
307  U.db(),
310  ),
311  U.mesh(),
312  dimensionedScalar("products", dimless, 0)
313  )
314  );
315 
316  volScalarField& products = tproducts.ref();
317 
318  forAll(productsIndex, j)
319  {
320  label specieI = productsIndex[j];
321  const volScalarField& Yp = this->thermoPtr_->composition().Y()[specieI];
322  products += Yp;
323  }
324 
326  (
327  max(scalar(1) - products/max(pc, scalar(1e-5)), scalar(0))
328  );
329 
330  pc = min(C_*c, scalar(1));
331 
332  const volScalarField fres(this->singleMixturePtr_->fres(fuelI));
333 
334  this->wFuel_ == mgft*pc*omegaFuelBar;
335 }
336 
337 
338 template<class CombThermoType, class ThermoType>
340 {
341  this->wFuel_ ==
343 
344  if (this->active())
345  {
346  calculateSourceNorm();
347  }
348 }
349 
350 
351 template<class CombThermoType, class ThermoType>
353 {
355  {
356  this->coeffs().lookup("Cv") >> Cv_ ;
357  this->coeffs().lookup("ftVarMin") >> ftVarMin_;
358  reactionRateFlameArea_->read(this->coeffs());
359  return true;
360  }
361  else
362  {
363  return false;
364  }
365 }
366 
367 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
368 
369 } // End namespace combustionModels
370 } // End namespace Foam
371 
372 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
scalar delta
Flame Surface Dennsity (FDS) combustion model.
Definition: FSD.H:80
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
U
Definition: pEqn.H:83
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
virtual bool read()
Update properties.
Definition: FSD.C:352
const double e
Elementary charge.
Definition: doubleFloat.H:78
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:644
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Macros for easy insertion into run-time selection tables.
const dimensionSet & dimensions() const
Return dimensions.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
dimensionedScalar exp(const dimensionedScalar &ds)
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Calculate the gradient of the given field.
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
LESModel< EddyDiffusivity< turbulenceModel > > LESModel
Typedefs for turbulence, RAS and LES models for compressible flow based on the standard laminar trans...
const Type & value() const
Return const reference to value.
word timeName
Definition: getTimeIndex.H:3
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
Calculate the divergence of the given field.
const Mesh & mesh() const
Return mesh.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Internal & ref()
Return a reference to the dimensioned internal field.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const Time & time() const
Return time.
Definition: IOobject.C:337
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const dimensionedScalar c
Speed of light in a vacuum.
virtual void correct()
Correct combustion rate.
Definition: FSD.C:339
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
dimensioned< scalar > mag(const dimensioned< Type > &)
Base class for combustion models using singleStepReactingMixture.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: PtrList.H:53
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:331
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
static autoPtr< reactionRateFlameArea > New(const dictionary &dict, const fvMesh &mesh, const combustionModel &combModel)
virtual ~FSD()
Destructor.
Definition: FSD.C:91
Namespace for OpenFOAM.