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