alphatWallBoilingWallFunctionFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2015-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
29 #include "phaseSystem.H"
35 
39 #include "mixedFvPatchFields.H"
40 
42 
43 using namespace Foam::constant::mathematical;
44 
45 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
46 
47 template<>
48 const char* Foam::NamedEnum
49 <
52  3
53 >::names[] = {"vapour", "vapor", "liquid"};
54 
55 const Foam::NamedEnum
56 <
59  3
60 >
63 
64 
65 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66 
67 namespace Foam
68 {
69 namespace compressible
70 {
71 
72 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
73 
75 {
76  // Data
77 
78  //- Wall function field
80 
81  //- Phase
82  const phaseModel& phase;
83 
84  //- Other phase
86 
87  //- Volume fraction
89 
90  //- Other volume fraction
92 
93  //- Interface
95 
96  //- Phase thermophysical transport model
98 
99  //- Phase convective turbulent thermal diffusivity
101 
102 
103  //- Constructor
105  (
107  const phaseModel& phase,
108  const phaseModel& otherPhase
109  )
110  :
111  field(field),
112  phase(phase),
113  otherPhase(otherPhase),
114  alphaw(phase.boundaryField()[patchi()]),
115  otherAlphaw(otherPhase.boundaryField()[patchi()]),
116  interface(phase, otherPhase),
117  ttm
118  (
119  field.db().lookupType<fluidThermophysicalTransportModel>
120  (
121  phase.name()
122  )
123  ),
124  alphatConv
125  (
127  (
128  ttm,
129  field.Prt_,
130  patchi()
131  )
132  )
133  {}
134 
135 
136  // Member Functions
137 
138  //- Patch
139  inline const fvPatch& patch() const
140  {
141  return field.patch();
142  }
143 
144  //- Patch index
145  inline label patchi() const
146  {
147  return patch().index();
148  }
149 };
150 
151 
153 :
155 {
156  // Data
157 
158  //- Name of the volatile specie
160 
161  //- Patch area by neighbouring cell volume ratio
163 
164  //- Liquid density
167 
168  //- Vapour density
171 
172  //- Liquid heat capacity
173  const scalarField& Cpw;
174 
175  //- Liquid laminar kinematic viscosity
177  const scalarField& nuw;
178 
179  //- Liquid laminar thermal diffusivity
181 
182  //- Liquid viscosity wall function
184 
185  //- Dimensionless wall distance
187 
188  //- Smoothing function
189  const scalarField P;
190 
191  //- Cell temperature
193 
194  //- Saturation temperature
196 
197  //- Latent heat
198  const scalarField L;
199 
200 
201  //- Constructor
203  (
205  const phaseModel& liquid,
206  const phaseModel& vapour
207  )
208  :
209  properties(field, liquid, vapour),
210  volatileSpecie
211  (
212  liquid.fluid().lookupOrDefault<word>("volatile", "none")
213  ),
214  AbyV
215  (
216  patch().magSf()
217  /scalarField
218  (
219  patch().boundaryMesh().mesh().V(),
220  patch().faceCells()
221  )
222  ),
223  trhoLiquidw(liquid.thermo().rho(patchi())),
224  rhoLiquidw(trhoLiquidw()),
225  trhoVapourw(vapour.thermo().rho(patchi())),
226  rhoVapourw(trhoVapourw()),
227  Cpw(liquid.thermo().Cp().boundaryField()[patchi()]),
228  tnuw(liquid.fluidThermo().nu(patchi())),
229  nuw(tnuw()),
230  kappaByCp(liquid.thermo().kappa().boundaryField()[patchi()]/Cpw),
231  nutw
232  (
234  (
235  ttm.momentumTransport(),
236  patchi()
237  )
238  ),
239  yPlus
240  (
241  pow025(nutw.Cmu())
242  *sqrt(ttm.momentumTransport().k()().boundaryField()[patchi()])
243  *ttm.momentumTransport().yb()[patchi()]
244  /nuw
245  ),
246  P
247  (
249  (
250  rhoLiquidw*nuw/kappaByCp/field.Prt_
251  )
252  ),
253  Tc
254  (
255  liquid.thermo().T().boundaryField()[patchi()].patchInternalField()
256  ),
257  Tsat
258  (
259  liquid.fluid().lookupInterfacialModel
260  <
262  >
263  (interface)
264  .Tsat(liquid.fluidThermo().p())()
265  .boundaryField()[patchi()]
266  ),
267  L
268  (
269  volatileSpecie != "none"
270  ? -refCast<const heatTransferPhaseSystem>(liquid.fluid())
271  .Li
272  (
273  interface,
274  volatileSpecie,
275  scalarField(patch().size(), +1),
276  Tsat,
277  patch().faceCells(),
278  heatTransferPhaseSystem::latentHeatScheme::upwind
279  )
280  : -refCast<const heatTransferPhaseSystem>(liquid.fluid())
281  .L
282  (
283  interface,
284  scalarField(patch().size(), +1),
285  Tsat,
286  patch().faceCells(),
287  heatTransferPhaseSystem::latentHeatScheme::upwind
288  )
289  )
290  {}
291 };
292 
293 
294 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
295 
297 alphatWallBoilingWallFunctionFvPatchScalarField::calcBoiling
298 (
299  const boilingLiquidProperties& props,
300  const scalarField& Tw,
301  scalarField& dDeparture,
302  scalarField& fDeparture,
303  scalarField& nucleationSiteDensity,
304  scalarField& qQuenching,
305  scalarField& qEvaporative,
306  scalarField& dmdtf
307 ) const
308 {
309  scalarField Tl;
310  if (!useLiquidTemperatureWallFunction_)
311  {
312  Tl = props.Tc;
313  }
314  else
315  {
316  // Liquid temperature at y+=250 is estimated from the
317  // logarithmic thermal wall function of Koncar, Krepper
318  // & Egorov (2005)
319  const scalarField TyPlus250
320  (
321  Prt_
322  *(
323  log(props.nutw.E()*250)
324  /props.nutw.kappa()
325  + props.P
326  )
327  );
328 
329  const scalarField TyPlus
330  (
331  Prt_
332  *(
333  log(props.nutw.E()*max(props.yPlus, scalar(11)))
334  /props.nutw.kappa()
335  + props.P
336  )
337  );
338 
339  Tl = Tw - (TyPlus250/TyPlus)*(Tw - props.Tc);
340  }
341 
342  // Bubble departure diameter
343  dDeparture =
344  departureDiameterModel_->dDeparture
345  (
346  props.phase,
347  props.otherPhase,
348  patch().index(),
349  Tl,
350  props.Tsat,
351  props.L
352  );
353 
354  // Bubble departure frequency
355  fDeparture =
356  departureFrequencyModel_->fDeparture
357  (
358  props.phase,
359  props.otherPhase,
360  patch().index(),
361  Tl,
362  props.Tsat,
363  props.L,
364  dDeparture
365  );
366 
367  // Nucleation site density
368  nucleationSiteDensity =
369  nucleationSiteModel_->nucleationSiteDensity
370  (
371  props.phase,
372  props.otherPhase,
373  patch().index(),
374  Tl,
375  props.Tsat,
376  props.L,
377  dDeparture,
378  fDeparture
379  );
380 
381  // Del Valle & Kenning (1985)
382  const scalarField Ja
383  (
384  props.rhoLiquidw
385  *props.Cpw
386  *max(props.Tsat - Tl, scalar(0))
387  /(props.rhoVapourw*props.L)
388  );
389 
390  const scalarField Al
391  (
392  wetFraction_*4.8*exp(min(-Ja/80, log(vGreat)))
393  );
394 
395  scalarField A2
396  (
397  min(pi*sqr(dDeparture)*nucleationSiteDensity*Al/4, scalar(1))
398  );
399  const scalarField A1(max(1 - A2, scalar(1e-4)));
400  scalarField A2E
401  (
402  min(pi*sqr(dDeparture)*nucleationSiteDensity*Al/4, scalar(5))
403  );
404 
405  if (props.volatileSpecie != "none" && !props.phase.pure())
406  {
407  const scalarField& Yvolatile =
408  props.phase
409  .Y(props.volatileSpecie)
410  .boundaryField()[patch().index()];
411  A2E *= Yvolatile;
412  A2 *= Yvolatile;
413  }
414 
415  // Volumetric mass source in the near wall cell due to the
416  // wall boiling
417  dmdtf = (1.0/6.0)*A2E*dDeparture*props.rhoVapourw*fDeparture*props.AbyV;
418 
419  // Quenching heat transfer coefficient
420  const scalarField hQ
421  (
422  2*props.kappaByCp*props.Cpw*fDeparture
423  *sqrt
424  (
425  tau_
426  /max(fDeparture, small)
427  /(pi*props.kappaByCp/props.rhoLiquidw)
428  )
429  );
430 
431  // Quenching heat flux
432  qQuenching = A2*hQ*max(Tw - Tl, scalar(0));
433 
434  // Evaporation heat flux
435  qEvaporative = dmdtf*props.L/props.AbyV;
436 
437  // Return total sum of convective, quenching and evaporative heat fluxes
438  const scalarField gradTw
439  (
440  patch().deltaCoeffs()*max(Tw - props.Tc, small*props.Tc)
441  );
442  return A1*props.alphatConv*props.Cpw*gradTw + qQuenching_ + qEvaporative_;
443 }
444 
445 
446 tmp<scalarField>
447 alphatWallBoilingWallFunctionFvPatchScalarField::calcBoiling
448 (
449  const boilingLiquidProperties& props,
450  const scalarField& Tw
451 ) const
452 {
453  scalarField dDeparture(dDeparture_);
454  scalarField fDeparture(fDeparture_);
455  scalarField nucleationSiteDensity(nucleationSiteDensity_);
456  scalarField qQuenching(qQuenching_);
457  scalarField qEvaporative(qEvaporative_);
458  scalarField dmdtf(dmdtf_);
459 
460  return
461  calcBoiling
462  (
463  props,
464  Tw,
465  dDeparture,
466  fDeparture,
467  nucleationSiteDensity,
468  qQuenching,
469  qEvaporative,
470  dmdtf
471  );
472 }
473 
474 
475 tmp<scalarField>
476 alphatWallBoilingWallFunctionFvPatchScalarField::evaluateBoiling
477 (
478  const boilingLiquidProperties& props,
479  const scalarField& Tw
480 )
481 {
482  return
483  calcBoiling
484  (
485  props,
486  Tw,
487  dDeparture_,
488  fDeparture_,
489  nucleationSiteDensity_,
490  qQuenching_,
491  qEvaporative_,
492  dmdtf_
493  );
494 }
495 
496 
497 const fvPatchScalarField&
498 alphatWallBoilingWallFunctionFvPatchScalarField::getTemperaturePatchField
499 (
500  const boilingLiquidProperties& props,
501  scalarField& isFixed,
502  scalarField& h,
503  scalarField& hTaPlusQa
504 ) const
505 {
506  isFixed.setSize(patch().size());
507  h.setSize(patch().size());
508  hTaPlusQa.setSize(patch().size());
509 
510  const fvPatchScalarField& Tw =
511  props.phase.thermo().T().boundaryField()[patch().index()];
512 
513  if (isA<fixedValueFvPatchScalarField>(Tw))
514  {
515  isFixed = 1;
516  h = rootVGreat;
517  hTaPlusQa = rootVGreat*Tw;
518  }
519  else if (isA<zeroGradientFvPatchScalarField>(Tw))
520  {
521  isFixed = 0;
522  h = 0;
523  hTaPlusQa = 0;
524  }
525  else if (isA<fixedGradientFvPatchScalarField>(Tw))
526  {
527  const fixedGradientFvPatchScalarField& Twm =
528  refCast<const fixedGradientFvPatchScalarField>(Tw);
529 
530  isFixed = 0;
531  h = 0;
532  hTaPlusQa = (*this)*props.Cpw*Twm.gradient();
533  }
534  else if (isA<mixedFvPatchScalarField>(Tw))
535  {
536  const mixedFvPatchScalarField& Twm =
537  refCast<const mixedFvPatchScalarField>(Tw);
538 
539  isFixed = pos(Twm.valueFraction() - 1 + rootSmall);
540  h =
541  Twm.valueFraction()
542  /max(1 - Twm.valueFraction(), rootVSmall)
543  *(*this)*props.Cpw*patch().deltaCoeffs();
544  hTaPlusQa =
545  h*Twm.refValue()
546  + (*this)*props.Cpw*Twm.refGrad();
547  }
548  else
549  {
551  << "Temperature boundary condition type not recognised"
552  << exit(FatalError);
553  }
554 
555  return Tw;
556 }
557 
558 
559 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
560 
561 alphatWallBoilingWallFunctionFvPatchScalarField::
562 alphatWallBoilingWallFunctionFvPatchScalarField
563 (
564  const fvPatch& p,
566  const dictionary& dict
567 )
568 :
569  fixedValueFvPatchScalarField(p, iF, dict),
571 
572  phaseType_(phaseTypeNames_.read(dict.lookup("phaseType"))),
573  useLiquidTemperatureWallFunction_
574  (
575  dict.lookupOrDefault<Switch>("useLiquidTemperatureWallFunction", true)
576  ),
577  tolerance_(dict.lookupOrDefault<scalar>("tolerance", rootSmall)),
578 
579  Prt_(dict.lookupOrDefault<scalar>("Prt", dimless, 0.85)),
580  tau_(dict.lookupOrDefault<scalar>("bubbleWaitingTimeRatio", dimless, 0.8)),
581 
582  partitioningModel_(nullptr),
583  nucleationSiteModel_(nullptr),
584  departureDiameterModel_(nullptr),
585  departureFrequencyModel_(nullptr),
586 
587  wetFraction_(p.size(), 0),
588  dDeparture_(p.size(), 1e-5),
589  fDeparture_(p.size(), 0),
590  nucleationSiteDensity_(p.size(), 0),
591  qQuenching_(p.size(), 0),
592  qEvaporative_(p.size(), 0),
593  dmdtf_(p.size(), 0)
594 {
595  // Sub-Models
596  partitioningModel_ =
598  (
599  dict.subDict("partitioningModel")
600  );
601 
602  if (phaseType_ == liquidPhase)
603  {
604  nucleationSiteModel_ =
606  (
607  dict.subDict("nucleationSiteModel")
608  );
609  departureDiameterModel_ =
611  (
612  dict.subDictBackwardsCompatible
613  (
614  {"departureDiameterModel", "departureDiamModel"}
615  )
616  );
617  departureFrequencyModel_ =
619  (
620  dict.subDictBackwardsCompatible
621  (
622  {"departureFrequencyModel", "departureFreqModel"}
623  )
624  );
625  }
626 
627  // Backwards compatible reading for old field keywords
628  auto readFieldBackwardsCompatible = [&p]
629  (
630  const dictionary& dict,
631  const unitConversion& units,
632  const wordList& keywords,
633  scalarField& field
634  )
635  {
636  forAll(keywords, i)
637  {
638  if (dict.found(keywords[i]))
639  {
640  field = scalarField(keywords[i], units, dict, p.size());
641  return;
642  }
643  }
644  };
645 
646  // State
647  readFieldBackwardsCompatible
648  (
649  dict,
650  unitFraction,
651  {"wetFraction", "wallLiquidFraction"},
652  wetFraction_
653  );
654 
655  if (phaseType_ == liquidPhase)
656  {
657  readFieldBackwardsCompatible
658  (
659  dict,
660  dimLength,
661  {"dDeparture"},
662  dDeparture_
663  );
664 
665  readFieldBackwardsCompatible
666  (
667  dict,
669  {"fDeparture", "depFrequency"},
670  fDeparture_
671  );
672 
673  readFieldBackwardsCompatible
674  (
675  dict,
677  {"nucleationSiteDensity", "nucSiteDensity"},
678  nucleationSiteDensity_
679  );
680 
681  readFieldBackwardsCompatible
682  (
683  dict,
685  {"qQuenching"},
686  qQuenching_
687  );
688 
689  readFieldBackwardsCompatible
690  (
691  dict,
693  {"qEvaporative"},
694  qEvaporative_
695  );
696 
697  readFieldBackwardsCompatible
698  (
699  dict,
701  {"dmdtf"},
702  dmdtf_
703  );
704  }
705 }
706 
707 
710 (
712  const fvPatch& p,
714  const fieldMapper& mapper
715 )
716 :
717  fixedValueFvPatchScalarField(psf, p, iF, mapper),
719 
720  phaseType_(psf.phaseType_),
721  useLiquidTemperatureWallFunction_(psf.useLiquidTemperatureWallFunction_),
722  tolerance_(psf.tolerance_),
723 
724  Prt_(psf.Prt_),
725  tau_(psf.tau_),
726 
727  partitioningModel_(psf.partitioningModel_, false),
728  nucleationSiteModel_(psf.nucleationSiteModel_, false),
729  departureDiameterModel_(psf.departureDiameterModel_, false),
730  departureFrequencyModel_(psf.departureFrequencyModel_, false),
731 
732  wetFraction_(mapper(psf.wetFraction_)),
733  dDeparture_(mapper(psf.dDeparture_)),
734  fDeparture_(mapper(psf.fDeparture_)),
735  nucleationSiteDensity_(mapper(psf.nucleationSiteDensity_)),
736  qQuenching_(mapper(psf.qQuenching_)),
737  qEvaporative_(mapper(psf.qEvaporative_)),
738  dmdtf_(mapper(psf.dmdtf_))
739 {}
740 
741 
744 (
747 )
748 :
749  fixedValueFvPatchScalarField(psf, iF),
751 
752  phaseType_(psf.phaseType_),
753  useLiquidTemperatureWallFunction_(psf.useLiquidTemperatureWallFunction_),
754  tolerance_(psf.tolerance_),
755 
756  Prt_(psf.Prt_),
757  tau_(psf.tau_),
758 
759  partitioningModel_(psf.partitioningModel_, false),
760  nucleationSiteModel_(psf.nucleationSiteModel_, false),
761  departureDiameterModel_(psf.departureDiameterModel_, false),
762  departureFrequencyModel_(psf.departureFrequencyModel_, false),
763 
764  wetFraction_(psf.wetFraction_),
765  dDeparture_(psf.dDeparture_),
766  fDeparture_(psf.fDeparture_),
767  nucleationSiteDensity_(psf.nucleationSiteDensity_),
768  qQuenching_(psf.qQuenching_),
769  qEvaporative_(psf.qEvaporative_),
770  dmdtf_(psf.dmdtf_)
771 {}
772 
773 
774 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
775 
777 (
778  const fvPatchScalarField& ptf,
779  const fieldMapper& mapper
780 )
781 {
782  fixedValueFvPatchScalarField::map(ptf, mapper);
783 
785  refCast<const alphatWallBoilingWallFunctionFvPatchScalarField>(ptf);
786 
787  mapper(wetFraction_, tiptf.wetFraction_);
788  mapper(dDeparture_, tiptf.dDeparture_);
789  mapper(fDeparture_, tiptf.fDeparture_);
790  mapper(nucleationSiteDensity_, tiptf.nucleationSiteDensity_);
791  mapper(qQuenching_, tiptf.qQuenching_);
792  mapper(qEvaporative_, tiptf.qEvaporative_);
793  mapper(dmdtf_, tiptf.dmdtf_);
794 }
795 
796 
798 (
799  const fvPatchScalarField& ptf
800 )
801 {
802  fixedValueFvPatchScalarField::reset(ptf);
803 
805  refCast<const alphatWallBoilingWallFunctionFvPatchScalarField>(ptf);
806 
807  wetFraction_.reset(tiptf.wetFraction_);
808  dDeparture_.reset(tiptf.dDeparture_);
809  fDeparture_.reset(tiptf.fDeparture_);
810  nucleationSiteDensity_.reset(tiptf.nucleationSiteDensity_);
811  qQuenching_.reset(tiptf.qQuenching_);
812  qEvaporative_.reset(tiptf.qEvaporative_);
813  dmdtf_.reset(tiptf.dmdtf_);
814 }
815 
816 
818 {
819  if (updated())
820  {
821  return;
822  }
823 
824  // Lookup the fluid model and the phases
825  const phaseSystem& fluid =
826  db().lookupObject<phaseSystem>(phaseSystem::propertiesName);
827 
828  switch (phaseType_)
829  {
830  case vapourPhase:
831  case vaporPhase:
832  {
833  const phaseModel& vapour = fluid.phases()[internalField().group()];
834  const phaseModel& liquid = fluid.phases()[otherPhaseName_];
835 
836  // Construct boiling properties
837  const properties props(*this, vapour, liquid);
838 
839  // Partitioning. Note: Assumes that there is only only one liquid
840  // phase and all other phases are vapour.
841  wetFraction_ = partitioningModel_->wetFraction(props.otherAlphaw);
842 
843  operator==
844  (
845  (1 - wetFraction_)
846  /max(1 - props.otherAlphaw, rootSmall)
847  *props.alphatConv
848  );
849 
850  break;
851  }
852 
853  case liquidPhase:
854  {
855  const phaseModel& liquid = fluid.phases()[internalField().group()];
856  const phaseModel& vapour = fluid.phases()[otherPhaseName_];
857 
858  // Boiling is enabled by the presence of saturation temperature
859  // modelling. This is consistent with interfacial thermal phase
860  // changes.
861  if
862  (
863  !fluid.foundInterfacialModel
864  <
866  >(phaseInterface(liquid, vapour))
867  )
868  {
869  // Construct non-boiling properties
870  const properties props(*this, liquid, vapour);
871 
872  Info<< "Saturation model for interface "
873  << props.interface.name()
874  << " not found. Wall boiling disabled." << endl;
875 
876  operator==(props.alphatConv);
877  }
878  else
879  {
880  // Construct boiling properties
881  const boilingLiquidProperties props(*this, liquid, vapour);
882 
883  // Partitioning. Note: Assumes that there is only only one
884  // liquid phase and all other phases are vapour.
885  wetFraction_ = partitioningModel_->wetFraction(props.alphaw);
886 
887  // Get the temperature boundary condition and extract its
888  // physical parameters
889  scalarField TwpfIsFixed, TwpfH, TwpfHTaPlusQa;
890  const fvPatchScalarField& Twpf =
891  getTemperaturePatchField
892  (
893  props,
894  TwpfIsFixed,
895  TwpfH,
896  TwpfHTaPlusQa
897  );
898 
899  // Define the residual. This should be monotonic in Tw.
900  auto R = [&](const scalarField& Tw)
901  {
902  return calcBoiling(props, Tw) - TwpfHTaPlusQa + TwpfH*Tw;
903  };
904 
905  // Solve using interval bisection. Boiling cannot occur below
906  // the saturation temperature, so that is taken to be the lower
907  // bound. The upper bound is harder to define. For now, take
908  // twice the current superheat as the maximum. The solution is
909  // likely to be below this value. If it is not, then the
910  // iteration will converge to this upper limit, creating a new
911  // superheat of twice the current value. Subsequent time-steps
912  // will then double this superheat again and again until it
913  // does eventually bound the solution.
914  const scalarField isBoiling(neg(R(props.Tsat)));
915  scalarField Tw0(props.Tsat);
916  scalarField Tw1
917  (
918  max
919  (
920  Twpf + (Twpf - props.Tsat),
921  props.Tsat*(1 + sqrt(tolerance_))
922  )
923  );
924  scalar e =
925  gMax((1 - TwpfIsFixed)*isBoiling*(Tw1 - Tw0)/(Tw0 + Tw1));
926  for (; e > tolerance_; e /= 2)
927  {
928  const scalarField TwM((Tw0 + Tw1)/2);
929  const scalarField rM(R(TwM));
930  Tw0 = pos(rM)*Tw0 + neg0(rM)*TwM;
931  Tw1 = pos(rM)*TwM + neg0(rM)*Tw1;
932  }
933  const scalarField Tw
934  (
935  TwpfIsFixed*Twpf + (1 - TwpfIsFixed)*(Tw0 + Tw1)/2
936  );
937 
938  // Use solution to re-evaluate the boiling and set the thermal
939  // diffusivity to recover the calculated heat flux
940  const scalarField gradTw
941  (
942  patch().deltaCoeffs()*max(Tw - props.Tc, small*props.Tc)
943  );
944 
945  const scalarField q(evaluateBoiling(props, Tw));
946 
947  operator==
948  (
949  isBoiling*q/props.Cpw/gradTw/max(props.alphaw, rootSmall)
950  + (1 - isBoiling)*props.alphatConv
951  );
952  }
953 
954  break;
955  }
956  }
957 
958  fixedValueFvPatchScalarField::updateCoeffs();
959 }
960 
961 
963 {
966 
967  // Controls
968  writeEntry(os, "phaseType", phaseTypeNames_[phaseType_]);
969  writeEntry
970  (
971  os,
972  "useLiquidTemperatureWallFunction",
973  useLiquidTemperatureWallFunction_
974  );
975  writeEntry(os, "tolerance", tolerance_);
976 
977  // Parameters
978  writeEntry(os, "Prt", Prt_);
979  writeEntry(os, "bubbleWaitingTimeRatio", tau_);
980 
981  // Sub-models
982  writeKeyword(os, "partitioningModel") << nl;
983  os << indent << token::BEGIN_BLOCK << incrIndent << nl;
984  partitioningModel_->write(os);
985  os << decrIndent << indent << token::END_BLOCK << nl;
986 
987  if (phaseType_ == liquidPhase)
988  {
989  writeKeyword(os, "nucleationSiteModel") << nl;
990  os << indent << token::BEGIN_BLOCK << incrIndent << nl;
991  nucleationSiteModel_->write(os);
992  os << decrIndent << indent << token::END_BLOCK << nl;
993 
994  writeKeyword(os, "departureDiameterModel") << nl;
995  os << indent << token::BEGIN_BLOCK << incrIndent << nl;
996  departureDiameterModel_->write(os);
997  os << decrIndent << indent << token::END_BLOCK << nl;
998 
999  writeKeyword(os, "departureFrequencyModel") << nl;
1000  os << indent << token::BEGIN_BLOCK << incrIndent << nl;
1001  departureFrequencyModel_->write(os);
1002  os << decrIndent << indent << token::END_BLOCK << nl;
1003  }
1004 
1005  // State
1006  writeEntry(os, "wetFraction", wetFraction_);
1007 
1008  if (phaseType_ == liquidPhase)
1009  {
1010  writeEntry(os, "dDeparture", dDeparture_);
1011  writeEntry(os, "fDeparture", fDeparture_);
1012  writeEntry(os, "nucleationSiteDensity", nucleationSiteDensity_);
1013  writeEntry(os, "qQuenching", qQuenching_);
1014  writeEntry(os, "qEvaporative", qEvaporative_);
1015  writeEntry(os, "dmdtf", dmdtf_);
1016  }
1017 }
1018 
1019 
1020 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1021 
1023 (
1026 );
1027 
1028 
1029 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1030 
1031 } // End namespace compressible
1032 } // End namespace Foam
1033 
1034 // ************************************************************************* //
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
label k
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
tmp< Field< Type > > T() const
Return the field transpose (only defined for second rank tensors)
Definition: Field.C:539
void reset(const Field< Type > &)
Reset the field values to the given field.
Definition: Field.C:456
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:54
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
This boundary condition provides a thermal wall function for turbulent thermal diffusivity (usuallyal...
Abstract base-class for all alphatWallFunctions supporting phase-change.
A thermal wall function for simulation of subcooled nucleate wall boiling with runtime selectable sub...
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void map(const fvPatchScalarField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
alphatWallBoilingWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:57
Abstract base class for fluid thermophysical transport models RAS, LES and laminar.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
Wrapper around saturationTemperatureModel to facilitate convenient construction on interfaces.
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition: liquid.H:53
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions,...
Class to represent an interface between phases. Derivations can further specify the configuration of ...
virtual word name() const
Name.
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:73
bool foundInterfacialModel(const phaseInterface &interface) const
Check availability of a sub model for a given interface.
static const word propertiesName
Default name of the phase properties dictionary.
Definition: phaseSystem.H:226
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:102
A class for managing temporary objects.
Definition: tmp.H:55
@ BEGIN_BLOCK
Definition: token.H:113
@ END_BLOCK
Definition: token.H:114
Unit conversion structure. Contains the associated dimensions and the multiplier with which to conver...
Upwind interpolation scheme class.
Definition: upwind.H:54
static autoPtr< departureDiameterModel > New(const dictionary &dict)
Select null constructed.
static autoPtr< departureFrequencyModel > New(const dictionary &dict)
Select null constructed.
static autoPtr< nucleationSiteModel > New(const dictionary &dict)
Select null constructed.
static autoPtr< partitioningModel > New(const dictionary &dict)
Select null constructed.
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const scalar yPlus
label patchi
makePatchTypeField(fvPatchScalarField, alphatJayatillekeWallFunctionFvPatchScalarField)
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
mathematical constants.
const dimensionedScalar h
Planck constant.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensionedScalar pos(const dimensionedScalar &ds)
const doubleScalar e
Definition: doubleScalar.H:106
dimensionedScalar exp(const dimensionedScalar &ds)
bool read(const char *, int32_t &)
Definition: int32IO.C:85
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:241
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:129
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimPower
const dimensionSet dimless
messageStream Info
const dimensionSet dimLength
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:234
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const HashTable< unitConversion > & units()
Get the table of unit conversions.
const dimensionSet dimTime
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionSet dimDensity
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
static scalar R(const scalar a, const scalar x)
Definition: invIncGamma.C:102
dimensionedScalar neg(const dimensionedScalar &ds)
error FatalError
const dimensionSet dimArea
dimensionedScalar neg0(const dimensionedScalar &ds)
Ostream & writeKeyword(Foam::Ostream &os, const keyType &kw)
Write the keyword to the Ostream with the current level of indentation.
Definition: keyType.C:155
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:227
static const char nl
Definition: Ostream.H:266
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Type gMax(const FieldField< Field, Type > &f)
fvPatchField< scalar > fvPatchScalarField
dimensionedScalar pow025(const dimensionedScalar &ds)
const unitConversion unitFraction
dictionary dict
volScalarField & p
fluidMulticomponentThermo & thermo
Definition: createFields.H:31
boilingLiquidProperties(const alphatWallBoilingWallFunctionFvPatchScalarField &field, const phaseModel &liquid, const phaseModel &vapour)
Constructor.
properties(const alphatWallBoilingWallFunctionFvPatchScalarField &field, const phaseModel &phase, const phaseModel &otherPhase)
Constructor.
const fluidThermophysicalTransportModel & ttm
Phase thermophysical transport model.
const alphatWallBoilingWallFunctionFvPatchScalarField & field
Wall function field.