heheuPsiThermo.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-2021 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 "heheuPsiThermo.H"
27 #include "fvMesh.H"
29 
30 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
31 
32 template<class BasicPsiThermo, class MixtureType>
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& psiCells = this->psi_.primitiveFieldRef();
42  scalarField& muCells = this->mu_.primitiveFieldRef();
43  scalarField& alphaCells = this->alpha_.primitiveFieldRef();
44 
45  forAll(TCells, celli)
46  {
47  const typename MixtureType::thermoMixtureType& thermoMixture =
48  this->cellThermoMixture(celli);
49 
50  const typename MixtureType::transportMixtureType& transportMixture =
51  this->cellTransportMixture(celli, thermoMixture);
52 
53  TCells[celli] = thermoMixture.THE
54  (
55  hCells[celli],
56  pCells[celli],
57  TCells[celli]
58  );
59 
60  psiCells[celli] = thermoMixture.psi(pCells[celli], TCells[celli]);
61 
62  muCells[celli] = transportMixture.mu(pCells[celli], TCells[celli]);
63  alphaCells[celli] =
64  transportMixture.kappa(pCells[celli], TCells[celli])
65  /thermoMixture.Cp(pCells[celli], TCells[celli]);
66 
67  TuCells[celli] = this->cellReactants(celli).THE
68  (
69  heuCells[celli],
70  pCells[celli],
71  TuCells[celli]
72  );
73  }
74 
75  volScalarField::Boundary& pBf =
76  this->p_.boundaryFieldRef();
77 
78  volScalarField::Boundary& TBf =
79  this->T_.boundaryFieldRef();
80 
81  volScalarField::Boundary& TuBf =
82  this->Tu_.boundaryFieldRef();
83 
84  volScalarField::Boundary& psiBf =
85  this->psi_.boundaryFieldRef();
86 
87  volScalarField::Boundary& heBf =
88  this->he().boundaryFieldRef();
89 
90  volScalarField::Boundary& heuBf =
91  this->heu().boundaryFieldRef();
92 
93  volScalarField::Boundary& muBf =
94  this->mu_.boundaryFieldRef();
95 
96  volScalarField::Boundary& alphaBf =
97  this->alpha_.boundaryFieldRef();
98 
99  forAll(this->T_.boundaryField(), patchi)
100  {
101  fvPatchScalarField& pp = pBf[patchi];
102  fvPatchScalarField& pT = TBf[patchi];
103  fvPatchScalarField& pTu = TuBf[patchi];
104  fvPatchScalarField& ppsi = psiBf[patchi];
105  fvPatchScalarField& phe = heBf[patchi];
106  fvPatchScalarField& pheu = heuBf[patchi];
107  fvPatchScalarField& pmu = muBf[patchi];
108  fvPatchScalarField& palpha = alphaBf[patchi];
109 
110  if (pT.fixesValue())
111  {
112  forAll(pT, facei)
113  {
114  const typename MixtureType::thermoMixtureType&
115  thermoMixture = this->patchFaceThermoMixture(patchi, facei);
116 
117  const typename MixtureType::transportMixtureType&
118  transportMixture =
119  this->patchFaceTransportMixture
120  (patchi, facei, thermoMixture);
121 
122  phe[facei] = thermoMixture.HE(pp[facei], pT[facei]);
123 
124  ppsi[facei] = thermoMixture.psi(pp[facei], pT[facei]);
125  pmu[facei] = transportMixture.mu(pp[facei], pT[facei]);
126  palpha[facei] =
127  transportMixture.kappa(pp[facei], pT[facei])
128  /thermoMixture.Cp(pp[facei], pT[facei]);
129  }
130  }
131  else
132  {
133  forAll(pT, facei)
134  {
135  const typename MixtureType::thermoMixtureType&
136  thermoMixture = this->patchFaceThermoMixture(patchi, facei);
137 
138  const typename MixtureType::transportMixtureType&
139  transportMixture =
140  this->patchFaceTransportMixture
141  (patchi, facei, thermoMixture);
142 
143  pT[facei] = thermoMixture.THE(phe[facei], pp[facei], pT[facei]);
144 
145  ppsi[facei] = thermoMixture.psi(pp[facei], pT[facei]);
146  pmu[facei] = transportMixture.mu(pp[facei], pT[facei]);
147  palpha[facei] =
148  transportMixture.kappa(pp[facei], pT[facei])
149  /thermoMixture.Cp(pp[facei], pT[facei]);
150 
151  pTu[facei] =
152  this->patchFaceReactants(patchi, facei)
153  .THE(pheu[facei], pp[facei], pTu[facei]);
154  }
155  }
156  }
157 }
158 
159 
160 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
161 
162 template<class BasicPsiThermo, class MixtureType>
164 (
165  const fvMesh& mesh,
166  const word& phaseName
167 )
168 :
170  Tu_
171  (
172  IOobject
173  (
174  "Tu",
175  mesh.time().timeName(),
176  mesh,
177  IOobject::MUST_READ,
178  IOobject::AUTO_WRITE
179  ),
180  mesh
181  ),
182  heu_
183  (
184  IOobject
185  (
186  MixtureType::thermoType::heName() + 'u',
187  mesh.time().timeName(),
188  mesh,
189  IOobject::NO_READ,
190  IOobject::NO_WRITE
191  ),
192  this->volScalarFieldProperty
193  (
194  MixtureType::thermoType::heName() + 'u',
196  &MixtureType::cellReactants,
197  &MixtureType::patchFaceReactants,
198  &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 BasicPsiThermo, class MixtureType>
217 {}
218 
219 
220 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
221 
222 template<class BasicPsiThermo, class MixtureType>
224 {
225  if (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 (debug)
236  {
237  Info<< " Finished" << endl;
238  }
239 }
240 
241 
242 template<class BasicPsiThermo, class MixtureType>
245 (
246  const scalarField& Tu,
247  const labelList& cells
248 ) const
249 {
250  return this->cellSetProperty
251  (
252  &MixtureType::cellReactants,
253  &MixtureType::thermoMixtureType::HE,
254  cells,
255  UIndirectList<scalar>(this->p_, cells),
256  Tu
257  );
258 }
259 
260 
261 template<class BasicPsiThermo, class MixtureType>
264 (
265  const scalarField& Tu,
266  const label patchi
267 ) const
268 {
269  return this->patchFieldProperty
270  (
271  &MixtureType::patchFaceReactants,
272  &MixtureType::thermoMixtureType::HE,
273  patchi,
274  this->p_.boundaryField()[patchi],
275  Tu
276  );
277 }
278 
279 
280 template<class BasicPsiThermo, class MixtureType>
283 {
284  return this->volScalarFieldProperty
285  (
286  "Tb",
288  &MixtureType::cellProducts,
289  &MixtureType::patchFaceProducts,
290  &MixtureType::thermoMixtureType::THE,
291  this->he_,
292  this->p_,
293  this->T_
294  );
295 }
296 
297 
298 template<class BasicPsiThermo, class MixtureType>
301 {
302  return this->volScalarFieldProperty
303  (
304  "psiu",
305  this->psi_.dimensions(),
306  &MixtureType::cellReactants,
307  &MixtureType::patchFaceReactants,
309  this->p_,
310  this->Tu_
311  );
312 }
313 
314 
315 template<class BasicPsiThermo, class MixtureType>
318 {
319  const volScalarField Tb(this->Tb());
320 
321  return this->volScalarFieldProperty
322  (
323  "psib",
324  this->psi_.dimensions(),
325  &MixtureType::cellProducts,
326  &MixtureType::patchFaceProducts,
328  this->p_,
329  Tb
330  );
331 }
332 
333 
334 template<class BasicPsiThermo, class MixtureType>
337 {
338  return this->volScalarFieldProperty
339  (
340  "muu",
342  &MixtureType::cellReactants,
343  &MixtureType::patchFaceReactants,
345  this->p_,
346  this->Tu_
347  );
348 }
349 
350 
351 template<class BasicPsiThermo, class MixtureType>
354 {
355  const volScalarField Tb(this->Tb());
356 
357  return this->volScalarFieldProperty
358  (
359  "mub",
361  &MixtureType::cellProducts,
362  &MixtureType::patchFaceProducts,
364  this->p_,
365  Tb
366  );
367 }
368 
369 
370 // ************************************************************************* //
virtual tmp< volScalarField > mub() const
Dynamic viscosity of burnt gas [kg/m/s].
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual tmp< volScalarField > psiu() const
Unburnt gas compressibility [s^2/m^2].
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
virtual tmp< volScalarField > Tb() const
Burnt gas temperature [K].
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual tmp< volScalarField > psib() const
Burnt gas compressibility [s^2/m^2].
virtual tmp< volScalarField > muu() const
Dynamic viscosity of unburnt gas [kg/m/s].
dynamicFvMesh & mesh
fvPatchField< scalar > fvPatchScalarField
A class for handling words, derived from string.
Definition: word.H:59
const dimensionSet dimDynamicViscosity
virtual void correct()
Update properties.
volScalarField scalarField(fieldObject, mesh)
virtual ~heheuPsiThermo()
Destructor.
const dimensionedScalar mu
Atomic mass unit.
thermo he()
const dimensionSet dimEnergy
const dimensionSet dimMass
label patchi
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
messageStream Info
const volScalarField & psi
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
virtual volScalarField & heu()
Unburnt gas enthalpy [J/kg].
heheuPsiThermo(const fvMesh &, const word &phaseName)
Construct from mesh and phase name.
const dimensionSet dimTemperature
#define InfoInFunction
Report an information message using Foam::Info.