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