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-2023 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 
54  volScalarField& psi = tPsi.ref();
55 
56  auto Yslicer = this->Yslicer();
57 
58  forAll(psi, celli)
59  {
60  auto composition = this->cellComposition(Yslicer, celli);
61 
62  psi[celli] =
63  ((this->*mixture)(composition).*psiMethod)(args[celli] ...);
64  }
65 
66  volScalarField::Boundary& psiBf = psi.boundaryFieldRef();
67 
68  forAll(psiBf, patchi)
69  {
70  forAll(psiBf[patchi], patchFacei)
71  {
72  auto composition =
73  this->patchFaceComposition(Yslicer, patchi, patchFacei);
74 
75  psiBf[patchi][patchFacei] =
76  ((this->*mixture)(composition).*psiMethod)
77  (
78  args.boundaryField()[patchi][patchFacei] ...
79  );
80  }
81  }
82 
83  return tPsi;
84 }
85 
86 
87 template<class MixtureType, class BasicThermoType>
88 template<class Mixture, class Method, class ... Args>
91 (
92  Mixture mixture,
93  Method psiMethod,
94  const labelList& cells,
95  const Args& ... args
96 ) const
97 {
98  // Note: Args are fields for the set, not for the mesh as a whole. The
99  // cells list is only used to get the mixture.
100 
101  tmp<scalarField> tPsi(new scalarField(cells.size()));
102  scalarField& psi = tPsi.ref();
103 
104  auto Yslicer = this->Yslicer();
105 
106  forAll(cells, i)
107  {
108  auto composition = this->cellComposition(Yslicer, cells[i]);
109 
110  psi[i] = ((this->*mixture)(composition).*psiMethod)(args[i] ...);
111  }
112 
113  return tPsi;
114 }
115 
116 
117 template<class MixtureType, class BasicThermoType>
118 template<class Mixture, class Method, class ... Args>
121 (
122  Mixture mixture,
123  Method psiMethod,
124  const label patchi,
125  const Args& ... args
126 ) const
127 {
128  tmp<scalarField> tPsi
129  (
130  new scalarField(this->T_.boundaryField()[patchi].size())
131  );
132  scalarField& psi = tPsi.ref();
133 
134  auto Yslicer = this->Yslicer();
135 
136  forAll(psi, patchFacei)
137  {
138  auto composition =
139  this->patchFaceComposition(Yslicer, patchi, patchFacei);
140 
141  psi[patchFacei] =
142  ((this->*mixture)(composition).*psiMethod)(args[patchFacei] ...);
143  }
144 
145  return tPsi;
146 }
147 
148 
149 template<class MixtureType, class BasicThermoType>
150 template<class Mixture, class Method, class ... Args>
153 (
154  Mixture mixture,
155  Method psiMethod,
156  const fvSource& source,
157  const Args& ... args
158 ) const
159 {
160  const labelUList cells = source.cells();
161 
162  tmp<scalarField> tPsi(new scalarField(cells.size()));
163  scalarField& psi = tPsi.ref();
164 
165  auto Yslicer = this->Yslicer(source);
166 
167  forAll(cells, i)
168  {
169  auto composition =
170  this->sourceCellComposition(Yslicer, i);
171 
172  psi[i] =
173  ((this->*mixture)(composition).*psiMethod)(args[i] ...);
174  }
175 
176  return tPsi;
177 }
178 
179 
180 template<class MixtureType, class BasicThermoType>
183 (
184  const volScalarField& psi,
185  const labelUList& cells
186 )
187 {
189 }
190 
191 
192 template<class MixtureType, class BasicThermoType>
195 (
197  const labelUList&
198 )
199 {
200  return psi.primitiveField();
201 }
202 
203 
204 template<class MixtureType, class BasicThermoType>
207 {
208  volScalarField::Boundary& hBf = h.boundaryFieldRef();
209 
210  forAll(hBf, patchi)
211  {
212  if (isA<gradientEnergyFvPatchScalarField>(hBf[patchi]))
213  {
214  refCast<gradientEnergyFvPatchScalarField>(hBf[patchi]).gradient()
215  = hBf[patchi].fvPatchField::snGrad();
216  }
217  else if (isA<mixedEnergyFvPatchScalarField>(hBf[patchi]))
218  {
219  refCast<mixedEnergyFvPatchScalarField>(hBf[patchi]).refGrad()
220  = hBf[patchi].fvPatchField::snGrad();
221  }
222  }
223 }
224 
225 
226 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
227 
228 template<class MixtureType, class BasicThermoType>
230 (
231  const fvMesh& mesh,
232  const word& phaseName
233 )
234 :
235  physicalProperties(mesh, phaseName),
236  MixtureType(properties()),
237  BasicThermoType(properties(), mixture(), mesh, phaseName),
238 
239  he_
240  (
241  IOobject
242  (
243  BasicThermoType::phasePropertyName
244  (
245  MixtureType::thermoType::heName(),
246  phaseName
247  ),
248  mesh.time().name(),
249  mesh,
250  IOobject::NO_READ,
251  IOobject::NO_WRITE
252  ),
253  volScalarFieldProperty
254  (
255  "he",
257  &MixtureType::thermoMixture,
258  &MixtureType::thermoMixtureType::he,
259  this->p_,
260  this->T_
261  ),
262  this->heBoundaryTypes(),
263  this->heBoundaryBaseTypes(),
264  this->heSourcesTypes()
265  ),
266 
267  Cp_
268  (
269  IOobject
270  (
271  BasicThermoType::phasePropertyName("Cp", phaseName),
272  mesh.time().name(),
273  mesh
274  ),
275  mesh,
277  ),
278 
279  Cv_
280  (
281  IOobject
282  (
283  BasicThermoType::phasePropertyName("Cv", phaseName),
284  mesh.time().name(),
285  mesh
286  ),
287  mesh,
289  )
290 {
292 }
293 
294 
295 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
296 
297 template<class MixtureType, class BasicThermoType>
299 {}
300 
301 
302 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
303 
304 template<class MixtureType, class BasicThermoType>
307 {
308  return volScalarFieldProperty
309  (
310  "W",
312  &MixtureType::thermoMixture,
314  );
315 }
316 
317 
318 template<class MixtureType, class BasicThermoType>
321 (
322  const label patchi
323 ) const
324 {
325  return patchFieldProperty
326  (
327  &MixtureType::thermoMixture,
329  patchi
330  );
331 }
332 
333 
334 template<class MixtureType, class BasicThermoType>
337 {
338  if (MixtureType::thermoType::enthalpy())
339  {
340  return Cp_;
341  }
342  else
343  {
344  return Cv_;
345  }
346 }
347 
348 
349 template<class MixtureType, class BasicThermoType>
352 (
353  const volScalarField& p,
354  const volScalarField& T
355 ) const
356 {
357  return volScalarFieldProperty
358  (
359  "he",
361  &MixtureType::thermoMixture,
363  p,
364  T
365  );
366 }
367 
368 
369 template<class MixtureType, class BasicThermoType>
372 (
373  const scalarField& T,
374  const labelList& cells
375 ) const
376 {
377  return cellSetProperty
378  (
379  &MixtureType::thermoMixture,
381  cells,
382  cellSetScalarList(this->p_, cells),
383  T
384  );
385 }
386 
387 
388 template<class MixtureType, class BasicThermoType>
391 (
392  const scalarField& T,
393  const label patchi
394 ) const
395 {
396  return patchFieldProperty
397  (
398  &MixtureType::thermoMixture,
400  patchi,
401  this->p_.boundaryField()[patchi],
402  T
403  );
404 }
405 
406 
407 template<class MixtureType, class BasicThermoType>
410 (
411  const scalarField& T,
412  const fvSource& source
413 ) const
414 {
415  return fieldSourceProperty
416  (
417  &MixtureType::thermoMixture,
419  source,
420  cellSetScalarList(this->p_, source.cells()),
421  T
422  );
423 }
424 
425 
426 template<class MixtureType, class BasicThermoType>
429 {
430  return volScalarFieldProperty
431  (
432  "hs",
434  &MixtureType::thermoMixture,
436  this->p_,
437  this->T_
438  );
439 }
440 
441 
442 template<class MixtureType, class BasicThermoType>
445 (
446  const volScalarField& p,
447  const volScalarField& T
448 ) const
449 {
450  return volScalarFieldProperty
451  (
452  "hs",
454  &MixtureType::thermoMixture,
456  p,
457  T
458  );
459 }
460 
461 
462 template<class MixtureType, class BasicThermoType>
465 (
466  const scalarField& T,
467  const labelList& cells
468 ) const
469 {
470  return cellSetProperty
471  (
472  &MixtureType::thermoMixture,
474  cells,
475  cellSetScalarList(this->p_, cells),
476  T
477  );
478 }
479 
480 
481 template<class MixtureType, class BasicThermoType>
484 (
485  const scalarField& T,
486  const label patchi
487 ) const
488 {
489  return patchFieldProperty
490  (
491  &MixtureType::thermoMixture,
493  patchi,
494  this->p_.boundaryField()[patchi],
495  T
496  );
497 }
498 
499 
500 template<class MixtureType, class BasicThermoType>
503 {
504  return volScalarFieldProperty
505  (
506  "ha",
508  &MixtureType::thermoMixture,
510  this->p_,
511  this->T_
512  );
513 }
514 
515 
516 template<class MixtureType, class BasicThermoType>
519 (
520  const volScalarField& p,
521  const volScalarField& T
522 ) const
523 {
524  return volScalarFieldProperty
525  (
526  "ha",
528  &MixtureType::thermoMixture,
530  p,
531  T
532  );
533 }
534 
535 
536 template<class MixtureType, class BasicThermoType>
539 (
540  const scalarField& T,
541  const labelList& cells
542 ) const
543 {
544  return cellSetProperty
545  (
546  &MixtureType::thermoMixture,
548  cells,
549  cellSetScalarList(this->p_, cells),
550  T
551  );
552 }
553 
554 
555 template<class MixtureType, class BasicThermoType>
558 (
559  const scalarField& T,
560  const label patchi
561 ) const
562 {
563  return patchFieldProperty
564  (
565  &MixtureType::thermoMixture,
567  patchi,
568  this->p_.boundaryField()[patchi],
569  T
570  );
571 }
572 
573 
574 template<class MixtureType, class BasicThermoType>
577 (
578  const scalarField& T,
579  const label patchi
580 ) const
581 {
582  return patchFieldProperty
583  (
584  &MixtureType::thermoMixture,
586  patchi,
587  this->p_.boundaryField()[patchi],
588  T
589  );
590 }
591 
592 
593 template<class MixtureType, class BasicThermoType>
596 (
597  const scalarField& T,
598  const label patchi
599 ) const
600 {
601  return patchFieldProperty
602  (
603  &MixtureType::thermoMixture,
605  patchi,
606  this->p_.boundaryField()[patchi],
607  T
608  );
609 }
610 
611 
612 template<class MixtureType, class BasicThermoType>
615 (
616  const scalarField& T,
617  const label patchi
618 ) const
619 {
620  if (MixtureType::thermoType::enthalpy())
621  {
622  return Cp(T, patchi);
623  }
624  else
625  {
626  return Cv(T, patchi);
627  }
628 }
629 
630 
631 template<class MixtureType, class BasicThermoType>
634 (
635  const volScalarField& h,
636  const volScalarField& p,
637  const volScalarField& T0
638 ) const
639 {
640  return volScalarFieldProperty
641  (
642  "T",
644  &MixtureType::thermoMixture,
645  &MixtureType::thermoMixtureType::The,
646  h,
647  p,
648  T0
649  );
650 }
651 
652 
653 template<class MixtureType, class BasicThermoType>
656 (
657  const scalarField& h,
658  const scalarField& T0,
659  const labelList& cells
660 ) const
661 {
662  return cellSetProperty
663  (
664  &MixtureType::thermoMixture,
665  &MixtureType::thermoMixtureType::The,
666  cells,
667  h,
668  cellSetScalarList(this->p_, cells),
669  T0
670  );
671 }
672 
673 
674 template<class MixtureType, class BasicThermoType>
677 (
678  const scalarField& h,
679  const scalarField& T0,
680  const label patchi
681 ) const
682 {
683  return patchFieldProperty
684  (
685  &MixtureType::thermoMixture,
686  &MixtureType::thermoMixtureType::The,
687  patchi,
688  h,
689  this->p_.boundaryField()[patchi],
690  T0
691  );
692 }
693 
694 
695 template<class MixtureType, class BasicThermoType>
697 {
699  {
700  MixtureType::read(*this);
701  BasicThermoType::read(*this);
702  return true;
703  }
704  else
705  {
706  return false;
707  }
708 }
709 
710 
711 // ************************************************************************* //
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:434
BasicThermo(const fvMesh &, const word &phaseName)
Construct from mesh and phase name.
Definition: BasicThermo.C:230
virtual const volScalarField & Cv() const
Heat capacity at constant volume [J/kg/K].
Definition: BasicThermo.H:242
virtual tmp< volScalarField > W() const
Molecular weight [kg/kmol].
Definition: BasicThermo.C:306
virtual const volScalarField & Cp() const
Heat capacity at constant pressure [J/kg/K].
Definition: BasicThermo.H:236
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 > fieldSourceProperty(Mixture mixture, Method psiMethod, const fvSource &source, const Args &... args) const
Return a scalarField of the given property for a source.
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:183
virtual tmp< volScalarField > ha() const
Absolute enthalpy [J/kg/K].
Definition: BasicThermo.C:502
virtual const volScalarField & he() const
Enthalpy/Internal energy [J/kg].
Definition: BasicThermo.H:223
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.
void heBoundaryCorrection(volScalarField &he)
Correct the enthalpy/internal energy field boundaries.
Definition: BasicThermo.C:206
virtual const volScalarField & Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
Definition: BasicThermo.C:336
virtual ~BasicThermo()
Destructor.
Definition: BasicThermo.C:298
virtual tmp< volScalarField > hs() const
Sensible enthalpy [J/kg/K].
Definition: BasicThermo.C:428
virtual tmp< volScalarField > The(const volScalarField &h, const volScalarField &p, const volScalarField &T0) const
Temperature from enthalpy/internal energy.
Definition: BasicThermo.C:634
virtual bool read()
Read thermophysical properties dictionary.
Definition: BasicThermo.C:696
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:99
Base class for finite volume sources.
Definition: fvSource.H:52
virtual labelUList cells() const =0
Return the cells that the source applies to.
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: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
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimTemperature
const dimensionSet dimMoles
const dimensionSet dimMass
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const scalarList W(::W(thermo))
thermo he()
Foam::argList args(argc, argv)
volScalarField & p
scalar T0
Definition: createFields.H:22