epsilonWallFunctionFvPatchScalarField.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016, 2019 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::epsilonWallFunctionFvPatchScalarField
29 
30 Group
31  grpWallFunctions
32 
33 Description
34  This boundary condition provides a wall constraint on turbulent kinetic
35  energy dissipation rate, i.e. \c epsilon, for low- and high-Reynolds number
36  turbulence models.
37 
38  The condition can be applied to wall boundaries for which it
39  - calculates \c epsilon and \c G
40  - specifies the near-wall \c epsilon value
41 
42  where
43  \vartable
44  epsilon | turbulent kinetic energy dissipation rate field
45  G | turbulent kinetic energy production field (divergence-free)
46  \endvartable
47 
48  The low-Re correction is activated by setting the entry
49  \c lowReCorrection to 'on'; in this mode the model switches between
50  viscous and turbulent functions based on the viscous-to-turbulent
51  \c y+ value derived from the \c kappa and \c E.
52 
53  When the \c lowReCorrection is inactive, the wall function operates
54  in high-Re mode.
55 
56 Usage
57  \table
58  Property | Description | Required | Default value
59  lowReCorrection | Low-Re correction active | no | off
60  \endtable
61 
62  Example of the boundary condition specification:
63  \verbatim
64  <patchName>
65  {
66  // Mandatory entries
67  type epsilonWallFunction;
68 
69  // Optional entries
70  }
71  \endverbatim
72 
73 Note
74  The coefficients \c Cmu, \c kappa, and \c E are obtained from
75  the specified \c nutWallFunction in order to ensure that each patch
76  possesses the same set of values for these coefficients.
77 
78 See also
79  Foam::fixedInternalValueFvPatchField
80 
81 SourceFiles
82  epsilonWallFunctionFvPatchScalarField.C
83 
84 \*---------------------------------------------------------------------------*/
85 
86 #ifndef epsilonWallFunctionFvPatchScalarField_H
87 #define epsilonWallFunctionFvPatchScalarField_H
88 
89 #include "fixedValueFvPatchField.H"
90 
91 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
92 
93 namespace Foam
94 {
95 
96 class turbulenceModel;
97 
98 /*---------------------------------------------------------------------------*\
99  Class epsilonWallFunctionFvPatchScalarField Declaration
100 \*---------------------------------------------------------------------------*/
101 
102 class epsilonWallFunctionFvPatchScalarField
103 :
104  public fixedValueFvPatchField<scalar>
105 {
106 protected:
107 
108  // Protected Data
109 
110  //- Tolerance used in weighted calculations
111  static scalar tolerance_;
112 
113  //- Apply low-Re correction term; default = no
114  bool lowReCorrection_;
115 
116  //- Initialised flag
117  bool initialised_;
118 
119  //- Master patch ID
120  label master_;
121 
122  //- Local copy of turbulence G field
123  scalarField G_;
124 
125  //- Local copy of turbulence epsilon field
127 
128  //- List of averaging corner weights
130 
131 
132  // Protected Member Functions
133 
134  //- Set the master patch - master is responsible for updating all
135  //- wall function patches
136  virtual void setMaster();
137 
138  //- Create the averaging weights for cells which are bounded by
139  //- multiple wall function faces
140  virtual void createAveragingWeights();
141 
142  //- Helper function to return non-const access to an epsilon patch
144  (
145  const label patchi
146  );
147 
148  //- Main driver to calculate the turbulence fields
149  virtual void calculateTurbulenceFields
150  (
152  scalarField& G0,
154  );
155 
156  //- Calculate the epsilon and G
157  virtual void calculate
158  (
160  const List<scalar>& cornerWeights,
161  const fvPatch& patch,
162  scalarField& G,
164  );
165 
166  //- Return non-const access to the master patch ID
167  virtual label& master()
168  {
169  return master_;
170  }
171 
172 
173 public:
174 
175  //- Runtime type information
176  TypeName("epsilonWallFunction");
177 
178 
179  // Constructors
180 
181  //- Construct from patch and internal field
183  (
184  const fvPatch&,
186  );
187 
188  //- Construct from patch, internal field and dictionary
190  (
191  const fvPatch&,
193  const dictionary&
194  );
195 
196  //- Construct by mapping given
197  //- epsilonWallFunctionFvPatchScalarField
198  //- onto a new patch
200  (
202  const fvPatch&,
204  const fvPatchFieldMapper&
205  );
206 
207  //- Construct as copy
209  (
211  );
212 
213  //- Construct and return a clone
214  virtual tmp<fvPatchScalarField> clone() const
215  {
217  (
219  );
220  }
221 
222  //- Construct as copy setting internal field reference
224  (
227  );
228 
229  //- Construct and return a clone setting internal field reference
231  (
233  ) const
234  {
236  (
238  );
239  }
240 
241 
242  //- Destructor
243  virtual ~epsilonWallFunctionFvPatchScalarField() = default;
244 
245 
246  // Member Functions
247 
248  // Access
249 
250  //- Return non-const access to the master's G field
251  scalarField& G(bool init = false);
252 
253  //- Return non-const access to the master's epsilon field
254  scalarField& epsilon(bool init = false);
255 
256 
257  // Evaluation functions
258 
259  //- Update the coefficients associated with the patch field
260  virtual void updateCoeffs();
261 
262  //- Update the coefficients associated with the patch field
263  virtual void updateWeightedCoeffs(const scalarField& weights);
264 
265  //- Manipulate matrix
266  virtual void manipulateMatrix(fvMatrix<scalar>& matrix);
267 
268  //- Manipulate matrix with given weights
269  virtual void manipulateMatrix
270  (
271  fvMatrix<scalar>& matrix,
272  const scalarField& weights
273  );
274 
275 
276  // I-O
277 
278  //- Write
279  virtual void write(Ostream&) const;
280 };
281 
282 
283 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
284 
285 } // End namespace Foam
286 
287 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
288 
289 #endif
290 
291 // ************************************************************************* //
Foam::epsilonWallFunctionFvPatchScalarField::initialised_
bool initialised_
Initialised flag.
Definition: epsilonWallFunctionFvPatchScalarField.H:134
Foam::epsilonWallFunctionFvPatchScalarField::master
virtual label & master()
Return non-const access to the master patch ID.
Definition: epsilonWallFunctionFvPatchScalarField.H:184
turbulence
Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\n"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::epsilonWallFunctionFvPatchScalarField::updateWeightedCoeffs
virtual void updateWeightedCoeffs(const scalarField &weights)
Update the coefficients associated with the patch field.
Definition: epsilonWallFunctionFvPatchScalarField.C:420
Foam::epsilonWallFunctionFvPatchScalarField::tolerance_
static scalar tolerance_
Tolerance used in weighted calculations.
Definition: epsilonWallFunctionFvPatchScalarField.H:128
Foam::epsilonWallFunctionFvPatchScalarField::createAveragingWeights
virtual void createAveragingWeights()
Definition: epsilonWallFunctionFvPatchScalarField.C:72
Foam::epsilonWallFunctionFvPatchScalarField::epsilonWallFunctionFvPatchScalarField
epsilonWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: epsilonWallFunctionFvPatchScalarField.C:246
Foam::fixedValueFvPatchField
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
Definition: fixedValueFvPatchField.H:80
Foam::constant::electromagnetic::epsilon0
const dimensionedScalar epsilon0
Electric constant: default SI units: [F/m].
Foam::epsilonWallFunctionFvPatchScalarField::epsilon
scalarField & epsilon(bool init=false)
Return non-const access to the master's epsilon field.
Definition: epsilonWallFunctionFvPatchScalarField.C:356
Foam::label
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:62
Foam::Field< scalar >
Foam::compressible::turbulenceModel
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Definition: turbulentFluidThermoModel.H:63
Foam::constant::electromagnetic::G0
const dimensionedScalar G0
Conductance quantum: default SI units: [S].
Foam::epsilonWallFunctionFvPatchScalarField::lowReCorrection_
bool lowReCorrection_
Apply low-Re correction term; default = no.
Definition: epsilonWallFunctionFvPatchScalarField.H:131
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Foam::epsilonWallFunctionFvPatchScalarField
This boundary condition provides a wall constraint on turbulent kinetic energy dissipation rate,...
Definition: epsilonWallFunctionFvPatchScalarField.H:119
Foam::epsilonWallFunctionFvPatchScalarField::calculate
virtual void calculate(const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &epsilon)
Calculate the epsilon and G.
Definition: epsilonWallFunctionFvPatchScalarField.C:183
Foam::epsilonWallFunctionFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: epsilonWallFunctionFvPatchScalarField.C:536
Foam::epsilonWallFunctionFvPatchScalarField::G
scalarField & G(bool init=false)
Return non-const access to the master's G field.
Definition: epsilonWallFunctionFvPatchScalarField.C:337
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:63
Foam::epsilonWallFunctionFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: epsilonWallFunctionFvPatchScalarField.C:374
Foam::epsilonWallFunctionFvPatchScalarField::G_
scalarField G_
Local copy of turbulence G field.
Definition: epsilonWallFunctionFvPatchScalarField.H:140
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::epsilonWallFunctionFvPatchScalarField::manipulateMatrix
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
Manipulate matrix.
Definition: epsilonWallFunctionFvPatchScalarField.C:477
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::epsilonWallFunctionFvPatchScalarField::calculateTurbulenceFields
virtual void calculateTurbulenceFields(const turbulenceModel &turbulence, scalarField &G0, scalarField &epsilon0)
Main driver to calculate the turbulence fields.
Definition: epsilonWallFunctionFvPatchScalarField.C:150
Foam::epsilonWallFunctionFvPatchScalarField::epsilonPatch
virtual epsilonWallFunctionFvPatchScalarField & epsilonPatch(const label patchi)
Helper function to return non-const access to an epsilon patch.
Definition: epsilonWallFunctionFvPatchScalarField.C:133
Foam::epsilonWallFunctionFvPatchScalarField::epsilon_
scalarField epsilon_
Local copy of turbulence epsilon field.
Definition: epsilonWallFunctionFvPatchScalarField.H:143
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:102
Foam::epsilonWallFunctionFvPatchScalarField::TypeName
TypeName("epsilonWallFunction")
Runtime type information.
Foam::epsilonWallFunctionFvPatchScalarField::master_
label master_
Master patch ID.
Definition: epsilonWallFunctionFvPatchScalarField.H:137
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:76
Foam::epsilonWallFunctionFvPatchScalarField::setMaster
virtual void setMaster()
Definition: epsilonWallFunctionFvPatchScalarField.C:42
fixedValueFvPatchField.H
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::epsilonWallFunctionFvPatchScalarField::clone
virtual tmp< fvPatchScalarField > clone() const
Construct and return a clone.
Definition: epsilonWallFunctionFvPatchScalarField.H:231
Foam::fvPatchField< scalar >::patch
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:343
Foam::epsilonWallFunctionFvPatchScalarField::~epsilonWallFunctionFvPatchScalarField
virtual ~epsilonWallFunctionFvPatchScalarField()=default
Destructor.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::epsilonWallFunctionFvPatchScalarField::cornerWeights_
List< List< scalar > > cornerWeights_
List of averaging corner weights.
Definition: epsilonWallFunctionFvPatchScalarField.H:146
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54