BasicThermo.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-2025 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 "BasicThermo.H"
29 
30 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
31 
32 template<class MixtureType, class BasicThermoType>
33 template<class Mixture, class Method, class ... Args>
36 (
37  const word& psiName,
38  const dimensionSet& psiDim,
39  Mixture mixture,
40  Method psiMethod,
41  const Args& ... args
42 ) const
43 {
45  (
47  (
48  IOobject::groupName(psiName, this->group()),
49  this->mesh(),
50  psiDim
51  )
52  );
53  volScalarField& psi = tPsi.ref();
54 
55  auto Yslicer = this->Yslicer();
56 
57  forAll(psi, celli)
58  {
59  auto composition = this->cellComposition(Yslicer, celli);
60 
61  psi[celli] =
62  ((this->*mixture)(composition).*psiMethod)(args[celli] ...);
63  }
64 
65  volScalarField::Boundary& psiBf = psi.boundaryFieldRef();
66 
67  forAll(psiBf, patchi)
68  {
69  forAll(psiBf[patchi], patchFacei)
70  {
71  auto composition =
72  this->patchFaceComposition(Yslicer, patchi, patchFacei);
73 
74  psiBf[patchi][patchFacei] =
75  ((this->*mixture)(composition).*psiMethod)
76  (
77  args.boundaryField()[patchi][patchFacei] ...
78  );
79  }
80  }
81 
82  return tPsi;
83 }
84 
85 
86 template<class MixtureType, class BasicThermoType>
87 template<class Mixture, class Method, class ... Args>
90 (
91  const word& psiName,
92  const dimensionSet& psiDim,
93  Mixture mixture,
94  Method psiMethod,
95  const Args& ... args
96 ) const
97 {
99  (
101  (
102  IOobject::groupName(psiName, this->group()),
103  this->mesh(),
104  psiDim
105  )
106  );
107  volScalarField::Internal& psi = tPsi.ref();
108 
109  auto Yslicer = this->Yslicer();
110 
111  forAll(psi, celli)
112  {
113  auto composition = this->cellComposition(Yslicer, celli);
114 
115  psi[celli] =
116  ((this->*mixture)(composition).*psiMethod)(args[celli] ...);
117  }
118 
119  return tPsi;
120 }
121 
122 
123 template<class MixtureType, class BasicThermoType>
124 template<class Mixture, class Method, class ... Args>
127 (
128  Mixture mixture,
129  Method psiMethod,
130  const labelList& cells,
131  const Args& ... args
132 ) const
133 {
134  // Note: Args are fields for the set, not for the mesh as a whole. The
135  // cells list is only used to get the mixture.
136 
137  tmp<scalarField> tPsi(new scalarField(cells.size()));
138  scalarField& psi = tPsi.ref();
139 
140  auto Yslicer = this->Yslicer();
141 
142  forAll(cells, i)
143  {
144  auto composition = this->cellComposition(Yslicer, cells[i]);
145 
146  psi[i] = ((this->*mixture)(composition).*psiMethod)(args[i] ...);
147  }
148 
149  return tPsi;
150 }
151 
152 
153 template<class MixtureType, class BasicThermoType>
154 template<class Mixture, class Method, class ... Args>
157 (
158  Mixture mixture,
159  Method psiMethod,
160  const label patchi,
161  const Args& ... args
162 ) const
163 {
164  tmp<scalarField> tPsi
165  (
166  new scalarField(this->T_.boundaryField()[patchi].size())
167  );
168  scalarField& psi = tPsi.ref();
169 
170  auto Yslicer = this->Yslicer();
171 
172  forAll(psi, patchFacei)
173  {
174  auto composition =
175  this->patchFaceComposition(Yslicer, patchi, patchFacei);
176 
177  psi[patchFacei] =
178  ((this->*mixture)(composition).*psiMethod)(args[patchFacei] ...);
179  }
180 
181  return tPsi;
182 }
183 
184 
185 template<class MixtureType, class BasicThermoType>
186 template<class Mixture, class Method, class ... Args>
189 (
190  const word& psiName,
191  const dimensionSet& psiDim,
192  Mixture mixture,
193  Method psiMethod,
194  const fvSource& model,
195  const volScalarField::Internal& source,
196  const Args& ... args
197 ) const
198 {
200  (
202  (
203  IOobject::groupName(psiName, this->group()),
204  this->mesh(),
205  psiDim
206  )
207  );
208  volScalarField::Internal& psi = tPsi.ref();
209 
210  auto Yslicer = this->Yslicer(model, source);
211 
212  forAll(psi, celli)
213  {
214  auto composition = this->sourceCellComposition(Yslicer, celli);
215 
216  psi[celli] =
217  ((this->*mixture)(composition).*psiMethod)(args[celli] ...);
218  }
219 
220  return tPsi;
221 }
222 
223 
224 template<class MixtureType, class BasicThermoType>
225 template<class Mixture, class Method, class ... Args>
228 (
229  Mixture mixture,
230  Method psiMethod,
231  const fvSource& model,
232  const scalarField& source,
233  const labelUList& cells,
234  const Args& ... args
235 ) const
236 {
237  tmp<scalarField> tPsi(new scalarField(cells.size()));
238  scalarField& psi = tPsi.ref();
239 
240  auto Yslicer = this->Yslicer(model, source, cells);
241 
242  forAll(cells, i)
243  {
244  auto composition =
245  this->sourceCellComposition(Yslicer, i);
246 
247  psi[i] =
248  ((this->*mixture)(composition).*psiMethod)(args[i] ...);
249  }
250 
251  return tPsi;
252 }
253 
254 
255 template<class MixtureType, class BasicThermoType>
258 (
259  const volScalarField& psi,
260  const labelUList& cells
261 )
262 {
264 }
265 
266 
267 template<class MixtureType, class BasicThermoType>
270 (
272  const labelUList&
273 )
274 {
275  return psi.primitiveField();
276 }
277 
278 
279 template<class MixtureType, class BasicThermoType>
282 {
283  volScalarField::Boundary& hBf = h.boundaryFieldRef();
284 
285  forAll(hBf, patchi)
286  {
287  if (isA<gradientEnergyFvPatchScalarField>(hBf[patchi]))
288  {
289  refCast<gradientEnergyFvPatchScalarField>(hBf[patchi]).gradient()
290  = hBf[patchi].fvPatchField::snGrad();
291  }
292  else if (isA<mixedEnergyFvPatchScalarField>(hBf[patchi]))
293  {
294  refCast<mixedEnergyFvPatchScalarField>(hBf[patchi]).refGrad()
295  = hBf[patchi].fvPatchField::snGrad();
296  }
297  }
298 }
299 
300 
301 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
302 
303 template<class MixtureType, class BasicThermoType>
305 (
306  const fvMesh& mesh,
307  const word& phaseName
308 )
309 :
310  physicalProperties(mesh, phaseName),
311  MixtureType(properties()),
312  BasicThermoType
313  (
314  properties(),
315  static_cast<const MixtureType&>(*this),
316  mesh,
317  phaseName
318  ),
319 
320  he_
321  (
322  IOobject
323  (
324  BasicThermoType::phasePropertyName
325  (
326  MixtureType::thermoType::heName(),
327  phaseName
328  ),
329  mesh.time().name(),
330  mesh,
331  IOobject::NO_READ,
332  IOobject::NO_WRITE
333  ),
334  volScalarFieldProperty
335  (
336  "he",
338  &MixtureType::thermoMixture,
339  &MixtureType::thermoMixtureType::he,
340  this->p_,
341  this->T_
342  ),
343  this->heBoundaryTypes(),
344  this->heBoundaryBaseTypes(),
345  this->heSourcesTypes()
346  ),
347 
348  Cp_
349  (
350  IOobject
351  (
352  BasicThermoType::phasePropertyName("Cp", phaseName),
353  mesh.time().name(),
354  mesh
355  ),
356  mesh,
358  ),
359 
360  Cv_
361  (
362  IOobject
363  (
364  BasicThermoType::phasePropertyName("Cv", phaseName),
365  mesh.time().name(),
366  mesh
367  ),
368  mesh,
370  )
371 {
373 }
374 
375 
376 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
377 
378 template<class MixtureType, class BasicThermoType>
380 {}
381 
382 
383 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
384 
385 template<class MixtureType, class BasicThermoType>
388 {
389  return volScalarFieldProperty
390  (
391  "W",
393  &MixtureType::thermoMixture,
395  );
396 }
397 
398 
399 template<class MixtureType, class BasicThermoType>
402 (
403  const label patchi
404 ) const
405 {
406  return patchFieldProperty
407  (
408  &MixtureType::thermoMixture,
410  patchi
411  );
412 }
413 
414 
415 template<class MixtureType, class BasicThermoType>
418 {
419  if (MixtureType::thermoType::enthalpy())
420  {
421  return Cp_;
422  }
423  else
424  {
425  return Cv_;
426  }
427 }
428 
429 
430 template<class MixtureType, class BasicThermoType>
433 (
434  const volScalarField& p,
435  const volScalarField& T
436 ) const
437 {
438  return volScalarFieldProperty
439  (
440  "he",
442  &MixtureType::thermoMixture,
444  p,
445  T
446  );
447 }
448 
449 
450 template<class MixtureType, class BasicThermoType>
453 (
456 ) const
457 {
458  return volInternalScalarFieldProperty
459  (
460  "he",
462  &MixtureType::thermoMixture,
464  p,
465  T
466  );
467 }
468 
469 
470 template<class MixtureType, class BasicThermoType>
473 (
474  const scalarField& T,
475  const labelList& cells
476 ) const
477 {
478  return cellSetProperty
479  (
480  &MixtureType::thermoMixture,
482  cells,
483  cellSetScalarList(this->p_, cells),
484  T
485  );
486 }
487 
488 
489 template<class MixtureType, class BasicThermoType>
492 (
493  const scalarField& T,
494  const label patchi
495 ) const
496 {
497  return patchFieldProperty
498  (
499  &MixtureType::thermoMixture,
501  patchi,
502  this->p_.boundaryField()[patchi],
503  T
504  );
505 }
506 
507 
508 template<class MixtureType, class BasicThermoType>
511 (
513  const fvSource& model,
514  const volScalarField::Internal& source
515 ) const
516 {
517  return fieldSourceProperty
518  (
519  "he",
521  &MixtureType::thermoMixture,
523  model,
524  source,
525  this->p_.internalField(),
526  T
527  );
528 }
529 
530 
531 template<class MixtureType, class BasicThermoType>
534 (
535  const scalarField& T,
536  const fvSource& model,
537  const scalarField& source,
538  const labelUList& cells
539 ) const
540 {
541  return fieldSourceProperty
542  (
543  &MixtureType::thermoMixture,
545  model,
546  source,
547  cells,
548  cellSetScalarList(this->p_, cells),
549  T
550  );
551 }
552 
553 
554 template<class MixtureType, class BasicThermoType>
557 {
558  return volScalarFieldProperty
559  (
560  "hs",
562  &MixtureType::thermoMixture,
564  this->p_,
565  this->T_
566  );
567 }
568 
569 
570 template<class MixtureType, class BasicThermoType>
573 (
574  const volScalarField& p,
575  const volScalarField& T
576 ) const
577 {
578  return volScalarFieldProperty
579  (
580  "hs",
582  &MixtureType::thermoMixture,
584  p,
585  T
586  );
587 }
588 
589 
590 template<class MixtureType, class BasicThermoType>
593 (
596 ) const
597 {
598  return volInternalScalarFieldProperty
599  (
600  "hs",
602  &MixtureType::thermoMixture,
604  p,
605  T
606  );
607 }
608 
609 
610 template<class MixtureType, class BasicThermoType>
613 (
614  const scalarField& T,
615  const labelList& cells
616 ) const
617 {
618  return cellSetProperty
619  (
620  &MixtureType::thermoMixture,
622  cells,
623  cellSetScalarList(this->p_, cells),
624  T
625  );
626 }
627 
628 
629 template<class MixtureType, class BasicThermoType>
632 (
633  const scalarField& T,
634  const label patchi
635 ) const
636 {
637  return patchFieldProperty
638  (
639  &MixtureType::thermoMixture,
641  patchi,
642  this->p_.boundaryField()[patchi],
643  T
644  );
645 }
646 
647 
648 template<class MixtureType, class BasicThermoType>
651 {
652  return volScalarFieldProperty
653  (
654  "ha",
656  &MixtureType::thermoMixture,
658  this->p_,
659  this->T_
660  );
661 }
662 
663 
664 template<class MixtureType, class BasicThermoType>
667 (
668  const volScalarField& p,
669  const volScalarField& T
670 ) const
671 {
672  return volScalarFieldProperty
673  (
674  "ha",
676  &MixtureType::thermoMixture,
678  p,
679  T
680  );
681 }
682 
683 
684 template<class MixtureType, class BasicThermoType>
687 (
690 ) const
691 {
692  return volInternalScalarFieldProperty
693  (
694  "ha",
696  &MixtureType::thermoMixture,
698  p,
699  T
700  );
701 }
702 
703 
704 template<class MixtureType, class BasicThermoType>
707 (
708  const scalarField& T,
709  const labelList& cells
710 ) const
711 {
712  return cellSetProperty
713  (
714  &MixtureType::thermoMixture,
716  cells,
717  cellSetScalarList(this->p_, cells),
718  T
719  );
720 }
721 
722 
723 template<class MixtureType, class BasicThermoType>
726 (
727  const scalarField& T,
728  const label patchi
729 ) const
730 {
731  return patchFieldProperty
732  (
733  &MixtureType::thermoMixture,
735  patchi,
736  this->p_.boundaryField()[patchi],
737  T
738  );
739 }
740 
741 
742 template<class MixtureType, class BasicThermoType>
745 (
746  const scalarField& T,
747  const label patchi
748 ) const
749 {
750  return patchFieldProperty
751  (
752  &MixtureType::thermoMixture,
754  patchi,
755  this->p_.boundaryField()[patchi],
756  T
757  );
758 }
759 
760 
761 template<class MixtureType, class BasicThermoType>
764 (
765  const scalarField& T,
766  const label patchi
767 ) const
768 {
769  return patchFieldProperty
770  (
771  &MixtureType::thermoMixture,
773  patchi,
774  this->p_.boundaryField()[patchi],
775  T
776  );
777 }
778 
779 
780 template<class MixtureType, class BasicThermoType>
783 (
784  const scalarField& T,
785  const label patchi
786 ) const
787 {
788  if (MixtureType::thermoType::enthalpy())
789  {
790  return Cp(T, patchi);
791  }
792  else
793  {
794  return Cv(T, patchi);
795  }
796 }
797 
798 
799 template<class MixtureType, class BasicThermoType>
802 (
803  const volScalarField& h,
804  const volScalarField& p,
805  const volScalarField& T0
806 ) const
807 {
808  return volScalarFieldProperty
809  (
810  "T",
812  &MixtureType::thermoMixture,
813  &MixtureType::thermoMixtureType::The,
814  h,
815  p,
816  T0
817  );
818 }
819 
820 
821 template<class MixtureType, class BasicThermoType>
824 (
825  const scalarField& h,
826  const scalarField& T0,
827  const labelList& cells
828 ) const
829 {
830  return cellSetProperty
831  (
832  &MixtureType::thermoMixture,
833  &MixtureType::thermoMixtureType::The,
834  cells,
835  h,
836  cellSetScalarList(this->p_, cells),
837  T0
838  );
839 }
840 
841 
842 template<class MixtureType, class BasicThermoType>
845 (
846  const scalarField& h,
847  const scalarField& T0,
848  const label patchi
849 ) const
850 {
851  return patchFieldProperty
852  (
853  &MixtureType::thermoMixture,
854  &MixtureType::thermoMixtureType::The,
855  patchi,
856  h,
857  this->p_.boundaryField()[patchi],
858  T0
859  );
860 }
861 
862 
863 template<class MixtureType, class BasicThermoType>
865 {
867  {
868  MixtureType::read(*this);
869  BasicThermoType::read(*this);
870  return true;
871  }
872  else
873  {
874  return false;
875  }
876 }
877 
878 
879 // ************************************************************************* //
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 ha(const scalar p, const scalar T) const
Definition: EtoHthermo.H:20
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:433
BasicThermo(const fvMesh &, const word &phaseName)
Construct from mesh and phase name.
Definition: BasicThermo.C:305
virtual const volScalarField & Cv() const
Heat capacity at constant volume [J/kg/K].
Definition: BasicThermo.H:262
virtual tmp< volScalarField > W() const
Molecular weight [kg/kmol].
Definition: BasicThermo.C:387
virtual const volScalarField & Cp() const
Heat capacity at constant pressure [J/kg/K].
Definition: BasicThermo.H:256
tmp< volScalarField > volScalarFieldProperty(const word &psiName, const dimensionSet &psiDim, Mixture mixture, Method psiMethod, const Args &... args) const
Return a volScalarField of the given property.
tmp< scalarField > cellSetProperty(Mixture mixture, Method psiMethod, const labelList &cells, const Args &... args) const
Return a scalarField of the given property on a cell set.
static UIndirectList< scalar > cellSetScalarList(const volScalarField &psi, const labelUList &cells)
Return an indirect list of a field for the given set of cells.
Definition: BasicThermo.C:258
virtual tmp< volScalarField > ha() const
Absolute enthalpy [J/kg/K].
Definition: BasicThermo.C:650
virtual const volScalarField & he() const
Enthalpy/Internal energy [J/kg].
Definition: BasicThermo.H:243
volScalarField he_
Energy field.
Definition: BasicThermo.H:72
tmp< scalarField > patchFieldProperty(Mixture mixture, Method psiMethod, const label patchi, const Args &... args) const
Return a scalarField of the given property on a patch.
tmp< volScalarField::Internal > fieldSourceProperty(const word &psiName, const dimensionSet &psiDim, Mixture mixture, Method psiMethod, const fvSource &model, const volScalarField::Internal &source, const Args &... args) const
Return a scalarField of the given property for a source.
void heBoundaryCorrection(volScalarField &he)
Correct the enthalpy/internal energy field boundaries.
Definition: BasicThermo.C:281
virtual const volScalarField & Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
Definition: BasicThermo.C:417
virtual ~BasicThermo()
Destructor.
Definition: BasicThermo.C:379
virtual tmp< volScalarField > hs() const
Sensible enthalpy [J/kg/K].
Definition: BasicThermo.C:556
virtual tmp< volScalarField > The(const volScalarField &h, const volScalarField &p, const volScalarField &T0) const
Temperature from enthalpy/internal energy.
Definition: BasicThermo.C:802
tmp< volScalarField::Internal > volInternalScalarFieldProperty(const word &psiName, const dimensionSet &psiDim, Mixture mixture, Method psiMethod, const Args &... args) const
Return a volScalarField::Internal of the given property.
virtual bool read()
Read thermophysical properties dictionary.
Definition: BasicThermo.C:864
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
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:125
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
Base class for finite volume sources.
Definition: fvSource.H:52
A base class for physical properties.
virtual bool read()
Read physicalProperties dictionary.
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:197
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
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
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimTemperature
const dimensionSet dimMoles
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimMass
const scalarList W(::W(thermo))
thermo he()
Foam::argList args(argc, argv)
volScalarField & p
scalar T0
Definition: createFields.H:22