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-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 
27 #include "fvMesh.H"
29 
30 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
31 
32 template<class BaseThermo>
34 {
35  const scalarField& hCells = this->he_;
36  const scalarField& heuCells = this->heu_;
37  const scalarField& pCells = this->p_;
38 
39  scalarField& TCells = this->T_.primitiveFieldRef();
40  scalarField& TuCells = this->Tu_.primitiveFieldRef();
41  scalarField& CpCells = this->Cp_.primitiveFieldRef();
42  scalarField& CvCells = this->Cv_.primitiveFieldRef();
43  scalarField& psiCells = this->psi_.primitiveFieldRef();
44  scalarField& muCells = this->mu_.primitiveFieldRef();
45  scalarField& kappaCells = this->kappa_.primitiveFieldRef();
46 
47  auto Yslicer = this->Yslicer();
48 
49  forAll(TCells, celli)
50  {
51  auto composition = this->cellComposition(Yslicer, celli);
52 
53  const typename BaseThermo::mixtureType::thermoMixtureType&
54  thermoMixture = this->thermoMixture(composition);
55 
56  const typename BaseThermo::mixtureType::transportMixtureType&
57  transportMixture =
58  this->transportMixture(composition, thermoMixture);
59 
60  TCells[celli] = thermoMixture.The
61  (
62  hCells[celli],
63  pCells[celli],
64  TCells[celli]
65  );
66 
67  CpCells[celli] = thermoMixture.Cp(pCells[celli], TCells[celli]);
68  CvCells[celli] = thermoMixture.Cv(pCells[celli], TCells[celli]);
69  psiCells[celli] = thermoMixture.psi(pCells[celli], TCells[celli]);
70 
71  muCells[celli] = transportMixture.mu(pCells[celli], TCells[celli]);
72  kappaCells[celli] =
73  transportMixture.kappa(pCells[celli], TCells[celli]);
74 
75  TuCells[celli] = this->reactants(composition).The
76  (
77  heuCells[celli],
78  pCells[celli],
79  TuCells[celli]
80  );
81  }
82 
83  volScalarField::Boundary& pBf =
84  this->p_.boundaryFieldRef();
85 
86  volScalarField::Boundary& TBf =
87  this->T_.boundaryFieldRef();
88 
89  volScalarField::Boundary& TuBf =
90  this->Tu_.boundaryFieldRef();
91 
92  volScalarField::Boundary& CpBf =
93  this->Cp_.boundaryFieldRef();
94 
95  volScalarField::Boundary& CvBf =
96  this->Cv_.boundaryFieldRef();
97 
98  volScalarField::Boundary& psiBf =
99  this->psi_.boundaryFieldRef();
100 
101  volScalarField::Boundary& heBf =
102  this->he().boundaryFieldRef();
103 
104  volScalarField::Boundary& heuBf =
105  this->heu().boundaryFieldRef();
106 
107  volScalarField::Boundary& muBf =
108  this->mu_.boundaryFieldRef();
109 
110  volScalarField::Boundary& kappaBf =
111  this->kappa_.boundaryFieldRef();
112 
113  forAll(this->T_.boundaryField(), patchi)
114  {
115  fvPatchScalarField& pp = pBf[patchi];
116  fvPatchScalarField& pT = TBf[patchi];
117  fvPatchScalarField& pTu = TuBf[patchi];
118  fvPatchScalarField& pCp = CpBf[patchi];
119  fvPatchScalarField& pCv = CvBf[patchi];
120  fvPatchScalarField& ppsi = psiBf[patchi];
121  fvPatchScalarField& phe = heBf[patchi];
122  fvPatchScalarField& pheu = heuBf[patchi];
123  fvPatchScalarField& pmu = muBf[patchi];
124  fvPatchScalarField& pkappa = kappaBf[patchi];
125 
126  if (pT.fixesValue())
127  {
128  forAll(pT, facei)
129  {
130  auto composition =
131  this->patchFaceComposition(Yslicer, patchi, facei);
132 
133  const typename BaseThermo::mixtureType::thermoMixtureType&
134  thermoMixture = this->thermoMixture(composition);
135 
136  const typename BaseThermo::mixtureType::transportMixtureType&
137  transportMixture =
138  this->transportMixture(composition, thermoMixture);
139 
140  phe[facei] = thermoMixture.he(pp[facei], pT[facei]);
141 
142  pCp[facei] = thermoMixture.Cp(pp[facei], pT[facei]);
143  pCv[facei] = thermoMixture.Cv(pp[facei], pT[facei]);
144  ppsi[facei] = thermoMixture.psi(pp[facei], pT[facei]);
145  pmu[facei] = transportMixture.mu(pp[facei], pT[facei]);
146  pkappa[facei] = transportMixture.kappa(pp[facei], pT[facei]);
147  }
148  }
149  else
150  {
151  forAll(pT, facei)
152  {
153  auto composition =
154  this->patchFaceComposition(Yslicer, patchi, facei);
155 
156  const typename BaseThermo::mixtureType::thermoMixtureType&
157  thermoMixture = this->thermoMixture(composition);
158 
159  const typename BaseThermo::mixtureType::transportMixtureType&
160  transportMixture =
161  this->transportMixture(composition, thermoMixture);
162 
163  pT[facei] = thermoMixture.The(phe[facei], pp[facei], pT[facei]);
164 
165  pCp[facei] = thermoMixture.Cp(pp[facei], pT[facei]);
166  pCv[facei] = thermoMixture.Cv(pp[facei], pT[facei]);
167  ppsi[facei] = thermoMixture.psi(pp[facei], pT[facei]);
168  pmu[facei] = transportMixture.mu(pp[facei], pT[facei]);
169  pkappa[facei] = transportMixture.kappa(pp[facei], pT[facei]);
170 
171  pTu[facei] =
172  this->reactants(composition)
173  .The(pheu[facei], pp[facei], pTu[facei]);
174  }
175  }
176  }
177 }
178 
179 
180 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
181 
182 template<class BaseThermo>
184 (
185  const fvMesh& mesh,
186  const word& phaseName
187 )
188 :
189  BaseThermo(mesh, phaseName),
190  Tu_
191  (
192  IOobject
193  (
194  "Tu",
195  mesh.time().name(),
196  mesh,
197  IOobject::MUST_READ,
198  IOobject::AUTO_WRITE
199  ),
200  mesh
201  ),
202  heu_
203  (
204  IOobject
205  (
206  BaseThermo::mixtureType::thermoType::heName() + 'u',
207  mesh.time().name(),
208  mesh,
209  IOobject::NO_READ,
210  IOobject::NO_WRITE
211  ),
212  this->volScalarFieldProperty
213  (
214  BaseThermo::mixtureType::thermoType::heName() + 'u',
216  &BaseThermo::mixtureType::reactants,
217  &BaseThermo::mixtureType::thermoMixtureType::he,
218  this->p_,
219  this->Tu_
220  ),
221  this->heuBoundaryTypes()
222  )
223 {
224  this->heuBoundaryCorrection(this->heu_);
225 
226  calculate();
227 
228  this->psi_.oldTime(); // Switch on saving old time
229 }
230 
231 
232 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
233 
234 template<class BaseThermo>
236 {}
237 
238 
239 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
240 
241 template<class BaseThermo>
243 {
244  if (BaseThermo::debug)
245  {
246  InfoInFunction << endl;
247  }
248 
249  // force the saving of the old-time values
250  this->psi_.oldTime();
251 
252  calculate();
253 
254  if (BaseThermo::debug)
255  {
256  Info<< " Finished" << endl;
257  }
258 }
259 
260 
261 template<class BaseThermo>
264 (
265  const scalarField& Tu,
266  const labelList& cells
267 ) const
268 {
269  return this->cellSetProperty
270  (
271  &BaseThermo::mixtureType::reactants,
273  cells,
274  UIndirectList<scalar>(this->p_, cells),
275  Tu
276  );
277 }
278 
279 
280 template<class BaseThermo>
283 (
284  const scalarField& Tu,
285  const label patchi
286 ) const
287 {
288  return this->patchFieldProperty
289  (
290  &BaseThermo::mixtureType::reactants,
292  patchi,
293  this->p_.boundaryField()[patchi],
294  Tu
295  );
296 }
297 
298 
299 template<class BaseThermo>
302 {
303  return this->volScalarFieldProperty
304  (
305  "Tb",
307  &BaseThermo::mixtureType::products,
308  &BaseThermo::mixtureType::thermoMixtureType::The,
309  this->he_,
310  this->p_,
311  this->T_
312  );
313 }
314 
315 
316 template<class BaseThermo>
319 {
320  return this->volScalarFieldProperty
321  (
322  "psiu",
323  this->psi_.dimensions(),
324  &BaseThermo::mixtureType::reactants,
326  this->p_,
327  this->Tu_
328  );
329 }
330 
331 
332 template<class BaseThermo>
335 {
336  const volScalarField Tb(this->Tb());
337 
338  return this->volScalarFieldProperty
339  (
340  "psib",
341  this->psi_.dimensions(),
342  &BaseThermo::mixtureType::products,
344  this->p_,
345  Tb
346  );
347 }
348 
349 
350 template<class BaseThermo>
353 {
354  return this->volScalarFieldProperty
355  (
356  "muu",
358  &BaseThermo::mixtureType::reactants,
360  this->p_,
361  this->Tu_
362  );
363 }
364 
365 
366 template<class BaseThermo>
369 {
370  const volScalarField Tb(this->Tb());
371 
372  return this->volScalarFieldProperty
373  (
374  "mub",
376  &BaseThermo::mixtureType::products,
378  this->p_,
379  Tb
380  );
381 }
382 
383 
384 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Generic GeometricField class.
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 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 ~PsiuMulticomponentThermo()
Destructor.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
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, PatchField, 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:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
messageStream Info
const dimensionSet dimTemperature
const dimensionSet dimMass
fvPatchField< scalar > fvPatchScalarField
thermo he()