epsilonWallFunctionFvPatchScalarField.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 
27 #include "turbulenceModel.H"
28 #include "fvPatchFieldMapper.H"
29 #include "fvMatrix.H"
30 #include "volFields.H"
31 #include "wallFvPatch.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
37 
38 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
39 
41 {
42  if (!isA<wallFvPatch>(patch()))
43  {
45  << "Invalid wall function specification" << nl
46  << " Patch type for patch " << patch().name()
47  << " must be wall" << nl
48  << " Current patch type is " << patch().type() << nl << endl
49  << abort(FatalError);
50  }
51 }
52 
53 
55 (
56  Ostream& os
57 ) const
58 {
59  os.writeKeyword("Cmu") << Cmu_ << token::END_STATEMENT << nl;
60  os.writeKeyword("kappa") << kappa_ << token::END_STATEMENT << nl;
61  os.writeKeyword("E") << E_ << token::END_STATEMENT << nl;
62 }
63 
64 
66 {
67  if (master_ != -1)
68  {
69  return;
70  }
71 
72  const volScalarField& epsilon =
73  static_cast<const volScalarField&>(this->internalField());
74 
75  const volScalarField::Boundary& bf = epsilon.boundaryField();
76 
77  label master = -1;
78  forAll(bf, patchi)
79  {
80  if (isA<epsilonWallFunctionFvPatchScalarField>(bf[patchi]))
81  {
83 
84  if (master == -1)
85  {
86  master = patchi;
87  }
88 
89  epf.master() = master;
90  }
91  }
92 }
93 
94 
96 {
97  const volScalarField& epsilon =
98  static_cast<const volScalarField&>(this->internalField());
99 
100  const volScalarField::Boundary& bf = epsilon.boundaryField();
101 
102  const fvMesh& mesh = epsilon.mesh();
103 
104  if (initialised_ && !mesh.changing())
105  {
106  return;
107  }
108 
109  volScalarField weights
110  (
111  IOobject
112  (
113  "weights",
114  mesh.time().timeName(),
115  mesh,
118  false // do not register
119  ),
120  mesh,
121  dimensionedScalar("zero", dimless, 0.0)
122  );
123 
124  DynamicList<label> epsilonPatches(bf.size());
125  forAll(bf, patchi)
126  {
127  if (isA<epsilonWallFunctionFvPatchScalarField>(bf[patchi]))
128  {
129  epsilonPatches.append(patchi);
130 
131  const labelUList& faceCells = bf[patchi].patch().faceCells();
132  forAll(faceCells, i)
133  {
134  weights[faceCells[i]]++;
135  }
136  }
137  }
138 
139  cornerWeights_.setSize(bf.size());
140  forAll(epsilonPatches, i)
141  {
142  label patchi = epsilonPatches[i];
143  const fvPatchScalarField& wf = weights.boundaryField()[patchi];
145  }
146 
147  G_.setSize(internalField().size(), 0.0);
149 
150  initialised_ = true;
151 }
152 
153 
156 {
157  const volScalarField& epsilon =
158  static_cast<const volScalarField&>(this->internalField());
159 
160  const volScalarField::Boundary& bf = epsilon.boundaryField();
161 
163  refCast<const epsilonWallFunctionFvPatchScalarField>(bf[patchi]);
164 
165  return const_cast<epsilonWallFunctionFvPatchScalarField&>(epf);
166 }
167 
168 
170 (
171  const turbulenceModel& turbulence,
172  scalarField& G0,
173  scalarField& epsilon0
174 )
175 {
176  // accumulate all of the G and epsilon contributions
178  {
179  if (!cornerWeights_[patchi].empty())
180  {
182 
183  const List<scalar>& w = cornerWeights_[patchi];
184 
185  epf.calculate(turbulence, w, epf.patch(), G0, epsilon0);
186  }
187  }
188 
189  // apply zero-gradient condition for epsilon
191  {
192  if (!cornerWeights_[patchi].empty())
193  {
195 
196  epf == scalarField(epsilon0, epf.patch().faceCells());
197  }
198  }
199 }
200 
201 
203 (
204  const turbulenceModel& turbulence,
205  const List<scalar>& cornerWeights,
206  const fvPatch& patch,
207  scalarField& G,
208  scalarField& epsilon
209 )
210 {
211  const label patchi = patch.index();
212 
213  const scalarField& y = turbulence.y()[patchi];
214 
215  const scalar Cmu25 = pow025(Cmu_);
216  const scalar Cmu75 = pow(Cmu_, 0.75);
217 
218  const tmp<volScalarField> tk = turbulence.k();
219  const volScalarField& k = tk();
220 
221  const tmp<scalarField> tnuw = turbulence.nu(patchi);
222  const scalarField& nuw = tnuw();
223 
224  const tmp<scalarField> tnutw = turbulence.nut(patchi);
225  const scalarField& nutw = tnutw();
226 
227  const fvPatchVectorField& Uw = turbulence.U().boundaryField()[patchi];
228 
229  const scalarField magGradUw(mag(Uw.snGrad()));
230 
231  // Set epsilon and G
232  forAll(nutw, facei)
233  {
234  label celli = patch.faceCells()[facei];
235 
236  scalar w = cornerWeights[facei];
237 
238  epsilon[celli] += w*Cmu75*pow(k[celli], 1.5)/(kappa_*y[facei]);
239 
240  G[celli] +=
241  w
242  *(nutw[facei] + nuw[facei])
243  *magGradUw[facei]
244  *Cmu25*sqrt(k[celli])
245  /(kappa_*y[facei]);
246  }
247 }
248 
249 
250 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
251 
254 (
255  const fvPatch& p,
257 )
258 :
260  Cmu_(0.09),
261  kappa_(0.41),
262  E_(9.8),
263  G_(),
264  epsilon_(),
265  initialised_(false),
266  master_(-1),
268 {
269  checkType();
270 }
271 
272 
275 (
277  const fvPatch& p,
279  const fvPatchFieldMapper& mapper
280 )
281 :
282  fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
283  Cmu_(ptf.Cmu_),
284  kappa_(ptf.kappa_),
285  E_(ptf.E_),
286  G_(),
287  epsilon_(),
288  initialised_(false),
289  master_(-1),
291 {
292  checkType();
293 }
294 
295 
298 (
299  const fvPatch& p,
301  const dictionary& dict
302 )
303 :
305  Cmu_(dict.lookupOrDefault<scalar>("Cmu", 0.09)),
306  kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
307  E_(dict.lookupOrDefault<scalar>("E", 9.8)),
308  G_(),
309  epsilon_(),
310  initialised_(false),
311  master_(-1),
313 {
314  checkType();
315 
316  // apply zero-gradient condition on start-up
318 }
319 
320 
323 (
325 )
326 :
328  Cmu_(ewfpsf.Cmu_),
329  kappa_(ewfpsf.kappa_),
330  E_(ewfpsf.E_),
331  G_(),
332  epsilon_(),
333  initialised_(false),
334  master_(-1),
336 {
337  checkType();
338 }
339 
340 
343 (
346 )
347 :
348  fixedValueFvPatchField<scalar>(ewfpsf, iF),
349  Cmu_(ewfpsf.Cmu_),
350  kappa_(ewfpsf.kappa_),
351  E_(ewfpsf.E_),
352  G_(),
353  epsilon_(),
354  initialised_(false),
355  master_(-1),
357 {
358  checkType();
359 }
360 
361 
362 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
363 
365 {
366  if (patch().index() == master_)
367  {
368  if (init)
369  {
370  G_ = 0.0;
371  }
372 
373  return G_;
374  }
375 
376  return epsilonPatch(master_).G();
377 }
378 
379 
381 (
382  bool init
383 )
384 {
385  if (patch().index() == master_)
386  {
387  if (init)
388  {
389  epsilon_ = 0.0;
390  }
391 
392  return epsilon_;
393  }
394 
395  return epsilonPatch(master_).epsilon(init);
396 }
397 
398 
400 {
401  if (updated())
402  {
403  return;
404  }
405 
406  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
407  (
409  (
411  internalField().group()
412  )
413  );
414 
415  setMaster();
416 
417  if (patch().index() == master_)
418  {
420  calculateTurbulenceFields(turbModel, G(true), epsilon(true));
421  }
422 
423  const scalarField& G0 = this->G();
424  const scalarField& epsilon0 = this->epsilon();
425 
426  typedef DimensionedField<scalar, volMesh> FieldType;
427 
428  FieldType& G =
429  const_cast<FieldType&>
430  (
431  db().lookupObject<FieldType>(turbModel.GName())
432  );
433 
434  FieldType& epsilon = const_cast<FieldType&>(internalField());
435 
436  forAll(*this, facei)
437  {
438  label celli = patch().faceCells()[facei];
439 
440  G[celli] = G0[celli];
441  epsilon[celli] = epsilon0[celli];
442  }
443 
445 }
446 
447 
449 (
450  const scalarField& weights
451 )
452 {
453  if (updated())
454  {
455  return;
456  }
457 
458  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
459  (
461  (
463  internalField().group()
464  )
465  );
466 
467  setMaster();
468 
469  if (patch().index() == master_)
470  {
472  calculateTurbulenceFields(turbModel, G(true), epsilon(true));
473  }
474 
475  const scalarField& G0 = this->G();
476  const scalarField& epsilon0 = this->epsilon();
477 
478  typedef DimensionedField<scalar, volMesh> FieldType;
479 
480  FieldType& G =
481  const_cast<FieldType&>
482  (
483  db().lookupObject<FieldType>(turbModel.GName())
484  );
485 
486  FieldType& epsilon = const_cast<FieldType&>(internalField());
487 
488  scalarField& epsilonf = *this;
489 
490  // only set the values if the weights are > tolerance
491  forAll(weights, facei)
492  {
493  scalar w = weights[facei];
494 
495  if (w > tolerance_)
496  {
497  label celli = patch().faceCells()[facei];
498 
499  G[celli] = (1.0 - w)*G[celli] + w*G0[celli];
500  epsilon[celli] = (1.0 - w)*epsilon[celli] + w*epsilon0[celli];
501  epsilonf[facei] = epsilon[celli];
502  }
503  }
504 
506 }
507 
508 
510 (
511  fvMatrix<scalar>& matrix
512 )
513 {
514  if (manipulatedMatrix())
515  {
516  return;
517  }
518 
519  matrix.setValues(patch().faceCells(), patchInternalField());
520 
522 }
523 
524 
526 (
527  fvMatrix<scalar>& matrix,
528  const Field<scalar>& weights
529 )
530 {
531  if (manipulatedMatrix())
532  {
533  return;
534  }
535 
536  DynamicList<label> constraintCells(weights.size());
537  DynamicList<scalar> constraintEpsilon(weights.size());
538  const labelUList& faceCells = patch().faceCells();
539 
540  const DimensionedField<scalar, volMesh>& epsilon
541  = internalField();
542 
543  label nConstrainedCells = 0;
544 
545 
546  forAll(weights, facei)
547  {
548  // only set the values if the weights are > tolerance
549  if (weights[facei] > tolerance_)
550  {
551  nConstrainedCells++;
552 
553  label celli = faceCells[facei];
554 
555  constraintCells.append(celli);
556  constraintEpsilon.append(epsilon[celli]);
557  }
558  }
559 
560  if (debug)
561  {
562  Pout<< "Patch: " << patch().name()
563  << ": number of constrained cells = " << nConstrainedCells
564  << " out of " << patch().size()
565  << endl;
566  }
567 
568  matrix.setValues
569  (
570  constraintCells,
571  scalarField(constraintEpsilon.xfer())
572  );
573 
575 }
576 
577 
579 {
580  writeLocalEntries(os);
582 }
583 
584 
585 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
586 
587 namespace Foam
588 {
590  (
593  );
594 }
595 
596 
597 // ************************************************************************* //
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:339
const char *const group
Group name for atomic constants.
word GName() const
Helper function to return the name of the turbulence G field.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
const dimensionedScalar G0
Conductance quantum: default SI units: [S].
virtual void calculate(const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &epsilon)
Calculate the epsilon and G.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:225
bool empty() const
Return true if the UList is empty (ie, size() is zero)
virtual void write(Ostream &) const
Write.
scalarField G_
Local copy of turbulence G field.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
const word & name() const
Return name.
Definition: fvPatch.H:149
const volVectorField & U() const
Access function to velocity 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)
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
Manipulate matrix.
static scalar tolerance_
Tolerance used in weighted calculations.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:715
label k
Boltzmann constant.
const dimensionedScalar epsilon0
Electric constant: default SI units: [F/m].
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
Abstract base class for turbulence models (RAS, LES and laminar).
Macros for easy insertion into run-time selection tables.
virtual void operator==(const fvPatchField< Type > &)
Definition: fvPatchField.C:562
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
scalar y
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:78
dynamicFvMesh & mesh
const Boundary & boundaryField() const
Return const-reference to the boundary field.
virtual label & master()
Return non-const access to the master patch ID.
static const word propertiesName
Default name of the turbulence properties dictionary.
virtual void calculateTurbulenceFields(const turbulenceModel &turbulence, scalarField &G0, scalarField &epsilon0)
Main driver to calculate the turbulence fields.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
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:474
Foam::fvPatchFieldMapper.
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:370
virtual label size() const
Return size.
Definition: fvPatch.H:161
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:292
bool changing() const
Is mesh changing (topology changing and/or moving)
Definition: polyMesh.H:521
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:71
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
Definition: fvPatchField.C:345
This boundary condition provides a turbulence dissipation wall function condition for high Reynolds n...
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:60
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 tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:217
virtual epsilonWallFunctionFvPatchScalarField & epsilonPatch(const label patchi)
Helper function to return non-const access to an epsilon patch.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:93
static const char nl
Definition: Ostream.H:262
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
label size() const
Return the number of elements in the UList.
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:198
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:179
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void setSize(const label)
Reset size of List.
Definition: List.C:295
label patchi
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
List< List< scalar > > cornerWeights_
List of averaging corner weights.
virtual void updateWeightedCoeffs(const scalarField &weights)
Update the coefficients associated with the patch field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const Mesh & mesh() const
Return mesh.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:312
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...
virtual void checkType()
Check the type of the patch.
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:54
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:91
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
epsilonWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
const nearWallDist & y() const
Return the near wall distances.
bool manipulatedMatrix() const
Return true if the matrix has already been manipulated.
Definition: fvPatchField.H:376
const DimensionedField< Type, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:346
Namespace for OpenFOAM.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
scalarField & epsilon(bool init=false)
Return non-const access to the master&#39;s epsilon field.