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-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
35template<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
57template<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 heBf[patchi].useImplicit(T.boundaryField()[patchi].useImplicit());
87 }
88
89 this->heBoundaryCorrection(he);
90
91 // Note: T does not have oldTime
92 if (p.nOldTimes() > 0)
93 {
94 init(p.oldTime(), T.oldTime(), he.oldTime());
95 }
96}
97
98
99// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
100
101template<class BasicThermo, class MixtureType>
103(
104 const fvMesh& mesh,
105 const word& phaseName
106)
107:
108 BasicThermo(mesh, phaseName),
109 MixtureType(*this, mesh, phaseName),
111 he_
112 (
114 (
115 BasicThermo::phasePropertyName
116 (
117 MixtureType::thermoType::heName()
118 ),
119 mesh.time().timeName(),
120 mesh,
121 IOobject::NO_READ,
122 IOobject::NO_WRITE
123 ),
124 mesh,
126 this->heBoundaryTypes(),
127 this->heBoundaryBaseTypes()
128 )
129{
130 init(this->p_, this->T_, he_);
131}
132
133
134template<class BasicThermo, class MixtureType>
136(
137 const fvMesh& mesh,
138 const dictionary& dict,
139 const word& phaseName
140)
141:
142 BasicThermo(mesh, dict, phaseName),
143 MixtureType(*this, mesh, phaseName),
144
145 he_
146 (
148 (
149 BasicThermo::phasePropertyName
150 (
151 MixtureType::thermoType::heName()
152 ),
153 mesh.time().timeName(),
154 mesh,
155 IOobject::NO_READ,
156 IOobject::NO_WRITE
157 ),
158 mesh,
160 this->heBoundaryTypes(),
161 this->heBoundaryBaseTypes()
162 )
163{
164 init(this->p_, this->T_, he_);
165}
166
167
168template<class BasicThermo, class MixtureType>
170(
171 const fvMesh& mesh,
172 const word& phaseName,
173 const word& dictionaryName
174)
175:
176 BasicThermo(mesh, phaseName, dictionaryName),
177 MixtureType(*this, mesh, phaseName),
178
179 he_
180 (
182 (
183 BasicThermo::phasePropertyName
184 (
185 MixtureType::thermoType::heName()
186 ),
187 mesh.time().timeName(),
188 mesh,
189 IOobject::NO_READ,
190 IOobject::NO_WRITE
191 ),
192 mesh,
194 this->heBoundaryTypes(),
195 this->heBoundaryBaseTypes()
196 )
197{
198 init(this->p_, this->T_, he_);
199}
200
201
203// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
204
205template<class BasicThermo, class MixtureType>
207{}
208
209
210// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
211
212template<class BasicThermo, class MixtureType>
215 const volScalarField& p,
216 const volScalarField& T
217) const
218{
219 const fvMesh& mesh = this->T_.mesh();
220
222 (
224 (
226 (
227 "he",
228 mesh.time().timeName(),
229 mesh,
232 false
233 ),
234 mesh,
235 he_.dimensions()
236 )
237 );
238
240 scalarField& heCells = he.primitiveFieldRef();
241 const scalarField& pCells = p;
242 const scalarField& TCells = T;
243
244 forAll(heCells, celli)
245 {
246 heCells[celli] =
247 this->cellMixture(celli).HE(pCells[celli], TCells[celli]);
248 }
249
251
252 forAll(heBf, patchi)
253 {
254 scalarField& hep = heBf[patchi];
255 const scalarField& pp = p.boundaryField()[patchi];
256 const scalarField& Tp = T.boundaryField()[patchi];
257
258 forAll(hep, facei)
259 {
260 hep[facei] =
261 this->patchFaceMixture(patchi, facei).HE(pp[facei], Tp[facei]);
262 }
263 }
265 return the;
266}
267
268
269template<class BasicThermo, class MixtureType>
271(
273 const scalarField& T,
274 const labelList& cells
275) const
276{
278 scalarField& he = the.ref();
279
280 forAll(T, celli)
281 {
282 he[celli] = this->cellMixture(cells[celli]).HE(p[celli], T[celli]);
284
285 return the;
287
288
289template<class BasicThermo, class MixtureType>
291(
292 const scalarField& p,
293 const scalarField& T,
294 const label patchi
295) const
296{
298 scalarField& he = the.ref();
299
300 forAll(T, facei)
301 {
302 he[facei] =
303 this->patchFaceMixture(patchi, facei).HE(p[facei], T[facei]);
305
306 return the;
307}
308
309
310template<class BasicThermo, class MixtureType>
314 const fvMesh& mesh = this->T_.mesh();
315
319 (
322 "hc",
323 mesh.time().timeName(),
324 mesh,
327 false
328 ),
329 mesh,
330 he_.dimensions()
331 )
332 );
333
334 volScalarField& hcf = thc.ref();
335 scalarField& hcCells = hcf.primitiveFieldRef();
336
337 forAll(hcCells, celli)
338 {
339 hcCells[celli] = this->cellMixture(celli).Hc();
340 }
341
343
344 forAll(hcfBf, patchi)
345 {
346 scalarField& hcp = hcfBf[patchi];
347
348 forAll(hcp, facei)
349 {
350 hcp[facei] = this->patchFaceMixture(patchi, facei).Hc();
351 }
352 }
353
354 return thc;
355}
356
357
358template<class BasicThermo, class MixtureType>
360(
361 const scalarField& p,
362 const scalarField& T,
363 const label patchi
364) const
365{
367 scalarField& cp = tCp.ref();
368
369 forAll(T, facei)
370 {
371 cp[facei] =
372 this->patchFaceMixture(patchi, facei).Cp(p[facei], T[facei]);
373 }
374
375 return tCp;
376}
377
378
379template<class BasicThermo, class MixtureType>
382(
383 const scalarField& p,
384 const scalarField& T,
385 const labelList& cells
386) const
387{
389 auto& Cp = tCp.ref();
390
391 forAll(cells, i)
392 {
393 const label celli = cells[i];
394 Cp[i] = this->cellMixture(celli).Cp(p[i], T[i]);
395 }
396
397 return tCp;
398}
399
400
401template<class BasicThermo, class MixtureType>
404{
405 const fvMesh& mesh = this->T_.mesh();
406
408 (
410 (
412 (
413 "Cp",
414 mesh.time().timeName(),
415 mesh,
418 false
419 ),
420 mesh,
422 )
423 );
424
425 volScalarField& cp = tCp.ref();
426
427 forAll(this->T_, celli)
428 {
429 cp[celli] =
430 this->cellMixture(celli).Cp(this->p_[celli], this->T_[celli]);
431 }
432
434
435 forAll(cpBf, patchi)
436 {
437 const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
438 const fvPatchScalarField& pT = this->T_.boundaryField()[patchi];
439 fvPatchScalarField& pCp = cpBf[patchi];
440
441 forAll(pT, facei)
442 {
443 pCp[facei] =
444 this->patchFaceMixture(patchi, facei).Cp(pp[facei], pT[facei]);
445 }
446 }
447
448 return tCp;
449}
450
451
452template<class BasicThermo, class MixtureType>
455(
456 const scalarField& p,
457 const scalarField& T,
458 const label patchi
459) const
460{
462 scalarField& cv = tCv.ref();
463
464 forAll(T, facei)
465 {
466 cv[facei] =
467 this->patchFaceMixture(patchi, facei).Cv(p[facei], T[facei]);
468 }
469
470 return tCv;
471}
472
473
474template<class BasicThermo, class MixtureType>
477(
478 const scalarField& p,
479 const scalarField& T,
480 const labelList& cells
481) const
482{
483 auto tRho = tmp<scalarField>::New(T.size());
484 auto& rho = tRho.ref();
485
486 forAll(cells, i)
487 {
488 const label celli = cells[i];
489 rho[i] = this->cellMixture(celli).rho(p[i], T[i]);
490 }
491
492 return tRho;
493}
494
495
496template<class BasicThermo, class MixtureType>
499{
500 const fvMesh& mesh = this->T_.mesh();
501
503 (
505 (
507 (
508 "Cv",
509 mesh.time().timeName(),
510 mesh,
513 false
514 ),
515 mesh,
517 )
518 );
519
520 volScalarField& cv = tCv.ref();
521
522 forAll(this->T_, celli)
523 {
524 cv[celli] =
525 this->cellMixture(celli).Cv(this->p_[celli], this->T_[celli]);
526 }
527
529
530 forAll(cvBf, patchi)
531 {
532 cvBf[patchi] = Cv
533 (
534 this->p_.boundaryField()[patchi],
535 this->T_.boundaryField()[patchi],
536 patchi
537 );
538 }
539
540 return tCv;
541}
542
543
544template<class BasicThermo, class MixtureType>
546(
547 const scalarField& p,
548 const scalarField& T,
549 const label patchi
550) const
551{
552 tmp<scalarField> tgamma(new scalarField(T.size()));
553 scalarField& gamma = tgamma.ref();
554
555 forAll(T, facei)
556 {
557 gamma[facei] =
558 this->patchFaceMixture(patchi, facei).gamma(p[facei], T[facei]);
559 }
560
561 return tgamma;
562}
563
564
565template<class BasicThermo, class MixtureType>
568{
569 const fvMesh& mesh = this->T_.mesh();
570
572 (
574 (
576 (
577 "gamma",
578 mesh.time().timeName(),
579 mesh,
582 false
583 ),
584 mesh,
585 dimless
586 )
587 );
588
589 volScalarField& gamma = tgamma.ref();
590
591 forAll(this->T_, celli)
592 {
593 gamma[celli] =
594 this->cellMixture(celli).gamma(this->p_[celli], this->T_[celli]);
595 }
596
597 volScalarField::Boundary& gammaBf = gamma.boundaryFieldRef();
598
599 forAll(gammaBf, patchi)
600 {
601 const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
602 const fvPatchScalarField& pT = this->T_.boundaryField()[patchi];
603 fvPatchScalarField& pgamma = gammaBf[patchi];
604
605 forAll(pT, facei)
606 {
607 pgamma[facei] = this->patchFaceMixture(patchi, facei).gamma
608 (
609 pp[facei],
610 pT[facei]
611 );
612 }
613 }
614
615 return tgamma;
616}
617
618
619template<class BasicThermo, class MixtureType>
621(
622 const scalarField& p,
623 const scalarField& T,
624 const label patchi
625) const
626{
627 tmp<scalarField> tCpv(new scalarField(T.size()));
628 scalarField& Cpv = tCpv.ref();
629
630 forAll(T, facei)
631 {
632 Cpv[facei] =
633 this->patchFaceMixture(patchi, facei).Cpv(p[facei], T[facei]);
634 }
635
636 return tCpv;
637}
638
639
640template<class BasicThermo, class MixtureType>
643{
644 const fvMesh& mesh = this->T_.mesh();
645
647 (
649 (
651 (
652 "Cpv",
653 mesh.time().timeName(),
654 mesh,
657 false
658 ),
659 mesh,
661 )
662 );
663
664 volScalarField& Cpv = tCpv.ref();
665
666 forAll(this->T_, celli)
667 {
668 Cpv[celli] =
669 this->cellMixture(celli).Cpv(this->p_[celli], this->T_[celli]);
670 }
671
673
674 forAll(CpvBf, patchi)
675 {
676 const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
677 const fvPatchScalarField& pT = this->T_.boundaryField()[patchi];
678 fvPatchScalarField& pCpv = CpvBf[patchi];
679
680 forAll(pT, facei)
681 {
682 pCpv[facei] =
683 this->patchFaceMixture(patchi, facei).Cpv(pp[facei], pT[facei]);
684 }
685 }
686
687 return tCpv;
688}
689
690
691template<class BasicThermo, class MixtureType>
693(
694 const scalarField& p,
695 const scalarField& T,
696 const label patchi
697) const
698{
699 tmp<scalarField> tCpByCpv(new scalarField(T.size()));
700 scalarField& CpByCpv = tCpByCpv.ref();
701
702 forAll(T, facei)
703 {
704 CpByCpv[facei] =
705 this->patchFaceMixture(patchi, facei).CpByCpv(p[facei], T[facei]);
706 }
707
708 return tCpByCpv;
709}
710
711
712template<class BasicThermo, class MixtureType>
715{
716 const fvMesh& mesh = this->T_.mesh();
717
718 tmp<volScalarField> tCpByCpv
719 (
721 (
723 (
724 "CpByCpv",
725 mesh.time().timeName(),
726 mesh,
729 false
730 ),
731 mesh,
732 dimless
733 )
734 );
735
736 volScalarField& CpByCpv = tCpByCpv.ref();
737
738 forAll(this->T_, celli)
739 {
740 CpByCpv[celli] = this->cellMixture(celli).CpByCpv
741 (
742 this->p_[celli],
743 this->T_[celli]
744 );
745 }
746
747 volScalarField::Boundary& CpByCpvBf =
748 CpByCpv.boundaryFieldRef();
749
750 forAll(CpByCpvBf, patchi)
751 {
752 const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
753 const fvPatchScalarField& pT = this->T_.boundaryField()[patchi];
754 fvPatchScalarField& pCpByCpv = CpByCpvBf[patchi];
755
756 forAll(pT, facei)
757 {
758 pCpByCpv[facei] = this->patchFaceMixture(patchi, facei).CpByCpv
759 (
760 pp[facei],
761 pT[facei]
762 );
763 }
764 }
765
766 return tCpByCpv;
767}
768
769
770template<class BasicThermo, class MixtureType>
772(
773 const scalarField& h,
774 const scalarField& p,
775 const scalarField& T0,
776 const labelList& cells
777) const
778{
780 scalarField& T = tT.ref();
781
782 forAll(h, celli)
783 {
784 T[celli] =
785 this->cellMixture(cells[celli]).THE(h[celli], p[celli], T0[celli]);
786 }
787
788 return tT;
789}
790
791
792template<class BasicThermo, class MixtureType>
794(
795 const scalarField& h,
796 const scalarField& p,
797 const scalarField& T0,
798 const label patchi
799) const
800{
801
803 scalarField& T = tT.ref();
804 forAll(h, facei)
805 {
806 T[facei] = this->patchFaceMixture
807 (
808 patchi,
809 facei
810 ).THE(h[facei], p[facei], T0[facei]);
811 }
812
813 return tT;
814}
815
816
817template<class BasicThermo, class MixtureType>
819(
820) const
821{
822 const fvMesh& mesh = this->T_.mesh();
823
825 (
827 (
829 (
830 "W",
831 mesh.time().timeName(),
832 mesh,
835 false
836 ),
837 mesh,
839 )
840 );
841
842 volScalarField& W = tW.ref();
843 scalarField& WCells = W.primitiveFieldRef();
844
845 forAll(WCells, celli)
846 {
847 WCells[celli] = this->cellMixture(celli).W();
848 }
849
851
852 forAll(WBf, patchi)
853 {
854 scalarField& Wp = WBf[patchi];
855 forAll(Wp, facei)
856 {
857 Wp[facei] = this->patchFaceMixture(patchi, facei).W();
858 }
859 }
860
861 return tW;
862}
863
864
865template<class BasicThermo, class MixtureType>
868{
869 tmp<Foam::volScalarField> kappa(Cp()*this->alpha_);
870 kappa.ref().rename("kappa");
871 return kappa;
872}
873
874
875template<class BasicThermo, class MixtureType>
877(
878 const label patchi
879) const
880{
881 return
882 Cp
883 (
884 this->p_.boundaryField()[patchi],
885 this->T_.boundaryField()[patchi],
886 patchi
887 )*this->alpha_.boundaryField()[patchi];
888}
889
890
891template<class BasicThermo, class MixtureType>
894{
895 tmp<Foam::volScalarField> alphaEff(this->CpByCpv()*this->alpha_);
896 alphaEff.ref().rename("alphahe");
897 return alphaEff;
898}
899
900
901template<class BasicThermo, class MixtureType>
904{
905 return
906 this->CpByCpv
907 (
908 this->p_.boundaryField()[patchi],
909 this->T_.boundaryField()[patchi],
910 patchi
911 )
912 *this->alpha_.boundaryField()[patchi];
913}
914
915
916template<class BasicThermo, class MixtureType>
919(
920 const volScalarField& alphat
921) const
922{
923 tmp<Foam::volScalarField> kappaEff(Cp()*(this->alpha_ + alphat));
924 kappaEff.ref().rename("kappaEff");
925 return kappaEff;
926}
927
928
929template<class BasicThermo, class MixtureType>
932(
933 const scalarField& alphat,
934 const label patchi
935) const
936{
937 return
938 Cp
939 (
940 this->p_.boundaryField()[patchi],
941 this->T_.boundaryField()[patchi],
942 patchi
943 )
944 *(
945 this->alpha_.boundaryField()[patchi]
946 + alphat
947 );
948}
949
950
951template<class BasicThermo, class MixtureType>
954(
955 const volScalarField& alphat
956) const
957{
958 tmp<Foam::volScalarField> alphaEff(this->CpByCpv()*(this->alpha_ + alphat));
959 alphaEff.ref().rename("alphaEff");
960 return alphaEff;
961}
962
963
964template<class BasicThermo, class MixtureType>
967(
968 const scalarField& alphat,
969 const label patchi
970) const
971{
972 return
973 this->CpByCpv
974 (
975 this->p_.boundaryField()[patchi],
976 this->T_.boundaryField()[patchi],
977 patchi
978 )
979 *(
980 this->alpha_.boundaryField()[patchi]
981 + alphat
982 );
983}
984
985
986template<class BasicThermo, class MixtureType>
988{
989 if (BasicThermo::read())
990 {
991 MixtureType::read(*this);
992 return true;
993 }
994
995 return false;
996}
997
998
999// ************************************************************************* //
volScalarField & he
Definition: YEEqn.H:52
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
Enthalpy/Internal energy for a mixture.
Definition: heThermo.H:56
virtual tmp< volScalarField > W() const
Molecular weight [kg/kmol].
Definition: heThermo.C:819
virtual tmp< volScalarField > Cp() const
Heat capacity at constant pressure [J/kg/K].
Definition: heThermo.C:403
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:772
virtual tmp< volScalarField > alphahe() const
Thermal diffusivity for energy of mixture [kg/m/s].
Definition: heThermo.C:893
volScalarField he_
Energy field.
Definition: heThermo.H:62
virtual tmp< volScalarField > CpByCpv() const
Heat capacity ratio [].
Definition: heThermo.C:714
virtual ~heThermo()
Destructor.
Definition: heThermo.C:206
virtual tmp< volScalarField > hc() const
Chemical enthalpy [J/kg].
Definition: heThermo.C:312
void heBoundaryCorrection(volScalarField &he)
Correct the enthalpy/internal energy field boundaries.
Definition: heThermo.C:37
virtual tmp< volScalarField > gamma() const
Gamma = Cp/Cv [].
Definition: heThermo.C:567
virtual volScalarField & he()
Enthalpy/Internal energy [J/kg].
Definition: heThermo.H:163
virtual tmp< volScalarField > kappa() const
Thermal diffusivity for temperature of mixture [J/m/s/K].
Definition: heThermo.C:867
virtual tmp< scalarField > rhoEoS(const scalarField &p, const scalarField &T, const labelList &cells) const
Density from pressure and temperature.
Definition: heThermo.C:477
virtual tmp< volScalarField > Cv() const
Heat capacity at constant volume [J/kg/K].
Definition: heThermo.C:498
virtual bool read()
Read thermophysical properties dictionary.
Definition: heThermo.C:987
virtual tmp< volScalarField > Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
Definition: heThermo.C:642
tmp< volScalarField > alphaEff() const
Effective thermal turbulent diffusivity of mixture [kg/m/s].
virtual void rename(const word &newName)
Rename.
Definition: regIOobject.C:417
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
const volScalarField & T
const tmp< volScalarField > & tCv
Definition: EEqn.H:5
const volScalarField & Cv
Definition: EEqn.H:8
const scalar gamma
Definition: EEqn.H:9
const tmp< volScalarField > & tCp
Definition: EEqn.H:4
const volScalarField & Cp
Definition: EEqn.H:7
dynamicFvMesh & mesh
const cellShapeList & cells
volScalarField alphaEff("alphaEff", turbulence->nu()/Pr+alphat)
word timeName
Definition: getTimeIndex.H:3
kappaEff
Definition: TEqn.H:10
const dimensionSet dimless
Dimensionless.
const dimensionSet dimEnergy
const dimensionSet dimMoles(0, 0, 0, 0, 1, 0, 0)
Definition: dimensionSets.H:55
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:54
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
dictionary dict
volScalarField & h
const volScalarField & cp
scalar T0
Definition: createFields.H:22
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333