omegaWallFunctionFvPatchScalarField.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) 2011-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 
28 #include "turbulenceModel.H"
29 #include "fvPatchFieldMapper.H"
30 #include "fvMatrix.H"
31 #include "volFields.H"
32 #include "wallFvPatch.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 scalar omegaWallFunctionFvPatchScalarField::tolerance_ = 1e-5;
43 
44 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
45 
47 {
48  if (!isA<wallFvPatch>(patch()))
49  {
51  << "Invalid wall function specification" << nl
52  << " Patch type for patch " << patch().name()
53  << " must be wall" << nl
54  << " Current patch type is " << patch().type() << nl << endl
55  << abort(FatalError);
56  }
57 }
58 
59 
61 {
62  os.writeKeyword("Cmu") << Cmu_ << token::END_STATEMENT << nl;
63  os.writeKeyword("kappa") << kappa_ << token::END_STATEMENT << nl;
64  os.writeKeyword("E") << E_ << token::END_STATEMENT << nl;
65  os.writeKeyword("beta1") << beta1_ << token::END_STATEMENT << nl;
66  os.writeKeyword("blended") << blended_ << token::END_STATEMENT << nl;
67 }
68 
69 
71 {
72  if (master_ != -1)
73  {
74  return;
75  }
76 
77  const volScalarField& omega =
78  static_cast<const volScalarField&>(this->internalField());
79 
80  const volScalarField::Boundary& bf = omega.boundaryField();
81 
82  label master = -1;
83  forAll(bf, patchi)
84  {
85  if (isA<omegaWallFunctionFvPatchScalarField>(bf[patchi]))
86  {
88 
89  if (master == -1)
90  {
91  master = patchi;
92  }
93 
94  opf.master() = master;
95  }
96  }
97 }
98 
99 
101 {
102  const volScalarField& omega =
103  static_cast<const volScalarField&>(this->internalField());
104 
105  const volScalarField::Boundary& bf = omega.boundaryField();
106 
107  const fvMesh& mesh = omega.mesh();
108 
109  if (initialised_ && !mesh.changing())
110  {
111  return;
112  }
113 
114  volScalarField weights
115  (
116  IOobject
117  (
118  "weights",
119  mesh.time().timeName(),
120  mesh,
123  false // do not register
124  ),
125  mesh,
126  dimensionedScalar("zero", dimless, 0.0)
127  );
128 
129  DynamicList<label> omegaPatches(bf.size());
130  forAll(bf, patchi)
131  {
132  if (isA<omegaWallFunctionFvPatchScalarField>(bf[patchi]))
133  {
134  omegaPatches.append(patchi);
135 
136  const labelUList& faceCells = bf[patchi].patch().faceCells();
137  forAll(faceCells, i)
138  {
139  label celli = faceCells[i];
140  weights[celli]++;
141  }
142  }
143  }
144 
145  cornerWeights_.setSize(bf.size());
146  forAll(omegaPatches, i)
147  {
148  label patchi = omegaPatches[i];
149  const fvPatchScalarField& wf = weights.boundaryField()[patchi];
151  }
152 
153  G_.setSize(internalField().size(), 0.0);
154  omega_.setSize(internalField().size(), 0.0);
155 
156  initialised_ = true;
157 }
158 
159 
162 {
163  const volScalarField& omega =
164  static_cast<const volScalarField&>(this->internalField());
165 
166  const volScalarField::Boundary& bf = omega.boundaryField();
167 
169  refCast<const omegaWallFunctionFvPatchScalarField>(bf[patchi]);
170 
171  return const_cast<omegaWallFunctionFvPatchScalarField&>(opf);
172 }
173 
174 
176 (
177  const turbulenceModel& turbModel,
178  scalarField& G0,
179  scalarField& omega0
180 )
181 {
182  // accumulate all of the G and omega contributions
184  {
185  if (!cornerWeights_[patchi].empty())
186  {
188 
189  const List<scalar>& w = cornerWeights_[patchi];
190 
191  opf.calculate(turbModel, w, opf.patch(), G0, omega0);
192  }
193  }
194 
195  // apply zero-gradient condition for omega
197  {
198  if (!cornerWeights_[patchi].empty())
199  {
201 
202  opf == scalarField(omega0, opf.patch().faceCells());
203  }
204  }
205 }
206 
207 
209 (
210  const turbulenceModel& turbModel,
211  const List<scalar>& cornerWeights,
212  const fvPatch& patch,
213  scalarField& G0,
214  scalarField& omega0
215 )
216 {
217  const label patchi = patch.index();
218 
219  const scalarField& y = turbModel.y()[patchi];
220 
221  const scalar Cmu25 = pow025(Cmu_);
222  const scalar Cmu5 = sqrt(Cmu_);
223 
224  const tmp<volScalarField> tk = turbModel.k();
225  const volScalarField& k = tk();
226 
227  const tmp<scalarField> tnuw = turbModel.nu(patchi);
228  const scalarField& nuw = tnuw();
229 
230  const tmp<scalarField> tnutw = turbModel.nut(patchi);
231  const scalarField& nutw = tnutw();
232 
233  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
234 
235  const scalarField magGradUw(mag(Uw.snGrad()));
236 
237  const FieldType& G =
238  db().lookupObject<FieldType>(turbModel.GName());
239 
240  // Set omega and G
241  forAll(nutw, facei)
242  {
243  const label celli = patch.faceCells()[facei];
244  const scalar w = cornerWeights[facei];
245 
246  const scalar Rey = y[facei]*sqrt(k[celli])/nuw[facei];
247  const scalar yPlus = Cmu25*Rey;
248  const scalar uPlus = (1/kappa_)*log(E_*yPlus);
249 
250  if (blended_)
251  {
252  const scalar lamFrac = exp(-Rey/11);
253  const scalar turbFrac = 1 - lamFrac;
254 
255  const scalar uStar = sqrt
256  (
257  lamFrac*nuw[facei]*magGradUw[facei] + turbFrac*Cmu5*k[celli]
258  );
259 
260  const scalar omegaVis = 6*nuw[facei]/(beta1_*sqr(y[facei]));
261  const scalar omegaLog = uStar/(Cmu5*kappa_*y[facei]);
262 
263  omega0[celli] += w*(lamFrac*omegaVis + turbFrac*omegaLog);
264 
265  G0[celli] +=
266  w
267  *(
268  lamFrac*G[celli]
269 
270  + turbFrac
271  *sqr(uStar*magGradUw[facei]*y[facei]/uPlus)
272  /(nuw[facei]*kappa_*yPlus)
273  );
274  }
275  else
276  {
277  if (yPlus < yPlusLam_)
278  {
279  const scalar omegaVis = 6*nuw[facei]/(beta1_*sqr(y[facei]));
280 
281  omega0[celli] += w*omegaVis;
282 
283  G0[celli] += w*G[celli];
284  }
285  else
286  {
287  const scalar uStar = sqrt(Cmu5*k[celli]);
288  const scalar omegaLog = uStar/(Cmu5*kappa_*y[facei]);
289 
290  omega0[celli] += w*omegaLog;
291 
292  G0[celli] +=
293  w*
294  sqr(uStar*magGradUw[facei]*y[facei]/uPlus)
295  /(nuw[facei]*kappa_*yPlus);
296  }
297  }
298  }
299 }
300 
301 
302 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
303 
305 (
306  const fvPatch& p,
308 )
309 :
311  Cmu_(0.09),
312  kappa_(0.41),
313  E_(9.8),
314  beta1_(0.075),
315  blended_(false),
317  G_(),
318  omega_(),
319  initialised_(false),
320  master_(-1),
322 {
323  checkType();
324 }
325 
326 
328 (
329  const fvPatch& p,
331  const dictionary& dict
332 )
333 :
335  Cmu_(dict.lookupOrDefault<scalar>("Cmu", 0.09)),
336  kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
337  E_(dict.lookupOrDefault<scalar>("E", 9.8)),
338  beta1_(dict.lookupOrDefault<scalar>("beta1", 0.075)),
339  blended_(dict.lookupOrDefault<Switch>("blended", false)),
341  G_(),
342  omega_(),
343  initialised_(false),
344  master_(-1),
346 {
347  checkType();
348 
349  // apply zero-gradient condition on start-up
351 }
352 
353 
355 (
357  const fvPatch& p,
359  const fvPatchFieldMapper& mapper
360 )
361 :
362  fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
363  Cmu_(ptf.Cmu_),
364  kappa_(ptf.kappa_),
365  E_(ptf.E_),
366  beta1_(ptf.beta1_),
367  blended_(ptf.blended_),
368  yPlusLam_(ptf.yPlusLam_),
369  G_(),
370  omega_(),
371  initialised_(false),
372  master_(-1),
374 {
375  checkType();
376 }
377 
378 
380 (
382 )
383 :
385  Cmu_(owfpsf.Cmu_),
386  kappa_(owfpsf.kappa_),
387  E_(owfpsf.E_),
388  beta1_(owfpsf.beta1_),
389  blended_(owfpsf.blended_),
390  yPlusLam_(owfpsf.yPlusLam_),
391  G_(),
392  omega_(),
393  initialised_(false),
394  master_(-1),
396 {
397  checkType();
398 }
399 
400 
402 (
405 )
406 :
407  fixedValueFvPatchField<scalar>(owfpsf, iF),
408  Cmu_(owfpsf.Cmu_),
409  kappa_(owfpsf.kappa_),
410  E_(owfpsf.E_),
411  beta1_(owfpsf.beta1_),
412  blended_(owfpsf.blended_),
413  yPlusLam_(owfpsf.yPlusLam_),
414  G_(),
415  omega_(),
416  initialised_(false),
417  master_(-1),
419 {
420  checkType();
421 }
422 
423 
424 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
425 
427 {
428  if (patch().index() == master_)
429  {
430  if (init)
431  {
432  G_ = 0.0;
433  }
434 
435  return G_;
436  }
437 
438  return omegaPatch(master_).G();
439 }
440 
441 
443 {
444  if (patch().index() == master_)
445  {
446  if (init)
447  {
448  omega_ = 0.0;
449  }
450 
451  return omega_;
452  }
453 
454  return omegaPatch(master_).omega(init);
455 }
456 
457 
459 {
460  if (updated())
461  {
462  return;
463  }
464 
465  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
466  (
468  (
470  internalField().group()
471  )
472  );
473 
474  setMaster();
475 
476  if (patch().index() == master_)
477  {
479  calculateTurbulenceFields(turbModel, G(true), omega(true));
480  }
481 
482  const scalarField& G0 = this->G();
483  const scalarField& omega0 = this->omega();
484 
485  FieldType& G =
486  const_cast<FieldType&>
487  (
488  db().lookupObject<FieldType>(turbModel.GName())
489  );
490 
491  FieldType& omega = const_cast<FieldType&>(internalField());
492 
493  forAll(*this, facei)
494  {
495  label celli = patch().faceCells()[facei];
496 
497  G[celli] = G0[celli];
498  omega[celli] = omega0[celli];
499  }
500 
502 }
503 
504 
506 (
507  const scalarField& weights
508 )
509 {
510  if (updated())
511  {
512  return;
513  }
514 
515  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
516  (
518  (
520  internalField().group()
521  )
522  );
523 
524  setMaster();
525 
526  if (patch().index() == master_)
527  {
529  calculateTurbulenceFields(turbModel, G(true), omega(true));
530  }
531 
532  const scalarField& G0 = this->G();
533  const scalarField& omega0 = this->omega();
534 
535  FieldType& G =
536  const_cast<FieldType&>
537  (
538  db().lookupObject<FieldType>(turbModel.GName())
539  );
540 
541  FieldType& omega = const_cast<FieldType&>(internalField());
542 
543  scalarField& omegaf = *this;
544 
545  // only set the values if the weights are > tolerance
546  forAll(weights, facei)
547  {
548  scalar w = weights[facei];
549 
550  if (w > tolerance_)
551  {
552  label celli = patch().faceCells()[facei];
553 
554  G[celli] = (1.0 - w)*G[celli] + w*G0[celli];
555  omega[celli] = (1.0 - w)*omega[celli] + w*omega0[celli];
556  omegaf[facei] = omega[celli];
557  }
558  }
559 
561 }
562 
563 
565 (
566  fvMatrix<scalar>& matrix
567 )
568 {
569  if (manipulatedMatrix())
570  {
571  return;
572  }
573 
574  matrix.setValues(patch().faceCells(), patchInternalField());
575 
577 }
578 
579 
581 (
582  fvMatrix<scalar>& matrix,
583  const Field<scalar>& weights
584 )
585 {
586  if (manipulatedMatrix())
587  {
588  return;
589  }
590 
591  DynamicList<label> constraintCells(weights.size());
592  DynamicList<scalar> constraintomega(weights.size());
593  const labelUList& faceCells = patch().faceCells();
594 
596  = internalField();
597 
598  label nConstrainedCells = 0;
599 
600 
601  forAll(weights, facei)
602  {
603  // only set the values if the weights are > tolerance
604  if (weights[facei] > tolerance_)
605  {
606  nConstrainedCells++;
607 
608  label celli = faceCells[facei];
609 
610  constraintCells.append(celli);
611  constraintomega.append(omega[celli]);
612  }
613  }
614 
615  if (debug)
616  {
617  Pout<< "Patch: " << patch().name()
618  << ": number of constrained cells = " << nConstrainedCells
619  << " out of " << patch().size()
620  << endl;
621  }
622 
623  matrix.setValues
624  (
625  constraintCells,
626  scalarField(constraintomega.xfer())
627  );
628 
630 }
631 
632 
634 {
635  writeLocalEntries(os);
637 }
638 
639 
640 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
641 
643 (
646 );
647 
648 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
649 
650 } // End namespace Foam
651 
652 // ************************************************************************* //
const char *const group
Group name for atomic constants.
scalar uPlus
bool changing() const
Is mesh changing (topology changing and/or moving)
Definition: polyMesh.H:530
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:313
Switch blended_
Blending switch (defaults to false)
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
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void write(Ostream &) const
Write.
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:216
scalar Rey
const dimensionedScalar G0
Conductance quantum: default SI units: [S].
dimensionedScalar log(const dimensionedScalar &ds)
scalar yPlusLam() const
Return the Y+ at the edge of the laminar sublayer.
const volVectorField & U() const
Access function to velocity field.
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:179
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
omegaWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:371
const Boundary & boundaryField() const
Return const-reference to the boundary field.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
virtual void updateWeightedCoeffs(const scalarField &weights)
Update the coefficients associated with the patch field.
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
dimensionedScalar pow025(const dimensionedScalar &ds)
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none.
Definition: Switch.H:60
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual void createAveragingWeights()
Create the averaging weights for cells which are bounded by.
virtual void checkType()
Check the type of the patch.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:626
label k
Boltzmann constant.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
bool manipulatedMatrix() const
Return true if the matrix has already been manipulated.
Definition: fvPatchField.H:377
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
word GName() const
Helper function to return the name of the turbulence G field.
Abstract base class for turbulence models (RAS, LES and laminar).
const word & name() const
Return name.
Definition: fvPatch.H:149
Macros for easy insertion into run-time selection tables.
virtual void operator==(const fvPatchField< Type > &)
Definition: fvPatchField.C:561
scalarField & G(bool init=false)
Return non-const access to the master&#39;s G field.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
scalar y
dynamicFvMesh & mesh
virtual omegaWallFunctionFvPatchScalarField & omegaPatch(const label patchi)
Helper function to return non-const access to an omega patch.
dimensionedScalar exp(const dimensionedScalar &ds)
static const word propertiesName
Default name of the turbulence properties dictionary.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:93
static word groupName(Name name, const word &group)
void setValues(const labelUList &cells, const UList< Type > &values)
Set solution in given cells to the specified values.
Definition: fvMatrix.C:486
Foam::fvPatchFieldMapper.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:292
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.
const nearWallDist & y() const
Return the near wall distances.
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:72
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
Definition: fvPatchField.C:344
errorManip< error > abort(error &err)
Definition: errorManip.H:131
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:61
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:340
virtual label size() const
Return size.
Definition: fvPatch.H:161
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:265
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:224
const Mesh & mesh() const
Return mesh.
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
List< List< scalar > > cornerWeights_
List of averaging corner weights.
scalarField G_
Local copy of turbulence G field.
scalarField & omega(bool init=false)
Return non-const access to the master&#39;s omega field.
void setSize(const label)
Reset size of List.
Definition: List.C:281
label patchi
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
This boundary condition provides a wall constraint on turbulnce specific dissipation, omega for both low and high Reynolds number turbulence models.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
scalar yPlus
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:311
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:197
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual void setMaster()
Set the master patch - master is responsible for updating all.
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:98
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
virtual label & master()
Return non-const access to the master patch ID.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
Manipulate matrix.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
label size() const
Return the number of elements in the UList.
Definition: ListI.H:170
scalar yPlusLam_
y+ at the edge of the laminar sublayer
scalarField omega_
Local copy of turbulence omega field.
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.
const DimensionedField< Type, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:347