epsilonWallFunctionFvPatchScalarField.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2019 OpenFOAM Foundation
9  Copyright (C) 2017-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
31 #include "turbulenceModel.H"
32 #include "fvMatrix.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
38 
39 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
40 
42 {
43  if (master_ != -1)
44  {
45  return;
46  }
47 
48  const auto& epsilon =
49  static_cast<const volScalarField&>(this->internalField());
50 
51  const volScalarField::Boundary& bf = epsilon.boundaryField();
52 
53  label master = -1;
54  forAll(bf, patchi)
55  {
56  if (isA<epsilonWallFunctionFvPatchScalarField>(bf[patchi]))
57  {
59 
60  if (master == -1)
61  {
62  master = patchi;
63  }
64 
65  epf.master() = master;
66  }
67  }
68 }
69 
70 
72 {
73  const auto& epsilon =
74  static_cast<const volScalarField&>(this->internalField());
75 
76  const volScalarField::Boundary& bf = epsilon.boundaryField();
77 
78  const fvMesh& mesh = epsilon.mesh();
79 
80  if (initialised_ && !mesh.changing())
81  {
82  return;
83  }
84 
85  volScalarField weights
86  (
87  IOobject
88  (
89  "weights",
90  mesh.time().timeName(),
91  mesh,
95  ),
96  mesh,
98  );
99 
100  DynamicList<label> epsilonPatches(bf.size());
101  forAll(bf, patchi)
102  {
103  if (isA<epsilonWallFunctionFvPatchScalarField>(bf[patchi]))
104  {
105  epsilonPatches.append(patchi);
106 
107  const labelUList& faceCells = bf[patchi].patch().faceCells();
108  for (const auto& faceCell : faceCells)
109  {
110  ++weights[faceCell];
111  }
112  }
113  }
114 
115  cornerWeights_.setSize(bf.size());
116 
117  for (const auto& patchi : epsilonPatches)
118  {
119  const fvPatchScalarField& wf = weights.boundaryField()[patchi];
120  cornerWeights_[patchi] = 1.0/wf.patchInternalField();
121  }
122 
123  G_.setSize(internalField().size(), Zero);
124  epsilon_.setSize(internalField().size(), Zero);
126  initialised_ = true;
127 }
128 
129 
132 (
133  const label patchi
134 )
135 {
136  const auto& epsilon =
137  static_cast<const volScalarField&>(this->internalField());
138 
139  const volScalarField::Boundary& bf = epsilon.boundaryField();
140 
141  const auto& epf =
142  refCast<const epsilonWallFunctionFvPatchScalarField>(bf[patchi]);
143 
144  return const_cast<epsilonWallFunctionFvPatchScalarField&>(epf);
145 }
146 
147 
149 (
151  scalarField& G0,
153 )
154 {
155  // Accumulate all of the G and epsilon contributions
156  forAll(cornerWeights_, patchi)
157  {
158  if (!cornerWeights_[patchi].empty())
159  {
160  epsilonWallFunctionFvPatchScalarField& epf = epsilonPatch(patchi);
161 
162  const List<scalar>& w = cornerWeights_[patchi];
163 
164  epf.calculate(turbulence, w, epf.patch(), G0, epsilon0);
165  }
166  }
167 
168  // Apply zero-gradient condition for epsilon
169  forAll(cornerWeights_, patchi)
170  {
171  if (!cornerWeights_[patchi].empty())
172  {
173  epsilonWallFunctionFvPatchScalarField& epf = epsilonPatch(patchi);
174 
175  epf == scalarField(epsilon0, epf.patch().faceCells());
176  }
177  }
178 }
179 
180 
182 (
183  const turbulenceModel& turbModel,
184  const List<scalar>& cornerWeights,
185  const fvPatch& patch,
186  scalarField& G0,
188 )
189 {
190  const label patchi = patch.index();
191 
192  const scalar Cmu25 = pow025(wallCoeffs_.Cmu());
193  const scalar Cmu75 = pow(wallCoeffs_.Cmu(), 0.75);
194  const scalar kappa = wallCoeffs_.kappa();
195  const scalar yPlusLam = wallCoeffs_.yPlusLam();
196 
197  const scalarField& y = turbModel.y()[patchi];
198 
199  const labelUList& faceCells = patch.faceCells();
200 
201  const tmp<scalarField> tnuw = turbModel.nu(patchi);
202  const scalarField& nuw = tnuw();
203 
204  const tmp<volScalarField> tk = turbModel.k();
205  const volScalarField& k = tk();
206 
207  // Calculate y-plus
208  const auto yPlus = [&](const label facei) -> scalar
209  {
210  return
211  (
212  Cmu25*y[facei]*sqrt(k[faceCells[facei]])/nuw[facei]
213  );
214  };
215 
216  // Contribution from the viscous sublayer
217  const auto epsilonVis = [&](const label facei) -> scalar
218  {
219  return
220  (
221  cornerWeights[facei]*2.0*k[faceCells[facei]]*nuw[facei]
222  / sqr(y[facei])
223  );
224  };
225 
226  // Contribution from the inertial sublayer
227  const auto epsilonLog = [&](const label facei) -> scalar
228  {
229  return
230  (
231  cornerWeights[facei]*Cmu75*pow(k[faceCells[facei]], 1.5)
232  / (kappa*y[facei])
233  );
234  };
235 
236  switch (blender_)
237  {
238  case blenderType::STEPWISE:
239  {
240  forAll(faceCells, facei)
241  {
242  if (lowReCorrection_ && yPlus(facei) < yPlusLam)
243  {
244  epsilon0[faceCells[facei]] += epsilonVis(facei);
245  }
246  else
247  {
248  epsilon0[faceCells[facei]] += epsilonLog(facei);
249  }
250  }
251  break;
252  }
253 
254  case blenderType::BINOMIAL:
255  {
256  forAll(faceCells, facei)
257  {
258  // (ME:Eqs. 15-16)
259  epsilon0[faceCells[facei]] +=
260  pow
261  (
262  pow(epsilonVis(facei), n_) + pow(epsilonLog(facei), n_),
263  scalar(1)/n_
264  );
265  }
266  break;
267  }
268 
269  case blenderType::MAX:
270  {
271  forAll(faceCells, facei)
272  {
273  // (PH:Eq. 27)
274  epsilon0[faceCells[facei]] +=
275  max(epsilonVis(facei), epsilonLog(facei));
276  }
277  break;
278  }
279 
280  case blenderType::EXPONENTIAL:
281  {
282  forAll(faceCells, facei)
283  {
284  // (PH:p. 193)
285  const scalar yPlusFace = yPlus(facei);
286  const scalar Gamma =
287  0.001*pow4(yPlusFace)/(scalar(1) + yPlusFace);
288  const scalar invGamma = scalar(1)/(Gamma + ROOTVSMALL);
289 
290  epsilon0[faceCells[facei]] +=
291  (
292  epsilonVis(facei)*exp(-Gamma)
293  + epsilonLog(facei)*exp(-invGamma)
294  );
295  }
296  break;
297  }
298 
299  case blenderType::TANH:
300  {
301  forAll(faceCells, facei)
302  {
303  // (KAS:Eqs. 33-34)
304  const scalar epsilonVisFace = epsilonVis(facei);
305  const scalar epsilonLogFace = epsilonLog(facei);
306  const scalar b1 = epsilonVisFace + epsilonLogFace;
307  const scalar b2 =
308  pow
309  (
310  pow(epsilonVisFace, 1.2) + pow(epsilonLogFace, 1.2),
311  1.0/1.2
312  );
313  const scalar phiTanh = tanh(pow4(0.1*yPlus(facei)));
314 
315  epsilon0[faceCells[facei]] += phiTanh*b1 + (1 - phiTanh)*b2;
316  }
317  break;
318  }
319  }
320 
321  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
322  const scalarField magGradUw(mag(Uw.snGrad()));
323 
324  const tmp<scalarField> tnutw = turbModel.nut(patchi);
325  const scalarField& nutw = tnutw();
326 
327  forAll(faceCells, facei)
328  {
329  if (!lowReCorrection_ || (yPlus(facei) > yPlusLam))
330  {
331  G0[faceCells[facei]] +=
332  cornerWeights[facei]
333  *(nutw[facei] + nuw[facei])
334  *magGradUw[facei]
335  *Cmu25*sqrt(k[faceCells[facei]])
336  /(kappa*y[facei]);
337  }
338  }
339 }
340 
341 
343 (
344  Ostream& os
345 ) const
346 {
348  os.writeEntryIfDifferent<bool>("lowReCorrection", false, lowReCorrection_);
349  wallCoeffs_.writeEntries(os);
350 }
351 
352 
353 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
354 
357 (
358  const fvPatch& p,
360 )
361 :
362  fixedValueFvPatchField<scalar>(p, iF),
364  lowReCorrection_(false),
365  initialised_(false),
366  master_(-1),
367  wallCoeffs_(),
368  G_(),
369  epsilon_(),
370  cornerWeights_()
371 {}
372 
373 
376 (
378  const fvPatch& p,
380  const fvPatchFieldMapper& mapper
381 )
382 :
383  fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
385  lowReCorrection_(ptf.lowReCorrection_),
386  initialised_(false),
387  master_(-1),
388  wallCoeffs_(ptf.wallCoeffs_),
389  G_(),
390  epsilon_(),
391  cornerWeights_()
392 {}
393 
394 
397 (
398  const fvPatch& p,
400  const dictionary& dict
401 )
402 :
403  fixedValueFvPatchField<scalar>(p, iF, dict),
404  wallFunctionBlenders(dict, blenderType::STEPWISE, scalar(2)),
405  lowReCorrection_(dict.getOrDefault("lowReCorrection", false)),
406  initialised_(false),
407  master_(-1),
408  wallCoeffs_(dict),
409  G_(),
410  epsilon_(),
411  cornerWeights_()
412 {
413  // Apply zero-gradient condition on start-up
414  this->extrapolateInternal();
415 }
416 
417 
420 (
422 )
423 :
424  fixedValueFvPatchField<scalar>(ewfpsf),
425  wallFunctionBlenders(ewfpsf),
426  lowReCorrection_(ewfpsf.lowReCorrection_),
427  initialised_(false),
428  master_(-1),
429  wallCoeffs_(ewfpsf.wallCoeffs_),
430  G_(),
431  epsilon_(),
432  cornerWeights_()
433 {}
434 
435 
438 (
441 )
442 :
443  fixedValueFvPatchField<scalar>(ewfpsf, iF),
444  wallFunctionBlenders(ewfpsf),
445  lowReCorrection_(ewfpsf.lowReCorrection_),
446  initialised_(false),
447  master_(-1),
448  wallCoeffs_(ewfpsf.wallCoeffs_),
449  G_(),
450  epsilon_(),
451  cornerWeights_()
452 {}
453 
454 
455 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
456 
458 (
459  bool init
460 )
461 {
462  if (patch().index() == master_)
463  {
464  if (init)
465  {
466  G_ = 0.0;
467  }
468 
469  return G_;
470  }
471 
472  return epsilonPatch(master_).G();
473 }
474 
475 
477 (
478  bool init
479 )
480 {
481  if (patch().index() == master_)
482  {
483  if (init)
484  {
485  epsilon_ = 0.0;
486  }
487 
488  return epsilon_;
489  }
490 
491  return epsilonPatch(master_).epsilon(init);
492 }
493 
494 
496 {
497  if (updated())
498  {
499  return;
500  }
501 
502  const auto& turbModel = db().lookupObject<turbulenceModel>
503  (
505  (
507  internalField().group()
508  )
509  );
510 
511  setMaster();
512 
513  if (patch().index() == master_)
514  {
515  createAveragingWeights();
516  calculateTurbulenceFields(turbModel, G(true), epsilon(true));
517  }
518 
519  const scalarField& G0 = this->G();
520  const scalarField& epsilon0 = this->epsilon();
521 
522  typedef DimensionedField<scalar, volMesh> FieldType;
523 
524  FieldType& G = db().lookupObjectRef<FieldType>(turbModel.GName());
525 
526  FieldType& epsilon = const_cast<FieldType&>(internalField());
527 
528  forAll(*this, facei)
529  {
530  const label celli = patch().faceCells()[facei];
531 
532  G[celli] = G0[celli];
533  epsilon[celli] = epsilon0[celli];
534  }
535 
537 }
538 
539 
541 (
542  const scalarField& weights
543 )
544 {
545  if (updated())
546  {
547  return;
548  }
549 
550  const auto& turbModel = db().lookupObject<turbulenceModel>
551  (
553  (
555  internalField().group()
556  )
557  );
558 
559  setMaster();
560 
561  if (patch().index() == master_)
562  {
563  createAveragingWeights();
564  calculateTurbulenceFields(turbModel, G(true), epsilon(true));
565  }
566 
567  const scalarField& G0 = this->G();
568  const scalarField& epsilon0 = this->epsilon();
569 
570  typedef DimensionedField<scalar, volMesh> FieldType;
571 
572  FieldType& G = db().lookupObjectRef<FieldType>(turbModel.GName());
573 
574  FieldType& epsilon = const_cast<FieldType&>(internalField());
575 
576  scalarField& epsilonf = *this;
577 
578  // Only set the values if the weights are > tolerance
579  forAll(weights, facei)
580  {
581  const scalar w = weights[facei];
582 
583  if (w > tolerance_)
584  {
585  const label celli = patch().faceCells()[facei];
586 
587  G[celli] = (1.0 - w)*G[celli] + w*G0[celli];
588  epsilon[celli] = (1.0 - w)*epsilon[celli] + w*epsilon0[celli];
589  epsilonf[facei] = epsilon[celli];
590  }
591  }
592 
594 }
595 
596 
598 (
599  fvMatrix<scalar>& matrix
600 )
601 {
602  if (manipulatedMatrix())
603  {
604  return;
605  }
606 
607  matrix.setValues(patch().faceCells(), patchInternalField());
608 
610 }
611 
612 
614 (
615  fvMatrix<scalar>& matrix,
616  const Field<scalar>& weights
617 )
618 {
619  if (manipulatedMatrix())
620  {
621  return;
622  }
623 
624  DynamicList<label> constraintCells(weights.size());
625  DynamicList<scalar> constraintValues(weights.size());
626  const labelUList& faceCells = patch().faceCells();
627 
628  const DimensionedField<scalar, volMesh>& fld = internalField();
629 
630  forAll(weights, facei)
631  {
632  // Only set the values if the weights are > tolerance
633  if (weights[facei] > tolerance_)
634  {
635  const label celli = faceCells[facei];
636 
637  constraintCells.append(celli);
638  constraintValues.append(fld[celli]);
639  }
640  }
641 
642  if (debug)
643  {
644  Pout<< "Patch: " << patch().name()
645  << ": number of constrained cells = " << constraintCells.size()
646  << " out of " << patch().size()
647  << endl;
648  }
649 
650  matrix.setValues(constraintCells, constraintValues);
651 
653 }
654 
655 
657 (
658  Ostream& os
659 ) const
660 {
662  writeLocalEntries(os);
664 }
665 
666 
667 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
668 
669 namespace Foam
670 {
672  (
674  epsilonWallFunctionFvPatchScalarField
675  );
676 }
677 
678 
679 // ************************************************************************* //
dimensionedScalar tanh(const dimensionedScalar &ds)
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
static bool initialised_(false)
fvPatchField< vector > fvPatchVectorField
void writeLocalEntries(Ostream &) const
Write local wall function variables.
const dimensionedScalar G0
Conductance quantum: default SI units: [S].
virtual void calculate(const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &epsilon)
Calculate the epsilon and G.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const fvPatch & patch() const noexcept
Return the patch.
Definition: fvPatchField.H:269
const dimensionedScalar G
Newtonian constant of gravitation.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
dimensionedScalar pow025(const dimensionedScalar &ds)
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
Manipulate matrix.
static scalar tolerance_
Tolerance used in weighted calculations.
void extrapolateInternal()
Assign the patch field from the internal field.
Definition: fvPatchField.C:62
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of boundary fields.
label k
Boltzmann constant.
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:52
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
const dimensionedScalar epsilon0
Electric constant: default SI units: [F/m].
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Abstract base class for turbulence models (RAS, LES and laminar).
The class wallFunctionBlenders is a base class that hosts common entries for various derived wall-fun...
Macros for easy insertion into run-time selection tables.
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
constexpr const char *const group
Group name for atomic constants.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
scalar y
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
void setSize(const label n)
Alias for resize()
Definition: List.H:316
dynamicFvMesh & mesh
dimensionedScalar exp(const dimensionedScalar &ds)
virtual label & master()
Return non-const access to the master patch ID.
fvPatchField< scalar > fvPatchScalarField
static const word propertiesName
Default name of the turbulence properties dictionary.
virtual void calculateTurbulenceFields(const turbulenceModel &turbulence, scalarField &G0, scalarField &epsilon0)
Main driver to calculate the turbulence fields.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
Definition: fvPatchField.C:342
This boundary condition provides wall functions for the turbulent kinetic energy dissipation rate (i...
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:336
scalarField & G(bool init=false)
Return non-const access to the master&#39;s G field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual epsilonWallFunctionFvPatchScalarField & epsilonPatch(const label patchi)
Helper function to return non-const access to an epsilon patch.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
int debug
Static debugging option.
OBJstream os(runTime.globalPath()/outputName)
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
bool changing() const noexcept
Is mesh changing (topology changing and/or moving)
Definition: polyMesh.H:767
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
virtual void updateWeightedCoeffs(const scalarField &weights)
Update the coefficients associated with the patch field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
scalar yPlus
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:309
void writeEntries(Ostream &) const
Write wall-function blending data as dictionary entries.
scalar epsilon
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
dimensionedScalar pow4(const dimensionedScalar &ds)
Nothing to be read.
const DimensionedField< scalar, volMesh > & internalField() const noexcept
Return const-reference to the dimensioned internal field.
Definition: fvPatchField.H:662
const std::string patch
OpenFOAM patch number as a std::string.
void setValues(const labelUList &cellLabels, const Type &value)
Set solution in given cells to the specified value and eliminate the corresponding equations from the...
Definition: fvMatrix.C:978
virtual void setMaster()
Set the master patch - master is responsible for updating all wall function patches.
volScalarField & p
virtual void createAveragingWeights()
Create the averaging weights for cells which are bounded by multiple wall function faces...
epsilonWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Do not request registration (bool: false)
Namespace for OpenFOAM.
scalarField & epsilon(bool init=false)
Return non-const access to the master&#39;s epsilon field.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127