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