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-2019 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::thermoType& mixture_ =
48  this->cellMixture(celli);
49 
50  TCells[celli] = mixture_.THE
51  (
52  hCells[celli],
53  pCells[celli],
54  TCells[celli]
55  );
56 
57  psiCells[celli] = mixture_.psi(pCells[celli], TCells[celli]);
58 
59  muCells[celli] = mixture_.mu(pCells[celli], TCells[celli]);
60  alphaCells[celli] = mixture_.alphah(pCells[celli], TCells[celli]);
61 
62  TuCells[celli] = this->cellReactants(celli).THE
63  (
64  heuCells[celli],
65  pCells[celli],
66  TuCells[celli]
67  );
68  }
69 
70  volScalarField::Boundary& pBf =
71  this->p_.boundaryFieldRef();
72 
73  volScalarField::Boundary& TBf =
74  this->T_.boundaryFieldRef();
75 
76  volScalarField::Boundary& TuBf =
77  this->Tu_.boundaryFieldRef();
78 
79  volScalarField::Boundary& psiBf =
80  this->psi_.boundaryFieldRef();
81 
82  volScalarField::Boundary& heBf =
83  this->he().boundaryFieldRef();
84 
85  volScalarField::Boundary& heuBf =
86  this->heu().boundaryFieldRef();
87 
88  volScalarField::Boundary& muBf =
89  this->mu_.boundaryFieldRef();
90 
91  volScalarField::Boundary& alphaBf =
92  this->alpha_.boundaryFieldRef();
93 
94  forAll(this->T_.boundaryField(), patchi)
95  {
96  fvPatchScalarField& pp = pBf[patchi];
97  fvPatchScalarField& pT = TBf[patchi];
98  fvPatchScalarField& pTu = TuBf[patchi];
99  fvPatchScalarField& ppsi = psiBf[patchi];
100  fvPatchScalarField& phe = heBf[patchi];
101  fvPatchScalarField& pheu = heuBf[patchi];
102  fvPatchScalarField& pmu = muBf[patchi];
103  fvPatchScalarField& palpha = alphaBf[patchi];
104 
105  if (pT.fixesValue())
106  {
107  forAll(pT, facei)
108  {
109  const typename MixtureType::thermoType& mixture_ =
110  this->patchFaceMixture(patchi, facei);
111 
112  phe[facei] = mixture_.HE(pp[facei], pT[facei]);
113 
114  ppsi[facei] = mixture_.psi(pp[facei], pT[facei]);
115  pmu[facei] = mixture_.mu(pp[facei], pT[facei]);
116  palpha[facei] = mixture_.alphah(pp[facei], pT[facei]);
117  }
118  }
119  else
120  {
121  forAll(pT, facei)
122  {
123  const typename MixtureType::thermoType& mixture_ =
124  this->patchFaceMixture(patchi, facei);
125 
126  pT[facei] = mixture_.THE(phe[facei], pp[facei], pT[facei]);
127 
128  ppsi[facei] = mixture_.psi(pp[facei], pT[facei]);
129  pmu[facei] = mixture_.mu(pp[facei], pT[facei]);
130  palpha[facei] = mixture_.alphah(pp[facei], pT[facei]);
131 
132  pTu[facei] =
133  this->patchFaceReactants(patchi, facei)
134  .THE(pheu[facei], pp[facei], pTu[facei]);
135  }
136  }
137  }
138 }
139 
140 
141 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
142 
143 template<class BasicPsiThermo, class MixtureType>
145 (
146  const fvMesh& mesh,
147  const word& phaseName
148 )
149 :
151  Tu_
152  (
153  IOobject
154  (
155  "Tu",
156  mesh.time().timeName(),
157  mesh,
158  IOobject::MUST_READ,
159  IOobject::AUTO_WRITE
160  ),
161  mesh
162  ),
163  heu_
164  (
165  IOobject
166  (
167  MixtureType::thermoType::heName() + 'u',
168  mesh.time().timeName(),
169  mesh,
170  IOobject::NO_READ,
171  IOobject::NO_WRITE
172  ),
173  this->volScalarFieldProperty
174  (
175  MixtureType::thermoType::heName() + 'u',
177  &MixtureType::cellReactants,
178  &MixtureType::patchFaceReactants,
179  &MixtureType::thermoType::HE,
180  this->p_,
181  this->Tu_
182  ),
183  this->heuBoundaryTypes()
184  )
185 {
186  this->heuBoundaryCorrection(this->heu_);
187 
188  calculate();
189 
190  this->psi_.oldTime(); // Switch on saving old time
191 }
192 
193 
194 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
195 
196 template<class BasicPsiThermo, class MixtureType>
198 {}
199 
200 
201 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
202 
203 template<class BasicPsiThermo, class MixtureType>
205 {
206  if (debug)
207  {
208  InfoInFunction << endl;
209  }
210 
211  // force the saving of the old-time values
212  this->psi_.oldTime();
213 
214  calculate();
215 
216  if (debug)
217  {
218  Info<< " Finished" << endl;
219  }
220 }
221 
222 
223 template<class BasicPsiThermo, class MixtureType>
226 (
227  const scalarField& Tu,
228  const labelList& cells
229 ) const
230 {
231  return this->cellSetProperty
232  (
233  &MixtureType::cellReactants,
234  &MixtureType::thermoType::HE,
235  cells,
236  UIndirectList<scalar>(this->p(), cells),
237  Tu
238  );
239 }
240 
241 
242 template<class BasicPsiThermo, class MixtureType>
245 (
246  const scalarField& Tu,
247  const label patchi
248 ) const
249 {
250  return this->patchFieldProperty
251  (
252  &MixtureType::patchFaceReactants,
253  &MixtureType::thermoType::HE,
254  patchi,
255  this->p().boundaryField()[patchi],
256  Tu
257  );
258 }
259 
260 
261 template<class BasicPsiThermo, class MixtureType>
264 {
265  return this->volScalarFieldProperty
266  (
267  "Tb",
269  &MixtureType::cellProducts,
270  &MixtureType::patchFaceProducts,
271  &MixtureType::thermoType::THE,
272  this->he_,
273  this->p_,
274  this->T_
275  );
276 }
277 
278 
279 template<class BasicPsiThermo, class MixtureType>
282 {
283  return this->volScalarFieldProperty
284  (
285  "psiu",
286  this->psi_.dimensions(),
287  &MixtureType::cellReactants,
288  &MixtureType::patchFaceReactants,
290  this->p_,
291  this->T_
292  );
293 }
294 
295 
296 template<class BasicPsiThermo, class MixtureType>
299 {
300  return this->volScalarFieldProperty
301  (
302  "psib",
303  this->psi_.dimensions(),
304  &MixtureType::cellProducts,
305  &MixtureType::patchFaceProducts,
307  this->p_,
308  this->T_
309  );
310 }
311 
312 
313 template<class BasicPsiThermo, class MixtureType>
316 {
317  return this->volScalarFieldProperty
318  (
319  "muu",
321  &MixtureType::cellReactants,
322  &MixtureType::patchFaceReactants,
324  this->p_,
325  this->T_
326  );
327 }
328 
329 
330 template<class BasicPsiThermo, class MixtureType>
333 {
334  return this->volScalarFieldProperty
335  (
336  "mub",
338  &MixtureType::cellProducts,
339  &MixtureType::patchFaceProducts,
341  this->p_,
342  this->T_
343  );
344 }
345 
346 
347 // ************************************************************************* //
volScalarField & he
Definition: YEEqn.H:50
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].
const dimensionSet dimDynamicViscosity
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 dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
virtual void correct()
Update properties.
volScalarField scalarField(fieldObject, mesh)
virtual ~heheuPsiThermo()
Destructor.
const dimensionSet dimEnergy
label patchi
const dimensionedScalar & mu
Atomic mass unit.
Enthalpy/Internal energy for a mixture.
Definition: heThermo.H:49
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
messageStream Info
const volScalarField & psi
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
volScalarField & p
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.
#define InfoInFunction
Report an information message using Foam::Info.