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-2022 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().timeName(),
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 
113  forAll(epsilonPatches, i)
114  {
115  label patchi = epsilonPatches[i];
116  const fvPatchScalarField& wf = weights.boundaryField()[patchi];
118  }
119 
120  G_.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,
146  scalarField& epsilon0
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
156  {
157  if (!cornerWeights_[patchi].empty())
158  {
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
169  {
170  if (!cornerWeights_[patchi].empty())
171  {
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,
186  scalarField& epsilon0
187 )
188 {
189  const label patchi = patch.index();
190 
192  nutWallFunctionFvPatchScalarField::nutw(turbModel, patchi);
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 )
246 :
248  G_(),
249  epsilon_(),
250  initialised_(false),
251  master_(-1),
253 {}
254 
255 
258 (
259  const fvPatch& p,
261  const dictionary& dict
262 )
263 :
265  G_(),
266  epsilon_(),
267  initialised_(false),
268  master_(-1),
270 {
271  // Apply zero-gradient condition on start-up
273 }
274 
275 
278 (
280  const fvPatch& p,
282  const fvPatchFieldMapper& mapper
283 )
284 :
285  fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
286  G_(),
287  epsilon_(),
288  initialised_(false),
289  master_(-1),
291 {}
292 
293 
296 (
299 )
300 :
301  fixedValueFvPatchField<scalar>(ewfpsf, iF),
302  G_(),
303  epsilon_(),
304  initialised_(false),
305  master_(-1),
307 {}
308 
309 
310 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
311 
313 {
314  if (patch().index() == master_)
315  {
316  if (init)
317  {
318  G_ = 0.0;
319  }
320 
321  return G_;
322  }
323 
324  return epsilonPatch(master_).G();
325 }
326 
327 
329 (
330  bool init
331 )
332 {
333  if (patch().index() == master_)
334  {
335  if (init)
336  {
337  epsilon_ = 0.0;
338  }
339 
340  return epsilon_;
341  }
342 
343  return epsilonPatch(master_).epsilon(init);
344 }
345 
346 
348 {
349  if (updated())
350  {
351  return;
352  }
353 
354  const momentumTransportModel& turbModel =
356  (
358  (
359  momentumTransportModel::typeName,
360  internalField().group()
361  )
362  );
363 
364  setMaster();
365 
366  if (patch().index() == master_)
367  {
369  calculateTurbulenceFields(turbModel, G(true), epsilon(true));
370  }
371 
372  const scalarField& G0 = this->G();
373  const scalarField& epsilon0 = this->epsilon();
374 
375  typedef DimensionedField<scalar, volMesh> FieldType;
376  FieldType& G =
377  const_cast<FieldType&>
378  (
379  db().lookupObject<FieldType>(turbModel.GName())
380  );
381  FieldType& epsilon =
382  const_cast<FieldType&>
383  (
384  internalField()
385  );
386 
387  scalarField weights(patch().magSf()/patch().patch().magFaceAreas());
388  forAll(weights, facei)
389  {
390  scalar& w = weights[facei];
391  w = w <= tolerance_ ? 0 : (w - tolerance_)/(1 - tolerance_);
392  }
393 
394  forAll(weights, facei)
395  {
396  const scalar w = weights[facei];
397  const label celli = patch().faceCells()[facei];
398 
399  G[celli] = (1 - w)*G[celli] + w*G0[celli];
400  epsilon[celli] = (1 - w)*epsilon[celli] + w*epsilon0[celli];
401  }
402 
403  this->operator==(scalarField(epsilon, patch().faceCells()));
404 
406 }
407 
408 
410 (
411  fvMatrix<scalar>& matrix
412 )
413 {
414  if (manipulatedMatrix())
415  {
416  return;
417  }
418 
419  const scalarField& epsilon0 = this->epsilon();
420 
421  scalarField weights(patch().magSf()/patch().patch().magFaceAreas());
422  forAll(weights, facei)
423  {
424  scalar& w = weights[facei];
425  w = w <= tolerance_ ? 0 : (w - tolerance_)/(1 - tolerance_);
426  }
427 
428  matrix.setValues
429  (
430  patch().faceCells(),
431  UIndirectList<scalar>(epsilon0, patch().faceCells()),
432  weights
433  );
434 
436 }
437 
438 
439 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
440 
441 namespace Foam
442 {
444  (
447  );
448 }
449 
450 
451 // ************************************************************************* //
const char *const group
Group name for atomic constants.
bool changing() const
Is mesh changing (topology changing and/or moving)
Definition: polyMesh.H:537
dictionary dict
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:325
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:166
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:175
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:369
const Boundary & boundaryField() const
Return const-reference to the boundary field.
scalarField G_
Local copy of turbulence G field.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
dimensionedScalar pow025(const dimensionedScalar &ds)
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
Manipulate matrix.
static scalar tolerance_
Tolerance used in weighted calculations.
virtual void calculate(const momentumTransportModel &turbModel, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &epsilon)
Calculate the epsilon and G.
virtual void calculateTurbulenceFields(const momentumTransportModel &turbModel, scalarField &G0, scalarField &epsilon0)
Main driver to calculate the turbulence fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
const dimensionSet dimless
label k
Boltzmann constant.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:372
bool manipulatedMatrix() const
Return true if the matrix has already been manipulated.
Definition: fvPatchField.H:375
fvMesh & mesh
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
const volScalarField::Boundary & y() const
Return the near wall distance.
Macros for easy insertion into run-time selection tables.
virtual void operator==(const fvPatchField< Type > &)
Definition: fvPatchField.C:424
static const nutWallFunctionFvPatchScalarField & nutw(const momentumTransportModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
scalar y
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions...
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
virtual label & master()
Return non-const access to the master patch ID.
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:666
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:99
static word groupName(Name name, const word &group)
Foam::fvPatchFieldMapper.
scalarField epsilon_
Local copy of turbulence epsilon field.
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:236
This boundary condition provides a turbulence dissipation wall constraint for low- and high-Reynolds ...
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:60
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:351
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.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:174
const Mesh & mesh() const
Return mesh.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Abstract base class for turbulence models (RAS, LES and laminar).
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
void setValues(const labelUList &cells, const ListType< Type > &values)
Set solution in given cells to the specified values.
Definition: fvMatrix.C:502
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void setSize(const label)
Reset size of List.
Definition: List.C:281
const volVectorField & U() const
Access function to velocity field.
label patchi
List< List< scalar > > cornerWeights_
List of averaging corner weights.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionedScalar G0
Conductance quantum: default SI units: [S].
scalar yPlus
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:216
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
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:147
dimensioned< scalar > mag(const dimensioned< Type > &)
const dimensionedScalar epsilon0
Electric constant: default SI units: [F/m].
word GName() const
Helper function to return the name of the turbulence G field.
virtual void setMaster()
Set the master patch - master is responsible for updating all.
static scalar yPlusLam(const scalar kappa, const scalar E)
Calculate the Y+ at the edge of the laminar sublayer.
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
virtual void createAveragingWeights()
Create the averaging weights for cells which are bounded by.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
epsilonWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
label size() const
Return the number of elements in the UList.
Definition: ListI.H:171
const dimensionedScalar e
Elementary charge.
Definition: doubleScalar.H:105
Namespace for OpenFOAM.
const DimensionedField< Type, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:357
scalarField & epsilon(bool init=false)
Return non-const access to the master&#39;s epsilon field.