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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2017-2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "heheuPsiThermo.H"
30 #include "fvMesh.H"
32 
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34 
35 template<class BasicPsiThermo, class MixtureType>
37 {
38  const scalarField& hCells = this->he_;
39  const scalarField& heuCells = this->heu_;
40  const scalarField& pCells = this->p_;
41 
42  scalarField& TCells = this->T_.primitiveFieldRef();
43  scalarField& TuCells = this->Tu_.primitiveFieldRef();
44  scalarField& psiCells = this->psi_.primitiveFieldRef();
45  scalarField& muCells = this->mu_.primitiveFieldRef();
46  scalarField& alphaCells = this->alpha_.primitiveFieldRef();
47 
48  forAll(TCells, celli)
49  {
50  const typename MixtureType::thermoType& mixture_ =
51  this->cellMixture(celli);
52 
53  if (this->updateT())
54  {
55  TCells[celli] = mixture_.THE
56  (
57  hCells[celli],
58  pCells[celli],
59  TCells[celli]
60  );
61  }
62 
63  psiCells[celli] = mixture_.psi(pCells[celli], TCells[celli]);
64 
65  muCells[celli] = mixture_.mu(pCells[celli], TCells[celli]);
66  alphaCells[celli] = mixture_.alphah(pCells[celli], TCells[celli]);
67 
68  TuCells[celli] = this->cellReactants(celli).THE
69  (
70  heuCells[celli],
71  pCells[celli],
72  TuCells[celli]
73  );
74  }
75 
76  volScalarField::Boundary& pBf =
77  this->p_.boundaryFieldRef();
78 
79  volScalarField::Boundary& TBf =
80  this->T_.boundaryFieldRef();
81 
82  volScalarField::Boundary& TuBf =
83  this->Tu_.boundaryFieldRef();
84 
85  volScalarField::Boundary& psiBf =
86  this->psi_.boundaryFieldRef();
87 
88  volScalarField::Boundary& heBf =
89  this->he().boundaryFieldRef();
90 
91  volScalarField::Boundary& heuBf =
92  this->heu().boundaryFieldRef();
93 
94  volScalarField::Boundary& muBf =
95  this->mu_.boundaryFieldRef();
96 
97  volScalarField::Boundary& alphaBf =
98  this->alpha_.boundaryFieldRef();
99 
100  forAll(this->T_.boundaryField(), patchi)
101  {
102  fvPatchScalarField& pp = pBf[patchi];
103  fvPatchScalarField& pT = TBf[patchi];
104  fvPatchScalarField& pTu = TuBf[patchi];
105  fvPatchScalarField& ppsi = psiBf[patchi];
106  fvPatchScalarField& phe = heBf[patchi];
107  fvPatchScalarField& pheu = heuBf[patchi];
108  fvPatchScalarField& pmu = muBf[patchi];
109  fvPatchScalarField& palpha = alphaBf[patchi];
110 
111  if (pT.fixesValue())
112  {
113  forAll(pT, facei)
114  {
115  const typename MixtureType::thermoType& mixture_ =
116  this->patchFaceMixture(patchi, facei);
117 
118  phe[facei] = mixture_.HE(pp[facei], pT[facei]);
119 
120  ppsi[facei] = mixture_.psi(pp[facei], pT[facei]);
121  pmu[facei] = mixture_.mu(pp[facei], pT[facei]);
122  palpha[facei] = mixture_.alphah(pp[facei], pT[facei]);
123  }
124  }
125  else
126  {
127  forAll(pT, facei)
128  {
129  const typename MixtureType::thermoType& mixture_ =
130  this->patchFaceMixture(patchi, facei);
131 
132  if (this->updateT())
133  {
134  pT[facei] = mixture_.THE(phe[facei], pp[facei], pT[facei]);
135  }
136 
137  ppsi[facei] = mixture_.psi(pp[facei], pT[facei]);
138  pmu[facei] = mixture_.mu(pp[facei], pT[facei]);
139  palpha[facei] = mixture_.alphah(pp[facei], pT[facei]);
140 
141  pTu[facei] =
142  this->patchFaceReactants(patchi, facei)
143  .THE(pheu[facei], pp[facei], pTu[facei]);
144  }
145  }
146  }
147 }
148 
149 
150 /*
151 template<class BasicPsiThermo, class MixtureType>
152 void Foam::heheuPsiThermo<BasicPsiThermo, MixtureType>::calculateT()
153 {
154  //const scalarField& hCells = this->he_.primitiveFieldRef();
155  const scalarField& heuCells = this->heu_.primitiveFieldRef();
156  const scalarField& pCells = this->p_.primitiveFieldRef();
157 
158  scalarField& TCells = this->T_.primitiveFieldRef();
159  scalarField& TuCells = this->Tu_.primitiveFieldRef();
160  scalarField& psiCells = this->psi_.primitiveFieldRef();
161  scalarField& muCells = this->mu_.primitiveFieldRef();
162  scalarField& alphaCells = this->alpha_.primitiveFieldRef();
163 
164  forAll(TCells, celli)
165  {
166  const typename MixtureType::thermoType& mixture_ =
167  this->cellMixture(celli);
168 
169 // TCells[celli] = mixture_.THE
170 // (
171 // hCells[celli],
172 // pCells[celli],
173 // TCells[celli]
174 // );
175 
176  psiCells[celli] = mixture_.psi(pCells[celli], TCells[celli]);
177 
178  muCells[celli] = mixture_.mu(pCells[celli], TCells[celli]);
179  alphaCells[celli] = mixture_.alphah(pCells[celli], TCells[celli]);
180 
181  TuCells[celli] = this->cellReactants(celli).THE
182  (
183  heuCells[celli],
184  pCells[celli],
185  TuCells[celli]
186  );
187  }
188 
189  forAll(this->T_.boundaryField(), patchi)
190  {
191  fvPatchScalarField& pp = this->p_.boundaryFieldRef()[patchi];
192  fvPatchScalarField& pT = this->T_.boundaryFieldRef()[patchi];
193  fvPatchScalarField& pTu = this->Tu_.boundaryFieldRef()[patchi];
194  fvPatchScalarField& ppsi = this->psi_.boundaryFieldRef()[patchi];
195 
196  fvPatchScalarField& ph = this->he_.boundaryFieldRef()[patchi];
197  fvPatchScalarField& pheu = this->heu_.boundaryFieldRef()[patchi];
198 
199  fvPatchScalarField& pmu_ = this->mu_.boundaryFieldRef()[patchi];
200  fvPatchScalarField& palpha_ = this->alpha_.boundaryFieldRef()[patchi];
201 
202  if (pT.fixesValue())
203  {
204  forAll(pT, facei)
205  {
206  const typename MixtureType::thermoType& mixture_ =
207  this->patchFaceMixture(patchi, facei);
208 
209  ph[facei] = mixture_.HE(pp[facei], pT[facei]);
210 
211  ppsi[facei] = mixture_.psi(pp[facei], pT[facei]);
212  pmu_[facei] = mixture_.mu(pp[facei], pT[facei]);
213  palpha_[facei] = mixture_.alphah(pp[facei], pT[facei]);
214  }
215  }
216  else
217  {
218  forAll(pT, facei)
219  {
220  const typename MixtureType::thermoType& mixture_ =
221  this->patchFaceMixture(patchi, facei);
222 
223  //pT[facei] = mixture_.THE(ph[facei], pp[facei], pT[facei]);
224 
225  ppsi[facei] = mixture_.psi(pp[facei], pT[facei]);
226  pmu_[facei] = mixture_.mu(pp[facei], pT[facei]);
227  palpha_[facei] = mixture_.alphah(pp[facei], pT[facei]);
228 
229  pTu[facei] =
230  this->patchFaceReactants(patchi, facei)
231  .THE(pheu[facei], pp[facei], pTu[facei]);
232  }
233  }
234  }
235 }
236 */
237 
238 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
239 
240 template<class BasicPsiThermo, class MixtureType>
242 (
243  const fvMesh& mesh,
244  const word& phaseName
245 )
246 :
248  Tu_
249  (
250  IOobject
251  (
252  "Tu",
253  mesh.time().timeName(),
254  mesh,
255  IOobject::MUST_READ,
256  IOobject::AUTO_WRITE
257  ),
258  mesh
259  ),
260 
261  heu_
262  (
263  IOobject
264  (
265  MixtureType::thermoType::heName() + 'u',
266  mesh.time().timeName(),
267  mesh,
268  IOobject::NO_READ,
269  IOobject::NO_WRITE
270  ),
271  mesh,
272  dimensionSet(0, 2, -2, 0, 0),
273  this->heuBoundaryTypes()
274  )
275 {
276  scalarField& heuCells = this->heu_.primitiveFieldRef();
277  const scalarField& pCells = this->p_;
278  const scalarField& TuCells = this->Tu_;
279 
280  forAll(heuCells, celli)
281  {
282  heuCells[celli] = this->cellReactants(celli).HE
283  (
284  pCells[celli],
285  TuCells[celli]
286  );
287  }
288 
289  volScalarField::Boundary& heuBf = heu_.boundaryFieldRef();
290 
291  forAll(heuBf, patchi)
292  {
293  fvPatchScalarField& pheu = heuBf[patchi];
294  const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
295  const fvPatchScalarField& pTu = this->Tu_.boundaryField()[patchi];
296 
297  forAll(pheu, facei)
298  {
299  pheu[facei] = this->patchFaceReactants(patchi, facei).HE
300  (
301  pp[facei],
302  pTu[facei]
303  );
304  }
305  }
306 
307  this->heuBoundaryCorrection(this->heu_);
308 
309  calculate();
310  this->psi_.oldTime(); // Switch on saving old time
311 }
312 
313 
314 template<class BasicPsiThermo, class MixtureType>
316 (
317  const fvMesh& mesh,
318  const word& phaseName,
319  const word& dictName
320 )
321 :
323  Tu_
324  (
325  IOobject
326  (
327  "Tu",
328  mesh.time().timeName(),
329  mesh,
330  IOobject::MUST_READ,
331  IOobject::AUTO_WRITE
332  ),
333  mesh
334  ),
335 
336  heu_
337  (
338  IOobject
339  (
340  MixtureType::thermoType::heName() + 'u',
341  mesh.time().timeName(),
342  mesh,
343  IOobject::NO_READ,
344  IOobject::NO_WRITE
345  ),
346  mesh,
347  dimensionSet(0, 2, -2, 0, 0),
348  this->heuBoundaryTypes()
349  )
350 {
351  scalarField& heuCells = this->heu_.primitiveFieldRef();
352  const scalarField& pCells = this->p_;
353  const scalarField& TuCells = this->Tu_;
354 
355  forAll(heuCells, celli)
356  {
357  heuCells[celli] = this->cellReactants(celli).HE
358  (
359  pCells[celli],
360  TuCells[celli]
361  );
362  }
363 
364  volScalarField::Boundary& heuBf = heu_.boundaryFieldRef();
365 
366  forAll(heuBf, patchi)
367  {
368  fvPatchScalarField& pheu = heuBf[patchi];
369  const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
370  const fvPatchScalarField& pTu = this->Tu_.boundaryField()[patchi];
371 
372  forAll(pheu, facei)
373  {
374  pheu[facei] = this->patchFaceReactants(patchi, facei).HE
375  (
376  pp[facei],
377  pTu[facei]
378  );
379  }
380  }
381 
382  this->heuBoundaryCorrection(this->heu_);
383 
384  calculate();
385  this->psi_.oldTime(); // Switch on saving old time
386 }
387 
388 
389 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
390 
391 template<class BasicPsiThermo, class MixtureType>
393 {}
394 
395 
396 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
397 
398 template<class BasicPsiThermo, class MixtureType>
400 {
401  if (debug)
402  {
403  InfoInFunction << endl;
404  }
405 
406  // force the saving of the old-time values
407  this->psi_.oldTime();
408 
409  calculate();
410 
411  if (debug)
412  {
413  Info<< " Finished" << endl;
414  }
415 }
416 
417 /*
418 template<class BasicPsiThermo, class MixtureType>
419 void Foam::heheuPsiThermo<BasicPsiThermo, MixtureType>::correctT()
420 {
421  // force the saving of the old-time values
422  this->psi_.oldTime();
423 
424  calculateT();
425 }
426 */
427 
428 template<class BasicPsiThermo, class MixtureType>
431 (
432  const scalarField& p,
433  const scalarField& Tu,
434  const labelList& cells
435 ) const
436 {
437  tmp<scalarField> theu(new scalarField(Tu.size()));
438  scalarField& heu = theu.ref();
439 
440  forAll(Tu, celli)
441  {
442  heu[celli] = this->cellReactants(cells[celli]).HE(p[celli], Tu[celli]);
443  }
444 
445  return theu;
446 }
447 
448 
449 template<class BasicPsiThermo, class MixtureType>
452 (
453  const scalarField& p,
454  const scalarField& Tu,
455  const label patchi
456 ) const
457 {
458  tmp<scalarField> theu(new scalarField(Tu.size()));
459  scalarField& heu = theu.ref();
460 
461  forAll(Tu, facei)
462  {
463  heu[facei] =
464  this->patchFaceReactants(patchi, facei).HE(p[facei], Tu[facei]);
465  }
466 
467  return theu;
468 }
469 
470 
471 template<class BasicPsiThermo, class MixtureType>
474 {
476  (
477  new volScalarField
478  (
479  IOobject
480  (
481  "Tb",
482  this->T_.time().timeName(),
483  this->T_.db(),
484  IOobject::NO_READ,
485  IOobject::NO_WRITE,
486  false
487  ),
488  this->T_
489  )
490  );
491 
492  volScalarField& Tb_ = tTb.ref();
493  scalarField& TbCells = Tb_.primitiveFieldRef();
494  const scalarField& pCells = this->p_;
495  const scalarField& TCells = this->T_;
496  const scalarField& hCells = this->he_;
497 
498  forAll(TbCells, celli)
499  {
500  TbCells[celli] = this->cellProducts(celli).THE
501  (
502  hCells[celli],
503  pCells[celli],
504  TCells[celli]
505  );
506  }
507 
508  volScalarField::Boundary& TbBf = Tb_.boundaryFieldRef();
509 
510  forAll(TbBf, patchi)
511  {
512  fvPatchScalarField& pTb = TbBf[patchi];
513 
514  const fvPatchScalarField& ph = this->he_.boundaryField()[patchi];
515  const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
516  const fvPatchScalarField& pT = this->T_.boundaryField()[patchi];
517 
518  forAll(pTb, facei)
519  {
520  pTb[facei] =
521  this->patchFaceProducts(patchi, facei)
522  .THE(ph[facei], pp[facei], pT[facei]);
523  }
524  }
525 
526  return tTb;
527 }
528 
529 
530 template<class BasicPsiThermo, class MixtureType>
533 {
534  tmp<volScalarField> tpsiu
535  (
536  new volScalarField
537  (
538  IOobject
539  (
540  "psiu",
541  this->psi_.time().timeName(),
542  this->psi_.db(),
543  IOobject::NO_READ,
544  IOobject::NO_WRITE,
545  false
546  ),
547  this->psi_.mesh(),
548  this->psi_.dimensions()
549  )
550  );
551 
552  volScalarField& psiu = tpsiu.ref();
553  scalarField& psiuCells = psiu.primitiveFieldRef();
554  const scalarField& TuCells = this->Tu_;
555  const scalarField& pCells = this->p_;
556 
557  forAll(psiuCells, celli)
558  {
559  psiuCells[celli] =
560  this->cellReactants(celli).psi(pCells[celli], TuCells[celli]);
561  }
562 
563  volScalarField::Boundary& psiuBf = psiu.boundaryFieldRef();
564 
565  forAll(psiuBf, patchi)
566  {
567  fvPatchScalarField& ppsiu = psiuBf[patchi];
568 
569  const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
570  const fvPatchScalarField& pTu = this->Tu_.boundaryField()[patchi];
571 
572  forAll(ppsiu, facei)
573  {
574  ppsiu[facei] =
575  this->
576  patchFaceReactants(patchi, facei).psi(pp[facei], pTu[facei]);
577  }
578  }
579 
580  return tpsiu;
581 }
582 
583 
584 template<class BasicPsiThermo, class MixtureType>
587 {
588  tmp<volScalarField> tpsib
589  (
590  new volScalarField
591  (
592  IOobject
593  (
594  "psib",
595  this->psi_.time().timeName(),
596  this->psi_.db(),
597  IOobject::NO_READ,
598  IOobject::NO_WRITE,
599  false
600  ),
601  this->psi_.mesh(),
602  this->psi_.dimensions()
603  )
604  );
605 
606  volScalarField& psib = tpsib.ref();
607  scalarField& psibCells = psib.primitiveFieldRef();
608  const volScalarField Tb_(Tb());
609  const scalarField& TbCells = Tb_;
610  const scalarField& pCells = this->p_;
611 
612  forAll(psibCells, celli)
613  {
614  psibCells[celli] =
615  this->cellProducts(celli).psi(pCells[celli], TbCells[celli]);
616  }
617 
618  volScalarField::Boundary& psibBf = psib.boundaryFieldRef();
619 
620  forAll(psibBf, patchi)
621  {
622  fvPatchScalarField& ppsib = psibBf[patchi];
623 
624  const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
625  const fvPatchScalarField& pTb = Tb_.boundaryField()[patchi];
626 
627  forAll(ppsib, facei)
628  {
629  ppsib[facei] =
630  this->patchFaceProducts
631  (patchi, facei).psi(pp[facei], pTb[facei]);
632  }
633  }
634 
635  return tpsib;
636 }
637 
638 
639 template<class BasicPsiThermo, class MixtureType>
642 {
644  (
645  new volScalarField
646  (
647  IOobject
648  (
649  "muu",
650  this->T_.time().timeName(),
651  this->T_.db(),
652  IOobject::NO_READ,
653  IOobject::NO_WRITE,
654  false
655  ),
656  this->T_.mesh(),
657  dimensionSet(1, -1, -1, 0, 0)
658  )
659  );
660 
661  volScalarField& muu_ = tmuu.ref();
662  scalarField& muuCells = muu_.primitiveFieldRef();
663  const scalarField& pCells = this->p_;
664  const scalarField& TuCells = this->Tu_;
665 
666  forAll(muuCells, celli)
667  {
668  muuCells[celli] = this->cellReactants(celli).mu
669  (
670  pCells[celli],
671  TuCells[celli]
672  );
673  }
674 
675  volScalarField::Boundary& muuBf = muu_.boundaryFieldRef();
676 
677  forAll(muuBf, patchi)
678  {
679  fvPatchScalarField& pMuu = muuBf[patchi];
680  const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
681  const fvPatchScalarField& pTu = this->Tu_.boundaryField()[patchi];
682 
683  forAll(pMuu, facei)
684  {
685  pMuu[facei] = this->patchFaceReactants(patchi, facei).mu
686  (
687  pp[facei],
688  pTu[facei]
689  );
690  }
691  }
692 
693  return tmuu;
694 }
695 
696 
697 template<class BasicPsiThermo, class MixtureType>
700 {
702  (
703  new volScalarField
704  (
705  IOobject
706  (
707  "mub",
708  this->T_.time().timeName(),
709  this->T_.db(),
710  IOobject::NO_READ,
711  IOobject::NO_WRITE,
712  false
713  ),
714  this->T_.mesh(),
715  dimensionSet(1, -1, -1, 0, 0)
716  )
717  );
718 
719  volScalarField& mub_ = tmub.ref();
720  scalarField& mubCells = mub_.primitiveFieldRef();
721  const volScalarField Tb_(Tb());
722  const scalarField& pCells = this->p_;
723  const scalarField& TbCells = Tb_;
724 
725  forAll(mubCells, celli)
726  {
727  mubCells[celli] = this->cellProducts(celli).mu
728  (
729  pCells[celli],
730  TbCells[celli]
731  );
732  }
733 
734  volScalarField::Boundary& mubBf = mub_.boundaryFieldRef();
735 
736  forAll(mubBf, patchi)
737  {
738  fvPatchScalarField& pMub = mubBf[patchi];
739  const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
740  const fvPatchScalarField& pTb = Tb_.boundaryField()[patchi];
741 
742  forAll(pMub, facei)
743  {
744  pMub[facei] = this->patchFaceProducts(patchi, facei).mu
745  (
746  pp[facei],
747  pTb[facei]
748  );
749  }
750  }
751 
752  return tmub;
753 }
754 
755 
756 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::fvPatchField< scalar >
Foam::fvPatchScalarField
fvPatchField< scalar > fvPatchScalarField
Definition: fvPatchFieldsFwd.H:40
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
InfoInFunction
#define InfoInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:316
Foam::heThermo
Enthalpy/Internal energy for a mixture.
Definition: heThermo.H:52
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
dictName
const word dictName("blockMeshDict")
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::heheuPsiThermo::correct
virtual void correct()
Update properties.
Definition: heheuPsiThermo.C:399
Foam::dimensionSet
Dimension set for the base types.
Definition: dimensionSet.H:65
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::heheuPsiThermo
Definition: heheuPsiThermo.H:52
Foam::heheuPsiThermo::mub
virtual tmp< volScalarField > mub() const
Dynamic viscosity of burnt gas [kg/ms].
Definition: heheuPsiThermo.C:699
Foam::heheuPsiThermo::psiu
virtual tmp< volScalarField > psiu() const
Unburnt gas compressibility [s^2/m^2].
Definition: heheuPsiThermo.C:532
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
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:62
Foam::Field< scalar >
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::heheuPsiThermo::Tb
virtual tmp< volScalarField > Tb() const
Burnt gas temperature [K].
Definition: heheuPsiThermo.C:473
heheuPsiThermo.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
fvMesh.H
Foam::heheuPsiThermo::heu
virtual volScalarField & heu()
Update properties based on T.
Definition: heheuPsiThermo.H:119
Foam::GeometricField::primitiveFieldRef
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Definition: GeometricField.C:735
Foam::heheuPsiThermo::muu
virtual tmp< volScalarField > muu() const
Dynamic viscosity of unburnt gas [kg/ms].
Definition: heheuPsiThermo.C:641
he
volScalarField & he
Definition: YEEqn.H:52
Foam::GeometricField::ref
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
Definition: GeometricField.C:718
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:752
Foam::List< label >
fixedValueFvPatchFields.H
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::heheuPsiThermo::~heheuPsiThermo
virtual ~heheuPsiThermo()
Destructor.
Definition: heheuPsiThermo.C:392
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::heheuPsiThermo::psib
virtual tmp< volScalarField > psib() const
Burnt gas compressibility [s^2/m^2].
Definition: heheuPsiThermo.C:586