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-2022 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->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 (
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().name(),
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  (
234  (
235  BasicThermo::phasePropertyName("Cp", phaseName),
236  mesh.time().name(),
237  mesh
238  ),
239  mesh,
241  ),
242 
243  Cv_
244  (
245  IOobject
246  (
247  BasicThermo::phasePropertyName("Cv", phaseName),
248  mesh.time().name(),
249  mesh
250  ),
251  mesh,
253  )
254 {
255  heBoundaryCorrection(he_);
256 }
257 
258 
259 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
260 
261 template<class BasicThermo, class MixtureType>
263 {}
264 
265 
266 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
267 
268 template<class BasicThermo, class MixtureType>
271 {
272  if (MixtureType::thermoType::enthalpy())
273  {
274  return Cp_;
275  }
276  else
277  {
278  return Cv_;
279  }
280 }
281 
282 
283 template<class BasicThermo, class MixtureType>
285 (
286  const volScalarField& p,
287  const volScalarField& T
288 ) const
289 {
290  return volScalarFieldProperty
291  (
292  "he",
294  &MixtureType::cellThermoMixture,
295  &MixtureType::patchFaceThermoMixture,
296  &MixtureType::thermoMixtureType::HE,
297  p,
298  T
299  );
300 }
301 
302 
303 template<class BasicThermo, class MixtureType>
305 (
306  const scalarField& T,
307  const labelList& cells
308 ) const
309 {
310  return cellSetProperty
311  (
312  &MixtureType::cellThermoMixture,
313  &MixtureType::thermoMixtureType::HE,
314  cells,
315  cellSetScalarList(this->p_, cells),
316  T
317  );
318 }
319 
320 
321 template<class BasicThermo, class MixtureType>
323 (
324  const scalarField& T,
325  const label patchi
326 ) const
327 {
328  return patchFieldProperty
329  (
330  &MixtureType::patchFaceThermoMixture,
331  &MixtureType::thermoMixtureType::HE,
332  patchi,
333  this->p_.boundaryField()[patchi],
335  );
336 }
337 
338 
339 template<class BasicThermo, class MixtureType>
342 {
343  return volScalarFieldProperty
344  (
345  "hs",
347  &MixtureType::cellThermoMixture,
348  &MixtureType::patchFaceThermoMixture,
350  this->p_,
351  this->T_
352  );
353 }
354 
355 
356 template<class BasicThermo, class MixtureType>
358 (
359  const volScalarField& p,
360  const volScalarField& T
361 ) const
362 {
363  return volScalarFieldProperty
364  (
365  "hs",
367  &MixtureType::cellThermoMixture,
368  &MixtureType::patchFaceThermoMixture,
370  p,
371  T
372  );
373 }
374 
375 
376 template<class BasicThermo, class MixtureType>
378 (
379  const scalarField& T,
380  const labelList& cells
381 ) const
382 {
383  return cellSetProperty
384  (
385  &MixtureType::cellThermoMixture,
387  cells,
388  cellSetScalarList(this->p_, cells),
389  T
390  );
391 }
392 
393 
394 template<class BasicThermo, class MixtureType>
396 (
397  const scalarField& T,
398  const label patchi
399 ) const
400 {
401  return patchFieldProperty
402  (
403  &MixtureType::patchFaceThermoMixture,
405  patchi,
406  this->p_.boundaryField()[patchi],
407  T
408  );
409 }
410 
411 
412 template<class BasicThermo, class MixtureType>
415 {
416  return volScalarFieldProperty
417  (
418  "ha",
420  &MixtureType::cellThermoMixture,
421  &MixtureType::patchFaceThermoMixture,
423  this->p_,
424  this->T_
425  );
426 }
427 
428 
429 template<class BasicThermo, class MixtureType>
431 (
432  const volScalarField& p,
433  const volScalarField& T
434 ) const
435 {
436  return volScalarFieldProperty
437  (
438  "ha",
440  &MixtureType::cellThermoMixture,
441  &MixtureType::patchFaceThermoMixture,
443  p,
444  T
445  );
446 }
447 
448 
449 template<class BasicThermo, class MixtureType>
451 (
452  const scalarField& T,
453  const labelList& cells
454 ) const
455 {
456  return cellSetProperty
457  (
458  &MixtureType::cellThermoMixture,
460  cells,
461  cellSetScalarList(this->p_, cells),
462  T
463  );
464 }
465 
466 
467 template<class BasicThermo, class MixtureType>
469 (
470  const scalarField& T,
471  const label patchi
472 ) const
473 {
474  return patchFieldProperty
475  (
476  &MixtureType::patchFaceThermoMixture,
478  patchi,
479  this->p_.boundaryField()[patchi],
480  T
481  );
482 }
483 
484 
485 template<class BasicThermo, class MixtureType>
488 {
489  return volScalarFieldProperty
490  (
491  "hc",
493  &MixtureType::cellThermoMixture,
494  &MixtureType::patchFaceThermoMixture,
495  &MixtureType::thermoMixtureType::Hf
496  );
497 }
498 
499 
500 template<class BasicThermo, class MixtureType>
502 (
503  const scalarField& T,
504  const label patchi
505 ) const
506 {
507  return patchFieldProperty
508  (
509  &MixtureType::patchFaceThermoMixture,
511  patchi,
512  this->p_.boundaryField()[patchi],
513  T
514  );
515 }
516 
517 
518 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>
539 (
540  const scalarField& T,
541  const label patchi
542 ) const
543 {
544  return patchFieldProperty
545  (
546  &MixtureType::patchFaceThermoMixture,
547  &MixtureType::thermoMixtureType::gamma,
548  patchi,
549  this->p_.boundaryField()[patchi],
550  T
551  );
552 }
553 
554 
555 template<class BasicThermo, class MixtureType>
558 {
559  return volScalarField::New("gamma", Cp_/Cv_);
560 }
561 
562 
563 template<class BasicThermo, class MixtureType>
565 (
566  const scalarField& T,
567  const label patchi
568 ) const
569 {
570  if (MixtureType::thermoType::enthalpy())
571  {
572  return Cp(T, patchi);
573  }
574  else
575  {
576  return Cv(T, patchi);
577  }
578 }
579 
580 
581 template<class BasicThermo, class MixtureType>
583 (
584  const volScalarField& h,
585  const volScalarField& p,
586  const volScalarField& T0
587 ) const
588 {
589  return volScalarFieldProperty
590  (
591  "T",
593  &MixtureType::cellThermoMixture,
594  &MixtureType::patchFaceThermoMixture,
595  &MixtureType::thermoMixtureType::THE,
596  h,
597  p,
598  T0
599  );
600 }
601 
602 
603 template<class BasicThermo, class MixtureType>
605 (
606  const scalarField& h,
607  const scalarField& T0,
608  const labelList& cells
609 ) const
610 {
611  return cellSetProperty
612  (
613  &MixtureType::cellThermoMixture,
614  &MixtureType::thermoMixtureType::THE,
615  cells,
616  h,
617  cellSetScalarList(this->p_, cells),
618  T0
619  );
620 }
621 
622 
623 template<class BasicThermo, class MixtureType>
625 (
626  const scalarField& h,
627  const scalarField& T0,
628  const label patchi
629 ) const
630 {
631  return patchFieldProperty
632  (
633  &MixtureType::patchFaceThermoMixture,
634  &MixtureType::thermoMixtureType::THE,
635  patchi,
636  h,
637  this->p_.boundaryField()[patchi],
638  T0
639  );
640 }
641 
642 
643 template<class BasicThermo, class MixtureType>
646 {
647  return volScalarFieldProperty
648  (
649  "W",
651  &MixtureType::cellThermoMixture,
652  &MixtureType::patchFaceThermoMixture,
654  );
655 }
656 
657 
658 template<class BasicThermo, class MixtureType>
660 (
661  const label patchi
662 ) const
663 {
664  return patchFieldProperty
665  (
666  &MixtureType::patchFaceThermoMixture,
668  patchi
669  );
670 }
671 
672 
673 template<class BasicThermo, class MixtureType>
675 {
676  if (BasicThermo::read())
677  {
678  MixtureType::read(*this);
679  return true;
680  }
681  else
682  {
683  return false;
684  }
685 }
686 
687 
688 // ************************************************************************* //
scalar Ha(const scalar p, const scalar T) const
Definition: EtoHthermo.H:20
scalar Hs(const scalar p, const scalar T) const
Definition: EtoHthermo.H:11
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
scalar Cv(const scalar p, const scalar T) const
Definition: HtoEthermo.H:2
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Generic GeometricBoundaryField class.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
A List with indirect addressing.
Definition: UIndirectList.H:60
A class representing the concept of a uniform field which stores only the single value and providing ...
Definition: UniformField.H:48
Dimension set for the base types.
Definition: dimensionSet.H:122
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
heThermo(const fvMesh &, const word &phaseName)
Construct from mesh.
Definition: heThermo.C:195
virtual const volScalarField & Cv() const
Heat capacity at constant volume [J/kg/K].
Definition: heThermo.H:206
virtual tmp< volScalarField > W() const
Molecular weight [kg/kmol].
Definition: heThermo.C:645
virtual const volScalarField & Cp() const
Heat capacity at constant pressure [J/kg/K].
Definition: heThermo.H:200
virtual tmp< volScalarField > ha() const
Absolute enthalpy [J/kg/K].
Definition: heThermo.C:414
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
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.
virtual ~heThermo()
Destructor.
Definition: heThermo.C:262
virtual tmp< volScalarField > hc() const
Enthalpy of formation [J/kg].
Definition: heThermo.C:487
void heBoundaryCorrection(volScalarField &he)
Correct the enthalpy/internal energy field boundaries.
Definition: heThermo.C:171
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 > gamma() const
Gamma = Cp/Cv [].
Definition: heThermo.C:557
virtual const volScalarField & Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
Definition: heThermo.C:270
virtual tmp< volScalarField > hs() const
Sensible enthalpy [J/kg/K].
Definition: heThermo.C:341
virtual volScalarField & he()
Enthalpy/Internal energy [J/kg].
Definition: heThermo.H:188
virtual tmp< volScalarField > THE(const volScalarField &h, const volScalarField &p, const volScalarField &T0) const
Temperature from enthalpy/internal energy.
Definition: heThermo.C:583
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 bool read()
Read thermophysical properties dictionary.
Definition: heThermo.C:674
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A class for handling words, derived from string.
Definition: word.H:62
volScalarField scalarField(fieldObject, mesh)
label patchi
const cellShapeList & cells
const volScalarField & psi
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
const char *const group
Group name for atomic constants.
const dimensionedScalar h
Planck constant.
static const zero Zero
Definition: zero.H:97
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
const dimensionSet dimEnergy
const dimensionSet dimTemperature
const dimensionSet dimMoles
const dimensionSet dimMass
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const scalarList W(::W(thermo))
Foam::argList args(argc, argv)
volScalarField & p
scalar T0
Definition: createFields.H:22