Loading...
Searching...
No Matches
alphatWallBoilingWallFunctionFvPatchScalarField.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) 2015-2018 OpenFOAM Foundation
9 Copyright (C) 2018-2021 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
30#include "fvPatchFieldMapper.H"
32
33#include "phaseSystem.H"
35#include "ThermalDiffusivity.H"
37#include "saturationModel.H"
38#include "wallFvPatch.H"
41
42using namespace Foam::constant::mathematical;
43
44// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45
46const Foam::Enum
47<
48 Foam::compressible::
49 alphatWallBoilingWallFunctionFvPatchScalarField::phaseType
50>
51Foam::compressible::
52alphatWallBoilingWallFunctionFvPatchScalarField::phaseTypeNames_
53{
54 { phaseType::vaporPhase, "vapor" },
55 { phaseType::liquidPhase, "liquid" },
56};
57
58
59// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60
61namespace Foam
62{
63namespace compressible
64{
65
66// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
67
70(
71 const fvPatch& p,
73)
74:
76 otherPhaseName_("vapor"),
77 phaseType_(liquidPhase),
78 relax_(),
79 AbyV_(p.size(), 0),
80 alphatConv_(p.size(), 0),
81 dDep_(p.size(), 1e-5),
82 qq_(p.size(), 0),
83 K_(4),
84 partitioningModel_(nullptr),
85 nucleationSiteModel_(nullptr),
86 departureDiamModel_(nullptr),
87 departureFreqModel_(nullptr),
88 nucleatingModel_(nullptr),
89 filmBoilingModel_(nullptr),
90 LeidenfrostModel_(nullptr),
91 CHFModel_(nullptr),
92 CHFSoobModel_(nullptr),
93 MHFModel_(nullptr),
94 TDNBModel_(nullptr),
95 wp_(1),
96 liquidTatYplus_(false),
97 regimeTypes_(p.size(), -1)
98{
99 AbyV_ = this->patch().magSf();
100 forAll(AbyV_, facei)
102 const label faceCelli = this->patch().faceCells()[facei];
103 AbyV_[facei] /= iF.mesh().V()[faceCelli];
104 }
105}
106
107
110(
111 const fvPatch& p,
113 const dictionary& dict
114)
115:
116 alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField(p, iF, dict),
117 otherPhaseName_(dict.get<word>("otherPhase")),
118 phaseType_(phaseTypeNames_.get("phaseType", dict)),
119 relax_(Function1<scalar>::New("relax", dict, &db())),
120 AbyV_(p.size(), 0),
121 alphatConv_(p.size(), 0),
122 dDep_(p.size(), 1e-5),
123 qq_(p.size(), 0),
124 K_(4),
125 partitioningModel_(nullptr),
126 nucleationSiteModel_(nullptr),
127 departureDiamModel_(nullptr),
128 departureFreqModel_(nullptr),
129 nucleatingModel_(nullptr),
130 filmBoilingModel_(nullptr),
131 LeidenfrostModel_(nullptr),
132 CHFModel_(nullptr),
133 CHFSoobModel_(nullptr),
134 MHFModel_(nullptr),
135 TDNBModel_(nullptr),
136 wp_(1),
137 liquidTatYplus_(dict.getOrDefault<bool>("liquidTatYplus", false)),
138 regimeTypes_(p.size(), -1)
139{
140 // Check that otherPhaseName != this phase
141 if (internalField().group() == otherPhaseName_)
142 {
144 << "otherPhase should be the name of the vapor phase that "
145 << "corresponds to the liquid base of vice versa" << nl
146 << "This phase: " << internalField().group() << nl
147 << "otherPhase: " << otherPhaseName_
148 << abort(FatalError);
149 }
150
151 partitioningModel_ =
153 (
154 dict.subDict("partitioningModel")
155 );
156
157 switch (phaseType_)
158 {
159 case vaporPhase:
160 {
161 dmdt_ = Zero;
162
163 break;
164 }
165 case liquidPhase:
166 {
167 partitioningModel_ =
169 (
170 dict.subDict("partitioningModel")
171 );
172
173 // If nucleating model is specifed use it. Otherwise use
174 // RPI wall boiling model
175 const dictionary* nucleatingDict
176 = dict.findDict("nucleateFluxModel");
177
178 if (nucleatingDict)
179 {
180 nucleatingModel_ =
182 }
183 else
184 {
185 nucleationSiteModel_ =
187 (
188 dict.subDict("nucleationSiteModel")
189 );
190
191 departureDiamModel_ =
193 (
194 dict.subDict("departureDiamModel")
195 );
196
197 departureFreqModel_ =
199 (
200 dict.subDict("departureFreqModel")
201 );
202 }
203
204
205 const dictionary* LeidenfrostDict =
206 dict.findDict("LeidenfrostModel");
207
208 if (LeidenfrostDict)
209 {
210 LeidenfrostModel_ =
212 }
213
214 const dictionary* CHFDict = dict.findDict("CHFModel");
215
216 if (CHFDict)
217 {
218 CHFModel_ =
220 }
221
222 const dictionary* HFSubCoolDict = dict.findDict("CHFSubCoolModel");
223
224 if (HFSubCoolDict)
225 {
226 CHFSoobModel_ =
228 }
229
230 const dictionary* MHFDict = dict.findDict("MHFModel");
231
232 if (MHFDict)
233 {
234 MHFModel_ =
236 }
237
238 const dictionary* TDNBDict = dict.findDict("TDNBModel");
239
240 if (TDNBDict)
241 {
242 TDNBModel_ =
244 }
245
246 const dictionary* filmDict = dict.findDict("filmBoilingModel");
247
248 if (filmDict)
249 {
250 filmBoilingModel_ =
252 }
253
254 if (dict.found("dDep"))
255 {
256 dDep_ = scalarField("dDep", dict, p.size());
257 }
258
259 dict.readIfPresent("K", K_);
260
261 dict.readIfPresent("wp", wp_);
262
263 if (dict.found("qQuenching"))
264 {
265 qq_ = scalarField("qQuenching", dict, p.size());
266 }
267
268 break;
269 }
270 }
271
272 if (dict.found("alphatConv"))
273 {
274 alphatConv_ = scalarField("alphatConv", dict, p.size());
275 }
276
277 AbyV_ = this->patch().magSf();
278 forAll(AbyV_, facei)
280 const label faceCelli = this->patch().faceCells()[facei];
281 AbyV_[facei] /= iF.mesh().V()[faceCelli];
282 }
283}
284
285
288(
289 const alphatWallBoilingWallFunctionFvPatchScalarField& psf,
290 const fvPatch& p,
292 const fvPatchFieldMapper& mapper
293)
294:
295 alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
296 (
297 psf,
298 p,
299 iF,
300 mapper
301 ),
302 otherPhaseName_(psf.otherPhaseName_),
303 phaseType_(psf.phaseType_),
304 relax_(psf.relax_.clone()),
305 AbyV_(psf.AbyV_),
306 alphatConv_(psf.alphatConv_, mapper),
307 dDep_(psf.dDep_, mapper),
308 qq_(psf.qq_, mapper),
309 K_(psf.K_),
310 partitioningModel_(psf.partitioningModel_),
311 nucleationSiteModel_(psf.nucleationSiteModel_),
312 departureDiamModel_(psf.departureDiamModel_),
313 nucleatingModel_(psf.nucleatingModel_),
314 filmBoilingModel_(psf.filmBoilingModel_),
315 LeidenfrostModel_(psf.LeidenfrostModel_),
316 CHFModel_(psf.CHFModel_),
317 CHFSoobModel_(psf.CHFSoobModel_),
318 MHFModel_(psf.MHFModel_),
319 TDNBModel_(psf.TDNBModel_),
320 wp_(psf.wp_),
321 liquidTatYplus_(psf.liquidTatYplus_),
322 regimeTypes_(psf.regimeTypes_)
323{}
324
325
328(
330)
331:
333 otherPhaseName_(psf.otherPhaseName_),
334 phaseType_(psf.phaseType_),
335 relax_(psf.relax_.clone()),
336 AbyV_(psf.AbyV_),
337 alphatConv_(psf.alphatConv_),
338 dDep_(psf.dDep_),
339 qq_(psf.qq_),
340 K_(psf.K_),
341 partitioningModel_(psf.partitioningModel_),
342 nucleationSiteModel_(psf.nucleationSiteModel_),
343 departureDiamModel_(psf.departureDiamModel_),
344 nucleatingModel_(psf.nucleatingModel_),
345 filmBoilingModel_(psf.filmBoilingModel_),
346 LeidenfrostModel_(psf.LeidenfrostModel_),
347 CHFModel_(psf.CHFModel_),
348 CHFSoobModel_(psf.CHFSoobModel_),
349 MHFModel_(psf.MHFModel_),
350 TDNBModel_(psf.TDNBModel_),
351 wp_(psf.wp_),
352 liquidTatYplus_(psf.liquidTatYplus_),
353 regimeTypes_(psf.regimeTypes_)
354{}
355
356
359(
362)
363:
365 otherPhaseName_(psf.otherPhaseName_),
366 phaseType_(psf.phaseType_),
367 relax_(psf.relax_.clone()),
368 AbyV_(psf.AbyV_),
369 alphatConv_(psf.alphatConv_),
370 dDep_(psf.dDep_),
371 qq_(psf.qq_),
372 K_(psf.K_),
373 partitioningModel_(psf.partitioningModel_),
374 nucleationSiteModel_(psf.nucleationSiteModel_),
375 departureDiamModel_(psf.departureDiamModel_),
376 nucleatingModel_(psf.nucleatingModel_),
377 filmBoilingModel_(psf.filmBoilingModel_),
378 LeidenfrostModel_(psf.LeidenfrostModel_),
379 CHFModel_(psf.CHFModel_),
380 CHFSoobModel_(psf.CHFSoobModel_),
381 MHFModel_(psf.MHFModel_),
382 TDNBModel_(psf.TDNBModel_),
383 wp_(psf.wp_),
384 liquidTatYplus_(psf.liquidTatYplus_),
385 regimeTypes_(psf.regimeTypes_)
386{}
387
388
389// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
390
393{
394 if (phasePair == phasePairKey(otherPhaseName_, internalField().group()))
396 return true;
397 }
398
399 return false;
400}
401
403dmdt(const phasePairKey& phasePair) const
404{
406 {
407 return dmdt_;
408 }
409
411 << " dmdt requested for invalid phasePair!"
412 << abort(FatalError);
413
414 return dmdt_;
415}
416
418mDotL(const phasePairKey& phasePair) const
419{
421 {
422 return mDotL_;
423 }
424
426 << " mDotL requested for invalid phasePair!"
427 << abort(FatalError);
428
429 return mDotL_;
430}
431
432
434{
435 if (updated())
436 {
437 return;
438 }
439
440 // Check that partitioningModel has been constructed
441 if (!partitioningModel_)
442 {
444 << "partitioningModel has not been constructed!"
445 << abort(FatalError);
446 }
447
448 // Lookup the fluid model
449 const phaseSystem& fluid =
451 (
452 db().lookupObject<phaseSystem>("phaseProperties")
453 );
454
455 const auto& satModel =
456 db().lookupObject<saturationModel>("saturationModel");
457
458 const label patchi = patch().index();
459
460 const scalar t = this->db().time().timeOutputValue();
461 const scalar relax = relax_->value(t);
462
463 switch (phaseType_)
464 {
465 case vaporPhase:
466 {
467 const phaseModel& vapor
468 (
469 fluid.phases()[internalField().group()]
470 );
471
472 const tmp<scalarField> talphaw = vapor.thermo().alpha(patchi);
473 const scalarField& alphaw = talphaw();
474
475 const fvPatchScalarField& hewv =
476 vapor.thermo().he().boundaryField()[patchi];
477
478 // Vapor Liquid phase fraction at the wall
479 const scalarField vaporw
480 (
481 max(vapor.boundaryField()[patchi], scalar(1e-16))
482 );
483 const scalarField liquidw(1.0 - vaporw);
484
485 // NOTE! Assumes 1-thisPhase for liquid fraction in
486 // multiphase simulations
487 const scalarField fLiquid
488 (
489 partitioningModel_->fLiquid(1-vaporw)
490 );
491
492 // Convective thermal diffusivity for single phase
493 const scalarField alphatv(calcAlphat(*this));
494
495 forAll(*this, i)
496 {
497 this->operator[](i) =
498 (
499 (1 - fLiquid[i])*(alphatv[i] + alphaw[i])
500 /max(vaporw[i], scalar(1e-8))
501 );
502 }
503
504 if (debug)
505 {
506 Info<< "alphat for vapour : " << nl << endl;
507
508 Info<< " alphatEffv: " << gMinMax(vaporw*(*this + alphaw))
509 << endl;
510
511 const scalarField qEff(vaporw*(*this + alphaw)*hewv.snGrad());
512
513 Info<< " qEffVap: " << gMinMax(qEff) << endl;
514
515 scalar Qeff = gWeightedSum(patch().magSf(), qEff);
516 Info<< " Effective heat transfer rate to vapor:" << Qeff
517 << nl << endl;
518 }
519 break;
520 }
521 case liquidPhase:
522 {
523 // Check that nucleationSiteModel has been constructed
524 if (!nucleatingModel_)
525 {
526 if (!nucleationSiteModel_)
527 {
529 << "nucleationSiteModel has not been constructed!"
530 << abort(FatalError);
531 }
532
533 // Check that departureDiameterModel has been constructed
534 if (!departureDiamModel_)
535 {
537 << "departureDiameterModel has not been constructed!"
538 << abort(FatalError);
539 }
540
541 // Check that nucleationSiteModel has been constructed
542 if (!departureFreqModel_)
543 {
545 << "departureFrequencyModel has not been constructed!"
546 << abort(FatalError);
547 }
548 }
549
550 const phaseModel& liquid
551 (
552 fluid.phases()[internalField().group()]
553 );
554
555 const phaseModel& vapor(fluid.phases()[otherPhaseName_]);
556
557 // Retrieve turbulence properties from models
558 const auto& turbModel =
559 db().lookupObject<phaseCompressibleTurbulenceModel>
560 (
562 (
564 liquid.name()
565 )
566 );
567 const auto& vaporTurbModel =
568 db().lookupObject<phaseCompressibleTurbulenceModel>
569 (
571 (
573 vapor.name()
574 )
575 );
576
577 const tmp<scalarField> tnutw = turbModel.nut(patchi);
578
579 const scalar Cmu25 = pow025(Cmu_);
580
581 const scalarField& y = turbModel.y()[patchi];
582
583 const tmp<scalarField> tmuw = turbModel.mu(patchi);
584 const scalarField& muw = tmuw();
585
586 const tmp<scalarField> talphaw = liquid.thermo().alphahe(patchi);
587 const scalarField& alphaw = talphaw();
588
589 const tmp<volScalarField> tk = turbModel.k();
590 const volScalarField& k = tk();
591 const fvPatchScalarField& kw = k.boundaryField()[patchi];
592
593 const fvPatchVectorField& Uw =
594 turbModel.U().boundaryField()[patchi];
595 const scalarField magUp(mag(Uw.patchInternalField() - Uw));
596 const scalarField magGradUw(mag(Uw.snGrad()));
597
598 const fvPatchScalarField& rhow =
599 turbModel.rho().boundaryField()[patchi];
600
601
602 const fvPatchScalarField& Tw =
603 liquid.thermo().T().boundaryField()[patchi];
604 const scalarField Tc(Tw.patchInternalField());
605
606 const scalarField uTau(Cmu25*sqrt(kw));
607
608 const scalarField yPlus(uTau*y/(muw/rhow));
609
610 const scalarField Pr(muw/alphaw);
611
612 // Molecular-to-turbulent Prandtl number ratio
613 const scalarField Prat(Pr/Prt_);
614
615 // Thermal sublayer thickness
616 const scalarField P(this->Psmooth(Prat));
617
618 const scalarField yPlusTherm(this->yPlusTherm(P, Prat));
619
620 const fvPatchScalarField& rhoVaporw =
621 vaporTurbModel.rho().boundaryField()[patchi];
622
623 tmp<volScalarField> tCp = liquid.thermo().Cp();
624 const volScalarField& Cp = tCp();
625 const fvPatchScalarField& Cpw = Cp.boundaryField()[patchi];
626
627 // Saturation temperature
628 const tmp<volScalarField> tTsat =
629 satModel.Tsat(liquid.thermo().p());
630
631 const volScalarField& Tsat = tTsat();
632 const fvPatchScalarField& Tsatw = Tsat.boundaryField()[patchi];
633 const scalarField Tsatc(Tsatw.patchInternalField());
634
635 const fvPatchScalarField& pw =
636 liquid.thermo().p().boundaryField()[patchi];
637
638 const fvPatchScalarField& hew =
639 liquid.thermo().he().boundaryField()[patchi];
640
641 const scalarField hwLiqSat
642 (
643 liquid.thermo().he().member() == "e"
644 ? liquid.thermo().he(pw, Tsatc, patchi)
645 + pw/rhow.patchInternalField()
646 : liquid.thermo().he(pw, Tsatc, patchi)
647 );
648
649 const scalarField L
650 (
651 vapor.thermo().he().member() == "e"
652 ? vapor.thermo().he(pw, Tsatc, patchi) + pw/rhoVaporw - hwLiqSat
653 : vapor.thermo().he(pw, Tsatc, patchi) - hwLiqSat
654 );
655
656 // Liquid phase fraction at the wall
657 const scalarField liquidw(liquid.boundaryField()[patchi]);
658
659 // Partition between phases
660 const scalarField fLiquid(partitioningModel_->fLiquid(liquidw));
661
662 scalarField Tl(Tc);
663 // Liquid temperature at y+=250 is estimated from logarithmic
664 // thermal wall function (Koncar, Krepper & Egorov, 2005)
665 if (liquidTatYplus_)
666 {
667 const scalarField Tplus_y250
668 (
669 Prt_*(log(E_*250)/kappa_ + P)
670 );
671 const scalarField Tplus(Prt_*(log(E_*yPlus)/kappa_ + P));
672 scalarField Tl(Tw - (Tplus_y250/Tplus)*(Tw - Tc));
673 Tl = max(Tc - 40, Tl);
674 }
675
676 // Film, transient boiling regimes
677 scalarField Qtb(this->size(), 0);
678 scalarField tDNB(this->size(), GREAT);
679 scalarField TLeiden(this->size(), GREAT);
680 scalarField htcFilmBoiling(this->size(), 0);
681
682 if
683 (
684 CHFModel_
685 && CHFSoobModel_
686 && TDNBModel_
687 && MHFModel_
688 && LeidenfrostModel_
689 && filmBoilingModel_
690 )
691 {
692
693 const scalarField CHF
694 (
695 CHFModel_->CHF
696 (
697 liquid,
698 vapor,
699 patchi,
700 Tl,
701 Tsatw,
702 L
703 )
704 );
705
706 if (debug == 2)
707 {
708 Info << "CHF : " << CHF << endl;
709 }
710
711 // Effect of sub-cooling to the CHF in saturated conditions
712 const scalarField CHFSubCool
713 (
714 CHFSoobModel_->CHFSubCool
715 (
716 liquid,
717 vapor,
718 patchi,
719 Tl,
720 Tsatw,
721 L
722 )
723 );
724
725 if (debug == 2)
726 {
727 Info << "CHF Sub Cool factor : " << CHFSubCool << endl;
728 }
729
730 const scalarField CHFtotal(CHF*CHFSubCool);
731
732 tDNB =
733 TDNBModel_->TDNB
734 (
735 liquid,
736 vapor,
737 patchi,
738 Tl,
739 Tsatw,
740 L
741 );
742
743 if (debug == 2)
744 {
745 Info<< "Temperature departure from biling : "
746 << tDNB << endl;
747 }
748
749 const scalarField MHF
750 (
751 MHFModel_->MHF
752 (
753 liquid,
754 vapor,
755 patchi,
756 Tl,
757 Tsatw,
758 L
759 )
760 );
761
762 if (debug == 2)
763 {
764 Info<< "MHF : " << MHF << endl;
765 }
766
767 TLeiden =
768 LeidenfrostModel_->TLeid
769 (
770 liquid,
771 vapor,
772 patchi,
773 Tl,
774 Tsatw,
775 L
776 );
777
778 if (debug == 2)
779 {
780 Info<< "Leidenfrost Temp : " << TLeiden << endl;
781 }
782
783 // htc for film boiling
784 htcFilmBoiling =
785 filmBoilingModel_->htcFilmBoil
786 (
787 liquid,
788 vapor,
789 patchi,
790 Tl,
791 Tsatw,
792 L
793 );
794
795 if (debug == 2)
796 {
797 Info<< "Htc film boiling : " << htcFilmBoiling << endl;
798 }
799
800 // htc for film transition boiling
801 // Indicator between CHF (phi = 0) and MHF (phi = 1)
802 const scalarField phi
803 (
804 min
805 (
806 max
807 (
808 wp_*(Tw - tDNB)/(TLeiden - tDNB),
809 scalar(0)
810 ),
811 scalar(1)
812 )
813 );
814
815 Qtb = CHFtotal*(1 - phi) + phi*MHF;
816
817 }
818
819 // Convective heat transfer area for Sub-cool boiling
820 scalarField A1(this->size(), 0);
821 qq_ = Zero;
822 scalarField dmdtSubCooling(this->size(), 0);
823
824 if (nucleatingModel_)
825 {
826 dmdtSubCooling =
827 nucleatingModel_->qNucleate
828 (
829 liquid,
830 vapor,
831 patchi,
832 Tl,
833 Tsatw,
834 L
835 )*AbyV_/L;
836
837
838 dmdtSubCooling *= fLiquid;
839 }
840 else
841 {
842 // Sub-cool boiling Nucleation
843 const scalarField N
844 (
845 nucleationSiteModel_->N
846 (
847 liquid,
848 vapor,
849 patchi,
850 Tl,
851 Tsatw,
852 L
853 )
854 );
855
856 // Bubble departure diameter:
857 dDep_ = departureDiamModel_->dDeparture
858 (
859 liquid,
860 vapor,
861 patchi,
862 Tl,
863 Tsatw,
864 L
865 );
866
867 // Bubble departure frequency:
868 const scalarField fDep
869 (
870 departureFreqModel_->fDeparture
871 (
872 liquid,
873 vapor,
874 patchi,
875 dDep_
876 )
877 );
878
879 scalarField Ja
880 (
881 rhow*Cpw*(Tsatw - Tl)/(rhoVaporw*L)
882 );
883
884 scalarField Al
885 (
886 fLiquid*4.8*exp(min(-Ja/80, log(VGREAT)))
887 );
888
889 scalarField A2
890 (
891 min(pi*sqr(dDep_)*N*Al/4, scalar(1))
892 );
893
894 A1 = max(1 - A2, scalar(1e-4));
895
896 // Following Bowring(1962)
897 scalarField A2E
898 (
899 min(pi*sqr(dDep_)*N*Al/4, scalar(5))
900 );
901
902 dmdtSubCooling =
903 (
904 (1.0/6.0)*A2E*dDep_*rhoVaporw*fDep*AbyV_
905 );
906
907 scalarField hQ
908 (
909 2*(alphaw*Cpw)*fDep
910 *sqrt
911 (
912 (0.8/max(fDep, SMALL))/(pi*alphaw/rhow)
913 )
914 );
915
916 // Quenching heat flux in Sub-cool boiling
917 qq_ =
918 (
919 (1 - relax)*qq_
920 + relax*A2*hQ*max(Tw - Tl, scalar(0))
921 );
922 }
923
924 // Convective thermal diffusivity for single phase
925 alphatConv_ = calcAlphat(alphatConv_);
926
927 const scalarField hewSn(hew.snGrad());
928
929 // AlphaEff for film regime
930 scalarField alphaFilm(this->size(), 0);
931
932 // Use to identify regimes per face
933 regimeTypes_ = -1;
934
935 forAll(*this, i)
936 {
937 if (Tw[i] > Tsatw[i])
938 {
939 // Sub-cool boiling
940 if (Tw[i] < tDNB[i])
941 {
942 // Sub-cool boiling
943 regimeTypes_[i] = regimeType::subcool;
944
945 dmdt_[i] =
946 (
947 (1 - relax)*dmdt_[i] + relax*dmdtSubCooling[i]
948 );
949
950 // Volumetric source in the near wall cell due to
951 // the wall boiling
952 mDotL_[i] = dmdt_[i]*L[i];
953
954 this->operator[](i) =
955 (
956 max
957 (
958 A1[i]*alphatConv_[i]
959 + (
960 (qq_[i] + mDotL_[i]/AbyV_[i])
961 / max(hewSn[i], scalar(1e-16))
962 )
963 /max(liquidw[i], scalar(1e-8)),
964 scalar(1e-8)
965 )
966 );
967
968 if (debug == 2)
969 {
970 Info<< "Sub-cool boiling: " << nl
971 << " fraction Liq: " << fLiquid[i] << nl
972 << " Heat flux: "
973 << (qq_[i] + mDotL_[i]/AbyV_[i]) << nl
974 << " delta Tsub: " << (Tw[i] - Tsatw[i])
975 << endl;
976 }
977 }
978 else if (Tw[i] > tDNB[i] && Tw[i] < TLeiden[i])
979 {
980 // transient boiling
981 regimeTypes_[i] = regimeType::transient;
982
983 // No convective heat transfer
984 alphatConv_[i] = 0.0;
985
986 // transient boiling
987 dmdt_[i] =
988 fLiquid[i]
989 *(
990 relax*Qtb[i]*AbyV_[i]/L[i]
991 + (1 - relax)*dmdt_[i]
992 );
993
994
995 mDotL_[i] = dmdt_[i]*L[i];
996
997 // No quenching flux
998 qq_[i] = 0.0;
999
1000 this->operator[](i) =
1001 max
1002 (
1003 (
1004 mDotL_[i]/AbyV_[i]
1005 /max(hewSn[i], scalar(1e-16))
1006 )/max(liquidw[i], scalar(1e-8)),
1007 scalar(1e-8)
1008 );
1009
1010 if (debug == 2)
1011 {
1012 Info<< "Transient boiling: " << nl
1013 << " fraction Liq: " << fLiquid[i] << nl
1014 << " Heat flux: " << Qtb[i] << nl
1015 << " delta Tsub: " << (Tw[i] - Tsatw[i])
1016 << endl;
1017 }
1018
1019 }
1020 else if (Tw[i] > TLeiden[i])
1021 {
1022 regimeTypes_[i] = regimeType::film; // film boiling
1023
1024 // No convective heat transfer
1025 alphatConv_[i] = 0.0;
1026
1027 // Film boiling
1028 dmdt_[i] =
1029 fLiquid[i]
1030 *(
1031 relax*htcFilmBoiling[i]
1032 *max(Tw[i] - Tsatw[i], scalar(0))
1033 *AbyV_[i]/L[i]
1034 + (1 - relax)*dmdt_[i]
1035 );
1036
1037
1038 mDotL_[i] = dmdt_[i]*L[i];
1039
1040 // No quenching flux
1041 qq_[i] = 0.0;
1042
1043 alphaFilm[i] =
1044 (
1045 mDotL_[i]/AbyV_[i]/max(hewSn[i], scalar(1e-16))
1046 );
1047
1048 // alphat is added alphal and multiplied by phase
1049 // alphaFilm in the coupled BC. We subtract
1050 // alpha and divide by phase to get a net alphaFilm
1051 this->operator[](i) =
1052 (
1053 alphaFilm[i]/max(liquidw[i], scalar(1e-8))
1054 - alphaw[i]
1055 );
1056
1057 if (debug == 2)
1058 {
1059 Info<< "Film boiling: " << nl
1060 << " fraction Liq: " << fLiquid[i] << nl
1061 << " Heat flux: "
1062 << htcFilmBoiling[i]*(Tw[i] - Tsatw[i]) << nl
1063 << " delta Tsub: " << (Tw[i] - Tsatw[i])
1064 << endl;
1065 }
1066 }
1067 }
1068 else
1069 {
1070 // Tw below Tsat. No boiling single phase convection
1071 // Single phase
1072 regimeTypes_[i] = regimeType::nonBoiling;
1073
1074 qq_[i] = 0.0;
1075 mDotL_[i] = 0.0;
1076 dmdt_[i] = 0.0;
1077
1078 // Turbulente thermal diffusivity for single phase.
1079 this->operator[](i) =
1080 (
1081 max
1082 (
1083 fLiquid[i]*(alphatConv_[i])
1084 /max(liquidw[i], scalar(1e-8)),
1085 scalar(1e-8)
1086 )
1087 );
1088 }
1089 }
1090
1091 if (debug)
1092 {
1093 const scalarField qEff
1094 (
1095 fLiquid*liquidw*(*this + alphaw)*hew.snGrad()
1096 );
1097
1098 Info<< "alphat for liquid: " << nl << nl;
1099 Info<< " qEffLiq: " << gMinMax(qEff) << nl;
1100 Info<< " alphatl: " << gMinMax(*this) << nl;
1101 Info<< " dmdt: " << gMinMax(dmdt_) << nl;
1102 Info<< " alphatlEff: "
1103 << gMinMax(liquidw*(*this + alphaw)) << nl;
1104
1105 Info<< " Effective heat transfer rate to liquid: "
1106 << gWeightedSum(patch().magSf(), qEff)
1107 << nl << endl;
1108
1109 if (debug == 2)
1110 {
1111 scalarField nSubCools(this->size(), 0);
1112 scalarField nTransients(this->size(), 0);
1113 scalarField nFilms(this->size(), 0);
1114 scalarField nNonBoilings(this->size(), 0);
1115
1116 forAll(*this, i)
1117 {
1118 //faceRegimes[i] = regimeTypes[i];
1119 switch (regimeTypes_[i])
1120 {
1121 case regimeType::subcool:
1122 nSubCools[i] = 1;
1123 break;
1124
1125 case regimeType::transient:
1126 nTransients[i] = 1;
1127 break;
1128
1129 case regimeType::film:
1130 nFilms[i] = 1;
1131 break;
1132
1133 case regimeType::nonBoiling:
1134 nNonBoilings[i] = 1;
1135 break;
1136 }
1137 }
1138
1139 scalar nSubCool(gSum(nSubCools));
1140 scalar nTransient(gSum(nTransients));
1141 scalar nFilm(gSum(nFilms));
1142 scalar nNonBoiling(gSum(nNonBoilings));
1143
1144 Info<< "Faces regime : " << nl << endl;
1145
1146 Info<< " sub Cool faces : " << nSubCool << endl;
1147 Info<< " transient faces : " << nTransient << endl;
1148 Info<< " film faces : " << nFilm << endl;
1149 Info<< " non-Boiling faces : " << nNonBoiling << endl;
1150 Info<< " total faces : "
1151 << nSubCool + nTransient + nFilm + nNonBoiling
1152 << endl << nl;
1153
1154 const scalarField qc
1155 (
1156 nNonBoilings*fLiquid*(alphatConv_ + alphaw)
1157 *hew.snGrad()
1158 );
1159
1160 scalar Qc = gWeightedSum(patch().magSf(), qc);
1161 Info<< " Convective heat transfer: " << Qc << endl;
1162
1163 const scalarField qFilm
1164 (
1165 relax*fLiquid*nFilms*htcFilmBoiling*(Tw - Tsatw)
1166 );
1167
1168 scalar QFilm = gWeightedSum(patch().magSf(), qFilm);
1169 Info<< " Film boiling heat transfer: " << QFilm << endl;
1170
1171 Info<< " Htc Film Boiling coeff: "
1172 << gMin(nFilms*htcFilmBoiling)
1173 << " - "
1174 << gMax(nFilms*htcFilmBoiling) << endl;
1175
1176 scalar Qtbtot =
1177 gSum(fLiquid*nTransients*Qtb*patch().magSf());
1178 Info<< " Transient boiling heat transfer:" << Qtbtot
1179 << endl;
1180
1181
1182 Info<< " TDNB: " << gMin(tDNB) << " - " << gMax(tDNB)
1183 << endl;
1184
1185 const scalarField qSubCool
1186 (
1187 fLiquid*nSubCools*
1188 (
1189 A1*alphatConv_*hew.snGrad() + qe() + qq()
1190 )
1191 );
1192
1193 scalar QsubCool = gWeightedSum(patch().magSf(), qSubCool);
1194
1195 Info<< " Sub Cool boiling heat transfer: " << QsubCool
1196 << endl;
1197
1198 Info<< nl;
1199 }
1200 }
1201 }
1202 break;
1203 default:
1204 {
1206 << "Unknown phase type. Valid types are: "
1207 << phaseTypeNames_ << nl << exit(FatalError);
1209 }
1210
1211 fixedValueFvPatchScalarField::updateCoeffs();
1212}
1213
1214
1216{
1218
1219 os.writeEntry("phaseType", phaseTypeNames_[phaseType_]);
1220
1221 relax_->writeData(os);
1222
1223 os.beginBlock("partitioningModel");
1224 partitioningModel_->write(os);
1225 os.endBlock();
1226
1227 switch (phaseType_)
1228 {
1229 case vaporPhase:
1230 break;
1231 case liquidPhase:
1232 {
1233 if (nucleationSiteModel_)
1234 {
1235 os.beginBlock("nucleationSiteModel");
1236 nucleationSiteModel_->write(os);
1237 os.endBlock();
1238 }
1239
1240 if (departureDiamModel_)
1241 {
1242 os.beginBlock("departureDiamModel");
1243 departureDiamModel_->write(os);
1244 os.endBlock();
1245 }
1246
1247 if (departureFreqModel_)
1248 {
1249 os.beginBlock("departureFreqModel");
1250 departureFreqModel_->write(os);
1251 os.endBlock();
1252 }
1253
1254 if (nucleatingModel_)
1255 {
1256 os.beginBlock("nucleateFluxModel");
1257 nucleatingModel_->write(os);
1258 os.endBlock();
1259 }
1260
1261 if (filmBoilingModel_)
1262 {
1263 os.beginBlock("filmBoilingModel");
1264 filmBoilingModel_->write(os);
1265 os.endBlock();
1266 }
1267
1268 if (LeidenfrostModel_)
1269 {
1270 os.beginBlock("LeidenfrostModel");
1271 LeidenfrostModel_->write(os);
1272 os.endBlock();
1273 }
1274
1275 if (CHFModel_)
1276 {
1277 os.beginBlock("CHFModel");
1278 CHFModel_->write(os);
1279 os.endBlock();
1280 }
1281
1282 if (CHFSoobModel_)
1283 {
1284 os.beginBlock("CHFSubCoolModel");
1285 CHFSoobModel_->write(os);
1286 os.endBlock();
1287 }
1288
1289 if (MHFModel_)
1290 {
1291 os.beginBlock("MHFModel");
1292 MHFModel_->write(os);
1293 os.endBlock();
1294 }
1295
1296 if (TDNBModel_)
1297 {
1298 os.beginBlock("TDNBModel");
1299 TDNBModel_->write(os);
1300 os.endBlock();
1301 }
1302
1303 os.writeEntry("K", K_);
1304 os.writeEntry("wp", wp_);
1305 os.writeEntry("liquidTatYplus", liquidTatYplus_);
1306 break;
1307 }
1308 }
1309
1310 os.writeEntry("otherPhase", otherPhaseName_);
1311
1312 dmdt_.writeEntry("dmdt", os);
1313 dDep_.writeEntry("dDep", os);
1314 qq_.writeEntry("qQuenching", os);
1315 alphatConv_.writeEntry("alphatConv", os);
1318}
1319
1320
1321// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1322
1324(
1326 alphatWallBoilingWallFunctionFvPatchScalarField
1327);
1328
1329// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1330
1331} // End namespace compressible
1332} // End namespace Foam
1333
1334// ************************************************************************* //
scalar y
label k
Macros for easy insertion into run-time selection tables.
twoPhaseSystem & fluid
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Mesh & mesh() const noexcept
Return const reference to mesh.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
tmp< Field< Type > > T() const
Return the field transpose (only defined for second rank tensors).
Definition Field.C:745
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition Function1.H:92
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
static word member(const word &name)
Return member (name without the extension).
Definition IOobject.C:313
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
virtual const volScalarField & alpha() const
Thermal diffusivity for enthalpy of mixture [kg/m/s].
const volScalarField & rho() const
Return the density field.
This boundary condition provides a thermal wall function for turbulent thermal diffusivity (usuallyal...
tmp< scalarField > calcAlphat(const scalarField &prevAlphat) const
Update turbulent thermal diffusivity.
tmp< scalarField > yPlusTherm(const scalarField &P, const scalarField &Prat) const
Calculate y+ at the edge of the thermal laminar sublayer.
alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual const scalarField & mDotL() const
Return the enthalpy source due to phase-change.
tmp< scalarField > qe() const
Return the evaporation surface heat flux [W/m2].
alphatWallBoilingWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual bool activePhasePair(const phasePairKey &) const
Is there phase change mass transfer for this phasePair.
const scalarField & qq() const
Return the quenching surface heat flux [W/m2].
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A FieldMapper for finite-volume patch fields.
virtual void write(Ostream &) const
Write.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition liquid.H:53
scalar Cp(scalar p, scalar T) const
Liquid heat capacity [J/(kg K)].
Definition liquidI.H:39
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phaseModel.H:56
const word & name() const
Definition phaseModel.H:166
virtual const rhoThermo & thermo() const =0
Return the thermophysical model.
An ordered or unorder pair of phase names. Typically specified as follows.
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition phasePair.H:52
Class to represent a system of phases and model interfacial transfers between them.
Definition phaseSystem.H:72
A class for managing temporary objects.
Definition tmp.H:75
static const word propertiesName
Default name of the turbulence properties dictionary.
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
static autoPtr< CHFModel > New(const dictionary &dict)
Select default constructed.
Definition CHFModel.C:39
static autoPtr< CHFSubCoolModel > New(const dictionary &dict)
Select default constructed.
static autoPtr< LeidenfrostModel > New(const dictionary &dict)
Select default constructed.
static autoPtr< MHFModel > New(const dictionary &dict)
Select default constructed.
Definition MHFModel.C:39
static autoPtr< TDNBModel > New(const dictionary &dict)
Select default constructed.
Definition TDNBModel.C:39
static autoPtr< departureDiameterModel > New(const dictionary &dict)
Select default constructed.
static autoPtr< departureFrequencyModel > New(const dictionary &dict)
Select default constructed.
static autoPtr< filmBoilingModel > New(const dictionary &dict)
Select default constructed.
static autoPtr< nucleateFluxModel > New(const dictionary &dict)
Select default constructed.
static autoPtr< nucleationSiteModel > New(const dictionary &dict)
Select default constructed.
static autoPtr< partitioningModel > New(const dictionary &dict)
Select default constructed.
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
UEqn relax()
const tmp< volScalarField > & tCp
Definition EEqn.H:4
const volScalarField & Cp
Definition EEqn.H:7
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
scalar yPlus
scalar magUp
scalar uTau
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
dictionary filmDict
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Mathematical constants.
constexpr scalar pi(M_PI)
constexpr const char *const group
Group name for mathematical constants.
Namespace for handling debugging switches.
Definition debug.C:45
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Type gSum(const FieldField< Field, Type > &f)
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
Type gWeightedSum(const UList< scalar > &weights, const UList< Type > &fld, const label comm)
The global weighted sum (integral) of a field, using the mag() of the weights.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
errorManip< error > abort(error &err)
Definition errorManip.H:139
ThermalDiffusivity< PhaseCompressibleTurbulenceModel< phaseModel > > phaseCompressibleTurbulenceModel
Typedef for phaseCompressibleTurbulenceModel.
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
Type gMin(const FieldField< Field, Type > &f)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
Type gMax(const FieldField< Field, Type > &f)
fvPatchField< vector > fvPatchVectorField
fvPatchField< scalar > fvPatchScalarField
dimensionedScalar pow025(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
volScalarField & e
dimensionedScalar Pr("Pr", dimless, laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Various UniformDimensionedField types.
const vector L(dict.get< vector >("L"))
const Vector< label > N(dict.get< Vector< label > >("N"))