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