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