omegaWallFunctionFvPatchScalarField.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-2016, 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& omega =
49  static_cast<const volScalarField&>(this->internalField());
50 
51  const volScalarField::Boundary& bf = omega.boundaryField();
52 
53  label master = -1;
54  forAll(bf, patchi)
55  {
56  if (isA<omegaWallFunctionFvPatchScalarField>(bf[patchi]))
57  {
59 
60  if (master == -1)
61  {
62  master = patchi;
63  }
64 
65  opf.master() = master;
66  }
67  }
68 }
69 
70 
72 {
73  const auto& omega =
74  static_cast<const volScalarField&>(this->internalField());
75 
76  const volScalarField::Boundary& bf = omega.boundaryField();
77 
78  const fvMesh& mesh = omega.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> omegaPatches(bf.size());
101  forAll(bf, patchi)
102  {
103  if (isA<omegaWallFunctionFvPatchScalarField>(bf[patchi]))
104  {
105  omegaPatches.append(patchi);
106 
107  const labelUList& faceCells = bf[patchi].patch().faceCells();
108  for (const auto& celli : faceCells)
109  {
110  ++weights[celli];
111  }
112  }
113  }
114 
115  cornerWeights_.setSize(bf.size());
116  for (const auto& patchi : omegaPatches)
117  {
118  const fvPatchScalarField& wf = weights.boundaryField()[patchi];
119  cornerWeights_[patchi] = 1.0/wf.patchInternalField();
120  }
121 
122  G_.setSize(internalField().size(), 0.0);
123  omega_.setSize(internalField().size(), 0.0);
125  initialised_ = true;
126 }
127 
128 
131 (
132  const label patchi
133 )
134 {
135  const auto& omega =
136  static_cast<const volScalarField&>(this->internalField());
137 
138  const volScalarField::Boundary& bf = omega.boundaryField();
139 
140  const auto& opf =
141  refCast<const omegaWallFunctionFvPatchScalarField>(bf[patchi]);
142 
143  return const_cast<omegaWallFunctionFvPatchScalarField&>(opf);
144 }
145 
146 
148 (
149  const turbulenceModel& turbModel,
150  scalarField& G0,
151  scalarField& omega0
152 )
153 {
154  // accumulate all of the G and omega contributions
155  forAll(cornerWeights_, patchi)
156  {
157  if (!cornerWeights_[patchi].empty())
158  {
159  omegaWallFunctionFvPatchScalarField& opf = omegaPatch(patchi);
160 
161  const List<scalar>& w = cornerWeights_[patchi];
162 
163  opf.calculate(turbModel, w, opf.patch(), G0, omega0);
164  }
165  }
166 
167  // apply zero-gradient condition for omega
168  forAll(cornerWeights_, patchi)
169  {
170  if (!cornerWeights_[patchi].empty())
171  {
172  omegaWallFunctionFvPatchScalarField& opf = omegaPatch(patchi);
173 
174  opf == scalarField(omega0, opf.patch().faceCells());
175  }
176  }
177 }
178 
179 
181 (
182  const turbulenceModel& turbModel,
183  const List<scalar>& cornerWeights,
184  const fvPatch& patch,
185  scalarField& G0,
186  scalarField& omega0
187 )
188 {
189  const label patchi = patch.index();
190 
191  const scalar Cmu25 = pow025(wallCoeffs_.Cmu());
192  const scalar kappa = wallCoeffs_.kappa();
193  const scalar yPlusLam = wallCoeffs_.yPlusLam();
194 
195  const scalarField& y = turbModel.y()[patchi];
196 
197  const labelUList& faceCells = patch.faceCells();
198 
199  const tmp<scalarField> tnutw = turbModel.nut(patchi);
200  const scalarField& nutw = tnutw();
201 
202  const tmp<scalarField> tnuw = turbModel.nu(patchi);
203  const scalarField& nuw = tnuw();
204 
205  const tmp<volScalarField> tk = turbModel.k();
206  const volScalarField& k = tk();
207 
208  // Calculate y-plus
209  const auto yPlus = [&](const label facei) -> scalar
210  {
211  return
212  (
213  Cmu25*y[facei]*sqrt(k[faceCells[facei]])/nuw[facei]
214  );
215  };
216 
217  // Contribution from the viscous sublayer
218  const auto omegaVis = [&](const label facei) -> scalar
219  {
220  return
221  (
222  cornerWeights[facei]*6.0*nuw[facei]/(beta1_*sqr(y[facei]))
223  );
224  };
225 
226  // Contribution from the inertial sublayer
227  const auto omegaLog = [&](const label facei) -> scalar
228  {
229  return
230  (
231  cornerWeights[facei]*sqrt(k[faceCells[facei]])
232  / (Cmu25*kappa*y[facei])
233  );
234  };
235 
236  switch (blender_)
237  {
238  case blenderType::STEPWISE:
239  {
240  forAll(faceCells, facei)
241  {
242  if (yPlus(facei) > yPlusLam)
243  {
244  omega0[faceCells[facei]] += omegaLog(facei);
245  }
246  else
247  {
248  omega0[faceCells[facei]] += omegaVis(facei);
249  }
250  }
251  break;
252  }
253 
254  case blenderType::BINOMIAL:
255  {
256  forAll(faceCells, facei)
257  {
258  omega0[faceCells[facei]] +=
259  pow
260  (
261  pow(omegaVis(facei), n_) + pow(omegaLog(facei), n_),
262  scalar(1)/n_
263  );
264  }
265  break;
266  }
267 
268  case blenderType::MAX:
269  {
270  forAll(faceCells, facei)
271  {
272  // (PH:Eq. 27)
273  omega0[faceCells[facei]] +=
274  max(omegaVis(facei), omegaLog(facei));
275  }
276  break;
277  }
278 
279  case blenderType::EXPONENTIAL:
280  {
281  forAll(faceCells, facei)
282  {
283  // (PH:Eq. 31)
284  const scalar yPlusFace = yPlus(facei);
285  const scalar Gamma = 0.01*pow4(yPlusFace)/(1 + 5*yPlusFace);
286  const scalar invGamma = scalar(1)/(Gamma + ROOTVSMALL);
287 
288  omega0[faceCells[facei]] +=
289  (
290  omegaVis(facei)*exp(-Gamma)
291  + omegaLog(facei)*exp(-invGamma)
292  );
293  }
294  break;
295  }
296 
297  case blenderType::TANH:
298  {
299  forAll(faceCells, facei)
300  {
301  // (KAS:Eqs. 33-34)
302  const scalar omegaVisFace = omegaVis(facei);
303  const scalar omegaLogFace = omegaLog(facei);
304  const scalar b1 = omegaVisFace + omegaLogFace;
305  const scalar b2 =
306  pow
307  (
308  pow(omegaVisFace, 1.2) + pow(omegaLogFace, 1.2),
309  1.0/1.2
310  );
311  const scalar phiTanh = tanh(pow4(0.1*yPlus(facei)));
312 
313  omega0[faceCells[facei]] += phiTanh*b1 + (1 - phiTanh)*b2;
314  }
315  break;
316  }
317  }
318 
319  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
320  const scalarField magGradUw(mag(Uw.snGrad()));
321 
322  forAll(faceCells, facei)
323  {
324  if (!(blender_ == blenderType::STEPWISE) || yPlus(facei) > yPlusLam)
325  {
326  G0[faceCells[facei]] +=
327  cornerWeights[facei]
328  *(nutw[facei] + nuw[facei])
329  *magGradUw[facei]
330  *Cmu25*sqrt(k[faceCells[facei]])
331  /(kappa*y[facei]);
332  }
333  }
334 }
335 
336 
338 (
339  Ostream& os
340 ) const
341 {
343  os.writeEntryIfDifferent<scalar>("beta1", 0.075, beta1_);
344  wallCoeffs_.writeEntries(os);
345 }
346 
347 
348 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
349 
351 (
352  const fvPatch& p,
354 )
355 :
356  fixedValueFvPatchField<scalar>(p, iF),
358  initialised_(false),
359  master_(-1),
360  beta1_(0.075),
361  wallCoeffs_(),
362  G_(),
363  omega_(),
364  cornerWeights_()
365 {}
366 
367 
369 (
371  const fvPatch& p,
373  const fvPatchFieldMapper& mapper
374 )
375 :
376  fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
378  initialised_(false),
379  master_(-1),
380  beta1_(ptf.beta1_),
381  wallCoeffs_(ptf.wallCoeffs_),
382  G_(),
383  omega_(),
384  cornerWeights_()
385 {}
386 
387 
389 (
390  const fvPatch& p,
392  const dictionary& dict
393 )
394 :
395  fixedValueFvPatchField<scalar>(p, iF, dict),
396  wallFunctionBlenders(dict, blenderType::BINOMIAL, scalar(2)),
397  initialised_(false),
398  master_(-1),
399  beta1_(dict.getOrDefault<scalar>("beta1", 0.075)),
400  wallCoeffs_(dict),
401  G_(),
402  omega_(),
403  cornerWeights_()
404 {
405  // Apply zero-gradient condition on start-up
406  this->extrapolateInternal();
407 }
408 
409 
411 (
413 )
414 :
415  fixedValueFvPatchField<scalar>(owfpsf),
416  wallFunctionBlenders(owfpsf),
417  initialised_(false),
418  master_(-1),
419  beta1_(owfpsf.beta1_),
420  wallCoeffs_(owfpsf.wallCoeffs_),
421  G_(),
422  omega_(),
423  cornerWeights_()
424 {}
425 
426 
428 (
431 )
432 :
433  fixedValueFvPatchField<scalar>(owfpsf, iF),
434  wallFunctionBlenders(owfpsf),
435  initialised_(false),
436  master_(-1),
437  beta1_(owfpsf.beta1_),
438  wallCoeffs_(owfpsf.wallCoeffs_),
439  G_(),
440  omega_(),
441  cornerWeights_()
442 {}
443 
444 
445 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
446 
448 (
449  bool init
450 )
451 {
452  if (patch().index() == master_)
453  {
454  if (init)
455  {
456  G_ = 0.0;
457  }
458 
459  return G_;
460  }
461 
462  return omegaPatch(master_).G();
463 }
464 
465 
467 (
468  bool init
469 )
470 {
471  if (patch().index() == master_)
472  {
473  if (init)
474  {
475  omega_ = 0.0;
476  }
477 
478  return omega_;
479  }
480 
481  return omegaPatch(master_).omega(init);
482 }
483 
484 
486 {
487  if (updated())
488  {
489  return;
490  }
491 
492  const auto& turbModel = db().lookupObject<turbulenceModel>
493  (
495  (
497  internalField().group()
498  )
499  );
500 
501  setMaster();
502 
503  if (patch().index() == master_)
504  {
505  createAveragingWeights();
506  calculateTurbulenceFields(turbModel, G(true), omega(true));
507  }
508 
509  const scalarField& G0 = this->G();
510  const scalarField& omega0 = this->omega();
511 
512  typedef DimensionedField<scalar, volMesh> FieldType;
513 
514  FieldType& G = db().lookupObjectRef<FieldType>(turbModel.GName());
515 
516  FieldType& omega = const_cast<FieldType&>(internalField());
517 
518  forAll(*this, facei)
519  {
520  const label celli = patch().faceCells()[facei];
521 
522  G[celli] = G0[celli];
523  omega[celli] = omega0[celli];
524  }
525 
527 }
528 
529 
531 (
532  const scalarField& weights
533 )
534 {
535  if (updated())
536  {
537  return;
538  }
539 
540  const auto& turbModel = db().lookupObject<turbulenceModel>
541  (
543  (
545  internalField().group()
546  )
547  );
548 
549  setMaster();
550 
551  if (patch().index() == master_)
552  {
553  createAveragingWeights();
554  calculateTurbulenceFields(turbModel, G(true), omega(true));
555  }
556 
557  const scalarField& G0 = this->G();
558  const scalarField& omega0 = this->omega();
559 
560  typedef DimensionedField<scalar, volMesh> FieldType;
561 
562  FieldType& G = db().lookupObjectRef<FieldType>(turbModel.GName());
563 
564  FieldType& omega = const_cast<FieldType&>(internalField());
565 
566  scalarField& omegaf = *this;
567 
568  // only set the values if the weights are > tolerance
569  forAll(weights, facei)
570  {
571  const scalar w = weights[facei];
572 
573  if (w > tolerance_)
574  {
575  const label celli = patch().faceCells()[facei];
576 
577  G[celli] = (1.0 - w)*G[celli] + w*G0[celli];
578  omega[celli] = (1.0 - w)*omega[celli] + w*omega0[celli];
579  omegaf[facei] = omega[celli];
580  }
581  }
582 
584 }
585 
586 
588 (
589  fvMatrix<scalar>& matrix
590 )
591 {
592  if (manipulatedMatrix())
593  {
594  return;
595  }
596 
597  matrix.setValues(patch().faceCells(), patchInternalField());
598 
600 }
601 
602 
604 (
605  fvMatrix<scalar>& matrix,
606  const Field<scalar>& weights
607 )
608 {
609  if (manipulatedMatrix())
610  {
611  return;
612  }
613 
614  DynamicList<label> constraintCells(weights.size());
615  DynamicList<scalar> constraintValues(weights.size());
616  const labelUList& faceCells = patch().faceCells();
617 
618  const DimensionedField<scalar, volMesh>& fld = internalField();
619 
620  forAll(weights, facei)
621  {
622  // only set the values if the weights are > tolerance
623  if (tolerance_ < weights[facei])
624  {
625  const label celli = faceCells[facei];
626 
627  constraintCells.append(celli);
628  constraintValues.append(fld[celli]);
629  }
630  }
631 
632  if (debug)
633  {
634  Pout<< "Patch: " << patch().name()
635  << ": number of constrained cells = " << constraintCells.size()
636  << " out of " << patch().size()
637  << endl;
638  }
639 
640  matrix.setValues(constraintCells, constraintValues);
641 
643 }
644 
645 
647 (
648  Ostream& os
649 ) const
650 {
652  writeLocalEntries(os);
654 }
655 
656 
657 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
658 
659 namespace Foam
660 {
662  (
664  omegaWallFunctionFvPatchScalarField
665  );
666 }
667 
668 
669 // ************************************************************************* //
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
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
static bool initialised_(false)
fvPatchField< vector > fvPatchVectorField
const dimensionedScalar G0
Conductance quantum: default SI units: [S].
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
omegaWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
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)
virtual void updateWeightedCoeffs(const scalarField &weights)
Update the coefficients associated with the patch field.
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 createAveragingWeights()
Create the averaging weights for cells which are bounded by multiple wall function faces...
scalarField & omega(bool init=false)
Return non-const access to the master&#39;s omega field.
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 kappa
Coulomb constant: default SI units: [N.m2/C2].
Abstract base class for turbulence models (RAS, LES and laminar).
virtual omegaWallFunctionFvPatchScalarField & omegaPatch(const label patchi)
Helper function to return non-const access to an omega patch.
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)
fvPatchField< scalar > fvPatchScalarField
static const word propertiesName
Default name of the turbulence properties dictionary.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
static scalar tolerance_
Tolerance used in weighted calculations.
virtual void calculateTurbulenceFields(const turbulenceModel &turbulence, scalarField &G0, scalarField &omega0)
Main driver to calculate the turbulence 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
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.
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)
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)
void writeLocalEntries(Ostream &) const
Write local wall function variables.
This boundary condition provides a wall function for the specific dissipation rate (i...
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.
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.
virtual void setMaster()
Set the master patch - master is responsible for updating all wall function patches.
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
volScalarField & p
virtual label & master()
Return non-const access to the master patch ID.
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
Manipulate matrix.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Do not request registration (bool: false)
Namespace for OpenFOAM.
virtual void calculate(const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &omega)
Calculate the omega and G.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127