heheuPsiThermo.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2017-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "heheuPsiThermo.H"
30 #include "fvMesh.H"
32 
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34 
35 template<class BasicPsiThermo, class MixtureType>
37 {
38  const scalarField& hCells = this->he_;
39  const scalarField& heuCells = this->heu_;
40  const scalarField& pCells = this->p_;
41 
42  scalarField& TCells = this->T_.primitiveFieldRef();
43  scalarField& TuCells = this->Tu_.primitiveFieldRef();
44  scalarField& psiCells = this->psi_.primitiveFieldRef();
45  scalarField& muCells = this->mu_.primitiveFieldRef();
46  scalarField& alphaCells = this->alpha_.primitiveFieldRef();
47 
48  forAll(TCells, celli)
49  {
50  const typename MixtureType::thermoType& mixture_ =
51  this->cellMixture(celli);
52 
53  if (this->updateT())
54  {
55  TCells[celli] = mixture_.THE
56  (
57  hCells[celli],
58  pCells[celli],
59  TCells[celli]
60  );
61  }
62 
63  psiCells[celli] = mixture_.psi(pCells[celli], TCells[celli]);
64 
65  muCells[celli] = mixture_.mu(pCells[celli], TCells[celli]);
66  alphaCells[celli] = mixture_.alphah(pCells[celli], TCells[celli]);
67 
68  TuCells[celli] = this->cellReactants(celli).THE
69  (
70  heuCells[celli],
71  pCells[celli],
72  TuCells[celli]
73  );
74  }
75 
76  volScalarField::Boundary& pBf =
77  this->p_.boundaryFieldRef();
78 
79  volScalarField::Boundary& TBf =
80  this->T_.boundaryFieldRef();
81 
82  volScalarField::Boundary& TuBf =
83  this->Tu_.boundaryFieldRef();
84 
85  volScalarField::Boundary& psiBf =
86  this->psi_.boundaryFieldRef();
87 
88  volScalarField::Boundary& heBf =
89  this->he().boundaryFieldRef();
90 
91  volScalarField::Boundary& heuBf =
92  this->heu().boundaryFieldRef();
93 
94  volScalarField::Boundary& muBf =
95  this->mu_.boundaryFieldRef();
96 
97  volScalarField::Boundary& alphaBf =
98  this->alpha_.boundaryFieldRef();
99 
100  forAll(this->T_.boundaryField(), patchi)
101  {
102  fvPatchScalarField& pp = pBf[patchi];
103  fvPatchScalarField& pT = TBf[patchi];
104  fvPatchScalarField& pTu = TuBf[patchi];
105  fvPatchScalarField& ppsi = psiBf[patchi];
106  fvPatchScalarField& phe = heBf[patchi];
107  fvPatchScalarField& pheu = heuBf[patchi];
108  fvPatchScalarField& pmu = muBf[patchi];
109  fvPatchScalarField& palpha = alphaBf[patchi];
110 
111  if (pT.fixesValue())
112  {
113  forAll(pT, facei)
114  {
115  const typename MixtureType::thermoType& mixture_ =
116  this->patchFaceMixture(patchi, facei);
117 
118  phe[facei] = mixture_.HE(pp[facei], pT[facei]);
119 
120  ppsi[facei] = mixture_.psi(pp[facei], pT[facei]);
121  pmu[facei] = mixture_.mu(pp[facei], pT[facei]);
122  palpha[facei] = mixture_.alphah(pp[facei], pT[facei]);
123  }
124  }
125  else
126  {
127  forAll(pT, facei)
128  {
129  const typename MixtureType::thermoType& mixture_ =
130  this->patchFaceMixture(patchi, facei);
131 
132  if (this->updateT())
133  {
134  pT[facei] = mixture_.THE(phe[facei], pp[facei], pT[facei]);
135  }
136 
137  ppsi[facei] = mixture_.psi(pp[facei], pT[facei]);
138  pmu[facei] = mixture_.mu(pp[facei], pT[facei]);
139  palpha[facei] = mixture_.alphah(pp[facei], pT[facei]);
140 
141  pTu[facei] =
142  this->patchFaceReactants(patchi, facei)
143  .THE(pheu[facei], pp[facei], pTu[facei]);
144  }
145  }
146  }
147 }
148 
149 
150 /*
151 template<class BasicPsiThermo, class MixtureType>
152 void Foam::heheuPsiThermo<BasicPsiThermo, MixtureType>::calculateT()
153 {
154  //const scalarField& hCells = this->he_.primitiveFieldRef();
155  const scalarField& heuCells = this->heu_.primitiveFieldRef();
156  const scalarField& pCells = this->p_.primitiveFieldRef();
157 
158  scalarField& TCells = this->T_.primitiveFieldRef();
159  scalarField& TuCells = this->Tu_.primitiveFieldRef();
160  scalarField& psiCells = this->psi_.primitiveFieldRef();
161  scalarField& muCells = this->mu_.primitiveFieldRef();
162  scalarField& alphaCells = this->alpha_.primitiveFieldRef();
163 
164  forAll(TCells, celli)
165  {
166  const typename MixtureType::thermoType& mixture_ =
167  this->cellMixture(celli);
168 
169 // TCells[celli] = mixture_.THE
170 // (
171 // hCells[celli],
172 // pCells[celli],
173 // TCells[celli]
174 // );
175 
176  psiCells[celli] = mixture_.psi(pCells[celli], TCells[celli]);
177 
178  muCells[celli] = mixture_.mu(pCells[celli], TCells[celli]);
179  alphaCells[celli] = mixture_.alphah(pCells[celli], TCells[celli]);
180 
181  TuCells[celli] = this->cellReactants(celli).THE
182  (
183  heuCells[celli],
184  pCells[celli],
185  TuCells[celli]
186  );
187  }
188 
189  forAll(this->T_.boundaryField(), patchi)
190  {
191  fvPatchScalarField& pp = this->p_.boundaryFieldRef()[patchi];
192  fvPatchScalarField& pT = this->T_.boundaryFieldRef()[patchi];
193  fvPatchScalarField& pTu = this->Tu_.boundaryFieldRef()[patchi];
194  fvPatchScalarField& ppsi = this->psi_.boundaryFieldRef()[patchi];
195 
196  fvPatchScalarField& ph = this->he_.boundaryFieldRef()[patchi];
197  fvPatchScalarField& pheu = this->heu_.boundaryFieldRef()[patchi];
198 
199  fvPatchScalarField& pmu_ = this->mu_.boundaryFieldRef()[patchi];
200  fvPatchScalarField& palpha_ = this->alpha_.boundaryFieldRef()[patchi];
201 
202  if (pT.fixesValue())
203  {
204  forAll(pT, facei)
205  {
206  const typename MixtureType::thermoType& mixture_ =
207  this->patchFaceMixture(patchi, facei);
208 
209  ph[facei] = mixture_.HE(pp[facei], pT[facei]);
210 
211  ppsi[facei] = mixture_.psi(pp[facei], pT[facei]);
212  pmu_[facei] = mixture_.mu(pp[facei], pT[facei]);
213  palpha_[facei] = mixture_.alphah(pp[facei], pT[facei]);
214  }
215  }
216  else
217  {
218  forAll(pT, facei)
219  {
220  const typename MixtureType::thermoType& mixture_ =
221  this->patchFaceMixture(patchi, facei);
222 
223  //pT[facei] = mixture_.THE(ph[facei], pp[facei], pT[facei]);
224 
225  ppsi[facei] = mixture_.psi(pp[facei], pT[facei]);
226  pmu_[facei] = mixture_.mu(pp[facei], pT[facei]);
227  palpha_[facei] = mixture_.alphah(pp[facei], pT[facei]);
228 
229  pTu[facei] =
230  this->patchFaceReactants(patchi, facei)
231  .THE(pheu[facei], pp[facei], pTu[facei]);
232  }
233  }
234  }
235 }
236 */
237 
238 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
239 
240 template<class BasicPsiThermo, class MixtureType>
242 (
243  const fvMesh& mesh,
244  const word& phaseName
245 )
246 :
247  heThermo<psiuReactionThermo, MixtureType>(mesh, phaseName),
248  Tu_
249  (
250  IOobject
251  (
252  "Tu",
253  mesh.time().timeName(),
254  mesh,
255  IOobject::MUST_READ,
256  IOobject::AUTO_WRITE,
257  IOobject::REGISTER
258  ),
259  mesh
260  ),
261 
262  heu_
263  (
264  IOobject
265  (
266  MixtureType::thermoType::heName() + 'u',
267  mesh.time().timeName(),
268  mesh,
269  IOobject::NO_READ,
270  IOobject::NO_WRITE
271  ),
272  mesh,
273  dimensionSet(0, 2, -2, 0, 0),
274  this->heuBoundaryTypes()
275  )
276 {
277  scalarField& heuCells = this->heu_.primitiveFieldRef();
278  const scalarField& pCells = this->p_;
279  const scalarField& TuCells = this->Tu_;
280 
281  forAll(heuCells, celli)
282  {
283  heuCells[celli] = this->cellReactants(celli).HE
284  (
285  pCells[celli],
286  TuCells[celli]
287  );
288  }
289 
291 
292  forAll(heuBf, patchi)
293  {
294  fvPatchScalarField& pheu = heuBf[patchi];
295  const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
296  const fvPatchScalarField& pTu = this->Tu_.boundaryField()[patchi];
297 
298  forAll(pheu, facei)
299  {
300  pheu[facei] = this->patchFaceReactants(patchi, facei).HE
301  (
302  pp[facei],
303  pTu[facei]
304  );
305  }
306  }
307 
308  this->heuBoundaryCorrection(this->heu_);
309 
310  calculate();
311  this->psi_.oldTime(); // Switch on saving old time
312 }
313 
314 
315 template<class BasicPsiThermo, class MixtureType>
317 (
318  const fvMesh& mesh,
319  const word& phaseName,
320  const word& dictName
321 )
322 :
323  heThermo<psiuReactionThermo, MixtureType>(mesh, phaseName, dictName),
324  Tu_
325  (
326  IOobject
327  (
328  "Tu",
329  mesh.time().timeName(),
330  mesh,
331  IOobject::MUST_READ,
332  IOobject::AUTO_WRITE,
333  IOobject::REGISTER
334  ),
335  mesh
336  ),
337 
338  heu_
339  (
340  IOobject
341  (
342  MixtureType::thermoType::heName() + 'u',
343  mesh.time().timeName(),
344  mesh,
345  IOobject::NO_READ,
346  IOobject::NO_WRITE
347  ),
348  mesh,
349  dimensionSet(0, 2, -2, 0, 0),
350  this->heuBoundaryTypes()
351  )
352 {
353  scalarField& heuCells = this->heu_.primitiveFieldRef();
354  const scalarField& pCells = this->p_;
355  const scalarField& TuCells = this->Tu_;
356 
357  forAll(heuCells, celli)
358  {
359  heuCells[celli] = this->cellReactants(celli).HE
360  (
361  pCells[celli],
362  TuCells[celli]
363  );
364  }
365 
367 
368  forAll(heuBf, patchi)
369  {
370  fvPatchScalarField& pheu = heuBf[patchi];
371  const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
372  const fvPatchScalarField& pTu = this->Tu_.boundaryField()[patchi];
373 
374  forAll(pheu, facei)
375  {
376  pheu[facei] = this->patchFaceReactants(patchi, facei).HE
377  (
378  pp[facei],
379  pTu[facei]
380  );
381  }
382  }
383 
384  this->heuBoundaryCorrection(this->heu_);
385 
386  calculate();
387  this->psi_.oldTime(); // Switch on saving old time
388 }
389 
390 
391 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
392 
393 template<class BasicPsiThermo, class MixtureType>
395 {}
396 
397 
398 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
399 
400 template<class BasicPsiThermo, class MixtureType>
402 {
404 
405  // force the saving of the old-time values
406  this->psi_.oldTime();
407 
408  calculate();
409 
410  DebugInfo << " Finished" << endl;
411 }
412 
413 
414 /*
415 template<class BasicPsiThermo, class MixtureType>
416 void Foam::heheuPsiThermo<BasicPsiThermo, MixtureType>::correctT()
417 {
418  // force the saving of the old-time values
419  this->psi_.oldTime();
420 
421  calculateT();
422 }
423 */
424 
425 template<class BasicPsiThermo, class MixtureType>
428 (
429  const scalarField& p,
430  const scalarField& Tu,
431  const labelList& cells
432 ) const
433 {
434  auto theu = tmp<scalarField>::New(Tu.size());
435  auto& heu = theu.ref();
436 
437  forAll(Tu, celli)
438  {
439  heu[celli] = this->cellReactants(cells[celli]).HE(p[celli], Tu[celli]);
440  }
441 
442  return theu;
443 }
444 
445 
446 template<class BasicPsiThermo, class MixtureType>
449 (
450  const scalarField& p,
451  const scalarField& Tu,
452  const label patchi
453 ) const
454 {
455  auto theu = tmp<scalarField>::New(Tu.size());
456  auto& heu = theu.ref();
457 
458  forAll(Tu, facei)
459  {
460  heu[facei] =
461  this->patchFaceReactants(patchi, facei).HE(p[facei], Tu[facei]);
462  }
464  return theu;
465 }
466 
467 
468 template<class BasicPsiThermo, class MixtureType>
471 {
472  auto tTb = volScalarField::New
473  (
474  "Tb",
476  this->T_
477  );
478  auto& Tb_ = tTb.ref();
479 
480  scalarField& TbCells = Tb_.primitiveFieldRef();
481  const scalarField& pCells = this->p_;
482  const scalarField& TCells = this->T_;
483  const scalarField& hCells = this->he_;
484 
485  forAll(TbCells, celli)
486  {
487  TbCells[celli] = this->cellProducts(celli).THE
488  (
489  hCells[celli],
490  pCells[celli],
491  TCells[celli]
492  );
493  }
494 
495  volScalarField::Boundary& TbBf = Tb_.boundaryFieldRef();
496 
497  forAll(TbBf, patchi)
498  {
499  fvPatchScalarField& pTb = TbBf[patchi];
500 
501  const fvPatchScalarField& ph = this->he_.boundaryField()[patchi];
502  const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
503  const fvPatchScalarField& pT = this->T_.boundaryField()[patchi];
504 
505  forAll(pTb, facei)
506  {
507  pTb[facei] =
508  this->patchFaceProducts(patchi, facei)
509  .THE(ph[facei], pp[facei], pT[facei]);
510  }
511  }
513  return tTb;
514 }
515 
516 
517 template<class BasicPsiThermo, class MixtureType>
520 {
521  auto tpsiu = volScalarField::New
522  (
523  "psiu",
525  this->psi_.mesh(),
526  this->psi_.dimensions()
527  );
528  auto& psiu = tpsiu.ref();
529 
530  scalarField& psiuCells = psiu.primitiveFieldRef();
531  const scalarField& TuCells = this->Tu_;
532  const scalarField& pCells = this->p_;
533 
534  forAll(psiuCells, celli)
535  {
536  psiuCells[celli] =
537  this->cellReactants(celli).psi(pCells[celli], TuCells[celli]);
538  }
539 
540  volScalarField::Boundary& psiuBf = psiu.boundaryFieldRef();
541 
542  forAll(psiuBf, patchi)
543  {
544  fvPatchScalarField& ppsiu = psiuBf[patchi];
545 
546  const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
547  const fvPatchScalarField& pTu = this->Tu_.boundaryField()[patchi];
548 
549  forAll(ppsiu, facei)
550  {
551  ppsiu[facei] =
552  this->
553  patchFaceReactants(patchi, facei).psi(pp[facei], pTu[facei]);
554  }
555  }
557  return tpsiu;
558 }
559 
560 
561 template<class BasicPsiThermo, class MixtureType>
564 {
565  auto tpsib = volScalarField::New
566  (
567  "psib",
569  this->psi_.mesh(),
570  this->psi_.dimensions()
571  );
572  auto& psib = tpsib.ref();
573 
574  scalarField& psibCells = psib.primitiveFieldRef();
575  const volScalarField Tb_(Tb());
576  const scalarField& TbCells = Tb_;
577  const scalarField& pCells = this->p_;
578 
579  forAll(psibCells, celli)
580  {
581  psibCells[celli] =
582  this->cellProducts(celli).psi(pCells[celli], TbCells[celli]);
583  }
584 
585  volScalarField::Boundary& psibBf = psib.boundaryFieldRef();
586 
587  forAll(psibBf, patchi)
588  {
589  fvPatchScalarField& ppsib = psibBf[patchi];
590 
591  const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
592  const fvPatchScalarField& pTb = Tb_.boundaryField()[patchi];
593 
594  forAll(ppsib, facei)
595  {
596  ppsib[facei] =
597  this->patchFaceProducts
598  (patchi, facei).psi(pp[facei], pTb[facei]);
599  }
600  }
602  return tpsib;
603 }
604 
605 
606 template<class BasicPsiThermo, class MixtureType>
609 {
610  auto tmuu = volScalarField::New
611  (
612  "muu",
614  this->T_.mesh(),
615  dimensionSet(1, -1, -1, 0, 0)
616  );
617  auto& muu_ = tmuu.ref();
618 
619  scalarField& muuCells = muu_.primitiveFieldRef();
620  const scalarField& pCells = this->p_;
621  const scalarField& TuCells = this->Tu_;
622 
623  forAll(muuCells, celli)
624  {
625  muuCells[celli] = this->cellReactants(celli).mu
626  (
627  pCells[celli],
628  TuCells[celli]
629  );
630  }
631 
632  volScalarField::Boundary& muuBf = muu_.boundaryFieldRef();
633 
634  forAll(muuBf, patchi)
635  {
636  fvPatchScalarField& pMuu = muuBf[patchi];
637  const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
638  const fvPatchScalarField& pTu = this->Tu_.boundaryField()[patchi];
639 
640  forAll(pMuu, facei)
641  {
642  pMuu[facei] = this->patchFaceReactants(patchi, facei).mu
643  (
644  pp[facei],
645  pTu[facei]
646  );
647  }
648  }
650  return tmuu;
651 }
652 
653 
654 template<class BasicPsiThermo, class MixtureType>
657 {
658  auto tmub = volScalarField::New
659  (
660  "mub",
662  this->T_.mesh(),
663  dimensionSet(1, -1, -1, 0, 0)
664  );
665  auto& mub_ = tmub.ref();
666 
667  scalarField& mubCells = mub_.primitiveFieldRef();
668  const volScalarField Tb_(Tb());
669  const scalarField& pCells = this->p_;
670  const scalarField& TbCells = Tb_;
671 
672  forAll(mubCells, celli)
673  {
674  mubCells[celli] = this->cellProducts(celli).mu
675  (
676  pCells[celli],
677  TbCells[celli]
678  );
679  }
680 
681  volScalarField::Boundary& mubBf = mub_.boundaryFieldRef();
682 
683  forAll(mubBf, patchi)
684  {
685  fvPatchScalarField& pMub = mubBf[patchi];
686  const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
687  const fvPatchScalarField& pTb = Tb_.boundaryField()[patchi];
688 
689  forAll(pMub, facei)
690  {
691  pMub[facei] = this->patchFaceProducts(patchi, facei).mu
692  (
693  pp[facei],
694  pTb[facei]
695  );
696  }
697  }
698 
699  return tmub;
700 }
701 
702 
703 // ************************************************************************* //
volScalarField & he
Definition: YEEqn.H:52
virtual tmp< volScalarField > mub() const
Dynamic viscosity of burnt gas [kg/ms].
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
virtual tmp< volScalarField > psiu() const
Unburnt gas compressibility [s^2/m^2].
virtual tmp< volScalarField > Tb() const
Burnt gas temperature [K].
const word dictName("faMeshDefinition")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
virtual tmp< volScalarField > psib() const
Burnt gas compressibility [s^2/m^2].
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of boundary fields.
virtual tmp< volScalarField > muu() const
Dynamic viscosity of unburnt gas [kg/ms].
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
word timeName
Definition: getTimeIndex.H:3
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:105
dynamicFvMesh & mesh
const cellShapeList & cells
fvPatchField< scalar > fvPatchScalarField
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
#define DebugInFunction
Report an information message using Foam::Info.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
#define DebugInfo
Report an information message using Foam::Info.
Foam::psiuReactionThermo.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
virtual void correct()
Update properties.
virtual ~heheuPsiThermo()
Destructor.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
Enthalpy/Internal energy for a mixture.
Definition: heThermo.H:48
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
virtual volScalarField & heu()
Update properties based on T.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
Do not request registration (bool: false)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())