Loading...
Searching...
No Matches
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-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
27Class
28 Foam::epsilonWallFunctionFvPatchScalarField
29
30Group
31 grpWallFunctions
32
33Description
34 This boundary condition provides wall functions for the turbulent kinetic
35 energy dissipation rate (i.e. \c epsilon) and the turbulent kinetic
36 energy production contribution (i.e. \c G) for low- and high-Reynolds
37 number simulations.
38
39Usage
40 Example of the boundary condition specification:
41 \verbatim
42 <patchName>
43 {
44 // Mandatory entries
45 type epsilonWallFunction;
46
47 // Optional entries
48 lowReCorrection false;
49
50 // Inherited entries
51 ...
52 }
53 \endverbatim
54
55 where the entries mean:
56 \table
57 Property | Description | Type | Reqd | Deflt
58 type | Type name: epsilonWallFunction | word | yes | -
59 lowReCorrection | Flag: apply low-Re correction | bool | no | false
60 \endtable
61
62 The inherited entries are elaborated in:
63 - \link wallFunctionCoefficients.H \endlink
64 - \link wallFunctionBlenders.H \endlink
65
66Note
67 - \c lowReCorrection operates with only \c stepwise blending treatment to
68 ensure the backward compatibility.
69 - If \c lowReCorrection is \c on, \c stepwise blending treatment is fully
70 active.
71 - If \c lowReCorrection is \c off, only the inertial sublayer prediction
72 is used in the wall function, hence high-Re mode operation.
73
74SourceFiles
75 epsilonWallFunctionFvPatchScalarField.C
76
77\*---------------------------------------------------------------------------*/
78
79#ifndef epsilonWallFunctionFvPatchScalarField_H
80#define epsilonWallFunctionFvPatchScalarField_H
81
85
86// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
87
88namespace Foam
89{
90
91class turbulenceModel;
92
93/*---------------------------------------------------------------------------*\
94 Class epsilonWallFunctionFvPatchScalarField Declaration
95\*---------------------------------------------------------------------------*/
96
98:
99 public fixedValueFvPatchField<scalar>,
101{
102protected:
103
104 // Protected Data
105
106 //- Tolerance used in weighted calculations
107 static scalar tolerance_;
108
109 //- Apply low-Re correction term (default = no)
111
112 //- Initialised flag
113 bool initialised_;
114
115 //- Master patch ID
116 label master_;
117
118 //- Wall-function coefficients
120
121 //- Local copy of turbulence G field
123
124 //- Local copy of turbulence epsilon field
126
127 //- List of averaging corner weights
129
130
131 // Protected Member Functions
133 //- Set the master patch - master is responsible for updating all
134 //- wall function patches
135 virtual void setMaster();
136
137 //- Create the averaging weights for cells which are bounded by
138 //- multiple wall function faces
139 virtual void createAveragingWeights();
140
141 //- Helper function to return non-const access to an epsilon patch
143 (
144 const label patchi
145 );
146
147 //- Main driver to calculate the turbulence fields
148 virtual void calculateTurbulenceFields
149 (
151 scalarField& G0,
152 scalarField& epsilon0
153 );
154
155 //- Calculate the epsilon and G
156 virtual void calculate
159 const List<scalar>& cornerWeights,
160 const fvPatch& patch,
161 scalarField& G,
163 );
164
165 //- Return non-const access to the master patch ID
166 virtual label& master()
167 {
168 return master_;
169 }
170
171 //- Write local wall function variables
172 void writeLocalEntries(Ostream&) const;
173
174
175public:
176
177 //- Runtime type information
178 TypeName("epsilonWallFunction");
179
180
181 // Constructors
182
183 //- Construct from patch and internal field
185 (
186 const fvPatch&,
188 );
189
190 //- Construct from patch, internal field and dictionary
192 (
193 const fvPatch&,
195 const dictionary&
196 );
197
198 //- Construct by mapping given
199 //- epsilonWallFunctionFvPatchScalarField
200 //- onto a new patch
202 (
204 const fvPatch&,
206 const fvPatchFieldMapper&
207 );
208
209 //- Construct as copy
211 (
213 );
214
215 //- Construct as copy setting internal field reference
217 (
220 );
221
222 //- Return a clone
224 {
225 return fvPatchField<scalar>::Clone(*this);
226 }
227
228 //- Clone with an internal field reference
230 (
232 ) const
233 {
234 return fvPatchField<scalar>::Clone(*this, iF);
235 }
236
237
238 //- Destructor
239 virtual ~epsilonWallFunctionFvPatchScalarField() = default;
240
241
242 // Member Functions
243
244 // Access
245
246 //- Return non-const access to the master's G field
247 scalarField& G(bool init = false);
248
249 //- Return non-const access to the master's epsilon field
250 scalarField& epsilon(bool init = false);
251
252
253 // Evaluation
254
255 //- Update the coefficients associated with the patch field
256 virtual void updateCoeffs();
257
258 //- Update the coefficients associated with the patch field
259 virtual void updateWeightedCoeffs(const scalarField& weights);
260
261 //- Manipulate matrix
262 virtual void manipulateMatrix(fvMatrix<scalar>& matrix);
263
264 //- Manipulate matrix with given weights
265 virtual void manipulateMatrix
266 (
267 fvMatrix<scalar>& matrix,
268 const scalarField& weights
269 );
270
271
272 // I-O
273
274 //- Write
275 virtual void write(Ostream&) const;
276};
277
278
279// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
281} // End namespace Foam
282
283// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
284
285#endif
286
287// ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
This boundary condition provides wall functions for the turbulent kinetic energy dissipation rate (i....
List< List< scalar > > cornerWeights_
List of averaging corner weights.
epsilonWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual tmp< fvPatchField< scalar > > clone(const DimensionedField< scalar, volMesh > &iF) const
Clone with an internal field reference.
virtual tmp< fvPatchField< scalar > > clone() const
Return a clone.
scalarField epsilon_
Local copy of turbulence epsilon field.
void writeLocalEntries(Ostream &) const
Write local wall function variables.
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
Manipulate matrix.
virtual void manipulateMatrix(fvMatrix< scalar > &matrix, const scalarField &weights)
Manipulate matrix with given weights.
virtual ~epsilonWallFunctionFvPatchScalarField()=default
Destructor.
static scalar tolerance_
Tolerance used in weighted calculations.
const bool lowReCorrection_
Apply low-Re correction term (default = no).
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void createAveragingWeights()
Create the averaging weights for cells which are bounded by multiple wall function faces.
wallFunctionCoefficients wallCoeffs_
Wall-function coefficients.
virtual void calculate(const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &epsilon)
Calculate the epsilon and G.
virtual label & master()
Return non-const access to the master patch ID.
TypeName("epsilonWallFunction")
Runtime type information.
scalarField & G(bool init=false)
Return non-const access to the master's G field.
virtual void updateWeightedCoeffs(const scalarField &weights)
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 void calculateTurbulenceFields(const turbulenceModel &turbulence, scalarField &G0, scalarField &epsilon0)
Main driver to calculate the turbulence fields.
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 wall function patches.
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
fixedValueFvPatchField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
const fvPatch & patch() const noexcept
Return the patch.
A FieldMapper for finite-volume patch fields.
static tmp< fvPatchField< Type > > Clone(const DerivedPatchField &pf, Args &&... args)
Clone a patch field, optionally with internal field reference etc.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
A class for managing temporary objects.
Definition tmp.H:75
Abstract base class for turbulence models (RAS, LES and laminar).
The class wallFunctionBlenders is a base class that hosts common entries for various derived wall-fun...
wallFunctionBlenders()
Default construct with default coefficients.
Class to host the wall-function coefficients being used in the wall function boundary conditions.
scalar epsilon
Namespace for OpenFOAM.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
runTime write()
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68