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-2020 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>
33 template
34 <
35  class CellMixture,
36  class PatchFaceMixture,
37  class Method,
38  class ... Args
39 >
42 (
43  const word& psiName,
44  const dimensionSet& psiDim,
45  CellMixture cellMixture,
46  PatchFaceMixture patchFaceMixture,
47  Method psiMethod,
48  const Args& ... args
49 ) const
50 {
52  (
54  (
55  IOobject::groupName(psiName, this->group()),
56  this->T_.mesh(),
57  psiDim
58  )
59  );
60 
61  volScalarField& psi = tPsi.ref();
62 
63  forAll(this->T_, celli)
64  {
65  psi[celli] = ((this->*cellMixture)(celli).*psiMethod)(args[celli] ...);
66  }
67 
68  volScalarField::Boundary& psiBf = psi.boundaryFieldRef();
69 
70  forAll(psiBf, patchi)
71  {
72  fvPatchScalarField& pPsi = psiBf[patchi];
73 
74  forAll(this->T_.boundaryField()[patchi], facei)
75  {
76  pPsi[facei] =
77  ((this->*patchFaceMixture)(patchi, facei).*psiMethod)
78  (
79  args.boundaryField()[patchi][facei] ...
80  );
81  }
82  }
83 
84  return tPsi;
85 }
86 
87 
88 template<class BasicThermo, class MixtureType>
89 template<class CellMixture, class Method, class ... Args>
92 (
93  CellMixture cellMixture,
94  Method psiMethod,
95  const labelList& cells,
96  const Args& ... args
97 ) const
98 {
99  // Note: Args are fields for the set, not for the mesh as a whole. The
100  // cells list is only used to get the mixture.
101 
102  tmp<scalarField> tPsi(new scalarField(cells.size()));
103  scalarField& psi = tPsi.ref();
104 
105  forAll(cells, celli)
106  {
107  psi[celli] =
108  ((this->*cellMixture)(cells[celli]).*psiMethod)(args[celli] ...);
109  }
110 
111  return tPsi;
112 }
113 
114 
115 template<class BasicThermo, class MixtureType>
116 template<class PatchFaceMixture, class Method, class ... Args>
119 (
120  PatchFaceMixture patchFaceMixture,
121  Method psiMethod,
122  const label patchi,
123  const Args& ... args
124 ) const
125 {
126  tmp<scalarField> tPsi
127  (
128  new scalarField(this->T_.boundaryField()[patchi].size())
129  );
130  scalarField& psi = tPsi.ref();
131 
132  forAll(this->T_.boundaryField()[patchi], facei)
133  {
134  psi[facei] =
135  ((this->*patchFaceMixture)(patchi, facei).*psiMethod)
136  (
137  args[facei] ...
138  );
139  }
140 
141  return tPsi;
142 }
143 
144 
145 template<class BasicThermo, class MixtureType>
148 (
149  const volScalarField& psi,
150  const labelList& cells
151 )
152 {
154 }
155 
156 
157 template<class BasicThermo, class MixtureType>
160 (
161  const uniformGeometricScalarField& psi,
162  const labelList&
163 )
164 {
165  return psi.primitiveField();
166 }
167 
168 
169 template<class BasicThermo, class MixtureType>
172 {
173  volScalarField::Boundary& hBf = h.boundaryFieldRef();
174 
175  forAll(hBf, patchi)
176  {
177  if (isA<gradientEnergyFvPatchScalarField>(hBf[patchi]))
178  {
179  refCast<gradientEnergyFvPatchScalarField>(hBf[patchi]).gradient()
180  = hBf[patchi].fvPatchField::snGrad();
181  }
182  else if (isA<mixedEnergyFvPatchScalarField>(hBf[patchi]))
183  {
184  refCast<mixedEnergyFvPatchScalarField>(hBf[patchi]).refGrad()
185  = hBf[patchi].fvPatchField::snGrad();
186  }
187  }
188 }
189 
190 
191 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
192 
193 template<class BasicThermo, class MixtureType>
195 (
196  const fvMesh& mesh,
197  const word& phaseName
198 )
199 :
200  BasicThermo(mesh, phaseName),
201  MixtureType(*this, mesh, phaseName),
202 
203  he_
204  (
205  IOobject
206  (
207  BasicThermo::phasePropertyName
208  (
209  MixtureType::thermoType::heName(),
210  phaseName
211  ),
212  mesh.time().timeName(),
213  mesh,
214  IOobject::NO_READ,
215  IOobject::NO_WRITE
216  ),
217  volScalarFieldProperty
218  (
219  "he",
221  &MixtureType::cellThermoMixture,
222  &MixtureType::patchFaceThermoMixture,
223  &MixtureType::thermoMixtureType::HE,
224  this->p_,
225  this->T_
226  ),
227  this->heBoundaryTypes(),
228  this->heBoundaryBaseTypes()
229  ),
230 
231  Cp_
232  (
233  IOobject
234  (
235  BasicThermo::phasePropertyName("Cp", phaseName),
236  mesh.time().timeName(),
237  mesh
238  ),
239  mesh,
241  ),
242 
243  Cv_
244  (
245  IOobject
246  (
247  BasicThermo::phasePropertyName("Cv", phaseName),
248  mesh.time().timeName(),
249  mesh
250  ),
251  mesh,
253  )
254 {
255  heBoundaryCorrection(he_);
256 }
257 
258 
259 template<class BasicThermo, class MixtureType>
261 (
262  const fvMesh& mesh,
263  const dictionary& dict,
264  const word& phaseName
265 )
266 :
267  BasicThermo(mesh, dict, phaseName),
268  MixtureType(*this, mesh, phaseName),
269 
270  he_
271  (
272  IOobject
273  (
274  BasicThermo::phasePropertyName
275  (
276  MixtureType::thermoType::heName(),
277  phaseName
278  ),
279  mesh.time().timeName(),
280  mesh,
281  IOobject::NO_READ,
282  IOobject::NO_WRITE
283  ),
284  volScalarFieldProperty
285  (
286  "he",
288  &MixtureType::cellThermoMixture,
289  &MixtureType::patchFaceThermoMixture,
290  &MixtureType::thermoMixtureType::HE,
291  this->p_,
292  this->T_
293  ),
294  this->heBoundaryTypes(),
295  this->heBoundaryBaseTypes()
296  ),
297 
298  Cp_
299  (
300  IOobject
301  (
302  BasicThermo::phasePropertyName("Cp", phaseName),
303  mesh.time().timeName(),
304  mesh
305  ),
306  mesh,
308  ),
309 
310  Cv_
311  (
312  IOobject
313  (
314  BasicThermo::phasePropertyName("Cv", phaseName),
315  mesh.time().timeName(),
316  mesh
317  ),
318  mesh,
320  )
321 {
322  heBoundaryCorrection(he_);
323 }
324 
325 
326 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
327 
328 template<class BasicThermo, class MixtureType>
330 {}
331 
332 
333 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
334 
335 template<class BasicThermo, class MixtureType>
337 (
338  const volScalarField& p,
339  const volScalarField& T
340 ) const
341 {
342  return volScalarFieldProperty
343  (
344  "he",
346  &MixtureType::cellThermoMixture,
347  &MixtureType::patchFaceThermoMixture,
348  &MixtureType::thermoMixtureType::HE,
349  p,
350  T
351  );
352 }
353 
354 
355 template<class BasicThermo, class MixtureType>
357 (
358  const scalarField& T,
359  const labelList& cells
360 ) const
361 {
362  return cellSetProperty
363  (
364  &MixtureType::cellThermoMixture,
365  &MixtureType::thermoMixtureType::HE,
366  cells,
367  cellSetScalarList(this->p_, cells),
368  T
369  );
370 }
371 
372 
373 template<class BasicThermo, class MixtureType>
375 (
376  const scalarField& T,
377  const label patchi
378 ) const
379 {
380  return patchFieldProperty
381  (
382  &MixtureType::patchFaceThermoMixture,
383  &MixtureType::thermoMixtureType::HE,
384  patchi,
385  this->p_.boundaryField()[patchi],
386  T
387  );
388 }
389 
390 
391 template<class BasicThermo, class MixtureType>
394 {
395  return volScalarFieldProperty
396  (
397  "hs",
399  &MixtureType::cellThermoMixture,
400  &MixtureType::patchFaceThermoMixture,
402  this->p_,
403  this->T_
404  );
405 }
406 
407 
408 template<class BasicThermo, class MixtureType>
410 (
411  const volScalarField& p,
412  const volScalarField& T
413 ) const
414 {
415  return volScalarFieldProperty
416  (
417  "hs",
419  &MixtureType::cellThermoMixture,
420  &MixtureType::patchFaceThermoMixture,
422  p,
423  T
424  );
425 }
426 
427 
428 template<class BasicThermo, class MixtureType>
430 (
431  const scalarField& T,
432  const labelList& cells
433 ) const
434 {
435  return cellSetProperty
436  (
437  &MixtureType::cellThermoMixture,
439  cells,
440  cellSetScalarList(this->p_, cells),
441  T
442  );
443 }
444 
445 
446 template<class BasicThermo, class MixtureType>
448 (
449  const scalarField& T,
450  const label patchi
451 ) const
452 {
453  return patchFieldProperty
454  (
455  &MixtureType::patchFaceThermoMixture,
457  patchi,
458  this->p_.boundaryField()[patchi],
459  T
460  );
461 }
462 
463 
464 template<class BasicThermo, class MixtureType>
467 {
468  return volScalarFieldProperty
469  (
470  "ha",
472  &MixtureType::cellThermoMixture,
473  &MixtureType::patchFaceThermoMixture,
475  this->p_,
476  this->T_
477  );
478 }
479 
480 
481 template<class BasicThermo, class MixtureType>
483 (
484  const volScalarField& p,
485  const volScalarField& T
486 ) const
487 {
488  return volScalarFieldProperty
489  (
490  "ha",
492  &MixtureType::cellThermoMixture,
493  &MixtureType::patchFaceThermoMixture,
495  p,
496  T
497  );
498 }
499 
500 
501 template<class BasicThermo, class MixtureType>
503 (
504  const scalarField& T,
505  const labelList& cells
506 ) const
507 {
508  return cellSetProperty
509  (
510  &MixtureType::cellThermoMixture,
512  cells,
513  cellSetScalarList(this->p_, cells),
514  T
515  );
516 }
517 
518 
519 template<class BasicThermo, class MixtureType>
521 (
522  const scalarField& T,
523  const label patchi
524 ) const
525 {
526  return patchFieldProperty
527  (
528  &MixtureType::patchFaceThermoMixture,
530  patchi,
531  this->p_.boundaryField()[patchi],
532  T
533  );
534 }
535 
536 
537 template<class BasicThermo, class MixtureType>
540 {
541  return volScalarFieldProperty
542  (
543  "hc",
545  &MixtureType::cellThermoMixture,
546  &MixtureType::patchFaceThermoMixture,
547  &MixtureType::thermoMixtureType::Hf
548  );
549 }
550 
551 
552 template<class BasicThermo, class MixtureType>
554 (
555  const scalarField& T,
556  const label patchi
557 ) const
558 {
559  return patchFieldProperty
560  (
561  &MixtureType::patchFaceThermoMixture,
563  patchi,
564  this->p_.boundaryField()[patchi],
565  T
566  );
567 }
568 
569 
570 template<class BasicThermo, class MixtureType>
573 (
574  const scalarField& T,
575  const label patchi
576 ) const
577 {
578  return patchFieldProperty
579  (
580  &MixtureType::patchFaceThermoMixture,
582  patchi,
583  this->p_.boundaryField()[patchi],
584  T
585  );
586 }
587 
588 
589 template<class BasicThermo, class MixtureType>
591 (
592  const scalarField& T,
593  const label patchi
594 ) const
595 {
596  return patchFieldProperty
597  (
598  &MixtureType::patchFaceThermoMixture,
599  &MixtureType::thermoMixtureType::gamma,
600  patchi,
601  this->p_.boundaryField()[patchi],
602  T
603  );
604 }
605 
606 
607 template<class BasicThermo, class MixtureType>
610 {
611  return volScalarField::New("gamma", Cp_/Cv_);
612 }
613 
614 
615 template<class BasicThermo, class MixtureType>
617 (
618  const scalarField& T,
619  const label patchi
620 ) const
621 {
622  if (MixtureType::thermoType::enthalpy())
623  {
624  return Cp(T, patchi);
625  }
626  else
627  {
628  return Cv(T, patchi);
629  }
630 }
631 
632 
633 template<class BasicThermo, class MixtureType>
636 {
637  if (MixtureType::thermoType::enthalpy())
638  {
639  return Cp_;
640  }
641  else
642  {
643  return Cv_;
644  }
645 }
646 
647 
648 template<class BasicThermo, class MixtureType>
650 (
651  const volScalarField& h,
652  const volScalarField& p,
653  const volScalarField& T0
654 ) const
655 {
656  return volScalarFieldProperty
657  (
658  "T",
660  &MixtureType::cellThermoMixture,
661  &MixtureType::patchFaceThermoMixture,
662  &MixtureType::thermoMixtureType::THE,
663  h,
664  p,
665  T0
666  );
667 }
668 
669 
670 template<class BasicThermo, class MixtureType>
672 (
673  const scalarField& h,
674  const scalarField& T0,
675  const labelList& cells
676 ) const
677 {
678  return cellSetProperty
679  (
680  &MixtureType::cellThermoMixture,
681  &MixtureType::thermoMixtureType::THE,
682  cells,
683  h,
684  cellSetScalarList(this->p_, cells),
685  T0
686  );
687 }
688 
689 
690 template<class BasicThermo, class MixtureType>
692 (
693  const scalarField& h,
694  const scalarField& T0,
695  const label patchi
696 ) const
697 {
698  return patchFieldProperty
699  (
700  &MixtureType::patchFaceThermoMixture,
701  &MixtureType::thermoMixtureType::THE,
702  patchi,
703  h,
704  this->p_.boundaryField()[patchi],
705  T0
706  );
707 }
708 
709 
710 template<class BasicThermo, class MixtureType>
713 {
714  return volScalarFieldProperty
715  (
716  "W",
718  &MixtureType::cellThermoMixture,
719  &MixtureType::patchFaceThermoMixture,
721  );
722 }
723 
724 
725 template<class BasicThermo, class MixtureType>
727 (
728  const label patchi
729 ) const
730 {
731  return patchFieldProperty
732  (
733  &MixtureType::patchFaceThermoMixture,
735  patchi
736  );
737 }
738 
739 
740 template<class BasicThermo, class MixtureType>
743 {
744  return volScalarField::New("kappa", Cp_*this->alpha_);
745 }
746 
747 
748 template<class BasicThermo, class MixtureType>
750 (
751  const label patchi
752 ) const
753 {
754  return
755  Cp
756  (
757  this->T_.boundaryField()[patchi],
758  patchi
759  )*this->alpha_.boundaryField()[patchi];
760 }
761 
762 
763 template<class BasicThermo, class MixtureType>
766 {
767  if (MixtureType::thermoType::enthalpy())
768  {
769  return volScalarField::New("alphahe", this->alpha_);
770  }
771  else
772  {
773  return volScalarField::New
774  (
775  "alphahe",
776  this->gamma()*this->alpha_
777  );
778  }
779 }
780 
781 
782 template<class BasicThermo, class MixtureType>
785 {
786  if (MixtureType::thermoType::enthalpy())
787  {
788  return this->alpha_.boundaryField()[patchi];
789  }
790  else
791  {
792  return
793  this->gamma
794  (
795  this->T_.boundaryField()[patchi],
796  patchi
797  )*this->alpha_.boundaryField()[patchi];
798  }
799 }
800 
801 
802 template<class BasicThermo, class MixtureType>
805 (
806  const volScalarField& alphat
807 ) const
808 {
809  return volScalarField::New("kappaEff", Cp_*(this->alpha_ + alphat));
810 }
811 
812 
813 template<class BasicThermo, class MixtureType>
816 (
817  const scalarField& alphat,
818  const label patchi
819 ) const
820 {
821  return
822  Cp
823  (
824  this->T_.boundaryField()[patchi],
825  patchi
826  )
827  *(this->alpha_.boundaryField()[patchi] + alphat);
828 }
829 
830 
831 template<class BasicThermo, class MixtureType>
834 (
835  const volScalarField& alphat
836 ) const
837 {
838  if (MixtureType::thermoType::enthalpy())
839  {
840  return volScalarField::New("alphaEff", this->alpha_ + alphat);
841  }
842  else
843  {
844  return volScalarField::New
845  (
846  "alphaEff",
847  this->gamma()*(this->alpha_ + alphat)
848  );
849  }
850 }
851 
852 
853 template<class BasicThermo, class MixtureType>
856 (
857  const scalarField& alphat,
858  const label patchi
859 ) const
860 {
861  if (MixtureType::thermoType::enthalpy())
862  {
863  return (this->alpha_.boundaryField()[patchi] + alphat);
864  }
865  else
866  {
867  return
868  this->gamma
869  (
870  this->T_.boundaryField()[patchi],
871  patchi
872  )*(this->alpha_.boundaryField()[patchi] + alphat);
873  }
874 }
875 
876 
877 template<class BasicThermo, class MixtureType>
879 {
880  if (BasicThermo::read())
881  {
882  MixtureType::read(*this);
883  return true;
884  }
885  else
886  {
887  return false;
888  }
889 }
890 
891 
892 // ************************************************************************* //
virtual tmp< volScalarField > kappa() const
Thermal diffusivity for temperature of mixture [W/m/K].
Definition: heThermo.C:742
const char *const group
Group name for atomic constants.
A class representing the concept of a uniform field which stores only the single value and providing ...
Definition: UniformField.H:47
scalar Cv(const scalar p, const scalar T) const
Definition: HtoEthermo.H:2
virtual tmp< volScalarField > Cp() const
Heat capacity at constant pressure [J/kg/K].
Definition: heThermo.H:208
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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:609
scalar Hs(const scalar p, const scalar T) const
Definition: EtoHthermo.H:11
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
virtual ~heThermo()
Destructor.
Definition: heThermo.C:329
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
virtual tmp< volScalarField > Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
Definition: heThermo.C:635
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:636
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
virtual volScalarField & he()
Enthalpy/Internal energy [J/kg].
Definition: heThermo.H:196
virtual tmp< volScalarField > hc() const
Enthalpy of formation [J/kg].
Definition: heThermo.C:539
Internal::FieldType primitiveField() const
tmp< scalarField > patchFieldProperty(PatchFaceMixture patchFaceMixture, Method psiMethod, const label patchi, const Args &... args) const
Return a scalarField of the given property on a patch.
Dimension set for the base types.
Definition: dimensionSet.H:120
dynamicFvMesh & mesh
const cellShapeList & cells
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
A class for handling words, derived from string.
Definition: word.H:59
static UIndirectList< scalar > cellSetScalarList(const volScalarField &psi, const labelList &cells)
Return an indirect list of a field for the given set of cells.
Definition: heThermo.C:148
static const zero Zero
Definition: zero.H:97
virtual tmp< volScalarField > Cv() const
Heat capacity at constant volume [J/kg/K].
Definition: heThermo.H:214
virtual tmp< volScalarField > THE(const volScalarField &h, const volScalarField &p, const volScalarField &T0) const
Temperature from enthalpy/internal energy.
Definition: heThermo.C:650
virtual tmp< volScalarField > alphahe() const
Thermal diffusivity for energy of mixture [kg/m/s].
Definition: heThermo.C:765
virtual tmp< volScalarField > hs() const
Sensible enthalpy [J/kg/K].
Definition: heThermo.C:393
virtual tmp< volScalarField > ha() const
Absolute enthalpy [J/kg/K].
Definition: heThermo.C:466
volScalarField scalarField(fieldObject, mesh)
virtual tmp< volScalarField > alphaEff(const volScalarField &alphat) const
Effective thermal turbulent diffusivity of mixture [kg/m/s].
Definition: heThermo.C:834
tmp< scalarField > cellSetProperty(CellMixture cellMixture, Method psiMethod, const labelList &cells, const Args &... args) const
Return a scalarField of the given property on a cell set.
const dimensionSet dimEnergy
const dimensionSet dimMass
const volScalarField & T
heThermo(const fvMesh &, const word &phaseName)
Construct from mesh.
Definition: heThermo.C:195
virtual bool read()
Read thermophysical properties dictionary.
Definition: heThermo.C:878
label patchi
virtual tmp< volScalarField > kappaEff(const volScalarField &) const
Effective thermal turbulent diffusivity for temperature.
Definition: heThermo.C:805
const dimensionSet dimMoles
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const scalarList W(::W(thermo))
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A List with indirect addressing.
Definition: fvMatrix.H:106
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
scalar Ha(const scalar p, const scalar T) const
Definition: EtoHthermo.H:20
const volScalarField & psi
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:171
tmp< volScalarField > volScalarFieldProperty(const word &psiName, const dimensionSet &psiDim, CellMixture cellMixture, PatchFaceMixture patchFaceMixture, Method psiMethod, const Args &... args) const
Return a volScalarField of the given property.
const dimensionSet dimTemperature
virtual tmp< volScalarField > W() const
Molecular weight [kg/kmol].
Definition: heThermo.C:712
scalar T0
Definition: createFields.H:22