heThermo.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-2018 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 "heThermo.H"
29 
30 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
31 
32 template<class BasicThermo, class MixtureType>
35 {
36  volScalarField::Boundary& hBf = h.boundaryFieldRef();
37 
38  forAll(hBf, patchi)
39  {
40  if (isA<gradientEnergyFvPatchScalarField>(hBf[patchi]))
41  {
42  refCast<gradientEnergyFvPatchScalarField>(hBf[patchi]).gradient()
43  = hBf[patchi].fvPatchField::snGrad();
44  }
45  else if (isA<mixedEnergyFvPatchScalarField>(hBf[patchi]))
46  {
47  refCast<mixedEnergyFvPatchScalarField>(hBf[patchi]).refGrad()
48  = hBf[patchi].fvPatchField::snGrad();
49  }
50  }
51 }
52 
53 
54 template<class BasicThermo, class MixtureType>
56 {
57  scalarField& heCells = he_.primitiveFieldRef();
58  const scalarField& pCells = this->p_;
59  const scalarField& TCells = this->T_;
60 
61  forAll(heCells, celli)
62  {
63  heCells[celli] =
64  this->cellMixture(celli).HE(pCells[celli], TCells[celli]);
65  }
66 
67  volScalarField::Boundary& heBf = he_.boundaryFieldRef();
68 
69  forAll(heBf, patchi)
70  {
71  heBf[patchi] == he
72  (
73  this->p_.boundaryField()[patchi],
74  this->T_.boundaryField()[patchi],
75  patchi
76  );
77  }
78 
79  this->heBoundaryCorrection(he_);
80 }
81 
82 
83 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
84 
85 
86 template<class BasicThermo, class MixtureType>
88 (
89  const fvMesh& mesh,
90  const word& phaseName
91 )
92 :
93  BasicThermo(mesh, phaseName),
94  MixtureType(*this, mesh, phaseName),
95 
96  he_
97  (
98  IOobject
99  (
100  BasicThermo::phasePropertyName
101  (
102  MixtureType::thermoType::heName()
103  ),
104  mesh.time().timeName(),
105  mesh,
106  IOobject::NO_READ,
107  IOobject::NO_WRITE
108  ),
109  mesh,
111  this->heBoundaryTypes(),
112  this->heBoundaryBaseTypes()
113  )
114 {
115  init();
116 }
117 
118 
119 template<class BasicThermo, class MixtureType>
121 (
122  const fvMesh& mesh,
123  const dictionary& dict,
124  const word& phaseName
125 )
126 :
127  BasicThermo(mesh, dict, phaseName),
128  MixtureType(*this, mesh, phaseName),
129 
130  he_
131  (
132  IOobject
133  (
134  BasicThermo::phasePropertyName
135  (
136  MixtureType::thermoType::heName()
137  ),
138  mesh.time().timeName(),
139  mesh,
140  IOobject::NO_READ,
141  IOobject::NO_WRITE
142  ),
143  mesh,
145  this->heBoundaryTypes(),
146  this->heBoundaryBaseTypes()
147  )
148 {
149  init();
150 }
151 
152 
153 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
154 
155 template<class BasicThermo, class MixtureType>
157 {}
158 
159 
160 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
161 
162 template<class BasicThermo, class MixtureType>
164 (
165  const volScalarField& p,
166  const volScalarField& T
167 ) const
168 {
169  const fvMesh& mesh = this->T_.mesh();
170 
172  (
173  new volScalarField
174  (
175  IOobject
176  (
177  "he",
178  mesh.time().timeName(),
179  mesh,
180  IOobject::NO_READ,
181  IOobject::NO_WRITE,
182  false
183  ),
184  mesh,
185  he_.dimensions()
186  )
187  );
188 
189  volScalarField& he = the.ref();
190  scalarField& heCells = he.primitiveFieldRef();
191  const scalarField& pCells = p;
192  const scalarField& TCells = T;
193 
194  forAll(heCells, celli)
195  {
196  heCells[celli] =
197  this->cellMixture(celli).HE(pCells[celli], TCells[celli]);
198  }
199 
200  volScalarField::Boundary& heBf = he.boundaryFieldRef();
201 
202  forAll(heBf, patchi)
203  {
204  scalarField& hep = heBf[patchi];
205  const scalarField& pp = p.boundaryField()[patchi];
206  const scalarField& Tp = T.boundaryField()[patchi];
207 
208  forAll(hep, facei)
209  {
210  hep[facei] =
211  this->patchFaceMixture(patchi, facei).HE(pp[facei], Tp[facei]);
212  }
213  }
214 
215  return the;
216 }
217 
218 
219 template<class BasicThermo, class MixtureType>
221 (
222  const scalarField& p,
223  const scalarField& T,
224  const labelList& cells
225 ) const
226 {
227  tmp<scalarField> the(new scalarField(T.size()));
228  scalarField& he = the.ref();
229 
230  forAll(T, celli)
231  {
232  he[celli] = this->cellMixture(cells[celli]).HE(p[celli], T[celli]);
233  }
234 
235  return the;
236 }
237 
238 
239 template<class BasicThermo, class MixtureType>
241 (
242  const scalarField& p,
243  const scalarField& T,
244  const label patchi
245 ) const
246 {
247  tmp<scalarField> the(new scalarField(T.size()));
248  scalarField& he = the.ref();
249 
250  forAll(T, facei)
251  {
252  he[facei] =
253  this->patchFaceMixture(patchi, facei).HE(p[facei], T[facei]);
254  }
255 
256  return the;
257 }
258 
259 
260 template<class BasicThermo, class MixtureType>
263 {
264  const fvMesh& mesh = this->T_.mesh();
265 
267  (
268  new volScalarField
269  (
270  IOobject
271  (
272  "hc",
273  mesh.time().timeName(),
274  mesh,
275  IOobject::NO_READ,
276  IOobject::NO_WRITE,
277  false
278  ),
279  mesh,
280  he_.dimensions()
281  )
282  );
283 
284  volScalarField& hcf = thc.ref();
285  scalarField& hcCells = hcf.primitiveFieldRef();
286 
287  forAll(hcCells, celli)
288  {
289  hcCells[celli] = this->cellMixture(celli).Hc();
290  }
291 
292  volScalarField::Boundary& hcfBf = hcf.boundaryFieldRef();
293 
294  forAll(hcfBf, patchi)
295  {
296  scalarField& hcp = hcfBf[patchi];
297 
298  forAll(hcp, facei)
299  {
300  hcp[facei] = this->patchFaceMixture(patchi, facei).Hc();
301  }
302  }
303 
304  return thc;
305 }
306 
307 
308 template<class BasicThermo, class MixtureType>
310 (
311  const scalarField& p,
312  const scalarField& T,
313  const label patchi
314 ) const
315 {
316  tmp<scalarField> tCp(new scalarField(T.size()));
317  scalarField& cp = tCp.ref();
318 
319  forAll(T, facei)
320  {
321  cp[facei] =
322  this->patchFaceMixture(patchi, facei).Cp(p[facei], T[facei]);
323  }
324 
325  return tCp;
326 }
327 
328 
329 template<class BasicThermo, class MixtureType>
332 {
333  const fvMesh& mesh = this->T_.mesh();
334 
336  (
337  new volScalarField
338  (
339  IOobject
340  (
341  "Cp",
342  mesh.time().timeName(),
343  mesh,
344  IOobject::NO_READ,
345  IOobject::NO_WRITE,
346  false
347  ),
348  mesh,
350  )
351  );
352 
353  volScalarField& cp = tCp.ref();
354 
355  forAll(this->T_, celli)
356  {
357  cp[celli] =
358  this->cellMixture(celli).Cp(this->p_[celli], this->T_[celli]);
359  }
360 
361  volScalarField::Boundary& cpBf = cp.boundaryFieldRef();
362 
363  forAll(cpBf, patchi)
364  {
365  const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
366  const fvPatchScalarField& pT = this->T_.boundaryField()[patchi];
367  fvPatchScalarField& pCp = cpBf[patchi];
368 
369  forAll(pT, facei)
370  {
371  pCp[facei] =
372  this->patchFaceMixture(patchi, facei).Cp(pp[facei], pT[facei]);
373  }
374  }
375 
376  return tCp;
377 }
378 
379 
380 template<class BasicThermo, class MixtureType>
383 (
384  const scalarField& p,
385  const scalarField& T,
386  const label patchi
387 ) const
388 {
389  tmp<scalarField> tCv(new scalarField(T.size()));
390  scalarField& cv = tCv.ref();
391 
392  forAll(T, facei)
393  {
394  cv[facei] =
395  this->patchFaceMixture(patchi, facei).Cv(p[facei], T[facei]);
396  }
397 
398  return tCv;
399 }
400 
401 
402 template<class BasicThermo, class MixtureType>
405 {
406  const fvMesh& mesh = this->T_.mesh();
407 
409  (
410  new volScalarField
411  (
412  IOobject
413  (
414  "Cv",
415  mesh.time().timeName(),
416  mesh,
417  IOobject::NO_READ,
418  IOobject::NO_WRITE,
419  false
420  ),
421  mesh,
423  )
424  );
425 
426  volScalarField& cv = tCv.ref();
427 
428  forAll(this->T_, celli)
429  {
430  cv[celli] =
431  this->cellMixture(celli).Cv(this->p_[celli], this->T_[celli]);
432  }
433 
434  volScalarField::Boundary& cvBf = cv.boundaryFieldRef();
435 
436  forAll(cvBf, patchi)
437  {
438  cvBf[patchi] = Cv
439  (
440  this->p_.boundaryField()[patchi],
441  this->T_.boundaryField()[patchi],
442  patchi
443  );
444  }
445 
446  return tCv;
447 }
448 
449 
450 template<class BasicThermo, class MixtureType>
452 (
453  const scalarField& p,
454  const scalarField& T,
455  const label patchi
456 ) const
457 {
458  tmp<scalarField> tgamma(new scalarField(T.size()));
459  scalarField& gamma = tgamma.ref();
460 
461  forAll(T, facei)
462  {
463  gamma[facei] =
464  this->patchFaceMixture(patchi, facei).gamma(p[facei], T[facei]);
465  }
466 
467  return tgamma;
468 }
469 
470 
471 template<class BasicThermo, class MixtureType>
474 {
475  const fvMesh& mesh = this->T_.mesh();
476 
477  tmp<volScalarField> tgamma
478  (
479  new volScalarField
480  (
481  IOobject
482  (
483  "gamma",
484  mesh.time().timeName(),
485  mesh,
486  IOobject::NO_READ,
487  IOobject::NO_WRITE,
488  false
489  ),
490  mesh,
491  dimless
492  )
493  );
494 
495  volScalarField& gamma = tgamma.ref();
496 
497  forAll(this->T_, celli)
498  {
499  gamma[celli] =
500  this->cellMixture(celli).gamma(this->p_[celli], this->T_[celli]);
501  }
502 
503  volScalarField::Boundary& gammaBf = gamma.boundaryFieldRef();
504 
505  forAll(gammaBf, patchi)
506  {
507  const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
508  const fvPatchScalarField& pT = this->T_.boundaryField()[patchi];
509  fvPatchScalarField& pgamma = gammaBf[patchi];
510 
511  forAll(pT, facei)
512  {
513  pgamma[facei] = this->patchFaceMixture(patchi, facei).gamma
514  (
515  pp[facei],
516  pT[facei]
517  );
518  }
519  }
520 
521  return tgamma;
522 }
523 
524 
525 template<class BasicThermo, class MixtureType>
527 (
528  const scalarField& p,
529  const scalarField& T,
530  const label patchi
531 ) const
532 {
533  tmp<scalarField> tCpv(new scalarField(T.size()));
534  scalarField& Cpv = tCpv.ref();
535 
536  forAll(T, facei)
537  {
538  Cpv[facei] =
539  this->patchFaceMixture(patchi, facei).Cpv(p[facei], T[facei]);
540  }
541 
542  return tCpv;
543 }
544 
545 
546 template<class BasicThermo, class MixtureType>
549 {
550  const fvMesh& mesh = this->T_.mesh();
551 
553  (
554  new volScalarField
555  (
556  IOobject
557  (
558  "Cpv",
559  mesh.time().timeName(),
560  mesh,
561  IOobject::NO_READ,
562  IOobject::NO_WRITE,
563  false
564  ),
565  mesh,
567  )
568  );
569 
570  volScalarField& Cpv = tCpv.ref();
571 
572  forAll(this->T_, celli)
573  {
574  Cpv[celli] =
575  this->cellMixture(celli).Cpv(this->p_[celli], this->T_[celli]);
576  }
577 
578  volScalarField::Boundary& CpvBf = Cpv.boundaryFieldRef();
579 
580  forAll(CpvBf, patchi)
581  {
582  const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
583  const fvPatchScalarField& pT = this->T_.boundaryField()[patchi];
584  fvPatchScalarField& pCpv = CpvBf[patchi];
585 
586  forAll(pT, facei)
587  {
588  pCpv[facei] =
589  this->patchFaceMixture(patchi, facei).Cpv(pp[facei], pT[facei]);
590  }
591  }
592 
593  return tCpv;
594 }
595 
596 
597 template<class BasicThermo, class MixtureType>
599 (
600  const scalarField& p,
601  const scalarField& T,
602  const label patchi
603 ) const
604 {
605  tmp<scalarField> tCpByCpv(new scalarField(T.size()));
606  scalarField& CpByCpv = tCpByCpv.ref();
607 
608  forAll(T, facei)
609  {
610  CpByCpv[facei] =
611  this->patchFaceMixture(patchi, facei).CpByCpv(p[facei], T[facei]);
612  }
613 
614  return tCpByCpv;
615 }
616 
617 
618 template<class BasicThermo, class MixtureType>
621 {
622  const fvMesh& mesh = this->T_.mesh();
623 
624  tmp<volScalarField> tCpByCpv
625  (
626  new volScalarField
627  (
628  IOobject
629  (
630  "CpByCpv",
631  mesh.time().timeName(),
632  mesh,
633  IOobject::NO_READ,
634  IOobject::NO_WRITE,
635  false
636  ),
637  mesh,
638  dimless
639  )
640  );
641 
642  volScalarField& CpByCpv = tCpByCpv.ref();
643 
644  forAll(this->T_, celli)
645  {
646  CpByCpv[celli] = this->cellMixture(celli).CpByCpv
647  (
648  this->p_[celli],
649  this->T_[celli]
650  );
651  }
652 
653  volScalarField::Boundary& CpByCpvBf =
654  CpByCpv.boundaryFieldRef();
655 
656  forAll(CpByCpvBf, patchi)
657  {
658  const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
659  const fvPatchScalarField& pT = this->T_.boundaryField()[patchi];
660  fvPatchScalarField& pCpByCpv = CpByCpvBf[patchi];
661 
662  forAll(pT, facei)
663  {
664  pCpByCpv[facei] = this->patchFaceMixture(patchi, facei).CpByCpv
665  (
666  pp[facei],
667  pT[facei]
668  );
669  }
670  }
671 
672  return tCpByCpv;
673 }
674 
675 
676 template<class BasicThermo, class MixtureType>
678 (
679  const scalarField& h,
680  const scalarField& p,
681  const scalarField& T0,
682  const labelList& cells
683 ) const
684 {
685  tmp<scalarField> tT(new scalarField(h.size()));
686  scalarField& T = tT.ref();
687 
688  forAll(h, celli)
689  {
690  T[celli] =
691  this->cellMixture(cells[celli]).THE(h[celli], p[celli], T0[celli]);
692  }
693 
694  return tT;
695 }
696 
697 
698 template<class BasicThermo, class MixtureType>
700 (
701  const scalarField& h,
702  const scalarField& p,
703  const scalarField& T0,
704  const label patchi
705 ) const
706 {
707 
708  tmp<scalarField> tT(new scalarField(h.size()));
709  scalarField& T = tT.ref();
710  forAll(h, facei)
711  {
712  T[facei] = this->patchFaceMixture
713  (
714  patchi,
715  facei
716  ).THE(h[facei], p[facei], T0[facei]);
717  }
718 
719  return tT;
720 }
721 
722 
723 template<class BasicThermo, class MixtureType>
725 (
726 ) const
727 {
728  const fvMesh& mesh = this->T_.mesh();
729 
731  (
732  new volScalarField
733  (
734  IOobject
735  (
736  "W",
737  mesh.time().timeName(),
738  mesh,
739  IOobject::NO_READ,
740  IOobject::NO_WRITE,
741  false
742  ),
743  mesh,
745  )
746  );
747 
748  volScalarField& W = tW.ref();
749  scalarField& WCells = W.primitiveFieldRef();
750 
751  forAll(WCells, celli)
752  {
753  WCells[celli] = this->cellMixture(celli).W();
754  }
755 
756  volScalarField::Boundary& WBf = W.boundaryFieldRef();
757 
758  forAll(WBf, patchi)
759  {
760  scalarField& Wp = WBf[patchi];
761  forAll(Wp, facei)
762  {
763  Wp[facei] = this->patchFaceMixture(patchi, facei).W();
764  }
765  }
766 
767  return tW;
768 }
769 
770 
771 template<class BasicThermo, class MixtureType>
774 {
775  tmp<Foam::volScalarField> kappa(Cp()*this->alpha_);
776  kappa.ref().rename("kappa");
777  return kappa;
778 }
779 
780 
781 template<class BasicThermo, class MixtureType>
783 (
784  const label patchi
785 ) const
786 {
787  return
788  Cp
789  (
790  this->p_.boundaryField()[patchi],
791  this->T_.boundaryField()[patchi],
792  patchi
793  )*this->alpha_.boundaryField()[patchi];
794 }
795 
796 
797 template<class BasicThermo, class MixtureType>
800 {
801  tmp<Foam::volScalarField> alphaEff(this->CpByCpv()*this->alpha_);
802  alphaEff.ref().rename("alphahe");
803  return alphaEff;
804 }
805 
806 
807 template<class BasicThermo, class MixtureType>
810 {
811  return
812  this->CpByCpv
813  (
814  this->p_.boundaryField()[patchi],
815  this->T_.boundaryField()[patchi],
816  patchi
817  )
818  *this->alpha_.boundaryField()[patchi];
819 }
820 
821 
822 template<class BasicThermo, class MixtureType>
825 (
826  const volScalarField& alphat
827 ) const
828 {
829  tmp<Foam::volScalarField> kappaEff(Cp()*(this->alpha_ + alphat));
830  kappaEff.ref().rename("kappaEff");
831  return kappaEff;
832 }
833 
834 
835 template<class BasicThermo, class MixtureType>
838 (
839  const scalarField& alphat,
840  const label patchi
841 ) const
842 {
843  return
844  Cp
845  (
846  this->p_.boundaryField()[patchi],
847  this->T_.boundaryField()[patchi],
848  patchi
849  )
850  *(
851  this->alpha_.boundaryField()[patchi]
852  + alphat
853  );
854 }
855 
856 
857 template<class BasicThermo, class MixtureType>
860 (
861  const volScalarField& alphat
862 ) const
863 {
864  tmp<Foam::volScalarField> alphaEff(this->CpByCpv()*(this->alpha_ + alphat));
865  alphaEff.ref().rename("alphaEff");
866  return alphaEff;
867 }
868 
869 
870 template<class BasicThermo, class MixtureType>
873 (
874  const scalarField& alphat,
875  const label patchi
876 ) const
877 {
878  return
879  this->CpByCpv
880  (
881  this->p_.boundaryField()[patchi],
882  this->T_.boundaryField()[patchi],
883  patchi
884  )
885  *(
886  this->alpha_.boundaryField()[patchi]
887  + alphat
888  );
889 }
890 
891 
892 template<class BasicThermo, class MixtureType>
894 {
895  if (BasicThermo::read())
896  {
897  MixtureType::read(*this);
898  return true;
899  }
900  else
901  {
902  return false;
903  }
904 }
905 
906 
907 // ************************************************************************* //
virtual tmp< volScalarField > kappa() const
Thermal diffusivity for temperature of mixture [J/m/s/K].
Definition: heThermo.C:773
volScalarField & he
Definition: YEEqn.H:51
virtual tmp< scalarField > THE(const scalarField &he, const scalarField &p, const scalarField &T0, const labelList &cells) const
Temperature from enthalpy/internal energy for cell-set.
Definition: heThermo.C:678
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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 tmp< volScalarField > gamma() const
Gamma = Cp/Cv [].
Definition: heThermo.C:473
virtual tmp< volScalarField > CpByCpv() const
Heat capacity ratio [].
Definition: heThermo.C:620
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
virtual ~heThermo()
Destructor.
Definition: heThermo.C:156
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
virtual tmp< volScalarField > Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
Definition: heThermo.C:548
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:626
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
virtual volScalarField & he()
Enthalpy/Internal energy [J/kg].
Definition: heThermo.H:147
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
virtual tmp< volScalarField > hc() const
Chemical enthalpy [J/kg].
Definition: heThermo.C:262
virtual tmp< volScalarField > Cp() const
Heat capacity at constant pressure [J/kg/K].
Definition: heThermo.C:331
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
virtual tmp< volScalarField > alphahe() const
Thermal diffusivity for energy of mixture [kg/m/s].
Definition: heThermo.C:799
volScalarField scalarField(fieldObject, mesh)
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
virtual tmp< volScalarField > alphaEff(const volScalarField &alphat) const
Effective thermal turbulent diffusivity of mixture [kg/m/s].
Definition: heThermo.C:860
const volScalarField & cp
const dimensionSet dimMoles(0, 0, 0, 0, 1, 0, 0)
Definition: dimensionSets.H:53
virtual void rename(const word &newName)
Rename.
Definition: regIOobject.C:409
const volScalarField & T
const dimensionSet dimEnergy
Internal & ref()
Return a reference to the dimensioned internal field.
virtual bool read()
Read thermophysical properties dictionary.
Definition: heThermo.C:893
label patchi
virtual tmp< volScalarField > kappaEff(const volScalarField &) const
Effective thermal turbulent diffusivity for temperature.
Definition: heThermo.C:825
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
virtual tmp< volScalarField > Cv() const
Heat capacity at constant volume [J/kg/K].
Definition: heThermo.C:404
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
Enthalpy/Internal energy for a mixture.
Definition: heThermo.H:49
const scalarList W(::W(thermo))
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
volScalarField alphaEff("alphaEff", turbulence->nu()/Pr+alphat)
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
void heBoundaryCorrection(volScalarField &he)
Correct the enthalpy/internal energy field boundaries.
Definition: heThermo.C:34
virtual tmp< volScalarField > W() const
Molecular weight [kg/kmol].
Definition: heThermo.C:725