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-------------------------------------------------------------------------------
11License
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
35template<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/*
151template<class BasicPsiThermo, class MixtureType>
152void 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
240template<class BasicPsiThermo, class MixtureType>
242(
243 const fvMesh& mesh,
244 const word& phaseName
245)
246:
247 heThermo<psiuReactionThermo, MixtureType>(mesh, phaseName),
248 Tu_
249 (
251 (
252 "Tu",
253 mesh.time().timeName(),
254 mesh,
255 IOobject::MUST_READ,
256 IOobject::AUTO_WRITE
257 ),
258 mesh
259 ),
260
261 heu_
262 (
264 (
265 MixtureType::thermoType::heName() + 'u',
266 mesh.time().timeName(),
267 mesh,
268 IOobject::NO_READ,
269 IOobject::NO_WRITE
270 ),
271 mesh,
272 dimensionSet(0, 2, -2, 0, 0),
273 this->heuBoundaryTypes()
274 )
275{
276 scalarField& heuCells = this->heu_.primitiveFieldRef();
277 const scalarField& pCells = this->p_;
278 const scalarField& TuCells = this->Tu_;
279
280 forAll(heuCells, celli)
281 {
282 heuCells[celli] = this->cellReactants(celli).HE
283 (
284 pCells[celli],
285 TuCells[celli]
286 );
287 }
288
290
291 forAll(heuBf, patchi)
292 {
293 fvPatchScalarField& pheu = heuBf[patchi];
294 const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
295 const fvPatchScalarField& pTu = this->Tu_.boundaryField()[patchi];
296
297 forAll(pheu, facei)
298 {
299 pheu[facei] = this->patchFaceReactants(patchi, facei).HE
300 (
301 pp[facei],
302 pTu[facei]
303 );
304 }
305 }
306
307 this->heuBoundaryCorrection(this->heu_);
308
309 calculate();
310 this->psi_.oldTime(); // Switch on saving old time
311}
312
313
314template<class BasicPsiThermo, class MixtureType>
316(
317 const fvMesh& mesh,
318 const word& phaseName,
319 const word& dictName
320)
321:
322 heThermo<psiuReactionThermo, MixtureType>(mesh, phaseName, dictName),
323 Tu_
324 (
326 (
327 "Tu",
328 mesh.time().timeName(),
329 mesh,
330 IOobject::MUST_READ,
331 IOobject::AUTO_WRITE
332 ),
333 mesh
334 ),
335
336 heu_
337 (
339 (
340 MixtureType::thermoType::heName() + 'u',
341 mesh.time().timeName(),
342 mesh,
343 IOobject::NO_READ,
344 IOobject::NO_WRITE
345 ),
346 mesh,
347 dimensionSet(0, 2, -2, 0, 0),
348 this->heuBoundaryTypes()
349 )
350{
351 scalarField& heuCells = this->heu_.primitiveFieldRef();
352 const scalarField& pCells = this->p_;
353 const scalarField& TuCells = this->Tu_;
354
355 forAll(heuCells, celli)
356 {
357 heuCells[celli] = this->cellReactants(celli).HE
358 (
359 pCells[celli],
360 TuCells[celli]
361 );
362 }
363
365
366 forAll(heuBf, patchi)
367 {
368 fvPatchScalarField& pheu = heuBf[patchi];
369 const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
370 const fvPatchScalarField& pTu = this->Tu_.boundaryField()[patchi];
371
372 forAll(pheu, facei)
373 {
374 pheu[facei] = this->patchFaceReactants(patchi, facei).HE
375 (
376 pp[facei],
377 pTu[facei]
378 );
379 }
380 }
381
382 this->heuBoundaryCorrection(this->heu_);
383
384 calculate();
385 this->psi_.oldTime(); // Switch on saving old time
386}
387
388
389// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
390
391template<class BasicPsiThermo, class MixtureType>
393{}
394
395
396// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
397
398template<class BasicPsiThermo, class MixtureType>
400{
402
403 // force the saving of the old-time values
404 this->psi_.oldTime();
405
406 calculate();
407
408 DebugInfo << " Finished" << endl;
409}
410
411
412/*
413template<class BasicPsiThermo, class MixtureType>
414void Foam::heheuPsiThermo<BasicPsiThermo, MixtureType>::correctT()
415{
416 // force the saving of the old-time values
417 this->psi_.oldTime();
418
419 calculateT();
420}
421*/
422
423template<class BasicPsiThermo, class MixtureType>
426(
427 const scalarField& p,
428 const scalarField& Tu,
429 const labelList& cells
430) const
431{
432 tmp<scalarField> theu(new scalarField(Tu.size()));
433 scalarField& heu = theu.ref();
434
435 forAll(Tu, celli)
436 {
437 heu[celli] = this->cellReactants(cells[celli]).HE(p[celli], Tu[celli]);
438 }
439
440 return theu;
441}
442
443
444template<class BasicPsiThermo, class MixtureType>
447(
448 const scalarField& p,
449 const scalarField& Tu,
450 const label patchi
451) const
452{
453 tmp<scalarField> theu(new scalarField(Tu.size()));
454 scalarField& heu = theu.ref();
455
456 forAll(Tu, facei)
457 {
458 heu[facei] =
459 this->patchFaceReactants(patchi, facei).HE(p[facei], Tu[facei]);
460 }
461
462 return theu;
463}
464
465
466template<class BasicPsiThermo, class MixtureType>
469{
471 (
473 (
475 (
476 "Tb",
477 this->T_.time().timeName(),
478 this->T_.db(),
481 false
482 ),
483 this->T_
484 )
485 );
486
487 volScalarField& Tb_ = tTb.ref();
488 scalarField& TbCells = Tb_.primitiveFieldRef();
489 const scalarField& pCells = this->p_;
490 const scalarField& TCells = this->T_;
491 const scalarField& hCells = this->he_;
492
493 forAll(TbCells, celli)
494 {
495 TbCells[celli] = this->cellProducts(celli).THE
496 (
497 hCells[celli],
498 pCells[celli],
499 TCells[celli]
500 );
501 }
502
504
505 forAll(TbBf, patchi)
506 {
507 fvPatchScalarField& pTb = TbBf[patchi];
508
509 const fvPatchScalarField& ph = this->he_.boundaryField()[patchi];
510 const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
511 const fvPatchScalarField& pT = this->T_.boundaryField()[patchi];
512
513 forAll(pTb, facei)
514 {
515 pTb[facei] =
516 this->patchFaceProducts(patchi, facei)
517 .THE(ph[facei], pp[facei], pT[facei]);
518 }
519 }
520
521 return tTb;
522}
523
524
525template<class BasicPsiThermo, class MixtureType>
528{
530 (
532 (
534 (
535 "psiu",
536 this->psi_.time().timeName(),
537 this->psi_.db(),
540 false
541 ),
542 this->psi_.mesh(),
543 this->psi_.dimensions()
544 )
545 );
546
547 volScalarField& psiu = tpsiu.ref();
548 scalarField& psiuCells = psiu.primitiveFieldRef();
549 const scalarField& TuCells = this->Tu_;
550 const scalarField& pCells = this->p_;
551
552 forAll(psiuCells, celli)
553 {
554 psiuCells[celli] =
555 this->cellReactants(celli).psi(pCells[celli], TuCells[celli]);
556 }
557
559
560 forAll(psiuBf, patchi)
561 {
562 fvPatchScalarField& ppsiu = psiuBf[patchi];
563
564 const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
565 const fvPatchScalarField& pTu = this->Tu_.boundaryField()[patchi];
566
567 forAll(ppsiu, facei)
568 {
569 ppsiu[facei] =
570 this->
571 patchFaceReactants(patchi, facei).psi(pp[facei], pTu[facei]);
572 }
573 }
574
575 return tpsiu;
576}
577
578
579template<class BasicPsiThermo, class MixtureType>
582{
584 (
586 (
588 (
589 "psib",
590 this->psi_.time().timeName(),
591 this->psi_.db(),
594 false
595 ),
596 this->psi_.mesh(),
597 this->psi_.dimensions()
598 )
599 );
600
601 volScalarField& psib = tpsib.ref();
602 scalarField& psibCells = psib.primitiveFieldRef();
603 const volScalarField Tb_(Tb());
604 const scalarField& TbCells = Tb_;
605 const scalarField& pCells = this->p_;
606
607 forAll(psibCells, celli)
608 {
609 psibCells[celli] =
610 this->cellProducts(celli).psi(pCells[celli], TbCells[celli]);
611 }
612
614
615 forAll(psibBf, patchi)
616 {
617 fvPatchScalarField& ppsib = psibBf[patchi];
618
619 const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
620 const fvPatchScalarField& pTb = Tb_.boundaryField()[patchi];
621
622 forAll(ppsib, facei)
623 {
624 ppsib[facei] =
625 this->patchFaceProducts
626 (patchi, facei).psi(pp[facei], pTb[facei]);
627 }
628 }
629
630 return tpsib;
631}
632
633
634template<class BasicPsiThermo, class MixtureType>
637{
639 (
641 (
643 (
644 "muu",
645 this->T_.time().timeName(),
646 this->T_.db(),
649 false
650 ),
651 this->T_.mesh(),
652 dimensionSet(1, -1, -1, 0, 0)
653 )
654 );
655
656 volScalarField& muu_ = tmuu.ref();
657 scalarField& muuCells = muu_.primitiveFieldRef();
658 const scalarField& pCells = this->p_;
659 const scalarField& TuCells = this->Tu_;
660
661 forAll(muuCells, celli)
662 {
663 muuCells[celli] = this->cellReactants(celli).mu
664 (
665 pCells[celli],
666 TuCells[celli]
667 );
668 }
669
671
672 forAll(muuBf, patchi)
673 {
674 fvPatchScalarField& pMuu = muuBf[patchi];
675 const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
676 const fvPatchScalarField& pTu = this->Tu_.boundaryField()[patchi];
677
678 forAll(pMuu, facei)
679 {
680 pMuu[facei] = this->patchFaceReactants(patchi, facei).mu
681 (
682 pp[facei],
683 pTu[facei]
684 );
685 }
686 }
687
688 return tmuu;
689}
690
691
692template<class BasicPsiThermo, class MixtureType>
695{
697 (
699 (
701 (
702 "mub",
703 this->T_.time().timeName(),
704 this->T_.db(),
707 false
708 ),
709 this->T_.mesh(),
710 dimensionSet(1, -1, -1, 0, 0)
711 )
712 );
713
714 volScalarField& mub_ = tmub.ref();
715 scalarField& mubCells = mub_.primitiveFieldRef();
716 const volScalarField Tb_(Tb());
717 const scalarField& pCells = this->p_;
718 const scalarField& TbCells = Tb_;
719
720 forAll(mubCells, celli)
721 {
722 mubCells[celli] = this->cellProducts(celli).mu
723 (
724 pCells[celli],
725 TbCells[celli]
726 );
727 }
728
730
731 forAll(mubBf, patchi)
732 {
733 fvPatchScalarField& pMub = mubBf[patchi];
734 const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
735 const fvPatchScalarField& pTb = Tb_.boundaryField()[patchi];
736
737 forAll(pMub, facei)
738 {
739 pMub[facei] = this->patchFaceProducts(patchi, facei).mu
740 (
741 pp[facei],
742 pTb[facei]
743 );
744 }
745 }
746
747 return tmub;
748}
749
750
751// ************************************************************************* //
volScalarField & he
Definition: YEEqn.H:52
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:109
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Enthalpy/Internal energy for a mixture.
Definition: heThermo.H:56
virtual tmp< volScalarField > mub() const
Dynamic viscosity of burnt gas [kg/ms].
virtual void correct()
Update properties.
virtual tmp< volScalarField > muu() const
Dynamic viscosity of unburnt gas [kg/ms].
virtual volScalarField & heu()
Update properties based on T.
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 ~heheuPsiThermo()
Destructor.
Foam::psiuReactionThermo.
ThermoType thermoType
The type of thermodynamics this mixture is instantiated for.
Definition: pureMixture.H:66
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
dynamicFvMesh & mesh
const word dictName("faMeshDefinition")
const cellShapeList & cells
word timeName
Definition: getTimeIndex.H:3
#define DebugInfo
Report an information message using Foam::Info.
#define DebugInFunction
Report an information message using Foam::Info.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
fvPatchField< scalar > fvPatchScalarField
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333