Loading...
Searching...
No Matches
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 {
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>
58void Foam::heThermo<BasicThermo, MixtureType>::init
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())
93 {
94 init(p.oldTime(), T.oldTime(), he.oldTime());
95 }
96}
97
98
99// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
100
101template<class BasicThermo, class MixtureType>
102Foam::heThermo<BasicThermo, MixtureType>::heThermo
103(
104 const fvMesh& mesh,
105 const word& phaseName
106)
107:
108 BasicThermo(mesh, phaseName),
109 MixtureType(*this, mesh, phaseName),
110
111 he_
112 (
113 IOobject
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,
125 dimEnergy/dimMass,
126 this->heBoundaryTypes(),
127 this->heBoundaryBaseTypes()
129{
130 init(this->p_, this->T_, he_);
132
133
134template<class BasicThermo, class MixtureType>
135Foam::heThermo<BasicThermo, MixtureType>::heThermo
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()
163{
164 init(this->p_, this->T_, he_);
165}
166
167
168template<class BasicThermo, class MixtureType>
169Foam::heThermo<BasicThermo, MixtureType>::heThermo
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_);
200
201
202
203// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
204
205template<class BasicThermo, class MixtureType>
207{}
208
209
210// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
211
212template<class BasicThermo, class MixtureType>
214(
215 const volScalarField& p,
216 const volScalarField& T
217) const
218{
219 const fvMesh& mesh = this->T_.mesh();
220
221 auto the = volScalarField::New
222 (
223 "he",
225 mesh,
226 he_.dimensions()
227 );
228 auto& he = the.ref();
229
230 scalarField& heCells = he.primitiveFieldRef();
231 const scalarField& pCells = p;
232 const scalarField& TCells = T;
233
234 forAll(heCells, celli)
235 {
236 heCells[celli] =
237 this->cellMixture(celli).HE(pCells[celli], TCells[celli]);
238 }
239
240 volScalarField::Boundary& heBf = he.boundaryFieldRef();
241
242 forAll(heBf, patchi)
243 {
244 scalarField& hep = heBf[patchi];
245 const scalarField& pp = p.boundaryField()[patchi];
246 const scalarField& Tp = T.boundaryField()[patchi];
247
248 forAll(hep, facei)
249 {
250 hep[facei] =
251 this->patchFaceMixture(patchi, facei).HE(pp[facei], Tp[facei]);
252 }
254
255 return the;
256}
257
258
259template<class BasicThermo, class MixtureType>
261(
262 const scalarField& p,
263 const scalarField& T,
264 const labelList& cells
265) const
266{
267 auto the = tmp<scalarField>::New(T.size());
268 auto& he = the.ref();
269
270 forAll(T, celli)
271 {
272 he[celli] = this->cellMixture(cells[celli]).HE(p[celli], T[celli]);
274
275 return the;
276}
277
278
279template<class BasicThermo, class MixtureType>
281(
282 const scalarField& p,
283 const scalarField& T,
284 const label patchi
285) const
286{
287 auto the = tmp<scalarField>::New(T.size());
288 auto& he = the.ref();
289
290 forAll(T, facei)
291 {
292 he[facei] =
293 this->patchFaceMixture(patchi, facei).HE(p[facei], T[facei]);
294 }
296 return the;
297}
298
299
300template<class BasicThermo, class MixtureType>
303{
304 const fvMesh& mesh = this->T_.mesh();
305
306 auto thc = volScalarField::New
307 (
308 "hc",
310 mesh,
311 he_.dimensions()
312 );
313 auto& hcf = thc.ref();
314
315 scalarField& hcCells = hcf.primitiveFieldRef();
316
317 forAll(hcCells, celli)
318 {
319 hcCells[celli] = this->cellMixture(celli).Hc();
320 }
321
322 volScalarField::Boundary& hcfBf = hcf.boundaryFieldRef();
323
324 forAll(hcfBf, patchi)
325 {
326 scalarField& hcp = hcfBf[patchi];
327
328 forAll(hcp, facei)
329 {
330 hcp[facei] = this->patchFaceMixture(patchi, facei).Hc();
331 }
333
334 return thc;
335}
336
337
338template<class BasicThermo, class MixtureType>
340(
341 const scalarField& p,
342 const scalarField& T,
343 const label patchi
344) const
345{
346 auto tCp = tmp<scalarField>::New(T.size());
347 auto& cp = tCp.ref();
348
349 forAll(T, facei)
350 {
351 cp[facei] =
352 this->patchFaceMixture(patchi, facei).Cp(p[facei], T[facei]);
353 }
355 return tCp;
356}
357
358
359template<class BasicThermo, class MixtureType>
363 const scalarField& p,
364 const scalarField& T,
365 const labelList& cells
366) const
367{
368 auto tCp = tmp<scalarField>::New(T.size());
369 auto& Cp = tCp.ref();
371 forAll(cells, i)
372 {
373 const label celli = cells[i];
374 Cp[i] = this->cellMixture(celli).Cp(p[i], T[i]);
375 }
377 return tCp;
378}
379
380
381template<class BasicThermo, class MixtureType>
384{
385 const fvMesh& mesh = this->T_.mesh();
386
388 (
389 "Cp",
391 mesh,
393 );
394 auto& cp = tCp.ref();
395
396 forAll(this->T_, celli)
398 cp[celli] =
399 this->cellMixture(celli).Cp(this->p_[celli], this->T_[celli]);
400 }
401
402 volScalarField::Boundary& cpBf = cp.boundaryFieldRef();
403
404 forAll(cpBf, patchi)
405 {
406 const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
407 const fvPatchScalarField& pT = this->T_.boundaryField()[patchi];
408 fvPatchScalarField& pCp = cpBf[patchi];
409
410 forAll(pT, facei)
411 {
412 pCp[facei] =
413 this->patchFaceMixture(patchi, facei).Cp(pp[facei], pT[facei]);
414 }
415 }
417 return tCp;
418}
419
420
421template<class BasicThermo, class MixtureType>
424(
425 const scalarField& p,
426 const scalarField& T,
427 const label patchi
428) const
429{
430 auto tCv = tmp<scalarField>::New(T.size());
431 auto& cv = tCv.ref();
432
433 forAll(T, facei)
434 {
435 cv[facei] =
436 this->patchFaceMixture(patchi, facei).Cv(p[facei], T[facei]);
437 }
439 return tCv;
440}
441
442
443template<class BasicThermo, class MixtureType>
446(
447 const scalarField& p,
448 const scalarField& T,
449 const labelList& cells
450) const
451{
452 auto tRho = tmp<scalarField>::New(T.size());
453 auto& rho = tRho.ref();
454
455 forAll(cells, i)
456 {
457 const label celli = cells[i];
458 rho[i] = this->cellMixture(celli).rho(p[i], T[i]);
459 }
461 return tRho;
462}
463
464
465template<class BasicThermo, class MixtureType>
468{
469 const fvMesh& mesh = this->T_.mesh();
470
472 (
473 "Cv",
475 mesh,
477 );
478 auto& cv = tCv.ref();
479
480 forAll(this->T_, celli)
481 {
482 cv[celli] =
483 this->cellMixture(celli).Cv(this->p_[celli], this->T_[celli]);
484 }
485
486 volScalarField::Boundary& cvBf = cv.boundaryFieldRef();
487
488 forAll(cvBf, patchi)
489 {
490 cvBf[patchi] = Cv
491 (
492 this->p_.boundaryField()[patchi],
493 this->T_.boundaryField()[patchi],
494 patchi
495 );
497
498 return tCv;
499}
500
501
502template<class BasicThermo, class MixtureType>
504(
505 const scalarField& p,
506 const scalarField& T,
507 const label patchi
508) const
509{
510 auto tgamma = tmp<scalarField>::New(T.size());
511 auto& gamma = tgamma.ref();
512
513 forAll(T, facei)
514 {
515 gamma[facei] =
516 this->patchFaceMixture(patchi, facei).gamma(p[facei], T[facei]);
517 }
519 return tgamma;
520}
521
522
523template<class BasicThermo, class MixtureType>
526{
527 const fvMesh& mesh = this->T_.mesh();
528
529 auto tgamma = volScalarField::New
530 (
531 "gamma",
533 mesh,
534 dimless
535 );
536 auto& gamma = tgamma.ref();
537
538 forAll(this->T_, celli)
539 {
540 gamma[celli] =
541 this->cellMixture(celli).gamma(this->p_[celli], this->T_[celli]);
542 }
543
544 volScalarField::Boundary& gammaBf = gamma.boundaryFieldRef();
545
546 forAll(gammaBf, patchi)
547 {
548 const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
549 const fvPatchScalarField& pT = this->T_.boundaryField()[patchi];
550 fvPatchScalarField& pgamma = gammaBf[patchi];
551
552 forAll(pT, facei)
553 {
554 pgamma[facei] = this->patchFaceMixture(patchi, facei).gamma
555 (
556 pp[facei],
557 pT[facei]
558 );
559 }
561
562 return tgamma;
563}
564
565
566template<class BasicThermo, class MixtureType>
568(
569 const scalarField& p,
570 const scalarField& T,
571 const label patchi
572) const
573{
574 auto tCpv = tmp<scalarField>::New(T.size());
575 auto& Cpv = tCpv.ref();
576
577 forAll(T, facei)
578 {
579 Cpv[facei] =
580 this->patchFaceMixture(patchi, facei).Cpv(p[facei], T[facei]);
581 }
583 return tCpv;
584}
585
586
587template<class BasicThermo, class MixtureType>
590{
591 const fvMesh& mesh = this->T_.mesh();
592
593 auto tCpv = volScalarField::New
594 (
595 "Cpv",
597 mesh,
599 );
600 auto& Cpv = tCpv.ref();
601
602 forAll(this->T_, celli)
603 {
604 Cpv[celli] =
605 this->cellMixture(celli).Cpv(this->p_[celli], this->T_[celli]);
606 }
607
608 volScalarField::Boundary& CpvBf = Cpv.boundaryFieldRef();
609
610 forAll(CpvBf, patchi)
611 {
612 const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
613 const fvPatchScalarField& pT = this->T_.boundaryField()[patchi];
614 fvPatchScalarField& pCpv = CpvBf[patchi];
615
616 forAll(pT, facei)
617 {
618 pCpv[facei] =
619 this->patchFaceMixture(patchi, facei).Cpv(pp[facei], pT[facei]);
620 }
622
623 return tCpv;
624}
625
626
627template<class BasicThermo, class MixtureType>
629(
630 const scalarField& p,
631 const scalarField& T,
632 const label patchi
633) const
634{
635 auto tCpByCpv = tmp<scalarField>::New(T.size());
636 auto& CpByCpv = tCpByCpv.ref();
637
638 forAll(T, facei)
639 {
640 CpByCpv[facei] =
641 this->patchFaceMixture(patchi, facei).CpByCpv(p[facei], T[facei]);
642 }
644 return tCpByCpv;
645}
646
647
648template<class BasicThermo, class MixtureType>
651{
652 const fvMesh& mesh = this->T_.mesh();
653
654 auto tCpByCpv = volScalarField::New
655 (
656 "CpByCpv",
658 mesh,
659 dimless
660 );
661 auto& CpByCpv = tCpByCpv.ref();
662
663 forAll(this->T_, celli)
664 {
665 CpByCpv[celli] = this->cellMixture(celli).CpByCpv
666 (
667 this->p_[celli],
668 this->T_[celli]
669 );
670 }
671
672 volScalarField::Boundary& CpByCpvBf =
673 CpByCpv.boundaryFieldRef();
674
675 forAll(CpByCpvBf, patchi)
676 {
677 const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
678 const fvPatchScalarField& pT = this->T_.boundaryField()[patchi];
679 fvPatchScalarField& pCpByCpv = CpByCpvBf[patchi];
680
681 forAll(pT, facei)
682 {
683 pCpByCpv[facei] = this->patchFaceMixture(patchi, facei).CpByCpv
684 (
685 pp[facei],
686 pT[facei]
687 );
688 }
690
691 return tCpByCpv;
692}
693
694
695template<class BasicThermo, class MixtureType>
697(
698 const scalarField& h,
699 const scalarField& p,
700 const scalarField& T0,
701 const labelList& cells
702) const
703{
704 auto tT = tmp<scalarField>::New(h.size());
705 auto& T = tT.ref();
706
707 forAll(h, celli)
708 {
709 T[celli] =
710 this->cellMixture(cells[celli]).THE(h[celli], p[celli], T0[celli]);
712
713 return tT;
714}
715
716
717template<class BasicThermo, class MixtureType>
719(
720 const scalarField& h,
721 const scalarField& p,
722 const scalarField& T0,
723 const label patchi
724) const
725{
726
727 auto tT = tmp<scalarField>::New(h.size());
728 auto& T = tT.ref();
729
730 forAll(h, facei)
731 {
732 T[facei] = this->patchFaceMixture
733 (
734 patchi,
735 facei
736 ).THE(h[facei], p[facei], T0[facei]);
738
739 return tT;
740}
741
742
743template<class BasicThermo, class MixtureType>
745(
746) const
747{
748 const fvMesh& mesh = this->T_.mesh();
749
750 auto tW = volScalarField::New
751 (
752 "W",
754 mesh,
756 );
757 auto& W = tW.ref();
758
759 scalarField& WCells = W.primitiveFieldRef();
760
761 forAll(WCells, celli)
762 {
763 WCells[celli] = this->cellMixture(celli).W();
764 }
765
766 auto& WBf = W.boundaryFieldRef();
767
768 forAll(WBf, patchi)
769 {
770 scalarField& Wp = WBf[patchi];
771 forAll(Wp, facei)
772 {
773 Wp[facei] = this->patchFaceMixture(patchi, facei).W();
774 }
775 }
777 return tW;
778}
779
780
781template<class BasicThermo, class MixtureType>
784{
786 kappa.ref().rename("kappa");
787 return kappa;
788}
789
790
791template<class BasicThermo, class MixtureType>
793(
794 const label patchi
795) const
796{
797 return
798 Cp
799 (
800 this->p_.boundaryField()[patchi],
801 this->T_.boundaryField()[patchi],
802 patchi
803 )*this->alpha_.boundaryField()[patchi];
804}
805
806
807template<class BasicThermo, class MixtureType>
810{
811 tmp<Foam::volScalarField> alphaEff(this->CpByCpv()*this->alpha_);
812 alphaEff.ref().rename("alphahe");
813 return alphaEff;
814}
815
816
817template<class BasicThermo, class MixtureType>
820{
821 return
822 this->CpByCpv
823 (
824 this->p_.boundaryField()[patchi],
825 this->T_.boundaryField()[patchi],
826 patchi
828 *this->alpha_.boundaryField()[patchi];
829}
830
831
832template<class BasicThermo, class MixtureType>
835(
836 const volScalarField& alphat
837) const
838{
839 tmp<Foam::volScalarField> kappaEff(Cp()*(this->alpha_ + alphat));
840 kappaEff.ref().rename("kappaEff");
841 return kappaEff;
842}
843
844
845template<class BasicThermo, class MixtureType>
848(
849 const scalarField& alphat,
850 const label patchi
851) const
852{
853 return
854 Cp
855 (
856 this->p_.boundaryField()[patchi],
857 this->T_.boundaryField()[patchi],
858 patchi
859 )
860 *(
861 this->alpha_.boundaryField()[patchi]
862 + alphat
863 );
864}
865
866
867template<class BasicThermo, class MixtureType>
870(
871 const volScalarField& alphat
872) const
873{
874 tmp<Foam::volScalarField> alphaEff(this->CpByCpv()*(this->alpha_ + alphat));
875 alphaEff.ref().rename("alphaEff");
876 return alphaEff;
877}
878
879
880template<class BasicThermo, class MixtureType>
883(
884 const scalarField& alphat,
885 const label patchi
886) const
887{
888 return
889 this->CpByCpv
890 (
891 this->p_.boundaryField()[patchi],
892 this->T_.boundaryField()[patchi],
893 patchi
894 )
895 *(
896 this->alpha_.boundaryField()[patchi]
897 + alphat
898 );
899}
900
901
902template<class BasicThermo, class MixtureType>
904{
905 if (BasicThermo::read())
906 {
907 MixtureType::read(*this);
908 return true;
909 }
910
911 return false;
912}
913
914
915// ************************************************************************* //
volScalarField & he
Definition YEEqn.H:52
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_REGISTER
Do not request registration (bool: false).
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
virtual tmp< scalarField > Cpv(const scalarField &p, const scalarField &T, const label patchi) const
Heat capacity at constant pressure/volume for patch [J/kg/K].
Definition heThermo.C:561
virtual tmp< volScalarField > W() const
Molecular weight [kg/kmol].
Definition heThermo.C:738
virtual tmp< volScalarField > Cp() const
Heat capacity at constant pressure [J/kg/K].
Definition heThermo.C:376
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:690
virtual tmp< volScalarField > kappaEff(const volScalarField &) const
Effective thermal diffusivity for temperature.
Definition heThermo.C:828
virtual tmp< volScalarField > alphahe() const
Thermal diffusivity for energy of mixture [kg/m/s].
Definition heThermo.C:802
virtual tmp< volScalarField > alphaEff(const volScalarField &alphat) const
Effective thermal diffusivity of mixture [kg/m/s].
Definition heThermo.C:863
volScalarField he_
Energy field.
Definition heThermo.H:60
virtual tmp< volScalarField > CpByCpv() const
Heat capacity ratio [].
Definition heThermo.C:643
virtual ~heThermo()
Destructor.
Definition heThermo.C:199
virtual tmp< volScalarField > hc() const
Chemical enthalpy [J/kg].
Definition heThermo.C:295
void heBoundaryCorrection(volScalarField &he)
Correct the enthalpy/internal energy field boundaries.
Definition heThermo.C:30
virtual tmp< volScalarField > gamma() const
Gamma = Cp/Cv [].
Definition heThermo.C:518
virtual volScalarField & he()
Enthalpy/Internal energy [J/kg].
Definition heThermo.H:190
virtual tmp< scalarField > Cp(const scalarField &p, const scalarField &T, const label patchi) const
Definition heThermo.C:333
virtual tmp< volScalarField > kappa() const
Thermal diffusivity for temperature of mixture [J/m/s/K].
Definition heThermo.C:776
virtual tmp< scalarField > rhoEoS(const scalarField &p, const scalarField &T, const labelList &cells) const
Density from pressure and temperature.
Definition heThermo.C:439
virtual tmp< volScalarField > Cv() const
Heat capacity at constant volume [J/kg/K].
Definition heThermo.C:460
virtual bool read()
Read thermophysical properties dictionary.
Definition heThermo.C:896
virtual tmp< volScalarField > Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
Definition heThermo.C:582
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
A class for handling words, derived from Foam::string.
Definition word.H:66
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
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
const dimensionSet dimless
Dimensionless.
List< label > labelList
A List of labels.
Definition List.H:62
const dimensionSet dimEnergy
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimMoles(0, 0, 0, 0, 1, 0, 0)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy the source to the destination (recursively if necessary).
Definition POSIX.C:1065
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
fvPatchField< scalar > fvPatchScalarField
dictionary dict
volScalarField & h
scalar T0
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299