omegaWallFunctionFvPatchScalarField.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::omegaWallFunctionFvPatchScalarField
29
30Group
31 grpWallFunctions
32
33Description
34 This boundary condition provides a wall function for the specific
35 dissipation rate (i.e. \c omega) and the turbulent kinetic energy
36 production contribution (i.e. \c G) for low- and high-Reynolds number
37 applications.
38
39Usage
40 Example of the boundary condition specification:
41 \verbatim
42 <patchName>
43 {
44 // Mandatory entries
45 type omegaWallFunction;
46
47 // Optional entries
48 beta1 0.075;
49
50 // Inherited entries
51 ...
52 }
53 \endverbatim
54
55 \table
56 Property | Description | Type | Reqd | Deflt
57 type | Type name: omegaWallFunction | word | yes | -
58 beta1 | Model coefficient | scalar | no | 0.075
59 \endtable
60
61 The inherited entries are elaborated in:
62 - \link wallFunctionCoefficients.H \endlink
63 - \link wallFunctionBlenders.H \endlink
64
65SourceFiles
66 omegaWallFunctionFvPatchScalarField.C
67
68\*---------------------------------------------------------------------------*/
69
70#ifndef omegaWallFunctionFvPatchScalarField_H
71#define omegaWallFunctionFvPatchScalarField_H
72
76
77// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
78
79namespace Foam
80{
81
82class turbulenceModel;
83
84/*---------------------------------------------------------------------------*\
85 Class omegaWallFunctionFvPatchScalarField Declaration
86\*---------------------------------------------------------------------------*/
87
88class omegaWallFunctionFvPatchScalarField
89:
90 public fixedValueFvPatchField<scalar>,
91 private wallFunctionBlenders
92{
93protected:
94
95 // Protected Data
96
97 //- Tolerance used in weighted calculations
98 static scalar tolerance_;
99
100 //- Initialised flag
101 bool initialised_;
102
103 //- Master patch ID
104 label master_;
106 //- beta1 coefficient
107 scalar beta1_;
108
109 //- Wall-function coefficients
111
112 //- Local copy of turbulence G field
114
115 //- Local copy of turbulence omega field
117
118 //- List of averaging corner weights
120
122 // Protected Member Functions
123
124 //- Set the master patch - master is responsible for updating all
125 //- wall function patches
126 virtual void setMaster();
128 //- Create the averaging weights for cells which are bounded by
129 //- multiple wall function faces
131
132 //- Helper function to return non-const access to an omega patch
134 (
135 const label patchi
136 );
137
138 //- Main driver to calculate the turbulence fields
139 virtual void calculateTurbulenceFields
140 (
141 const turbulenceModel& turbulence,
142 scalarField& G0,
143 scalarField& omega0
144 );
145
146 //- Calculate the omega and G
147 virtual void calculate
148 (
149 const turbulenceModel& turbulence,
150 const List<scalar>& cornerWeights,
151 const fvPatch& patch,
152 scalarField& G,
154 );
155
156 //- Return non-const access to the master patch ID
157 virtual label& master()
158 {
159 return master_;
160 }
161
162 //- Write local wall function variables
163 void writeLocalEntries(Ostream&) const;
164
165
166public:
167
168 //- Runtime type information
169 TypeName("omegaWallFunction");
170
171
172 // Constructors
173
174 //- Construct from patch and internal field
176 (
177 const fvPatch&,
179 );
180
181 //- Construct from patch, internal field and dictionary
183 (
184 const fvPatch&,
187 );
188
189 //- Construct by mapping given
190 //- omegaWallFunctionFvPatchScalarField
191 //- onto a new patch
193 (
195 const fvPatch&,
197 const fvPatchFieldMapper&
198 );
199
200 //- Construct as copy
202 (
204 );
205
206 //- Construct and return a clone
207 virtual tmp<fvPatchScalarField> clone() const
208 {
210 (
212 );
213 }
214
215 //- Construct as copy setting internal field reference
217 (
220 );
221
222 //- Construct and return a clone setting internal field reference
226 ) const
227 {
229 (
231 );
232 }
233
234 //- Destructor
235 virtual ~omegaWallFunctionFvPatchScalarField() = default;
236
237
238 // Member Functions
239
240 // Access
241
242 //- Return non-const access to the master's G field
243 scalarField& G(bool init = false);
244
245 //- Return non-const access to the master's omega field
246 scalarField& omega(bool init = false);
247
248
249 // Evaluation
250
251 //- Update the coefficients associated with the patch field
252 virtual void updateCoeffs();
253
254 //- Update the coefficients associated with the patch field
255 virtual void updateWeightedCoeffs(const scalarField& weights);
256
257 //- Manipulate matrix
258 virtual void manipulateMatrix(fvMatrix<scalar>& matrix);
259
260 //- Manipulate matrix with given weights
261 virtual void manipulateMatrix
262 (
263 fvMatrix<scalar>& matrix,
264 const scalarField& weights
265 );
266
267
268 // I-O
269
270 //- Write
271 virtual void write(Ostream&) const;
272};
273
274
275// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
276
277} // End namespace Foam
278
279// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
280
281#endif
282
283// ************************************************************************* //
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:77
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
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:121
A FieldMapper for finite-volume patch fields.
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:362
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
This boundary condition provides a wall function for the specific dissipation rate (i....
virtual tmp< fvPatchScalarField > clone(const DimensionedField< scalar, volMesh > &iF) const
Construct and return a clone setting internal field reference.
scalarField G_
Local copy of turbulence G field.
virtual void calculate(const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &omega)
Calculate the omega and G.
List< List< scalar > > cornerWeights_
List of averaging corner weights.
void writeLocalEntries(Ostream &) const
Write local wall function variables.
TypeName("omegaWallFunction")
Runtime type information.
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
Manipulate matrix.
static scalar tolerance_
Tolerance used in weighted calculations.
virtual tmp< fvPatchScalarField > clone() const
Construct and return a clone.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual ~omegaWallFunctionFvPatchScalarField()=default
Destructor.
wallFunctionCoefficients wallCoeffs_
Wall-function coefficients.
virtual label & master()
Return non-const access to the master patch ID.
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.
scalarField & omega(bool init=false)
Return non-const access to the master's omega field.
virtual void calculateTurbulenceFields(const turbulenceModel &turbulence, scalarField &G0, scalarField &omega0)
Main driver to calculate the turbulence fields.
scalarField omega_
Local copy of turbulence omega field.
virtual omegaWallFunctionFvPatchScalarField & omegaPatch(const label patchi)
Helper function to return non-const access to an omega patch.
A class for managing temporary objects.
Definition: tmp.H:65
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...
Class to host the wall-function coefficients being used in the wall function boundary conditions.
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Namespace for OpenFOAM.
runTime write()
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73