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  volScalarField::Boundary& hBf = h.boundaryFieldRef();
150 
151  forAll(hBf, patchi)
152  {
153  if (isA<gradientEnergyFvPatchScalarField>(hBf[patchi]))
154  {
155  refCast<gradientEnergyFvPatchScalarField>(hBf[patchi]).gradient()
156  = hBf[patchi].fvPatchField::snGrad();
157  }
158  else if (isA<mixedEnergyFvPatchScalarField>(hBf[patchi]))
159  {
160  refCast<mixedEnergyFvPatchScalarField>(hBf[patchi]).refGrad()
161  = hBf[patchi].fvPatchField::snGrad();
162  }
163  }
164 }
165 
166 
167 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
168 
169 template<class BasicThermo, class MixtureType>
171 (
172  const fvMesh& mesh,
173  const word& phaseName
174 )
175 :
176  BasicThermo(mesh, phaseName),
177  MixtureType(*this, mesh, phaseName),
178 
179  he_
180  (
181  IOobject
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  he(this->p(), this->T_),
193  this->heBoundaryTypes(),
194  this->heBoundaryBaseTypes()
195  )
196 {
197  heBoundaryCorrection(he_);
198 }
199 
200 
201 template<class BasicThermo, class MixtureType>
203 (
204  const fvMesh& mesh,
205  const dictionary& dict,
206  const word& phaseName
207 )
208 :
209  BasicThermo(mesh, dict, phaseName),
210  MixtureType(*this, mesh, phaseName),
211 
212  he_
213  (
214  IOobject
215  (
216  BasicThermo::phasePropertyName
217  (
218  MixtureType::thermoType::heName()
219  ),
220  mesh.time().timeName(),
221  mesh,
222  IOobject::NO_READ,
223  IOobject::NO_WRITE
224  ),
225  he(this->p(), this->T_),
226  this->heBoundaryTypes(),
227  this->heBoundaryBaseTypes()
228  )
229 {
230  heBoundaryCorrection(he_);
231 }
232 
233 
234 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
235 
236 template<class BasicThermo, class MixtureType>
238 {}
239 
240 
241 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
242 
243 template<class BasicThermo, class MixtureType>
245 (
246  const volScalarField& p,
247  const volScalarField& T
248 ) const
249 {
250  return volScalarFieldProperty
251  (
252  "he",
254  &MixtureType::cellMixture,
255  &MixtureType::patchFaceMixture,
256  &MixtureType::thermoType::HE,
257  p,
258  T
259  );
260 }
261 
262 
263 template<class BasicThermo, class MixtureType>
265 (
266  const scalarField& T,
267  const labelList& cells
268 ) const
269 {
270  return cellSetProperty
271  (
272  &MixtureType::cellMixture,
273  &MixtureType::thermoType::HE,
274  cells,
275  UIndirectList<scalar>(this->p(), cells),
276  T
277  );
278 }
279 
280 
281 template<class BasicThermo, class MixtureType>
283 (
284  const scalarField& T,
285  const label patchi
286 ) const
287 {
288  return patchFieldProperty
289  (
290  &MixtureType::patchFaceMixture,
291  &MixtureType::thermoType::HE,
292  patchi,
293  this->p().boundaryField()[patchi],
294  T
295  );
296 }
297 
298 
299 template<class BasicThermo, class MixtureType>
302 {
303  return volScalarFieldProperty
304  (
305  "hs",
307  &MixtureType::cellMixture,
308  &MixtureType::patchFaceMixture,
310  this->p(),
311  this->T_
312  );
313 }
314 
315 
316 template<class BasicThermo, class MixtureType>
318 (
319  const volScalarField& p,
320  const volScalarField& T
321 ) const
322 {
323  return volScalarFieldProperty
324  (
325  "hs",
327  &MixtureType::cellMixture,
328  &MixtureType::patchFaceMixture,
330  p,
331  T
332  );
333 }
334 
335 
336 template<class BasicThermo, class MixtureType>
338 (
339  const scalarField& T,
340  const labelList& cells
341 ) const
342 {
343  return cellSetProperty
344  (
345  &MixtureType::cellMixture,
347  cells,
348  UIndirectList<scalar>(this->p(), cells),
349  T
350  );
351 }
352 
353 
354 template<class BasicThermo, class MixtureType>
356 (
357  const scalarField& T,
358  const label patchi
359 ) const
360 {
361  return patchFieldProperty
362  (
363  &MixtureType::patchFaceMixture,
365  patchi,
366  this->p().boundaryField()[patchi],
367  T
368  );
369 }
370 
371 
372 template<class BasicThermo, class MixtureType>
375 {
376  return volScalarFieldProperty
377  (
378  "ha",
380  &MixtureType::cellMixture,
381  &MixtureType::patchFaceMixture,
383  this->p(),
384  this->T_
385  );
386 }
387 
388 
389 template<class BasicThermo, class MixtureType>
391 (
392  const volScalarField& p,
393  const volScalarField& T
394 ) const
395 {
396  return volScalarFieldProperty
397  (
398  "ha",
400  &MixtureType::cellMixture,
401  &MixtureType::patchFaceMixture,
403  p,
404  T
405  );
406 }
407 
408 
409 template<class BasicThermo, class MixtureType>
411 (
412  const scalarField& T,
413  const labelList& cells
414 ) const
415 {
416  return cellSetProperty
417  (
418  &MixtureType::cellMixture,
420  cells,
421  UIndirectList<scalar>(this->p(), cells),
422  T
423  );
424 }
425 
426 
427 template<class BasicThermo, class MixtureType>
429 (
430  const scalarField& T,
431  const label patchi
432 ) const
433 {
434  return patchFieldProperty
435  (
436  &MixtureType::patchFaceMixture,
438  patchi,
439  this->p().boundaryField()[patchi],
440  T
441  );
442 }
443 
444 
445 template<class BasicThermo, class MixtureType>
448 {
449  return volScalarFieldProperty
450  (
451  "hc",
453  &MixtureType::cellMixture,
454  &MixtureType::patchFaceMixture,
455  &MixtureType::thermoType::Hf
456  );
457 }
458 
459 
460 template<class BasicThermo, class MixtureType>
462 (
463  const scalarField& T,
464  const label patchi
465 ) const
466 {
467  return patchFieldProperty
468  (
469  &MixtureType::patchFaceMixture,
471  patchi,
472  this->p().boundaryField()[patchi],
473  T
474  );
475 }
476 
477 
478 template<class BasicThermo, class MixtureType>
481 {
482  return volScalarFieldProperty
483  (
484  "Cp",
486  &MixtureType::cellMixture,
487  &MixtureType::patchFaceMixture,
489  this->p(),
490  this->T_
491  );
492 }
493 
494 
495 template<class BasicThermo, class MixtureType>
498 (
499  const scalarField& T,
500  const label patchi
501 ) const
502 {
503  return patchFieldProperty
504  (
505  &MixtureType::patchFaceMixture,
507  patchi,
508  this->p().boundaryField()[patchi],
509  T
510  );
511 }
512 
513 
514 template<class BasicThermo, class MixtureType>
517 {
518  return volScalarFieldProperty
519  (
520  "Cv",
522  &MixtureType::cellMixture,
523  &MixtureType::patchFaceMixture,
525  this->p(),
526  this->T_
527  );
528 }
529 
530 
531 template<class BasicThermo, class MixtureType>
533 (
534  const scalarField& T,
535  const label patchi
536 ) const
537 {
538  return patchFieldProperty
539  (
540  &MixtureType::patchFaceMixture,
541  &MixtureType::thermoType::gamma,
542  patchi,
543  this->p().boundaryField()[patchi],
544  T
545  );
546 }
547 
548 
549 template<class BasicThermo, class MixtureType>
552 {
553  return volScalarFieldProperty
554  (
555  "gamma",
556  dimless,
557  &MixtureType::cellMixture,
558  &MixtureType::patchFaceMixture,
559  &MixtureType::thermoType::gamma,
560  this->p(),
561  this->T_
562  );
563 }
564 
565 
566 template<class BasicThermo, class MixtureType>
568 (
569  const scalarField& T,
570  const label patchi
571 ) const
572 {
573  return patchFieldProperty
574  (
575  &MixtureType::patchFaceMixture,
576  &MixtureType::thermoType::Cpv,
577  patchi,
578  this->p().boundaryField()[patchi],
579  T
580  );
581 }
582 
583 
584 template<class BasicThermo, class MixtureType>
587 {
588  return volScalarFieldProperty
589  (
590  "Cpv",
592  &MixtureType::cellMixture,
593  &MixtureType::patchFaceMixture,
594  &MixtureType::thermoType::Cpv,
595  this->p(),
596  this->T_
597  );
598 }
599 
600 
601 template<class BasicThermo, class MixtureType>
603 (
604  const scalarField& T,
605  const label patchi
606 ) const
607 {
608  return patchFieldProperty
609  (
610  &MixtureType::patchFaceMixture,
611  &MixtureType::thermoType::CpByCpv,
612  patchi,
613  this->p().boundaryField()[patchi],
614  T
615  );
616 }
617 
618 
619 template<class BasicThermo, class MixtureType>
622 {
623  return volScalarFieldProperty
624  (
625  "CpByCpv",
626  dimless,
627  &MixtureType::cellMixture,
628  &MixtureType::patchFaceMixture,
629  &MixtureType::thermoType::CpByCpv,
630  this->p(),
631  this->T_
632  );
633 }
634 
635 
636 template<class BasicThermo, class MixtureType>
638 (
639  const scalarField& h,
640  const scalarField& T0,
641  const labelList& cells
642 ) const
643 {
644  return cellSetProperty
645  (
646  &MixtureType::cellMixture,
647  &MixtureType::thermoType::THE,
648  cells,
649  h,
650  UIndirectList<scalar>(this->p(), cells),
651  T0
652  );
653 }
654 
655 
656 template<class BasicThermo, class MixtureType>
658 (
659  const scalarField& h,
660  const scalarField& T0,
661  const label patchi
662 ) const
663 {
664  return patchFieldProperty
665  (
666  &MixtureType::patchFaceMixture,
667  &MixtureType::thermoType::THE,
668  patchi,
669  h,
670  this->p().boundaryField()[patchi],
671  T0
672  );
673 }
674 
675 
676 template<class BasicThermo, class MixtureType>
679 {
680  return volScalarFieldProperty
681  (
682  "W",
684  &MixtureType::cellMixture,
685  &MixtureType::patchFaceMixture,
687  );
688 }
689 
690 
691 template<class BasicThermo, class MixtureType>
693 (
694  const label patchi
695 ) const
696 {
697  return patchFieldProperty
698  (
699  &MixtureType::patchFaceMixture,
701  patchi
702  );
703 }
704 
705 
706 template<class BasicThermo, class MixtureType>
709 {
710  return volScalarField::New
711  (
712  "kappa",
713  Cp()*this->alpha_
714  );
715 }
716 
717 
718 template<class BasicThermo, class MixtureType>
720 (
721  const label patchi
722 ) const
723 {
724  return
725  Cp
726  (
727  this->T_.boundaryField()[patchi],
728  patchi
729  )*this->alpha_.boundaryField()[patchi];
730 }
731 
732 
733 template<class BasicThermo, class MixtureType>
736 {
737  return volScalarField::New
738  (
739  "alphahe",
740  this->CpByCpv()*this->alpha_
741  );
742 }
743 
744 
745 template<class BasicThermo, class MixtureType>
748 {
749  return
750  this->CpByCpv
751  (
752  this->T_.boundaryField()[patchi],
753  patchi
754  )
755  *this->alpha_.boundaryField()[patchi];
756 }
757 
758 
759 template<class BasicThermo, class MixtureType>
762 (
763  const volScalarField& alphat
764 ) const
765 {
766  return volScalarField::New
767  (
768  "kappaEff",
769  Cp()*(this->alpha_ + alphat)
770  );
771 }
772 
773 
774 template<class BasicThermo, class MixtureType>
777 (
778  const scalarField& alphat,
779  const label patchi
780 ) const
781 {
782  return
783  Cp
784  (
785  this->T_.boundaryField()[patchi],
786  patchi
787  )
788  *(this->alpha_.boundaryField()[patchi] + alphat);
789 }
790 
791 
792 template<class BasicThermo, class MixtureType>
795 (
796  const volScalarField& alphat
797 ) const
798 {
799  return volScalarField::New
800  (
801  "alphaEff",
802  this->CpByCpv()*(this->alpha_ + alphat)
803  );
804 }
805 
806 
807 template<class BasicThermo, class MixtureType>
810 (
811  const scalarField& alphat,
812  const label patchi
813 ) const
814 {
815  return
816  this->CpByCpv
817  (
818  this->T_.boundaryField()[patchi],
819  patchi
820  )
821  *(this->alpha_.boundaryField()[patchi] + alphat);
822 }
823 
824 
825 template<class BasicThermo, class MixtureType>
827 {
828  if (BasicThermo::read())
829  {
830  MixtureType::read(*this);
831  return true;
832  }
833  else
834  {
835  return false;
836  }
837 }
838 
839 
840 // ************************************************************************* //
virtual tmp< volScalarField > kappa() const
Thermal diffusivity for temperature of mixture [W/m/K].
Definition: heThermo.C:708
volScalarField & he
Definition: YEEqn.H:50
const char *const group
Group name for atomic constants.
#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:551
virtual tmp< volScalarField > CpByCpv() const
Heat capacity ratio [].
Definition: heThermo.C:621
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:158
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
const volScalarField & Cv
virtual ~heThermo()
Destructor.
Definition: heThermo.C:237
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:586
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:622
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:174
virtual tmp< volScalarField > hc() const
Enthalpy of formation [J/kg].
Definition: heThermo.C:447
tmp< scalarField > patchFieldProperty(PatchFaceMixture patchFaceMixture, Method psiMethod, const label patchi, const Args &... args) const
Return a scalarField of the given property on a patch.
virtual tmp< volScalarField > Cp() const
Heat capacity at constant pressure [J/kg/K].
Definition: heThermo.C:480
Dimension set for the base types.
Definition: dimensionSet.H:120
dynamicFvMesh & mesh
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
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:735
virtual tmp< volScalarField > hs() const
Sensible enthalpy [J/kg/K].
Definition: heThermo.C:301
virtual tmp< volScalarField > ha() const
Absolute enthalpy [J/kg/K].
Definition: heThermo.C:374
volScalarField scalarField(fieldObject, mesh)
virtual tmp< volScalarField > alphaEff(const volScalarField &alphat) const
Effective thermal turbulent diffusivity of mixture [kg/m/s].
Definition: heThermo.C:795
const dimensionSet dimMoles(0, 0, 0, 0, 1, 0, 0)
Definition: dimensionSets.H:53
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
heThermo(const fvMesh &, const word &phaseName)
Construct from mesh.
Definition: heThermo.C:171
virtual bool read()
Read thermophysical properties dictionary.
Definition: heThermo.C:826
label patchi
virtual tmp< volScalarField > kappaEff(const volScalarField &) const
Effective thermal turbulent diffusivity for temperature.
Definition: heThermo.C:762
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:516
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
const scalarList W(::W(thermo))
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual tmp< scalarField > THE(const scalarField &he, const scalarField &T0, const labelList &cells) const
Temperature from enthalpy/internal energy for cell-set.
Definition: heThermo.C:638
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
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
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:147
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.
virtual tmp< volScalarField > W() const
Molecular weight [kg/kmol].
Definition: heThermo.C:678