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 | Copyright (C) 2011-2016 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 
223  const tmp<volScalarField> tk = turbModel.k();
224  const volScalarField& k = tk();
225 
226  const tmp<scalarField> tnuw = turbModel.nu(patchi);
227  const scalarField& nuw = tnuw();
228 
229  const tmp<scalarField> tnutw = turbModel.nut(patchi);
230  const scalarField& nutw = tnutw();
231 
232  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
233 
234  const scalarField magGradUw(mag(Uw.snGrad()));
235 
236  // Set omega and G
237  forAll(nutw, facei)
238  {
239  const label celli = patch.faceCells()[facei];
240 
241  const scalar yPlus = Cmu25*y[facei]*sqrt(k[celli])/nuw[facei];
242 
243  const scalar w = cornerWeights[facei];
244 
245  const scalar omegaVis = 6*nuw[facei]/(beta1_*sqr(y[facei]));
246  const scalar omegaLog = sqrt(k[celli])/(Cmu25*kappa_*y[facei]);
247 
248  // Switching between the laminar sub-layer and the log-region rather
249  // than blending has been found to provide more accurate results over a
250  // range of near-wall y+.
251  //
252  // For backward-compatibility the blending method is provided as an
253  // option
254 
255  if (blended_)
256  {
257  omega0[celli] += w*sqrt(sqr(omegaVis) + sqr(omegaLog));
258  }
259 
260  if (yPlus > yPlusLam_)
261  {
262  if (!blended_)
263  {
264  omega0[celli] += w*omegaLog;
265  }
266 
267  G0[celli] +=
268  w
269  *(nutw[facei] + nuw[facei])
270  *magGradUw[facei]
271  *Cmu25*sqrt(k[celli])
272  /(kappa_*y[facei]);
273  }
274  else
275  {
276  if (!blended_)
277  {
278  omega0[celli] += w*omegaVis;
279  }
280  }
281  }
282 }
283 
284 
285 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
286 
288 (
289  const fvPatch& p,
291 )
292 :
294  Cmu_(0.09),
295  kappa_(0.41),
296  E_(9.8),
297  beta1_(0.075),
298  blended_(false),
300  G_(),
301  omega_(),
302  initialised_(false),
303  master_(-1),
305 {
306  checkType();
307 }
308 
309 
311 (
313  const fvPatch& p,
315  const fvPatchFieldMapper& mapper
316 )
317 :
318  fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
319  Cmu_(ptf.Cmu_),
320  kappa_(ptf.kappa_),
321  E_(ptf.E_),
322  beta1_(ptf.beta1_),
323  blended_(ptf.blended_),
324  yPlusLam_(ptf.yPlusLam_),
325  G_(),
326  omega_(),
327  initialised_(false),
328  master_(-1),
330 {
331  checkType();
332 }
333 
334 
336 (
337  const fvPatch& p,
339  const dictionary& dict
340 )
341 :
343  Cmu_(dict.lookupOrDefault<scalar>("Cmu", 0.09)),
344  kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
345  E_(dict.lookupOrDefault<scalar>("E", 9.8)),
346  beta1_(dict.lookupOrDefault<scalar>("beta1", 0.075)),
347  blended_(dict.lookupOrDefault<Switch>("blended", false)),
349  G_(),
350  omega_(),
351  initialised_(false),
352  master_(-1),
354 {
355  checkType();
356 
357  // apply zero-gradient condition on start-up
359 }
360 
361 
363 (
365 )
366 :
368  Cmu_(owfpsf.Cmu_),
369  kappa_(owfpsf.kappa_),
370  E_(owfpsf.E_),
371  beta1_(owfpsf.beta1_),
372  blended_(owfpsf.blended_),
373  yPlusLam_(owfpsf.yPlusLam_),
374  G_(),
375  omega_(),
376  initialised_(false),
377  master_(-1),
379 {
380  checkType();
381 }
382 
383 
385 (
388 )
389 :
390  fixedValueFvPatchField<scalar>(owfpsf, iF),
391  Cmu_(owfpsf.Cmu_),
392  kappa_(owfpsf.kappa_),
393  E_(owfpsf.E_),
394  beta1_(owfpsf.beta1_),
395  blended_(owfpsf.blended_),
396  yPlusLam_(owfpsf.yPlusLam_),
397  G_(),
398  omega_(),
399  initialised_(false),
400  master_(-1),
402 {
403  checkType();
404 }
405 
406 
407 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
408 
410 {
411  if (patch().index() == master_)
412  {
413  if (init)
414  {
415  G_ = 0.0;
416  }
417 
418  return G_;
419  }
420 
421  return omegaPatch(master_).G();
422 }
423 
424 
426 {
427  if (patch().index() == master_)
428  {
429  if (init)
430  {
431  omega_ = 0.0;
432  }
433 
434  return omega_;
435  }
436 
437  return omegaPatch(master_).omega(init);
438 }
439 
440 
442 {
443  if (updated())
444  {
445  return;
446  }
447 
448  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
449  (
451  (
453  internalField().group()
454  )
455  );
456 
457  setMaster();
458 
459  if (patch().index() == master_)
460  {
462  calculateTurbulenceFields(turbModel, G(true), omega(true));
463  }
464 
465  const scalarField& G0 = this->G();
466  const scalarField& omega0 = this->omega();
467 
468  typedef DimensionedField<scalar, volMesh> FieldType;
469 
470  FieldType& G =
471  const_cast<FieldType&>
472  (
473  db().lookupObject<FieldType>(turbModel.GName())
474  );
475 
476  FieldType& omega = const_cast<FieldType&>(internalField());
477 
478  forAll(*this, facei)
479  {
480  label celli = patch().faceCells()[facei];
481 
482  G[celli] = G0[celli];
483  omega[celli] = omega0[celli];
484  }
485 
487 }
488 
489 
491 (
492  const scalarField& weights
493 )
494 {
495  if (updated())
496  {
497  return;
498  }
499 
500  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
501  (
503  (
505  internalField().group()
506  )
507  );
508 
509  setMaster();
510 
511  if (patch().index() == master_)
512  {
514  calculateTurbulenceFields(turbModel, G(true), omega(true));
515  }
516 
517  const scalarField& G0 = this->G();
518  const scalarField& omega0 = this->omega();
519 
520  typedef DimensionedField<scalar, volMesh> FieldType;
521 
522  FieldType& G =
523  const_cast<FieldType&>
524  (
525  db().lookupObject<FieldType>(turbModel.GName())
526  );
527 
528  FieldType& omega = const_cast<FieldType&>(internalField());
529 
530  scalarField& omegaf = *this;
531 
532  // only set the values if the weights are > tolerance
533  forAll(weights, facei)
534  {
535  scalar w = weights[facei];
536 
537  if (w > tolerance_)
538  {
539  label celli = patch().faceCells()[facei];
540 
541  G[celli] = (1.0 - w)*G[celli] + w*G0[celli];
542  omega[celli] = (1.0 - w)*omega[celli] + w*omega0[celli];
543  omegaf[facei] = omega[celli];
544  }
545  }
546 
548 }
549 
550 
552 (
553  fvMatrix<scalar>& matrix
554 )
555 {
556  if (manipulatedMatrix())
557  {
558  return;
559  }
560 
561  matrix.setValues(patch().faceCells(), patchInternalField());
562 
564 }
565 
566 
568 (
569  fvMatrix<scalar>& matrix,
570  const Field<scalar>& weights
571 )
572 {
573  if (manipulatedMatrix())
574  {
575  return;
576  }
577 
578  DynamicList<label> constraintCells(weights.size());
579  DynamicList<scalar> constraintomega(weights.size());
580  const labelUList& faceCells = patch().faceCells();
581 
583  = internalField();
584 
585  label nConstrainedCells = 0;
586 
587 
588  forAll(weights, facei)
589  {
590  // only set the values if the weights are > tolerance
591  if (weights[facei] > tolerance_)
592  {
593  nConstrainedCells++;
594 
595  label celli = faceCells[facei];
596 
597  constraintCells.append(celli);
598  constraintomega.append(omega[celli]);
599  }
600  }
601 
602  if (debug)
603  {
604  Pout<< "Patch: " << patch().name()
605  << ": number of constrained cells = " << nConstrainedCells
606  << " out of " << patch().size()
607  << endl;
608  }
609 
610  matrix.setValues
611  (
612  constraintCells,
613  scalarField(constraintomega.xfer())
614  );
615 
617 }
618 
619 
621 {
622  writeLocalEntries(os);
624 }
625 
626 
627 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
628 
630 (
633 );
634 
635 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
636 
637 } // End namespace Foam
638 
639 // ************************************************************************* //
const char *const group
Group name for atomic constants.
bool changing() const
Is mesh changing (topology changing and/or moving)
Definition: polyMesh.H:524
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
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
const dimensionedScalar G0
Conductance quantum: default SI units: [S].
scalar yPlusLam() const
Return the Y+ at the edge of the laminar sublayer.
const double e
Elementary charge.
Definition: doubleFloat.H:78
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:253
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:644
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.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
scalar y
dynamicFvMesh & mesh
virtual omegaWallFunctionFvPatchScalarField & omegaPatch(const label patchi)
Helper function to return non-const access to an omega patch.
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:262
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.
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.
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