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-2022 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,
94  false // do not register
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 tmp<scalarField> tnutw = turbModel.nut(patchi);
193  const scalarField& nutw = tnutw();
194 
195  const scalarField& y = turbModel.y()[patchi];
196 
197  const tmp<scalarField> tnuw = turbModel.nu(patchi);
198  const scalarField& nuw = tnuw();
199 
200  const tmp<volScalarField> tk = turbModel.k();
201  const volScalarField& k = tk();
202 
203  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
204 
205  const scalarField magGradUw(mag(Uw.snGrad()));
206 
207  const scalar Cmu25 = pow025(wallCoeffs_.Cmu());
208  const scalar Cmu75 = pow(wallCoeffs_.Cmu(), 0.75);
209  const scalar kappa = wallCoeffs_.kappa();
210  const scalar yPlusLam = wallCoeffs_.yPlusLam();
211 
212  // Set epsilon and G
213  forAll(nutw, facei)
214  {
215  const label celli = patch.faceCells()[facei];
216 
217  const scalar yPlus = Cmu25*y[facei]*sqrt(k[celli])/nuw[facei];
218 
219  const scalar w = cornerWeights[facei];
220 
221  // Contribution from the viscous sublayer
222  const scalar epsilonVis = w*2.0*k[celli]*nuw[facei]/sqr(y[facei]);
223 
224  // Contribution from the inertial sublayer
225  const scalar epsilonLog = w*Cmu75*pow(k[celli], 1.5)/(kappa*y[facei]);
226 
227  switch (blender_)
228  {
229  case blenderType::STEPWISE:
230  {
231  if (lowReCorrection_ && yPlus < yPlusLam)
232  {
233  epsilon0[celli] += epsilonVis;
234  }
235  else
236  {
237  epsilon0[celli] += epsilonLog;
238  }
239  break;
240  }
241 
242  case blenderType::BINOMIAL:
243  {
244  // (ME:Eqs. 15-16)
245  epsilon0[celli] +=
246  pow
247  (
248  pow(epsilonVis, n_) + pow(epsilonLog, n_),
249  scalar(1)/n_
250  );
251  break;
252  }
253 
254  case blenderType::MAX:
255  {
256  // (PH:Eq. 27)
257  epsilon0[celli] += max(epsilonVis, epsilonLog);
258  break;
259  }
260 
261  case blenderType::EXPONENTIAL:
262  {
263  // (PH:p. 193)
264  const scalar Gamma = 0.001*pow4(yPlus)/(scalar(1) + yPlus);
265  const scalar invGamma = scalar(1)/(Gamma + ROOTVSMALL);
266  epsilon0[celli] +=
267  epsilonVis*exp(-Gamma) + epsilonLog*exp(-invGamma);
268  break;
269  }
270 
271  case blenderType::TANH:
272  {
273  // (KAS:Eqs. 33-34)
274  const scalar phiTanh = tanh(pow4(0.1*yPlus));
275  const scalar b1 = epsilonVis + epsilonLog;
276  const scalar b2 =
277  pow(pow(epsilonVis, 1.2) + pow(epsilonLog, 1.2), 1.0/1.2);
278 
279  epsilon0[celli] += phiTanh*b1 + (1 - phiTanh)*b2;
280  break;
281  }
282  }
283 
284  if (!lowReCorrection_ || (yPlus > yPlusLam))
285  {
286  G0[celli] +=
287  w
288  *(nutw[facei] + nuw[facei])
289  *magGradUw[facei]
290  *Cmu25*sqrt(k[celli])
291  /(kappa*y[facei]);
292  }
293  }
294 }
295 
296 
298 (
299  Ostream& os
300 ) const
301 {
303  os.writeEntryIfDifferent<bool>("lowReCorrection", false, lowReCorrection_);
304  wallCoeffs_.writeEntries(os);
305 }
306 
307 
308 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
309 
312 (
313  const fvPatch& p,
315 )
316 :
317  fixedValueFvPatchField<scalar>(p, iF),
319  lowReCorrection_(false),
320  initialised_(false),
321  master_(-1),
322  wallCoeffs_(),
323  G_(),
324  epsilon_(),
325  cornerWeights_()
326 {}
327 
328 
331 (
333  const fvPatch& p,
335  const fvPatchFieldMapper& mapper
336 )
337 :
338  fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
340  lowReCorrection_(ptf.lowReCorrection_),
341  initialised_(false),
342  master_(-1),
343  wallCoeffs_(ptf.wallCoeffs_),
344  G_(),
345  epsilon_(),
346  cornerWeights_()
347 {}
348 
349 
352 (
353  const fvPatch& p,
355  const dictionary& dict
356 )
357 :
358  fixedValueFvPatchField<scalar>(p, iF, dict),
359  wallFunctionBlenders(dict, blenderType::STEPWISE, scalar(2)),
360  lowReCorrection_(dict.getOrDefault("lowReCorrection", false)),
361  initialised_(false),
362  master_(-1),
363  wallCoeffs_(dict),
364  G_(),
365  epsilon_(),
366  cornerWeights_()
367 {
368  // Apply zero-gradient condition on start-up
370 }
371 
372 
375 (
377 )
378 :
379  fixedValueFvPatchField<scalar>(ewfpsf),
380  wallFunctionBlenders(ewfpsf),
381  lowReCorrection_(ewfpsf.lowReCorrection_),
382  initialised_(false),
383  master_(-1),
384  wallCoeffs_(ewfpsf.wallCoeffs_),
385  G_(),
386  epsilon_(),
387  cornerWeights_()
388 {}
389 
390 
393 (
396 )
397 :
398  fixedValueFvPatchField<scalar>(ewfpsf, iF),
399  wallFunctionBlenders(ewfpsf),
400  lowReCorrection_(ewfpsf.lowReCorrection_),
401  initialised_(false),
402  master_(-1),
403  wallCoeffs_(ewfpsf.wallCoeffs_),
404  G_(),
405  epsilon_(),
406  cornerWeights_()
407 {}
408 
409 
410 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
411 
413 (
414  bool init
415 )
416 {
417  if (patch().index() == master_)
418  {
419  if (init)
420  {
421  G_ = 0.0;
422  }
423 
424  return G_;
425  }
426 
427  return epsilonPatch(master_).G();
428 }
429 
430 
432 (
433  bool init
434 )
435 {
436  if (patch().index() == master_)
437  {
438  if (init)
439  {
440  epsilon_ = 0.0;
441  }
442 
443  return epsilon_;
444  }
445 
446  return epsilonPatch(master_).epsilon(init);
447 }
448 
449 
451 {
452  if (updated())
453  {
454  return;
455  }
456 
457  const auto& turbModel = db().lookupObject<turbulenceModel>
458  (
460  (
462  internalField().group()
463  )
464  );
465 
466  setMaster();
467 
468  if (patch().index() == master_)
469  {
470  createAveragingWeights();
471  calculateTurbulenceFields(turbModel, G(true), epsilon(true));
472  }
473 
474  const scalarField& G0 = this->G();
475  const scalarField& epsilon0 = this->epsilon();
476 
477  typedef DimensionedField<scalar, volMesh> FieldType;
478 
479  FieldType& G = db().lookupObjectRef<FieldType>(turbModel.GName());
480 
481  FieldType& epsilon = const_cast<FieldType&>(internalField());
482 
483  forAll(*this, facei)
484  {
485  const label celli = patch().faceCells()[facei];
486 
487  G[celli] = G0[celli];
488  epsilon[celli] = epsilon0[celli];
489  }
490 
492 }
493 
494 
496 (
497  const scalarField& weights
498 )
499 {
500  if (updated())
501  {
502  return;
503  }
504 
505  const auto& turbModel = db().lookupObject<turbulenceModel>
506  (
508  (
510  internalField().group()
511  )
512  );
513 
514  setMaster();
515 
516  if (patch().index() == master_)
517  {
518  createAveragingWeights();
519  calculateTurbulenceFields(turbModel, G(true), epsilon(true));
520  }
521 
522  const scalarField& G0 = this->G();
523  const scalarField& epsilon0 = this->epsilon();
524 
525  typedef DimensionedField<scalar, volMesh> FieldType;
526 
527  FieldType& G = db().lookupObjectRef<FieldType>(turbModel.GName());
528 
529  FieldType& epsilon = const_cast<FieldType&>(internalField());
530 
531  scalarField& epsilonf = *this;
532 
533  // Only set the values if the weights are > tolerance
534  forAll(weights, facei)
535  {
536  const scalar w = weights[facei];
537 
538  if (w > tolerance_)
539  {
540  const label celli = patch().faceCells()[facei];
541 
542  G[celli] = (1.0 - w)*G[celli] + w*G0[celli];
543  epsilon[celli] = (1.0 - w)*epsilon[celli] + w*epsilon0[celli];
544  epsilonf[facei] = epsilon[celli];
545  }
546  }
547 
549 }
550 
551 
553 (
554  fvMatrix<scalar>& matrix
555 )
556 {
557  if (manipulatedMatrix())
558  {
559  return;
560  }
561 
562  matrix.setValues(patch().faceCells(), patchInternalField());
563 
565 }
566 
567 
569 (
570  fvMatrix<scalar>& matrix,
571  const Field<scalar>& weights
572 )
573 {
574  if (manipulatedMatrix())
575  {
576  return;
577  }
578 
579  DynamicList<label> constraintCells(weights.size());
580  DynamicList<scalar> constraintValues(weights.size());
581  const labelUList& faceCells = patch().faceCells();
582 
583  const DimensionedField<scalar, volMesh>& fld = internalField();
584 
585  forAll(weights, facei)
586  {
587  // Only set the values if the weights are > tolerance
588  if (weights[facei] > tolerance_)
589  {
590  const label celli = faceCells[facei];
591 
592  constraintCells.append(celli);
593  constraintValues.append(fld[celli]);
594  }
595  }
596 
597  if (debug)
598  {
599  Pout<< "Patch: " << patch().name()
600  << ": number of constrained cells = " << constraintCells.size()
601  << " out of " << patch().size()
602  << endl;
603  }
604 
605  matrix.setValues(constraintCells, constraintValues);
606 
608 }
609 
610 
612 (
613  Ostream& os
614 ) const
615 {
617  writeLocalEntries(os);
618  writeEntry("value", os);
619 }
620 
621 
622 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
623 
624 namespace Foam
625 {
627  (
629  epsilonWallFunctionFvPatchScalarField
630  );
631 }
632 
633 
634 // ************************************************************************* //
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:118
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:120
const fvPatch & patch() const noexcept
Return the patch.
Definition: fvPatchField.H:199
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:487
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:68
dimensionedScalar pow025(const dimensionedScalar &ds)
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
Manipulate matrix.
static scalar tolerance_
Tolerance used in weighted calculations.
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:361
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:333
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:80
virtual void operator==(const fvPatchField< scalar > &)
Definition: fvPatchField.C:507
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
constexpr const char *const group
Group name for atomic constants.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:84
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:289
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:303
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:327
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:55
static word timeName(const scalar t, const int precision=precision_)
Return time name of given scalar time formatted with the given precision.
Definition: Time.C:760
virtual tmp< Field< scalar > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:182
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:732
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:270
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 dimensioned internal field reference.
Definition: fvPatchField.H:576
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:977
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.
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:157