MulticomponentThermo.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 "MulticomponentThermo.H"
27 #include "fvMesh.H"
28 
29 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
30 
31 template<class BaseThermo>
32 template<class Method, class ... Args>
35 (
36  const word& psiName,
37  const dimensionSet& psiDim,
38  Method psiMethod,
39  const label speciei,
40  const Args& ... args
41 ) const
42 {
43  const typename BaseThermo::mixtureType::thermoType& thermo =
44  this->specieThermo(speciei);
45 
47  (
49  (
50  IOobject::groupName(psiName, this->T_.group()),
51  this->T_.mesh(),
52  psiDim
53  )
54  );
55 
56  volScalarField& psi = tPsi.ref();
57 
58  forAll(psi, celli)
59  {
60  psi[celli] = (thermo.*psiMethod)(args[celli] ...);
61  }
62 
63  volScalarField::Boundary& psiBf = psi.boundaryFieldRef();
64 
65  forAll(psiBf, patchi)
66  {
67  forAll(psiBf[patchi], patchFacei)
68  {
69  psiBf[patchi][patchFacei] =
70  (thermo.*psiMethod)
71  (
72  args.boundaryField()[patchi][patchFacei] ...
73  );
74  }
75  }
76 
77  return tPsi;
78 }
79 
80 
81 template<class BaseThermo>
82 template<class Method, class ... Args>
85 (
86  const word& psiName,
87  const dimensionSet& psiDim,
88  Method psiMethod,
89  const label speciei,
90  const Args& ... args
91 ) const
92 {
93  const typename BaseThermo::mixtureType::thermoType& thermo =
94  this->specieThermo(speciei);
95 
97  (
99  (
100  IOobject::groupName(psiName, this->T_.group()),
101  this->T_.mesh(),
102  psiDim
103  )
104  );
105 
106  volScalarField::Internal& psi = tPsi.ref();
107 
108  forAll(psi, celli)
109  {
110  psi[celli] = (thermo.*psiMethod)(args[celli] ...);
111  }
112 
113  return tPsi;
114 }
115 
116 
117 template<class BaseThermo>
118 template<class Method, class Arg, class ... Args>
121 (
122  Method psiMethod,
123  const label speciei,
124  const Arg& arg,
125  const Args& ... args
126 ) const
127 {
128  const typename BaseThermo::mixtureType::thermoType& thermo =
129  this->specieThermo(speciei);
130 
131  tmp<scalarField> tPsi(new scalarField(arg.size()));
132 
133  scalarField& psi = tPsi.ref();
134 
135  forAll(psi, i)
136  {
137  psi[i] = (thermo.*psiMethod)(arg[i], args[i] ...);
138  }
139 
140  return tPsi;
141 }
142 
143 
144 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
145 
146 template<class BaseThermo>
148 (
149  const fvMesh& mesh,
150  const word& phaseName
151 )
152 :
153  BaseThermo(mesh, phaseName)
154 {}
155 
156 
157 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
158 
159 template<class BaseThermo>
161 {}
162 
163 
164 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
165 
166 template<class BaseThermo>
168 (
169  const label speciei
170 ) const
171 {
172  return this->specieThermo(speciei).W();
173 }
174 
175 
176 template<class BaseThermo>
178 (
179  const label speciei
180 ) const
181 {
182  return
184  (
185  "W",
187  this->specieThermo(speciei).W()
188  );
189 }
190 
191 
192 template<class BaseThermo>
194 (
195  const label speciei,
196  const scalar p,
197  const scalar T
198 ) const
199 {
200  return this->specieThermo(speciei).rho(p, T);
201 }
202 
203 
204 template<class BaseThermo>
207 (
208  const label speciei,
209  const volScalarField& p,
210  const volScalarField& T
211 ) const
212 {
213  return volScalarFieldPropertyi
214  (
215  "rho",
216  dimDensity,
218  speciei,
219  p,
220  T
221  );
222 }
223 
224 
225 template<class BaseThermo>
227 (
228  const label speciei,
229  const scalar p,
230  const scalar T
231 ) const
232 {
233  return this->specieThermo(speciei).Cp(p, T);
234 }
235 
236 
237 template<class BaseThermo>
240 (
241  const label speciei,
242  const volScalarField& p,
243  const volScalarField& T
244 ) const
245 {
246  return volScalarFieldPropertyi
247  (
248  "Cp",
251  speciei,
252  p,
253  T
254  );
255 }
256 
257 
258 template<class BaseThermo>
260 (
261  const label speciei,
262  const scalar p,
263  const scalar T
264 ) const
265 {
266  return this->specieThermo(speciei).he(p, T);
267 }
268 
269 
270 template<class BaseThermo>
272 (
273  const label speciei,
274  const scalarField& p,
275  const scalarField& T
276 ) const
277 {
278  return scalarFieldPropertyi
279  (
281  speciei,
282  p,
283  T
284  );
285 }
286 
287 
288 template<class BaseThermo>
291 (
292  const label speciei,
295 ) const
296 {
297  return volInternalScalarFieldPropertyi
298  (
299  "he",
302  speciei,
303  p,
304  T
305  );
306 }
307 
308 
309 template<class BaseThermo>
312 (
313  const label speciei,
314  const volScalarField& p,
315  const volScalarField& T
316 ) const
317 {
318  return volScalarFieldPropertyi
319  (
320  "he",
323  speciei,
324  p,
325  T
326  );
327 }
328 
329 
330 template<class BaseThermo>
332 (
333  const label speciei,
334  const scalar p,
335  const scalar T
336 ) const
337 {
338  return this->specieThermo(speciei).hs(p, T);
339 }
340 
341 
342 template<class BaseThermo>
344 (
345  const label speciei,
346  const scalarField& p,
347  const scalarField& T
348 ) const
349 {
350  return scalarFieldPropertyi
351  (
353  speciei,
354  p,
355  T
356  );
357 }
358 
359 
360 template<class BaseThermo>
363 (
364  const label speciei,
367 ) const
368 {
369  return volInternalScalarFieldPropertyi
370  (
371  "hs",
374  speciei,
375  p,
376  T
377  );
378 }
379 
380 
381 template<class BaseThermo>
384 (
385  const label speciei,
386  const volScalarField& p,
387  const volScalarField& T
388 ) const
389 {
390  return volScalarFieldPropertyi
391  (
392  "hs",
395  speciei,
396  p,
397  T
398  );
399 }
400 
401 
402 template<class BaseThermo>
404 (
405  const label speciei,
406  const scalar p,
407  const scalar T
408 ) const
409 {
410  return this->specieThermo(speciei).ha(p, T);
411 }
412 
413 
414 template<class BaseThermo>
416 (
417  const label speciei,
418  const scalarField& p,
419  const scalarField& T
420 ) const
421 {
422  return scalarFieldPropertyi
423  (
425  speciei,
426  p,
427  T
428  );
429 }
430 
431 
432 template<class BaseThermo>
435 (
436  const label speciei,
439 ) const
440 {
441  return volInternalScalarFieldPropertyi
442  (
443  "ha",
446  speciei,
447  p,
448  T
449  );
450 }
451 
452 
453 template<class BaseThermo>
456 (
457  const label speciei,
458  const volScalarField& p,
459  const volScalarField& T
460 ) const
461 {
462  return volScalarFieldPropertyi
463  (
464  "ha",
467  speciei,
468  p,
469  T
470  );
471 }
472 
473 
474 template<class BaseThermo>
476 (
477  const label speciei
478 ) const
479 {
480  return
482  (
483  "hf",
485  this->specieThermo(speciei).hf()
486  );
487 }
488 
489 
490 template<class BaseThermo>
492 (
493  const label speciei
494 ) const
495 {
496  return this->specieThermo(speciei).hf();
497 }
498 
499 
500 template<class BaseThermo>
502 (
503  const label speciei,
504  const scalar p,
505  const scalar T
506 ) const
507 {
508  return this->specieThermo(speciei).kappa(p, T);
509 }
510 
511 
512 template<class BaseThermo>
515 (
516  const label speciei,
517  const volScalarField& p,
518  const volScalarField& T
519 ) const
520 {
521  return volScalarFieldPropertyi
522  (
523  "kappa",
526  speciei,
527  p,
528  T
529  );
530 }
531 
532 
533 // ************************************************************************* //
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
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricBoundaryField class.
Generic GeometricField class.
virtual ~MulticomponentThermo()
Destructor.
virtual scalar WiValue(const label speciei) const
Molecular weight [kg/kmol].
virtual scalar kappai(const label speciei, const scalar p, const scalar T) const
Thermal conductivity [W/m/K].
virtual scalar Cpi(const label speciei, const scalar p, const scalar T) const
Heat capacity at constant pressure [J/kg/K].
virtual scalar rhoi(const label speciei, const scalar p, const scalar T) const
Density [kg/m^3].
virtual dimensionedScalar hfi(const label speciei) const
Enthalpy of formation [J/kg].
tmp< volScalarField::Internal > volInternalScalarFieldPropertyi(const word &psiName, const dimensionSet &psiDim, Method psiMethod, const label speciei, const Args &... args) const
Return a volScalarField::Internal of the given property.
tmp< volScalarField > volScalarFieldPropertyi(const word &psiName, const dimensionSet &psiDim, Method psiMethod, const label speciei, const Args &... args) const
Return a volScalarField of the given property.
MulticomponentThermo(const fvMesh &, const word &phaseName)
Construct from mesh and phase name.
virtual scalar hfiValue(const label speciei) const
Enthalpy of formation [J/kg].
virtual dimensionedScalar Wi(const label speciei) const
Molecular weight [kg/kmol].
tmp< scalarField > scalarFieldPropertyi(Method psiMethod, const label speciei, const Arg &arg, const Args &... args) const
Return a scalarField of the given property.
virtual scalar hai(const label speciei, const scalar p, const scalar T) const
Absolute enthalpy [J/kg].
virtual scalar hei(const label speciei, const scalar p, const scalar T) const
Enthalpy/Internal energy [J/kg].
virtual scalar hsi(const label speciei, const scalar p, const scalar T) const
Sensible enthalpy [J/kg].
Dimension set for the base types.
Definition: dimensionSet.H:125
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
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 volScalarField & psi
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
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
scalarList W(const fluidMulticomponentThermo &thermo)
const dimensionSet dimDensity
const dimensionSet dimMoles
const dimensionSet dimMass
const dimensionSet dimThermalConductivity
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
thermo he()
Foam::argList args(argc, argv)
volScalarField & p
fluidMulticomponentThermo & thermo
Definition: createFields.H:31