epsilonWallFunctionFvPatchScalarField.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-2023 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 "momentumTransportModel.H"
29 #include "fvMatrix.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
35 
36 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
37 
39 {
40  if (master_ != -1)
41  {
42  return;
43  }
44 
45  const volScalarField& epsilon =
46  static_cast<const volScalarField&>(this->internalField());
47 
48  const volScalarField::Boundary& bf = epsilon.boundaryField();
49 
50  label master = -1;
51  forAll(bf, patchi)
52  {
53  if (isA<epsilonWallFunctionFvPatchScalarField>(bf[patchi]))
54  {
56 
57  if (master == -1)
58  {
59  master = patchi;
60  }
61 
62  epf.master() = master;
63  }
64  }
65 }
66 
67 
69 {
70  const volScalarField& epsilon =
71  static_cast<const volScalarField&>(this->internalField());
72 
73  const volScalarField::Boundary& bf = epsilon.boundaryField();
74 
75  const fvMesh& mesh = epsilon.mesh();
76 
77  if (initialised_ && !mesh.changing())
78  {
79  return;
80  }
81 
82  volScalarField weights
83  (
84  IOobject
85  (
86  "weights",
87  mesh.time().name(),
88  mesh,
91  false // do not register
92  ),
93  mesh,
95  );
96 
97  DynamicList<label> epsilonPatches(bf.size());
98  forAll(bf, patchi)
99  {
100  if (isA<epsilonWallFunctionFvPatchScalarField>(bf[patchi]))
101  {
102  epsilonPatches.append(patchi);
103 
104  const labelUList& faceCells = bf[patchi].patch().faceCells();
105  forAll(faceCells, i)
106  {
107  weights[faceCells[i]]++;
108  }
109  }
110  }
111 
112  cornerWeights_.setSize(bf.size());
113  forAll(epsilonPatches, i)
114  {
115  label patchi = epsilonPatches[i];
116  const fvPatchScalarField& wf = weights.boundaryField()[patchi];
117  cornerWeights_[patchi] = 1.0/wf.patchInternalField();
118  }
119 
120  G_.setSize(internalField().size(), 0.0);
121  epsilon_.setSize(internalField().size(), 0.0);
122 
123  initialised_ = true;
124 }
125 
126 
129 {
130  const volScalarField& epsilon =
131  static_cast<const volScalarField&>(this->internalField());
132 
133  const volScalarField::Boundary& bf = epsilon.boundaryField();
134 
136  refCast<const epsilonWallFunctionFvPatchScalarField>(bf[patchi]);
137 
138  return const_cast<epsilonWallFunctionFvPatchScalarField&>(epf);
139 }
140 
141 
143 (
144  const momentumTransportModel& turbModel,
145  scalarField& G0,
147 )
148 {
149  // Ask for the wall distance in advance. Some processes may not have any
150  // corner weights, so we cannot allow functions inside the if statements
151  // below to trigger the wall distance calculation.
152  turbModel.y();
153 
154  // Accumulate all of the G and epsilon contributions
155  forAll(cornerWeights_, patchi)
156  {
157  if (!cornerWeights_[patchi].empty())
158  {
159  epsilonWallFunctionFvPatchScalarField& epf = epsilonPatch(patchi);
160 
161  const List<scalar>& w = cornerWeights_[patchi];
162 
163  epf.calculate(turbModel, w, epf.patch(), G0, epsilon0);
164  }
165  }
166 
167  // Apply zero-gradient condition for epsilon
168  forAll(cornerWeights_, patchi)
169  {
170  if (!cornerWeights_[patchi].empty())
171  {
172  epsilonWallFunctionFvPatchScalarField& epf = epsilonPatch(patchi);
173 
174  epf == scalarField(epsilon0, epf.patch().faceCells());
175  }
176  }
177 }
178 
179 
181 (
182  const momentumTransportModel& turbModel,
183  const List<scalar>& cornerWeights,
184  const fvPatch& patch,
185  scalarField& G0,
187 )
188 {
189  const label patchi = patch.index();
190 
193 
194  const scalarField& y = turbModel.y()[patchi];
195 
196  const tmp<scalarField> tnuw = turbModel.nu(patchi);
197  const scalarField& nuw = tnuw();
198 
199  const tmp<volScalarField> tk = turbModel.k();
200  const volScalarField& k = tk();
201 
202  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
203 
204  const scalarField magGradUw(mag(Uw.snGrad()));
205 
206  const scalar Cmu25 = pow025(nutw.Cmu());
207  const scalar Cmu75 = pow(nutw.Cmu(), 0.75);
208 
209  // Set epsilon and G
210  forAll(nutw, facei)
211  {
212  const label celli = patch.faceCells()[facei];
213 
214  const scalar yPlus = Cmu25*y[facei]*sqrt(k[celli])/nuw[facei];
215 
216  const scalar w = cornerWeights[facei];
217 
218  if (yPlus > nutw.yPlusLam())
219  {
220  epsilon0[celli] +=
221  w*Cmu75*pow(k[celli], 1.5)/(nutw.kappa()*y[facei]);
222 
223  G0[celli] +=
224  w
225  *(nutw[facei] + nuw[facei])
226  *magGradUw[facei]
227  *Cmu25*sqrt(k[celli])
228  /(nutw.kappa()*y[facei]);
229  }
230  else
231  {
232  epsilon0[celli] += w*2.0*k[celli]*nuw[facei]/sqr(y[facei]);
233  }
234  }
235 }
236 
237 
238 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
239 
242 (
243  const fvPatch& p,
245  const dictionary& dict
246 )
247 :
248  fixedValueFvPatchField<scalar>(p, iF, dict),
249  G_(),
250  epsilon_(),
251  initialised_(false),
252  master_(-1),
253  cornerWeights_()
254 {
255  // Apply zero-gradient condition on start-up
257 }
258 
259 
262 (
264  const fvPatch& p,
266  const fvPatchFieldMapper& mapper
267 )
268 :
269  fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
270  G_(),
271  epsilon_(),
272  initialised_(false),
273  master_(-1),
274  cornerWeights_()
275 {}
276 
277 
280 (
283 )
284 :
285  fixedValueFvPatchField<scalar>(ewfpsf, iF),
286  G_(),
287  epsilon_(),
288  initialised_(false),
289  master_(-1),
290  cornerWeights_()
291 {}
292 
293 
294 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
295 
297 {
298  if (patch().index() == master_)
299  {
300  if (init)
301  {
302  G_ = 0.0;
303  }
304 
305  return G_;
306  }
307 
308  return epsilonPatch(master_).G();
309 }
310 
311 
313 (
314  bool init
315 )
316 {
317  if (patch().index() == master_)
318  {
319  if (init)
320  {
321  epsilon_ = 0.0;
322  }
323 
324  return epsilon_;
325  }
326 
327  return epsilonPatch(master_).epsilon(init);
328 }
329 
330 
332 {
333  if (updated())
334  {
335  return;
336  }
337 
338  const momentumTransportModel& turbModel =
339  db().lookupType<momentumTransportModel>(internalField().group());
340 
341  setMaster();
342 
343  if (patch().index() == master_)
344  {
345  createAveragingWeights();
346  calculateTurbulenceFields(turbModel, G(true), epsilon(true));
347  }
348 
349  const scalarField& G0 = this->G();
350  const scalarField& epsilon0 = this->epsilon();
351 
352  typedef DimensionedField<scalar, volMesh> FieldType;
353  FieldType& G =
354  const_cast<FieldType&>
355  (
356  db().lookupObject<FieldType>(turbModel.GName())
357  );
358  FieldType& epsilon =
359  const_cast<FieldType&>
360  (
361  internalField()
362  );
363 
364  scalarField weights(patch().magSf()/patch().patch().magFaceAreas());
365  forAll(weights, facei)
366  {
367  scalar& w = weights[facei];
368  w = w <= tolerance_ ? 0 : (w - tolerance_)/(1 - tolerance_);
369  }
370 
371  forAll(weights, facei)
372  {
373  const scalar w = weights[facei];
374  const label celli = patch().faceCells()[facei];
375 
376  G[celli] = (1 - w)*G[celli] + w*G0[celli];
377  epsilon[celli] = (1 - w)*epsilon[celli] + w*epsilon0[celli];
378  }
379 
380  this->operator==(scalarField(epsilon, patch().faceCells()));
381 
383 }
384 
385 
387 (
388  fvMatrix<scalar>& matrix
389 )
390 {
391  if (manipulatedMatrix())
392  {
393  return;
394  }
395 
396  const scalarField& epsilon0 = this->epsilon();
397 
398  scalarField weights(patch().magSf()/patch().patch().magFaceAreas());
399  forAll(weights, facei)
400  {
401  scalar& w = weights[facei];
402  w = w <= tolerance_ ? 0 : (w - tolerance_)/(1 - tolerance_);
403  }
404 
405  matrix.setValues
406  (
407  patch().faceCells(),
408  UIndirectList<scalar>(epsilon0, patch().faceCells()),
409  weights
410  );
411 
413 }
414 
415 
416 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
417 
418 namespace Foam
419 {
421  (
424  );
425 }
426 
427 
428 // ************************************************************************* //
scalar y
label k
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
Generic GeometricBoundaryField class.
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
void setSize(const label)
Reset size of List.
Definition: List.C:281
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
const word & name() const
Return const reference to name.
This boundary condition provides a turbulence dissipation wall constraint for low- and high-Reynolds ...
virtual label & master()
Return non-const access to the master patch ID.
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
Manipulate matrix.
static scalar tolerance_
Tolerance used in weighted calculations.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void createAveragingWeights()
Create the averaging weights for cells which are bounded by.
epsilonWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
scalarField & G(bool init=false)
Return non-const access to the master's G field.
virtual void calculateTurbulenceFields(const momentumTransportModel &turbModel, scalarField &G0, scalarField &epsilon0)
Main driver to calculate the turbulence fields.
virtual epsilonWallFunctionFvPatchScalarField & epsilonPatch(const label patchi)
Helper function to return non-const access to an epsilon patch.
scalarField & epsilon(bool init=false)
Return non-const access to the master's epsilon field.
virtual void setMaster()
Set the master patch - master is responsible for updating all.
virtual void calculate(const momentumTransportModel &turbModel, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &epsilon)
Calculate the epsilon and G.
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
void setValues(const labelUList &cells, const ListType< Type > &values)
Set solution in given cells to the specified values.
Definition: fvMatrix.C:502
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:402
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
const DimensionedField< Type, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:361
virtual void operator==(const fvPatchField< Type > &)
Definition: fvPatchField.C:417
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:172
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:164
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:204
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
Definition: fvPatchField.C:224
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:355
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:176
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:99
Abstract base class for turbulence models (RAS, LES and laminar).
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
word GName() const
Helper function to return the name of the turbulence G field.
const volVectorField & U() const
Access function to velocity field.
const volScalarField::Boundary & y() const
Return the near wall distance.
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions,...
static scalar yPlusLam(const scalar kappa, const scalar E)
Calculate the Y+ at the edge of the laminar sublayer.
static const nutWallFunctionFvPatchScalarField & nutw(const momentumTransportModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
bool changing() const
Is mesh changing.
Definition: polyMesh.H:490
A class for managing temporary objects.
Definition: tmp.H:55
const scalar yPlus
const scalar epsilon
label patchi
const dimensionedScalar G0
Conductance quantum: default SI units: [S].
const dimensionedScalar e
Elementary charge.
const dimensionedScalar epsilon0
Electric constant: default SI units: [F/m].
const dimensionedScalar G
Newtonian constant of gravitation.
Namespace for OpenFOAM.
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
const dimensionSet dimless
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dimensionedScalar pow025(const dimensionedScalar &ds)
dictionary dict
volScalarField & p