omegaWallFunctionFvPatchScalarField.H
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-2018 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 Class
25  Foam::omegaWallFunctionFvPatchScalarField
26 
27 Description
28  This boundary condition provides a wall constraint on turbulnce specific
29  dissipation, omega for both low and high Reynolds number turbulence models.
30 
31  The near-wall omega may be either blended between the viscous region and
32  logarithmic region values using:
33 
34  \f[
35  \omega = sqrt(\omega_{vis}^2 + \omega_{log}^2)
36  \f]
37 
38  where
39 
40  \vartable
41  \omega_{vis} | omega in viscous region
42  \omega_{log} | omega in logarithmic region
43  \endvartable
44 
45  see eq.(15) of:
46  \verbatim
47  Menter, F., Esch, T.
48  "Elements of Industrial Heat Transfer Prediction"
49  16th Brazilian Congress of Mechanical Engineering (COBEM),
50  Nov. 2001
51  \endverbatim
52 
53  or switched between these values based on the laminar-to-turbulent y+ value
54  derived from kappa and E. Recent tests have shown that the standard
55  switching method provides more accurate results for 10 < y+ < 30 when used
56  with high Reynolds number wall-functions and both methods provide accurate
57  results when used with continuous wall-functions. Based on this the
58  standard switching method is used by default.
59 
60 Usage
61  \table
62  Property | Description | Required | Default value
63  Cmu | Model coefficient | no | 0.09
64  kappa | von Karman constant | no | 0.41
65  E | Model coefficient | no | 9.8
66  beta1 | Model coefficient | no | 0.075
67  blended | Blending switch | no | false
68  \endtable
69 
70  Example of the boundary condition specification:
71  \verbatim
72  <patchName>
73  {
74  type omegaWallFunction;
75  }
76  \endverbatim
77 
78 See also
79  Foam::fixedInternalValueFvPatchField
80  Foam::epsilonWallFunctionFvPatchScalarField
81 
82 SourceFiles
83  omegaWallFunctionFvPatchScalarField.C
84 
85 \*---------------------------------------------------------------------------*/
86 
87 #ifndef omegaWallFunctionFvPatchScalarField_H
88 #define omegaWallFunctionFvPatchScalarField_H
89 
90 #include "fixedValueFvPatchField.H"
91 
92 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
93 
94 namespace Foam
95 {
96 
97 class turbulenceModel;
98 
99 /*---------------------------------------------------------------------------*\
100  Class omegaWallFunctionFvPatchScalarField Declaration
101 \*---------------------------------------------------------------------------*/
102 
103 class omegaWallFunctionFvPatchScalarField
104 :
105  public fixedValueFvPatchField<scalar>
106 {
107 protected:
108 
109  // Protected data
110 
111  //- Tolerance used in weighted calculations
112  static scalar tolerance_;
113 
114  //- Cmu coefficient
115  scalar Cmu_;
116 
117  //- Von Karman constant
118  scalar kappa_;
119 
120  //- E coefficient
121  scalar E_;
122 
123  //- beta1 coefficient
124  scalar beta1_;
125 
126  //- Blending switch (defaults to false)
127  Switch blended_;
128 
129  //- y+ at the edge of the laminar sublayer
130  scalar yPlusLam_;
131 
132  //- Local copy of turbulence G field
133  scalarField G_;
134 
135  //- Local copy of turbulence omega field
137 
138  //- Initialised flag
139  bool initialised_;
141  //- Master patch ID
142  label master_;
143 
144  //- List of averaging corner weights
146 
147 
148  // Protected Member Functions
150  //- Check the type of the patch
151  virtual void checkType();
153  //- Write local wall function variables
154  virtual void writeLocalEntries(Ostream&) const;
156  //- Set the master patch - master is responsible for updating all
157  // wall function patches
158  virtual void setMaster();
159 
160  //- Create the averaging weights for cells which are bounded by
161  // multiple wall function faces
162  virtual void createAveragingWeights();
163 
164  //- Helper function to return non-const access to an omega patch
166  (
167  const label patchi
168  );
169 
170  //- Main driver to calculate the turbulence fields
171  virtual void calculateTurbulenceFields
172  (
174  scalarField& G0,
175  scalarField& omega0
176  );
177 
178  //- Calculate the omega and G
179  virtual void calculate
180  (
182  const List<scalar>& cornerWeights,
183  const fvPatch& patch,
184  scalarField& G,
186  );
187 
188  //- Return non-const access to the master patch ID
189  virtual label& master()
190  {
191  return master_;
192  }
193 
195 
196 
197 public:
198 
199  //- Runtime type information
200  TypeName("omegaWallFunction");
201 
202 
203  // Constructors
204 
205  //- Construct from patch and internal field
207  (
208  const fvPatch&,
210  );
211 
212  //- Construct from patch, internal field and dictionary
214  (
215  const fvPatch&,
217  const dictionary&
218  );
219 
220  //- Construct by mapping given
221  // omegaWallFunctionFvPatchScalarField
222  // onto a new patch
224  (
226  const fvPatch&,
228  const fvPatchFieldMapper&
229  );
230 
231  //- Construct as copy
233  (
235  );
236 
237  //- Construct and return a clone
238  virtual tmp<fvPatchScalarField> clone() const
239  {
241  (
243  );
244  }
245 
246  //- Construct as copy setting internal field reference
248  (
251  );
252 
253  //- Construct and return a clone setting internal field reference
255  (
257  ) const
258  {
260  (
262  );
263  }
264 
265 
266  // Member functions
267 
268  // Access
269 
270  //- Return non-const access to the master's G field
271  scalarField& G(bool init = false);
272 
273  //- Return non-const access to the master's omega field
274  scalarField& omega(bool init = false);
276 
277  // Evaluation functions
278 
279  //- Update the coefficients associated with the patch field
280  virtual void updateCoeffs();
281 
282  //- Update the coefficients associated with the patch field
283  virtual void updateWeightedCoeffs(const scalarField& weights);
284 
285  //- Manipulate matrix
286  virtual void manipulateMatrix(fvMatrix<scalar>& matrix);
287 
288  //- Manipulate matrix with given weights
289  virtual void manipulateMatrix
290  (
291  fvMatrix<scalar>& matrix,
292  const scalarField& weights
293  );
294 
295 
296  // I-O
297 
298  //- Write
299  virtual void write(Ostream&) const;
300 };
301 
302 
303 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
304 
305 } // End namespace Foam
306 
307 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
308 
309 #endif
310 
311 // ************************************************************************* //
Switch blended_
Blending switch (defaults to false)
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
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
const dimensionedScalar G0
Conductance quantum: default SI units: [S].
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
omegaWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:60
virtual tmp< fvPatchScalarField > clone() const
Construct and return a clone.
virtual void updateWeightedCoeffs(const scalarField &weights)
Update the coefficients associated with the patch field.
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
virtual void createAveragingWeights()
Create the averaging weights for cells which are bounded by.
virtual void checkType()
Check the type of the patch.
Abstract base class for turbulence models (RAS, LES and laminar).
scalarField & G(bool init=false)
Return non-const access to the master&#39;s G field.
virtual omegaWallFunctionFvPatchScalarField & omegaPatch(const label patchi)
Helper function to return non-const access to an omega patch.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Foam::fvPatchFieldMapper.
static scalar tolerance_
Tolerance used in weighted calculations.
virtual void calculateTurbulenceFields(const turbulenceModel &turbulence, scalarField &G0, scalarField &omega0)
Main driver to calculate the turbulence fields.
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
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:340
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
Info<< "Reading field U\"<< 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\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
List< List< scalar > > cornerWeights_
List of averaging corner weights.
scalarField G_
Local copy of turbulence G field.
scalarField & omega(bool init=false)
Return non-const access to the master&#39;s omega field.
label patchi
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
This boundary condition provides a wall constraint on turbulnce specific dissipation, omega for both low and high Reynolds number turbulence models.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual void setMaster()
Set the master patch - master is responsible for updating all.
A class for managing temporary objects.
Definition: PtrList.H:53
TypeName("omegaWallFunction")
Runtime type information.
virtual label & master()
Return non-const access to the master patch ID.
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
Manipulate matrix.
scalar yPlusLam_
y+ at the edge of the laminar sublayer
scalarField omega_
Local copy of turbulence omega field.
Namespace for OpenFOAM.
virtual void calculate(const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &omega)
Calculate the omega and G.