PsiuMulticomponentThermo.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 
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 template<class BaseThermo>
33 {
34  const scalarField& hCells = this->he_;
35  const scalarField& heuCells = this->heu_;
36  const scalarField& pCells = this->p_;
37 
38  scalarField& TCells = this->T_.primitiveFieldRef();
39  scalarField& TuCells = this->Tu_.primitiveFieldRef();
40  scalarField& CpCells = this->Cp_.primitiveFieldRef();
41  scalarField& CvCells = this->Cv_.primitiveFieldRef();
42  scalarField& psiCells = this->psi_.primitiveFieldRef();
43  scalarField& muCells = this->mu_.primitiveFieldRef();
44  scalarField& kappaCells = this->kappa_.primitiveFieldRef();
45 
46  auto Yslicer = this->Yslicer();
47 
48  forAll(TCells, celli)
49  {
50  auto composition = this->cellComposition(Yslicer, celli);
51 
52  const typename BaseThermo::mixtureType::thermoMixtureType&
53  thermoMixture = this->thermoMixture(composition);
54 
55  const typename BaseThermo::mixtureType::transportMixtureType&
56  transportMixture =
57  this->transportMixture(composition, thermoMixture);
58 
59  TCells[celli] = thermoMixture.The
60  (
61  hCells[celli],
62  pCells[celli],
63  TCells[celli]
64  );
65 
66  CpCells[celli] = thermoMixture.Cp(pCells[celli], TCells[celli]);
67  CvCells[celli] = thermoMixture.Cv(pCells[celli], TCells[celli]);
68  psiCells[celli] = thermoMixture.psi(pCells[celli], TCells[celli]);
69 
70  muCells[celli] = transportMixture.mu(pCells[celli], TCells[celli]);
71  kappaCells[celli] =
72  transportMixture.kappa(pCells[celli], TCells[celli]);
73 
74  TuCells[celli] = this->reactants(composition).The
75  (
76  heuCells[celli],
77  pCells[celli],
78  TuCells[celli]
79  );
80  }
81 
82  volScalarField::Boundary& pBf = this->p_.boundaryFieldRef();
83  volScalarField::Boundary& TBf = this->T_.boundaryFieldRef();
84  volScalarField::Boundary& TuBf = this->Tu_.boundaryFieldRef();
85  volScalarField::Boundary& CpBf = this->Cp_.boundaryFieldRef();
86  volScalarField::Boundary& CvBf = this->Cv_.boundaryFieldRef();
87  volScalarField::Boundary& psiBf = this->psi_.boundaryFieldRef();
88  volScalarField::Boundary& heBf = this->he().boundaryFieldRef();
89  volScalarField::Boundary& heuBf = this->heu().boundaryFieldRef();
90  volScalarField::Boundary& muBf = this->mu_.boundaryFieldRef();
91  volScalarField::Boundary& kappaBf = this->kappa_.boundaryFieldRef();
92 
93  forAll(TBf, patchi)
94  {
95  fvPatchScalarField& pPf = pBf[patchi];
96  fvPatchScalarField& TPf = TBf[patchi];
97  fvPatchScalarField& TuPf = TuBf[patchi];
98  fvPatchScalarField& CpPf = CpBf[patchi];
99  fvPatchScalarField& CvPf = CvBf[patchi];
100  fvPatchScalarField& psiPf = psiBf[patchi];
101  fvPatchScalarField& hePf = heBf[patchi];
102  fvPatchScalarField& heuPf = heuBf[patchi];
103  fvPatchScalarField& muPf = muBf[patchi];
104  fvPatchScalarField& kappaPf = kappaBf[patchi];
105 
106  if (TPf.fixesValue())
107  {
108  forAll(TPf, facei)
109  {
110  auto composition =
111  this->patchFaceComposition(Yslicer, patchi, facei);
112 
113  const typename BaseThermo::mixtureType::thermoMixtureType&
114  thermoMixture = this->thermoMixture(composition);
115 
116  const typename BaseThermo::mixtureType::transportMixtureType&
117  transportMixture =
118  this->transportMixture(composition, thermoMixture);
119 
120  hePf[facei] = thermoMixture.he(pPf[facei], TPf[facei]);
121 
122  CpPf[facei] = thermoMixture.Cp(pPf[facei], TPf[facei]);
123  CvPf[facei] = thermoMixture.Cv(pPf[facei], TPf[facei]);
124  psiPf[facei] = thermoMixture.psi(pPf[facei], TPf[facei]);
125  muPf[facei] = transportMixture.mu(pPf[facei], TPf[facei]);
126  kappaPf[facei] = transportMixture.kappa(pPf[facei], TPf[facei]);
127  }
128  }
129  else
130  {
131  forAll(TPf, facei)
132  {
133  auto composition =
134  this->patchFaceComposition(Yslicer, patchi, facei);
135 
136  const typename BaseThermo::mixtureType::thermoMixtureType&
137  thermoMixture = this->thermoMixture(composition);
138 
139  const typename BaseThermo::mixtureType::transportMixtureType&
140  transportMixture =
141  this->transportMixture(composition, thermoMixture);
142 
143  TPf[facei] =
144  thermoMixture.The(hePf[facei], pPf[facei], TPf[facei]);
145 
146  CpPf[facei] = thermoMixture.Cp(pPf[facei], TPf[facei]);
147  CvPf[facei] = thermoMixture.Cv(pPf[facei], TPf[facei]);
148  psiPf[facei] = thermoMixture.psi(pPf[facei], TPf[facei]);
149  muPf[facei] = transportMixture.mu(pPf[facei], TPf[facei]);
150  kappaPf[facei] = transportMixture.kappa(pPf[facei], TPf[facei]);
151 
152  TuPf[facei] =
153  this->reactants(composition)
154  .The(heuPf[facei], pPf[facei], TuPf[facei]);
155  }
156  }
157  }
158 }
159 
160 
161 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
162 
163 template<class BaseThermo>
165 (
166  const fvMesh& mesh,
167  const word& phaseName
168 )
169 :
170  BaseThermo(mesh, phaseName),
171  Tu_
172  (
173  IOobject
174  (
175  "Tu",
176  mesh.time().name(),
177  mesh,
178  IOobject::MUST_READ,
179  IOobject::AUTO_WRITE
180  ),
181  mesh
182  ),
183  heu_
184  (
185  IOobject
186  (
187  BaseThermo::mixtureType::thermoType::heName() + 'u',
188  mesh.time().name(),
189  mesh,
190  IOobject::NO_READ,
191  IOobject::NO_WRITE
192  ),
193  this->volScalarFieldProperty
194  (
195  BaseThermo::mixtureType::thermoType::heName() + 'u',
197  &BaseThermo::mixtureType::reactants,
198  &BaseThermo::mixtureType::thermoMixtureType::he,
199  this->p_,
200  this->Tu_
201  ),
202  this->heuBoundaryTypes()
203  )
204 {
205  this->heuBoundaryCorrection(this->heu_);
206 
207  calculate();
208 
209  this->psi_.oldTime(); // Switch on saving old time
210 }
211 
212 
213 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
214 
215 template<class BaseThermo>
217 {}
218 
219 
220 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
221 
222 template<class BaseThermo>
224 {
225  if (BaseThermo::debug)
226  {
227  InfoInFunction << endl;
228  }
229 
230  // force the saving of the old-time values
231  this->psi_.oldTime();
232 
233  calculate();
234 
235  if (BaseThermo::debug)
236  {
237  Info<< " Finished" << endl;
238  }
239 }
240 
241 
242 template<class BaseThermo>
245 {
246  tmp<volScalarField> tfres
247  (
249  (
250  "fres",
251  this->mesh(),
252  dimless
253  )
254  );
255 
256  auto Yslicer = this->Yslicer();
257 
258  scalarField& fresCells = tfres.ref().primitiveFieldRef();
259 
260  forAll(fresCells, celli)
261  {
262  fresCells[celli] = BaseThermo::mixtureType::fres
263  (
264  this->cellComposition(Yslicer, celli)
265  );
266  }
267 
268  volScalarField::Boundary& fresBf = tfres.ref().boundaryFieldRef();
269 
270  forAll(fresBf, patchi)
271  {
272  fvPatchScalarField& fresPf = fresBf[patchi];
273 
274  forAll(fresPf, facei)
275  {
276  fresPf[facei] = BaseThermo::mixtureType::fres
277  (
278  this->patchFaceComposition(Yslicer, patchi, facei)
279  );
280  }
281  }
282 
283  return tfres;
284 }
285 
286 
287 template<class BaseThermo>
289 {
290  BaseThermo::mixtureType::reset(this->Y());
291 }
292 
293 
294 template<class BaseThermo>
297 (
298  const scalarField& Tu,
299  const labelList& cells
300 ) const
301 {
302  return this->cellSetProperty
303  (
304  &BaseThermo::mixtureType::reactants,
306  cells,
307  UIndirectList<scalar>(this->p_, cells),
308  Tu
309  );
310 }
311 
312 
313 template<class BaseThermo>
316 (
317  const scalarField& Tu,
318  const label patchi
319 ) const
320 {
321  return this->patchFieldProperty
322  (
323  &BaseThermo::mixtureType::reactants,
325  patchi,
326  this->p_.boundaryField()[patchi],
327  Tu
328  );
329 }
330 
331 
332 template<class BaseThermo>
335 {
336  return this->volScalarFieldProperty
337  (
338  "Tb",
340  &BaseThermo::mixtureType::products,
341  &BaseThermo::mixtureType::thermoMixtureType::The,
342  this->he_,
343  this->p_,
344  this->T_
345  );
346 }
347 
348 
349 template<class BaseThermo>
352 {
353  return
354  this->volScalarFieldProperty
355  (
356  "hf",
358  &BaseThermo::mixtureType::reactants,
359  &BaseThermo::mixtureType::thermoMixtureType::hf
360  )
361  - this->volScalarFieldProperty
362  (
363  "hf",
365  &BaseThermo::mixtureType::products,
366  &BaseThermo::mixtureType::thermoMixtureType::hf
367  );
368 }
369 
370 
371 template<class BaseThermo>
374 {
375  return this->volScalarFieldProperty
376  (
377  "psiu",
378  this->psi_.dimensions(),
379  &BaseThermo::mixtureType::reactants,
381  this->p_,
382  this->Tu_
383  );
384 }
385 
386 
387 template<class BaseThermo>
390 {
391  const volScalarField Tb(this->Tb());
392 
393  return this->volScalarFieldProperty
394  (
395  "psib",
396  this->psi_.dimensions(),
397  &BaseThermo::mixtureType::products,
399  this->p_,
400  Tb
401  );
402 }
403 
404 
405 template<class BaseThermo>
408 {
409  return this->volScalarFieldProperty
410  (
411  "muu",
413  &BaseThermo::mixtureType::reactants,
415  this->p_,
416  this->Tu_
417  );
418 }
419 
420 
421 template<class BaseThermo>
424 {
425  const volScalarField Tb(this->Tb());
426 
427  return this->volScalarFieldProperty
428  (
429  "mub",
431  &BaseThermo::mixtureType::products,
433  this->p_,
434  Tb
435  );
436 }
437 
438 
439 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Generic GeometricBoundaryField class.
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Thermo implementation based on compressibility with additional unburnt thermodynamic state.
virtual tmp< volScalarField > mub() const
Dynamic viscosity of burnt gas [kg/m/s].
PsiuMulticomponentThermo(const fvMesh &, const word &phaseName)
Construct from mesh and phase name.
virtual void correct()
Update properties.
virtual tmp< volScalarField > muu() const
Dynamic viscosity of unburnt gas [kg/m/s].
virtual tmp< volScalarField > hr() const
Standard enthalpy of reaction [J/kg].
virtual const volScalarField & heu() const
Unburnt gas enthalpy [J/kg].
virtual tmp< volScalarField > Tb() const
Burnt gas temperature [K].
virtual tmp< volScalarField > psib() const
Burnt gas compressibility [s^2/m^2].
virtual tmp< volScalarField > psiu() const
Unburnt gas compressibility [s^2/m^2].
virtual void reset()
Reset the mixture to an unburnt state and update EGR.
virtual ~PsiuMulticomponentThermo()
Destructor.
virtual tmp< volScalarField > fres() const
Return the residual fraction of fuel in the burnt mixture.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:91
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
#define InfoInFunction
Report an information message using Foam::Info.
const dimensionedScalar mu
Atomic mass unit.
label calculate(const fvMesh &mesh, const labelHashSet &patchIDs, const scalar minFaceFraction, GeometricField< scalar, GeoMesh > &distance)
Calculate distance data from patches.
const dimensionSet dimDynamicViscosity
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
const dimensionSet dimless
messageStream Info
const dimensionSet dimTemperature
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimMass
fvPatchField< scalar > fvPatchScalarField
thermo he()
PtrList< volScalarField > & Y