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-2024 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 Arg, class ... Args>
85 (
86  Method psiMethod,
87  const label speciei,
88  const Arg& arg,
89  const Args& ... args
90 ) const
91 {
92  const typename BaseThermo::mixtureType::thermoType& thermo =
93  this->specieThermo(speciei);
94 
95  tmp<scalarField> tPsi(new scalarField(arg.size()));
96 
97  scalarField& psi = tPsi.ref();
98 
99  forAll(psi, i)
100  {
101  psi[i] = (thermo.*psiMethod)(arg[i], args[i] ...);
102  }
103 
104  return tPsi;
105 }
106 
107 
108 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
109 
110 template<class BaseThermo>
112 (
113  const fvMesh& mesh,
114  const word& phaseName
115 )
116 :
117  BaseThermo(mesh, phaseName)
118 {}
119 
120 
121 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
122 
123 template<class BaseThermo>
125 {}
126 
127 
128 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
129 
130 template<class BaseThermo>
132 (
133  const label speciei
134 ) const
135 {
136  return this->specieThermo(speciei).W();
137 }
138 
139 
140 template<class BaseThermo>
142 (
143  const label speciei
144 ) const
145 {
146  return
148  (
149  "W",
151  this->specieThermo(speciei).W()
152  );
153 }
154 
155 
156 template<class BaseThermo>
158 (
159  const label speciei,
160  const scalar p,
161  const scalar T
162 ) const
163 {
164  return this->specieThermo(speciei).rho(p, T);
165 }
166 
167 
168 template<class BaseThermo>
171 (
172  const label speciei,
173  const volScalarField& p,
174  const volScalarField& T
175 ) const
176 {
177  return volScalarFieldPropertyi
178  (
179  "rho",
180  dimDensity,
182  speciei,
183  p,
184  T
185  );
186 }
187 
188 
189 template<class BaseThermo>
191 (
192  const label speciei,
193  const scalar p,
194  const scalar T
195 ) const
196 {
197  return this->specieThermo(speciei).Cp(p, T);
198 }
199 
200 
201 template<class BaseThermo>
204 (
205  const label speciei,
206  const volScalarField& p,
207  const volScalarField& T
208 ) const
209 {
210  return volScalarFieldPropertyi
211  (
212  "Cp",
215  speciei,
216  p,
217  T
218  );
219 }
220 
221 
222 template<class BaseThermo>
224 (
225  const label speciei,
226  const scalar p,
227  const scalar T
228 ) const
229 {
230  return this->specieThermo(speciei).he(p, T);
231 }
232 
233 
234 template<class BaseThermo>
236 (
237  const label speciei,
238  const scalarField& p,
239  const scalarField& T
240 ) const
241 {
242  return scalarFieldPropertyi
243  (
245  speciei,
246  p,
247  T
248  );
249 }
250 
251 
252 template<class BaseThermo>
255 (
256  const label speciei,
257  const volScalarField& p,
258  const volScalarField& T
259 ) const
260 {
261  return volScalarFieldPropertyi
262  (
263  "he",
266  speciei,
267  p,
268  T
269  );
270 }
271 
272 
273 template<class BaseThermo>
275 (
276  const label speciei,
277  const scalar p,
278  const scalar T
279 ) const
280 {
281  return this->specieThermo(speciei).hs(p, T);
282 }
283 
284 
285 template<class BaseThermo>
287 (
288  const label speciei,
289  const scalarField& p,
290  const scalarField& T
291 ) const
292 {
293  return scalarFieldPropertyi
294  (
296  speciei,
297  p,
298  T
299  );
300 }
301 
302 
303 template<class BaseThermo>
306 (
307  const label speciei,
308  const volScalarField& p,
309  const volScalarField& T
310 ) const
311 {
312  return volScalarFieldPropertyi
313  (
314  "hs",
317  speciei,
318  p,
319  T
320  );
321 }
322 
323 
324 template<class BaseThermo>
326 (
327  const label speciei,
328  const scalar p,
329  const scalar T
330 ) const
331 {
332  return this->specieThermo(speciei).ha(p, T);
333 }
334 
335 
336 template<class BaseThermo>
338 (
339  const label speciei,
340  const scalarField& p,
341  const scalarField& T
342 ) const
343 {
344  return scalarFieldPropertyi
345  (
347  speciei,
348  p,
349  T
350  );
351 }
352 
353 
354 template<class BaseThermo>
357 (
358  const label speciei,
359  const volScalarField& p,
360  const volScalarField& T
361 ) const
362 {
363  return volScalarFieldPropertyi
364  (
365  "ha",
368  speciei,
369  p,
370  T
371  );
372 }
373 
374 
375 template<class BaseThermo>
377 (
378  const label speciei
379 ) const
380 {
381  return
383  (
384  "hf",
386  this->specieThermo(speciei).hf()
387  );
388 }
389 
390 
391 template<class BaseThermo>
393 (
394  const label speciei
395 ) const
396 {
397  return this->specieThermo(speciei).hf();
398 }
399 
400 
401 template<class BaseThermo>
403 (
404  const label speciei,
405  const scalar p,
406  const scalar T
407 ) const
408 {
409  return this->specieThermo(speciei).kappa(p, T);
410 }
411 
412 
413 template<class BaseThermo>
416 (
417  const label speciei,
418  const volScalarField& p,
419  const volScalarField& T
420 ) const
421 {
422  return volScalarFieldPropertyi
423  (
424  "kappa",
427  speciei,
428  p,
429  T
430  );
431 }
432 
433 
434 // ************************************************************************* //
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:434
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 > 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:99
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 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
const dimensionSet dimTemperature
scalarList W(const fluidMulticomponentThermo &thermo)
const dimensionSet dimDensity
const dimensionSet dimMoles
const dimensionSet dimMass
const dimensionSet dimThermalConductivity
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
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