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-2020 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 
27 #include "phaseSystem.H"
30 #include "saturationModel.H"
31 #include "rhoReactionThermo.H"
33 
34 using namespace Foam::constant::mathematical;
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 template<>
39 const char* Foam::NamedEnum
40 <
43  2
44 >::names[] =
45 {
46  "vapor",
47  "liquid"
48 };
49 
50 const Foam::NamedEnum
51 <
54  2
55 >
58 
59 
60 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 
62 namespace Foam
63 {
64 namespace compressible
65 {
66 
67 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
68 
69 alphatWallBoilingWallFunctionFvPatchScalarField::
70 alphatWallBoilingWallFunctionFvPatchScalarField
71 (
72  const fvPatch& p,
73  const DimensionedField<scalar, volMesh>& iF
74 )
75 :
76  alphatPhaseChangeWallFunctionFvPatchScalarField(p, iF),
77  phaseType_(liquidPhase),
78  AbyV_(p.size(), 0),
79  alphatConv_(p.size(), 0),
80  dDep_(p.size(), 1e-5),
81  qq_(p.size(), 0),
82  partitioningModel_(nullptr),
83  nucleationSiteModel_(nullptr),
84  departureDiamModel_(nullptr),
85  departureFreqModel_(nullptr)
86 {
87  AbyV_ = this->patch().magSf();
88  forAll(AbyV_, facei)
89  {
90  const label faceCelli = this->patch().faceCells()[facei];
91  AbyV_[facei] /= iF.mesh().V()[faceCelli];
92  }
93 }
94 
95 
98 (
99  const fvPatch& p,
100  const DimensionedField<scalar, volMesh>& iF,
101  const dictionary& dict
102 )
103 :
104  alphatPhaseChangeWallFunctionFvPatchScalarField(p, iF, dict),
105  phaseType_(phaseTypeNames_.read(dict.lookup("phaseType"))),
106  AbyV_(p.size(), 0),
107  alphatConv_(p.size(), 0),
108  dDep_(p.size(), 1e-5),
109  qq_(p.size(), 0),
110  partitioningModel_(nullptr),
111  nucleationSiteModel_(nullptr),
112  departureDiamModel_(nullptr),
113  departureFreqModel_(nullptr)
114 {
115  // Check that otherPhaseName != this phase
116  if (internalField().group() == otherPhaseName_)
117  {
119  << "otherPhase should be the name of the vapor phase that "
120  << "corresponds to the liquid base or vice versa" << nl
121  << "This phase: " << internalField().group() << nl
122  << "otherPhase: " << otherPhaseName_
123  << abort(FatalError);
124  }
125 
126  switch (phaseType_)
127  {
128  case vaporPhase:
129  {
130  partitioningModel_ =
132  (
133  dict.subDict("partitioningModel")
134  );
135 
136  dmdtf_ = 0;
137 
138  break;
139  }
140  case liquidPhase:
141  {
142  partitioningModel_ =
144  (
145  dict.subDict("partitioningModel")
146  );
147 
148  nucleationSiteModel_ =
150  (
151  dict.subDict("nucleationSiteModel")
152  );
153 
154  departureDiamModel_ =
156  (
157  dict.subDict("departureDiamModel")
158  );
159 
160  departureFreqModel_ =
162  (
163  dict.subDict("departureFreqModel")
164  );
165 
166  if (dict.found("dDep"))
167  {
168  dDep_ = scalarField("dDep", dict, p.size());
169  }
170 
171  if (dict.found("qQuenching"))
172  {
173  qq_ = scalarField("qQuenching", dict, p.size());
174  }
175 
176  break;
177  }
178  }
179 
180  if (dict.found("alphatConv"))
181  {
182  alphatConv_ = scalarField("alphatConv", dict, p.size());
183  }
184 
185  AbyV_ = this->patch().magSf();
186  forAll(AbyV_, facei)
187  {
188  const label faceCelli = this->patch().faceCells()[facei];
189  AbyV_[facei] /= iF.mesh().V()[faceCelli];
190  }
191 }
192 
193 
196 (
197  const alphatWallBoilingWallFunctionFvPatchScalarField& psf,
198  const fvPatch& p,
199  const DimensionedField<scalar, volMesh>& iF,
200  const fvPatchFieldMapper& mapper
201 )
202 :
203  alphatPhaseChangeWallFunctionFvPatchScalarField
204  (
205  psf,
206  p,
207  iF,
208  mapper
209  ),
210  phaseType_(psf.phaseType_),
211  AbyV_(mapper(psf.AbyV_)),
212  alphatConv_(mapper(psf.alphatConv_)),
213  dDep_(mapper(psf.dDep_)),
214  qq_(mapper(psf.qq_)),
215  partitioningModel_(psf.partitioningModel_, false),
216  nucleationSiteModel_(psf.nucleationSiteModel_, false),
217  departureDiamModel_(psf.departureDiamModel_, false),
218  departureFreqModel_(psf.departureFreqModel_, false)
219 {}
220 
221 
224 (
225  const alphatWallBoilingWallFunctionFvPatchScalarField& psf
226 )
227 :
228  alphatPhaseChangeWallFunctionFvPatchScalarField(psf),
229  phaseType_(psf.phaseType_),
230  AbyV_(psf.AbyV_),
231  alphatConv_(psf.alphatConv_),
232  dDep_(psf.dDep_),
233  qq_(psf.qq_),
234  partitioningModel_(psf.partitioningModel_, false),
235  nucleationSiteModel_(psf.nucleationSiteModel_, false),
236  departureDiamModel_(psf.departureDiamModel_, false),
237  departureFreqModel_(psf.departureFreqModel_, false)
238 {}
239 
240 
243 (
244  const alphatWallBoilingWallFunctionFvPatchScalarField& psf,
245  const DimensionedField<scalar, volMesh>& iF
246 )
247 :
248  alphatPhaseChangeWallFunctionFvPatchScalarField(psf, iF),
249  phaseType_(psf.phaseType_),
250  AbyV_(psf.AbyV_),
251  alphatConv_(psf.alphatConv_),
252  dDep_(psf.dDep_),
253  qq_(psf.qq_),
254  partitioningModel_(psf.partitioningModel_, false),
255  nucleationSiteModel_(psf.nucleationSiteModel_, false),
256  departureDiamModel_(psf.departureDiamModel_, false),
257  departureFreqModel_(psf.departureFreqModel_, false)
258 {}
259 
260 
261 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
262 
264 (
265  const fvPatchFieldMapper& m
266 )
267 {
269 
270  m(AbyV_, AbyV_);
271  m(alphatConv_, alphatConv_);
272  m(dDep_, dDep_);
273  m(qq_, qq_);
274 }
275 
276 
278 (
279  const fvPatchScalarField& ptf,
280  const labelList& addr
281 )
282 {
284 
285  const alphatWallBoilingWallFunctionFvPatchScalarField& tiptf =
286  refCast<const alphatWallBoilingWallFunctionFvPatchScalarField>(ptf);
287 
288  AbyV_.rmap(tiptf.AbyV_, addr);
289  alphatConv_.rmap(tiptf.alphatConv_, addr);
290  dDep_.rmap(tiptf.dDep_, addr);
291  qq_.rmap(tiptf.qq_, addr);
292 }
293 
294 
296 {
297  if (updated())
298  {
299  return;
300  }
301 
302  // Lookup the fluid model
303  const phaseSystem& fluid =
304  db().lookupObject<phaseSystem>("phaseProperties");
305 
306  const word volatileSpecie(fluid.lookupOrDefault<word>("volatile", "none"));
307 
308  const label patchi = patch().index();
309 
310  switch (phaseType_)
311  {
312  case vaporPhase:
313  {
314  const phaseModel& vapor = fluid.phases()[internalField().group()];
315 
316  // Vapor phase fraction at the wall
317  const scalarField& vaporw = vapor.boundaryField()[patchi];
318 
319  // Partitioning
320  // NOTE! Assumes 1-thisPhase for liquid fraction in
321  // multiphase simulations
322  const scalarField fLiquid(partitioningModel_->fLiquid(1 - vaporw));
323 
324  operator==
325  (
326  calcAlphat(*this)*(1 - fLiquid)/max(vaporw, scalar(1e-8))
327  );
328  break;
329  }
330  case liquidPhase:
331  {
332  const phaseModel& liquid = fluid.phases()[internalField().group()];
333  const phaseModel& vapor = fluid.phases()[otherPhaseName_];
334 
335  const phasePair pair(vapor, liquid);
336 
337  if
338  (
339  db().foundObject<saturationModel>
340  (
341  IOobject::groupName("saturationModel", pair.name())
342  )
343  )
344  {
345  // Retrieve turbulence properties from models
346  const phaseCompressibleMomentumTransportModel& turbModel =
347  db().lookupObject<phaseCompressibleMomentumTransportModel>
348  (
350  (
351  momentumTransportModel::typeName,
352  liquid.name()
353  )
354  );
355  const phaseCompressibleMomentumTransportModel& vaporTurbModel =
356  db().lookupObject<phaseCompressibleMomentumTransportModel>
357  (
359  (
360  momentumTransportModel::typeName,
361  vapor.name()
362  )
363  );
364 
365  const nutWallFunctionFvPatchScalarField& nutw =
366  nutWallFunctionFvPatchScalarField::nutw(turbModel, patchi);
367 
368  const scalar Cmu25(pow025(nutw.Cmu()));
369 
370  const scalarField& y = turbModel.y()[patchi];
371 
372  const tmp<scalarField> tmuw = turbModel.mu(patchi);
373  const scalarField& muw = tmuw();
374 
375  const rhoThermo& lThermo = liquid.thermo();
376  const rhoThermo& vThermo = vapor.thermo();
377 
378  const scalarField& alphaw = lThermo.alpha(patchi);
379 
380  const tmp<volScalarField> tk = turbModel.k();
381  const volScalarField& k = tk();
382  const fvPatchScalarField& kw = k.boundaryField()[patchi];
383 
384  const fvPatchVectorField& Uw =
385  turbModel.U().boundaryField()[patchi];
386  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
387  const scalarField magGradUw(mag(Uw.snGrad()));
388 
389  const fvPatchScalarField& rhoLiquidw =
390  turbModel.rho().boundaryField()[patchi];
391 
392  const fvPatchScalarField& rhoVaporw =
393  vaporTurbModel.rho().boundaryField()[patchi];
394 
395  const fvPatchScalarField& hew =
396  lThermo.he().boundaryField()[patchi];
397 
398  const fvPatchScalarField& pw =
399  lThermo.p().boundaryField()[patchi];
400 
401  const fvPatchScalarField& Tw =
402  lThermo.T().boundaryField()[patchi];
403 
404  const scalarField Tc(Tw.patchInternalField());
405 
406  const scalarField uTau(Cmu25*sqrt(kw));
407 
408  const scalarField yPlus(uTau*y/(muw/rhoLiquidw));
409 
410  const scalarField Pr(muw/alphaw);
411 
412  // Molecular-to-turbulent Prandtl number ratio
413  const scalarField Prat(Pr/Prt_);
414 
415  // Thermal sublayer thickness
416  const scalarField P(this->Psmooth(Prat));
417 
418  const scalarField yPlusTherm(this->yPlusTherm(nutw, P, Prat));
419 
420  const scalarField Cpw(lThermo.Cp(Tw, patchi));
421 
422  // Saturation temperature
423  const saturationModel& satModel =
424  db().lookupObject<saturationModel>
425  (
426  IOobject::groupName("saturationModel", pair.name())
427  );
428  const tmp<volScalarField> tTsat = satModel.Tsat(lThermo.p());
429  const volScalarField& Tsat = tTsat();
430  const fvPatchScalarField& Tsatw(Tsat.boundaryField()[patchi]);
431 
432  // Latent heat
433  scalarField liquidHaw(lThermo.ha(Tc, patchi));
434 
435  scalarField vaporHaw(vThermo.ha(Tsatw, patchi));
436 
437  if (volatileSpecie != "none" && isA<rhoReactionThermo>(lThermo))
438  {
439  const basicSpecieMixture& composition =
440  refCast<const rhoReactionThermo>(lThermo).composition();
441 
442  liquidHaw =
443  composition.Ha
444  (
445  composition.species()[volatileSpecie],
446  pw,
447  Tc
448  );
449  }
450 
451  if (volatileSpecie != "none" && isA<rhoReactionThermo>(vThermo))
452  {
453  const basicSpecieMixture& composition =
454  refCast<const rhoReactionThermo>(vThermo).composition();
455 
456  vaporHaw =
457  composition.Ha
458  (
459  composition.species()[volatileSpecie],
460  pw,
461  Tsatw
462  );
463  }
464 
465  const scalarField L(vaporHaw - liquidHaw);
466 
467  // Liquid phase fraction at the wall
468  const scalarField liquidw(liquid.boundaryField()[patchi]);
469 
470  // Partitioning
471  const scalarField fLiquid(partitioningModel_->fLiquid(liquidw));
472 
473  // Convective thermal diffusivity
474  alphatConv_ = calcAlphat(alphatConv_);
475 
476  label maxIter(10);
477  for (label i=0; i<maxIter; i++)
478  {
479  // Liquid temperature at y+=250 is estimated from
480  // logarithmic thermal wall function (Koncar, Krepper &
481  // Egorov, 2005)
482  const scalarField Tplus_y250
483  (
484  Prt_*(log(nutw.E()*250)/nutw.kappa() + P)
485  );
486 
487  const scalarField Tplus
488  (
489  Prt_*(log(nutw.E()*yPlus)/nutw.kappa() + P)
490  );
491 
492  const scalarField Tl
493  (
494  max
495  (
496  Tc - 40,
497  Tw - (Tplus_y250/Tplus)*(Tw - Tc)
498  )
499  );
500 
501  // Bubble departure diameter:
502  dDep_ = departureDiamModel_->dDeparture
503  (
504  liquid,
505  vapor,
506  patchi,
507  Tl,
508  Tsatw,
509  L
510  );
511 
512  // Bubble departure frequency:
513  const scalarField fDep
514  (
515  departureFreqModel_->fDeparture
516  (
517  liquid,
518  vapor,
519  patchi,
520  Tl,
521  Tsatw,
522  L,
523  dDep_
524  )
525  );
526 
527  // Nucleation site density:
528  const scalarField N
529  (
530  nucleationSiteModel_->N
531  (
532  liquid,
533  vapor,
534  patchi,
535  Tl,
536  Tsatw,
537  L,
538  dDep_,
539  fDep
540  )
541  );
542 
543  // Area fractions:
544 
545  // Del Valle & Kenning (1985)
546  const scalarField Ja
547  (
548  rhoLiquidw*Cpw*(Tsatw - Tl)/(rhoVaporw*L)
549  );
550 
551  const scalarField Al
552  (
553  fLiquid*4.8*exp(min(-Ja/80, log(vGreat)))
554  );
555 
556  scalarField A2(min(pi*sqr(dDep_)*N*Al/4, scalar(1)));
557  const scalarField A1(max(1 - A2, scalar(1e-4)));
558  scalarField A2E(min(pi*sqr(dDep_)*N*Al/4, scalar(5)));
559 
560  if (volatileSpecie != "none" && !liquid.pure())
561  {
562  const volScalarField& Yvolatile =
563  liquid.Y(volatileSpecie);
564  A2E *= Yvolatile.boundaryField()[patchi];
565  A2 *= Yvolatile.boundaryField()[patchi];
566  }
567 
568  // Volumetric mass source in the near wall cell due to the
569  // wall boiling
570  dmdtf_ =
571  (1 - relax_)*dmdtf_
572  + relax_*(1.0/6.0)*A2E*dDep_*rhoVaporw*fDep*AbyV_;
573 
574  // Volumetric source in the near wall cell due to the wall
575  // boiling
576  dmdtLf_ = dmdtf_*L;
577 
578  // Quenching heat transfer coefficient
579  const scalarField hQ
580  (
581  2*(alphaw*Cpw)*fDep
582  *sqrt
583  (
584  (0.8/max(fDep, small))/(pi*alphaw/rhoLiquidw)
585  )
586  );
587 
588  // Quenching heat flux
589  qq_ =
590  (1 - relax_)*qq_
591  + relax_*(A2*hQ*max(Tw - Tl, scalar(0)));
592 
593  // Effective thermal diffusivity that corresponds to the
594  // calculated convective, quenching and evaporative heat
595  // fluxes
596 
597  operator==
598  (
599  (
600  A1*alphatConv_
601  + (qq_ + qe())/max(hew.snGrad(), scalar(1e-16))
602  )
603  /max(liquidw, scalar(1e-8))
604  );
605 
606  scalarField TsupPrev(max((Tw - Tsatw), scalar(0)));
607  const_cast<fvPatchScalarField&>(Tw).evaluate();
608  scalarField TsupNew(max((Tw - Tsatw), scalar(0)));
609 
610  scalar maxErr(max(mag(TsupPrev - TsupNew)));
611 
612  if (debug)
613  {
614  const scalarField qc
615  (
616  fLiquid*A1*(alphatConv_ + alphaw)*hew.snGrad()
617  );
618 
619  const scalarField qEff
620  (
621  liquidw*(*this + alphaw)*hew.snGrad()
622  );
623 
624  Info<< " L: " << gMin(L) << " - " << gMax(L) << endl;
625  Info<< " Tl: " << gMin(Tl) << " - " << gMax(Tl)
626  << endl;
627  Info<< " N: " << gMin(N) << " - " << gMax(N) << endl;
628  Info<< " dDep_: " << gMin(dDep_) << " - "
629  << gMax(dDep_) << endl;
630  Info<< " fDep: " << gMin(fDep) << " - "
631  << gMax(fDep) << endl;
632  Info<< " Al: " << gMin(Al) << " - " << gMax(Al)
633  << endl;
634  Info<< " A1: " << gMin(A1) << " - " << gMax(A1)
635  << endl;
636  Info<< " A2: " << gMin(A2) << " - " << gMax(A2)
637  << endl;
638  Info<< " A2E: " << gMin(A2E) << " - "
639  << gMax(A2E) << endl;
640  Info<< " dmdtW: " << gMin(dmdtf_) << " - "
641  << gMax(dmdtf_) << endl;
642  Info<< " qc: " << gMin(qc) << " - " << gMax(qc)
643  << endl;
644  Info<< " qq: " << gMin(fLiquid*qq()) << " - "
645  << gMax(fLiquid*qq()) << endl;
646  Info<< " qe: " << gMin(fLiquid*qe()) << " - "
647  << gMax(fLiquid*qe()) << endl;
648  Info<< " qEff: " << gMin(qEff) << " - "
649  << gMax(qEff) << endl;
650  Info<< " alphat: " << gMin(*this) << " - "
651  << gMax(*this) << endl;
652  Info<< " alphatConv: " << gMin(alphatConv_)
653  << " - " << gMax(alphatConv_) << endl;
654  }
655 
656  if (maxErr < 1e-1)
657  {
658  if (i > 0)
659  {
660  Info<< "Wall boiling wall function iterations: "
661  << i + 1 << endl;
662  }
663  break;
664  }
665  else if (i == (maxIter - 1))
666  {
667  Info<< "Maximum number of wall boiling wall function "
668  << "iterations (" << maxIter << ") reached." << endl
669  << "Maximum change in wall temperature on last "
670  << "iteration: " << maxErr << endl;
671  }
672 
673  }
674  break;
675  }
676  else
677  {
678  Info<< "Saturation model for phase pair " << pair.name()
679  << " not found. Wall boiling disabled." << endl;
680 
681  operator== (alphatConv_);
682  }
683  break;
684  }
685  default:
686  {
688  << "Unknown phase type. Valid types are: "
689  << phaseTypeNames_ << nl << exit(FatalError);
690  }
691  }
692 
693  fixedValueFvPatchScalarField::updateCoeffs();
694 }
695 
696 
698 {
700 
701  writeEntry(os, "phaseType", phaseTypeNames_[phaseType_]);
702  writeEntry(os, "alphatConv", alphatConv_);
703  writeEntry(os, "dDep", dDep_);
704  writeEntry(os, "qQuenching", qq_);
705 
706  switch (phaseType_)
707  {
708  case vaporPhase:
709  {
710  writeKeyword(os, "partitioningModel") << nl;
711  os << indent << token::BEGIN_BLOCK << incrIndent << nl;
712  partitioningModel_->write(os);
713  os << decrIndent << indent << token::END_BLOCK << nl;
714  break;
715  }
716  case liquidPhase:
717  {
718  writeKeyword(os, "partitioningModel") << nl;
719  os << indent << token::BEGIN_BLOCK << incrIndent << nl;
720  partitioningModel_->write(os);
721  os << decrIndent << indent << token::END_BLOCK << nl;
722 
723  writeKeyword(os, "nucleationSiteModel") << nl;
724  os << indent << token::BEGIN_BLOCK << incrIndent << nl;
725  nucleationSiteModel_->write(os);
726  os << decrIndent << indent << token::END_BLOCK << nl;
727 
728  writeKeyword(os, "departureDiamModel") << nl;
729  os << indent << token::BEGIN_BLOCK << incrIndent << nl;
730  departureDiamModel_->write(os);
731  os << decrIndent << indent << token::END_BLOCK << nl;
732 
733  writeKeyword(os, "departureFreqModel") << nl;
734  os << indent << token::BEGIN_BLOCK << incrIndent << nl;
735  departureFreqModel_->write(os);
736  os << decrIndent << indent << token::END_BLOCK << nl;
737 
738  break;
739  }
740  }
741 }
742 
743 
744 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
745 
747 (
749  alphatWallBoilingWallFunctionFvPatchScalarField
750 );
751 
752 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
753 
754 } // End namespace compressible
755 } // End namespace Foam
756 
757 // ************************************************************************* //
const char *const group
Group name for atomic constants.
static autoPtr< departureDiameterModel > New(const dictionary &dict)
Select null constructed.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
static autoPtr< partitioningModel > New(const dictionary &dict)
Select null constructed.
PhaseCompressibleMomentumTransportModel< phaseModel > phaseCompressibleMomentumTransportModel
Typedef for phaseCompressibleMomentumTransportModel.
fvPatchField< vector > fvPatchVectorField
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
Ostream & writeKeyword(Foam::Ostream &os, const keyType &kw)
Write the keyword to the Ostream with the current level of indentation.
Definition: keyType.C:155
dimensionedScalar log(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
tmp< scalarField > calcAlphat(const scalarField &prevAlphat) const
Evaluate the turbulent thermal diffusivity.
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Type gMin(const FieldField< Field, Type > &f)
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
basicSpecieMixture & composition
const Boundary & boundaryField() const
Return const-reference to the boundary field.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
dimensionedScalar pow025(const dimensionedScalar &ds)
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
static autoPtr< nucleationSiteModel > New(const dictionary &dict)
Select null constructed.
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:51
Macros for easy insertion into run-time selection tables.
scalar magUp
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
scalar uTau
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
static const nutWallFunctionFvPatchScalarField & nutw(const momentumTransportModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
stressControl lookup("compactNormalStress") >> compactNormalStress
dimensionedScalar exp(const dimensionedScalar &ds)
fvPatchField< scalar > fvPatchScalarField
mathematical constants.
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
List< label > labelList
A List of labels.
Definition: labelList.H:56
tmp< scalarField > yPlusTherm(const nutWallFunctionFvPatchScalarField &nutw, const scalarField &P, const scalarField &Prat) const
Calculate y+ at the edge of the thermal laminar sublayer.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
alphatWallBoilingWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
static const char nl
Definition: Ostream.H:260
Type gMax(const FieldField< Field, Type > &f)
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
tmp< scalarField > Psmooth(const scalarField &Prat) const
&#39;P&#39; function
label patchi
static autoPtr< departureFrequencyModel > New(const dictionary &dict)
Select null constructed.
scalar yPlus
messageStream Info
label N
Definition: createFields.H:22
dimensioned< scalar > mag(const dimensioned< Type > &)
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:228
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
makePatchTypeField(fvPatchScalarField, thermalBaffleFvPatchScalarField)
Namespace for OpenFOAM.