heheuPsiThermo.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 "heheuPsiThermo.H"
27 #include "fvMesh.H"
29 
30 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
31 
32 template<class BasicPsiThermo, class MixtureType>
34 {
35  const scalarField& hCells = this->he_.internalField();
36  const scalarField& heuCells = this->heu_.internalField();
37  const scalarField& pCells = this->p_.internalField();
38 
39  scalarField& TCells = this->T_.internalField();
40  scalarField& TuCells = this->Tu_.internalField();
41  scalarField& psiCells = this->psi_.internalField();
42  scalarField& muCells = this->mu_.internalField();
43  scalarField& alphaCells = this->alpha_.internalField();
44 
45  forAll(TCells, celli)
46  {
47  const typename MixtureType::thermoType& mixture_ =
48  this->cellMixture(celli);
49 
50  TCells[celli] = mixture_.THE
51  (
52  hCells[celli],
53  pCells[celli],
54  TCells[celli]
55  );
56 
57  psiCells[celli] = mixture_.psi(pCells[celli], TCells[celli]);
58 
59  muCells[celli] = mixture_.mu(pCells[celli], TCells[celli]);
60  alphaCells[celli] = mixture_.alphah(pCells[celli], TCells[celli]);
61 
62  TuCells[celli] = this->cellReactants(celli).THE
63  (
64  heuCells[celli],
65  pCells[celli],
66  TuCells[celli]
67  );
68  }
69 
70  forAll(this->T_.boundaryField(), patchi)
71  {
72  fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
73  fvPatchScalarField& pT = this->T_.boundaryField()[patchi];
74  fvPatchScalarField& pTu = this->Tu_.boundaryField()[patchi];
75  fvPatchScalarField& ppsi = this->psi_.boundaryField()[patchi];
76 
77  fvPatchScalarField& ph = this->he_.boundaryField()[patchi];
78  fvPatchScalarField& pheu = this->heu_.boundaryField()[patchi];
79 
80  fvPatchScalarField& pmu_ = this->mu_.boundaryField()[patchi];
81  fvPatchScalarField& palpha_ = this->alpha_.boundaryField()[patchi];
82 
83  if (pT.fixesValue())
84  {
85  forAll(pT, facei)
86  {
87  const typename MixtureType::thermoType& mixture_ =
88  this->patchFaceMixture(patchi, facei);
89 
90  ph[facei] = mixture_.HE(pp[facei], pT[facei]);
91 
92  ppsi[facei] = mixture_.psi(pp[facei], pT[facei]);
93  pmu_[facei] = mixture_.mu(pp[facei], pT[facei]);
94  palpha_[facei] = mixture_.alphah(pp[facei], pT[facei]);
95  }
96  }
97  else
98  {
99  forAll(pT, facei)
100  {
101  const typename MixtureType::thermoType& mixture_ =
102  this->patchFaceMixture(patchi, facei);
103 
104  pT[facei] = mixture_.THE(ph[facei], pp[facei], pT[facei]);
105 
106  ppsi[facei] = mixture_.psi(pp[facei], pT[facei]);
107  pmu_[facei] = mixture_.mu(pp[facei], pT[facei]);
108  palpha_[facei] = mixture_.alphah(pp[facei], pT[facei]);
109 
110  pTu[facei] =
111  this->patchFaceReactants(patchi, facei)
112  .THE(pheu[facei], pp[facei], pTu[facei]);
113  }
114  }
115  }
116 }
117 
118 
119 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
120 
121 template<class BasicPsiThermo, class MixtureType>
123 (
124  const fvMesh& mesh,
125  const word& phaseName
126 )
127 :
129  Tu_
130  (
131  IOobject
132  (
133  "Tu",
134  mesh.time().timeName(),
135  mesh,
136  IOobject::MUST_READ,
137  IOobject::AUTO_WRITE
138  ),
139  mesh
140  ),
141 
142  heu_
143  (
144  IOobject
145  (
146  MixtureType::thermoType::heName() + 'u',
147  mesh.time().timeName(),
148  mesh,
149  IOobject::NO_READ,
150  IOobject::NO_WRITE
151  ),
152  mesh,
153  dimensionSet(0, 2, -2, 0, 0),
154  this->heuBoundaryTypes()
155  )
156 {
157  scalarField& heuCells = this->heu_.internalField();
158  const scalarField& pCells = this->p_.internalField();
159  const scalarField& TuCells = this->Tu_.internalField();
160 
161  forAll(heuCells, celli)
162  {
163  heuCells[celli] = this->cellReactants(celli).HE
164  (
165  pCells[celli],
166  TuCells[celli]
167  );
168  }
169 
170  forAll(this->heu_.boundaryField(), patchi)
171  {
172  fvPatchScalarField& pheu = this->heu_.boundaryField()[patchi];
173  const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
174  const fvPatchScalarField& pTu = this->Tu_.boundaryField()[patchi];
175 
176  forAll(pheu, facei)
177  {
178  pheu[facei] = this->patchFaceReactants(patchi, facei).HE
179  (
180  pp[facei],
181  pTu[facei]
182  );
183  }
184  }
185 
186  this->heuBoundaryCorrection(this->heu_);
187 
188  calculate();
189  this->psi_.oldTime(); // Switch on saving old time
190 }
191 
192 
193 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
194 
195 template<class BasicPsiThermo, class MixtureType>
197 {}
198 
199 
200 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
201 
202 template<class BasicPsiThermo, class MixtureType>
204 {
205  if (debug)
206  {
207  Info<< "entering heheuPsiThermo<BasicPsiThermo, MixtureType>::correct()"
208  << endl;
209  }
210 
211  // force the saving of the old-time values
212  this->psi_.oldTime();
213 
214  calculate();
215 
216  if (debug)
217  {
218  Info<< "exiting heheuPsiThermo<BasicPsiThermo, MixtureType>::correct()"
219  << endl;
220  }
221 }
222 
223 
224 template<class BasicPsiThermo, class MixtureType>
227 (
228  const scalarField& p,
229  const scalarField& Tu,
230  const labelList& cells
231 ) const
232 {
233  tmp<scalarField> theu(new scalarField(Tu.size()));
234  scalarField& heu = theu();
235 
236  forAll(Tu, celli)
237  {
238  heu[celli] = this->cellReactants(cells[celli]).HE(p[celli], Tu[celli]);
239  }
240 
241  return theu;
242 }
243 
244 
245 template<class BasicPsiThermo, class MixtureType>
248 (
249  const scalarField& p,
250  const scalarField& Tu,
251  const label patchi
252 ) const
253 {
254  tmp<scalarField> theu(new scalarField(Tu.size()));
255  scalarField& heu = theu();
256 
257  forAll(Tu, facei)
258  {
259  heu[facei] =
260  this->patchFaceReactants(patchi, facei).HE(p[facei], Tu[facei]);
261  }
262 
263  return theu;
264 }
265 
266 
267 template<class BasicPsiThermo, class MixtureType>
270 {
272  (
273  new volScalarField
274  (
275  IOobject
276  (
277  "Tb",
278  this->T_.time().timeName(),
279  this->T_.db(),
280  IOobject::NO_READ,
281  IOobject::NO_WRITE,
282  false
283  ),
284  this->T_
285  )
286  );
287 
288  volScalarField& Tb_ = tTb();
289  scalarField& TbCells = Tb_.internalField();
290  const scalarField& pCells = this->p_.internalField();
291  const scalarField& TCells = this->T_.internalField();
292  const scalarField& hCells = this->he_.internalField();
293 
294  forAll(TbCells, celli)
295  {
296  TbCells[celli] = this->cellProducts(celli).THE
297  (
298  hCells[celli],
299  pCells[celli],
300  TCells[celli]
301  );
302  }
303 
304  forAll(Tb_.boundaryField(), patchi)
305  {
307 
308  const fvPatchScalarField& ph = this->he_.boundaryField()[patchi];
309  const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
310  const fvPatchScalarField& pT = this->T_.boundaryField()[patchi];
311 
312  forAll(pTb, facei)
313  {
314  pTb[facei] =
315  this->patchFaceProducts(patchi, facei)
316  .THE(ph[facei], pp[facei], pT[facei]);
317  }
318  }
319 
320  return tTb;
321 }
322 
323 
324 template<class BasicPsiThermo, class MixtureType>
327 {
328  tmp<volScalarField> tpsiu
329  (
330  new volScalarField
331  (
332  IOobject
333  (
334  "psiu",
335  this->psi_.time().timeName(),
336  this->psi_.db(),
337  IOobject::NO_READ,
338  IOobject::NO_WRITE,
339  false
340  ),
341  this->psi_.mesh(),
342  this->psi_.dimensions()
343  )
344  );
345 
346  volScalarField& psiu = tpsiu();
347  scalarField& psiuCells = psiu.internalField();
348  const scalarField& TuCells = this->Tu_.internalField();
349  const scalarField& pCells = this->p_.internalField();
350 
351  forAll(psiuCells, celli)
352  {
353  psiuCells[celli] =
354  this->cellReactants(celli).psi(pCells[celli], TuCells[celli]);
355  }
356 
357  forAll(psiu.boundaryField(), patchi)
358  {
359  fvPatchScalarField& ppsiu = psiu.boundaryField()[patchi];
360 
361  const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
362  const fvPatchScalarField& pTu = this->Tu_.boundaryField()[patchi];
363 
364  forAll(ppsiu, facei)
365  {
366  ppsiu[facei] =
367  this->
368  patchFaceReactants(patchi, facei).psi(pp[facei], pTu[facei]);
369  }
370  }
371 
372  return tpsiu;
373 }
374 
375 
376 template<class BasicPsiThermo, class MixtureType>
379 {
380  tmp<volScalarField> tpsib
381  (
382  new volScalarField
383  (
384  IOobject
385  (
386  "psib",
387  this->psi_.time().timeName(),
388  this->psi_.db(),
389  IOobject::NO_READ,
390  IOobject::NO_WRITE,
391  false
392  ),
393  this->psi_.mesh(),
394  this->psi_.dimensions()
395  )
396  );
397 
398  volScalarField& psib = tpsib();
399  scalarField& psibCells = psib.internalField();
400  const volScalarField Tb_(Tb());
401  const scalarField& TbCells = Tb_.internalField();
402  const scalarField& pCells = this->p_.internalField();
403 
404  forAll(psibCells, celli)
405  {
406  psibCells[celli] =
407  this->cellReactants(celli).psi(pCells[celli], TbCells[celli]);
408  }
409 
410  forAll(psib.boundaryField(), patchi)
411  {
412  fvPatchScalarField& ppsib = psib.boundaryField()[patchi];
413 
414  const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
415  const fvPatchScalarField& pTb = Tb_.boundaryField()[patchi];
416 
417  forAll(ppsib, facei)
418  {
419  ppsib[facei] =
420  this->patchFaceReactants
421  (patchi, facei).psi(pp[facei], pTb[facei]);
422  }
423  }
424 
425  return tpsib;
426 }
427 
428 
429 template<class BasicPsiThermo, class MixtureType>
432 {
434  (
435  new volScalarField
436  (
437  IOobject
438  (
439  "muu",
440  this->T_.time().timeName(),
441  this->T_.db(),
442  IOobject::NO_READ,
443  IOobject::NO_WRITE,
444  false
445  ),
446  this->T_.mesh(),
447  dimensionSet(1, -1, -1, 0, 0)
448  )
449  );
450 
451  volScalarField& muu_ = tmuu();
452  scalarField& muuCells = muu_.internalField();
453  const scalarField& pCells = this->p_.internalField();
454  const scalarField& TuCells = this->Tu_.internalField();
455 
456  forAll(muuCells, celli)
457  {
458  muuCells[celli] = this->cellReactants(celli).mu
459  (
460  pCells[celli],
461  TuCells[celli]
462  );
463  }
464 
465  forAll(muu_.boundaryField(), patchi)
466  {
467  fvPatchScalarField& pMuu = muu_.boundaryField()[patchi];
468  const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
469  const fvPatchScalarField& pTu = this->Tu_.boundaryField()[patchi];
470 
471  forAll(pMuu, facei)
472  {
473  pMuu[facei] = this->patchFaceReactants(patchi, facei).mu
474  (
475  pp[facei],
476  pTu[facei]
477  );
478  }
479  }
480 
481  return tmuu;
482 }
483 
484 
485 template<class BasicPsiThermo, class MixtureType>
488 {
490  (
491  new volScalarField
492  (
493  IOobject
494  (
495  "mub",
496  this->T_.time().timeName(),
497  this->T_.db(),
498  IOobject::NO_READ,
499  IOobject::NO_WRITE,
500  false
501  ),
502  this->T_.mesh(),
503  dimensionSet(1, -1, -1, 0, 0)
504  )
505  );
506 
507  volScalarField& mub_ = tmub();
508  scalarField& mubCells = mub_.internalField();
509  const volScalarField Tb_(Tb());
510  const scalarField& pCells = this->p_.internalField();
511  const scalarField& TbCells = Tb_.internalField();
512 
513  forAll(mubCells, celli)
514  {
515  mubCells[celli] = this->cellProducts(celli).mu
516  (
517  pCells[celli],
518  TbCells[celli]
519  );
520  }
521 
522  forAll(mub_.boundaryField(), patchi)
523  {
524  fvPatchScalarField& pMub = mub_.boundaryField()[patchi];
525  const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
526  const fvPatchScalarField& pTb = Tb_.boundaryField()[patchi];
527 
528  forAll(pMub, facei)
529  {
530  pMub[facei] = this->patchFaceProducts(patchi, facei).mu
531  (
532  pp[facei],
533  pTb[facei]
534  );
535  }
536  }
537 
538  return tmub;
539 }
540 
541 
542 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::heThermo
Enthalpy/Internal energy for a mixture.
Definition: heThermo.H:49
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::heheuPsiThermo::correct
virtual void correct()
Update properties.
Definition: heheuPsiThermo.C:203
Foam::dimensionSet
Dimension set for the base types.
Definition: dimensionSet.H:116
Foam::fvPatchField::fixesValue
virtual bool fixesValue() const
Return true if this patch field fixes a value.
Definition: fvPatchField.H:336
Foam::heheuPsiThermo::mub
virtual tmp< volScalarField > mub() const
Dynamic viscosity of burnt gas [kg/ms].
Definition: heheuPsiThermo.C:487
Foam::heheuPsiThermo::psiu
virtual tmp< volScalarField > psiu() const
Unburnt gas compressibility [s^2/m^2].
Definition: heheuPsiThermo.C:326
Foam::label
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
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::Info
messageStream Info
Foam::heheuPsiThermo::Tb
virtual tmp< volScalarField > Tb() const
Burnt gas temperature [K].
Definition: heheuPsiThermo.C:269
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
heheuPsiThermo.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fvMesh.H
Foam::heheuPsiThermo::heu
virtual volScalarField & heu()
Unburnt gas enthalpy [J/kg].
Definition: heheuPsiThermo.H:99
Foam::heheuPsiThermo::muu
virtual tmp< volScalarField > muu() const
Dynamic viscosity of unburnt gas [kg/ms].
Definition: heheuPsiThermo.C:431
Foam::heheuPsiThermo::calculate
void calculate()
Definition: heheuPsiThermo.C:33
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
scalarField
volScalarField scalarField(fieldObject, mesh)
fixedValueFvPatchFields.H
patchi
label patchi
Definition: getPatchFieldScalar.H:1
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::heheuPsiThermo::~heheuPsiThermo
virtual ~heheuPsiThermo()
Destructor.
Definition: heheuPsiThermo.C:196
Foam::heheuPsiThermo::heheuPsiThermo
heheuPsiThermo(const heheuPsiThermo< BasicPsiThermo, MixtureType > &)
Construct as copy (not implemented)
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::heheuPsiThermo::psib
virtual tmp< volScalarField > psib() const
Burnt gas compressibility [s^2/m^2].
Definition: heheuPsiThermo.C:378