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