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-2019 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 specified in the corresponding nutWallFunction.
55  Recent tests have shown that the standard switching method provides more
56  accurate results for 10 < y+ < 30 when used with high Reynolds number
57  wall-functions and both methods provide accurate results when used with
58  continuous wall-functions. Based on this the standard switching method is
59  used by default.
60 
61 Usage
62  \table
63  Property | Description | Required | Default value
64  beta1 | Model coefficient | no | 0.075
65  blended | Blending switch | no | false
66  \endtable
67 
68  Example of the boundary condition specification:
69  \verbatim
70  <patchName>
71  {
72  type omegaWallFunction;
73  }
74  \endverbatim
75 
76 See also
77  Foam::fixedInternalValueFvPatchField
78  Foam::epsilonWallFunctionFvPatchScalarField
79 
80 SourceFiles
81  omegaWallFunctionFvPatchScalarField.C
82 
83 \*---------------------------------------------------------------------------*/
84 
85 #ifndef omegaWallFunctionFvPatchScalarField_H
86 #define omegaWallFunctionFvPatchScalarField_H
87 
88 #include "fixedValueFvPatchField.H"
89 
90 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
91 
92 namespace Foam
93 {
94 
95 class turbulenceModel;
96 
97 /*---------------------------------------------------------------------------*\
98  Class omegaWallFunctionFvPatchScalarField Declaration
99 \*---------------------------------------------------------------------------*/
100 
101 class omegaWallFunctionFvPatchScalarField
102 :
103  public fixedValueFvPatchField<scalar>
104 {
105 protected:
106 
107  // Protected data
108 
109  //- Tolerance used in weighted calculations
110  static scalar tolerance_;
111 
112  //- beta1 coefficient
113  scalar beta1_;
114 
115  //- Blending switch (defaults to false)
116  Switch blended_;
117 
118  //- Local copy of turbulence G field
119  scalarField G_;
120 
121  //- Local copy of turbulence omega field
124  //- Initialised flag
125  bool initialised_;
126 
127  //- Master patch ID
128  label master_;
129 
130  //- List of averaging corner weights
133 
134  // Protected Member Functions
136  //- Set the master patch - master is responsible for updating all
137  // wall function patches
138  virtual void setMaster();
139 
140  //- Create the averaging weights for cells which are bounded by
141  // multiple wall function faces
142  virtual void createAveragingWeights();
143 
144  //- Helper function to return non-const access to an omega patch
146  (
147  const label patchi
148  );
149 
150  //- Main driver to calculate the turbulence fields
151  virtual void calculateTurbulenceFields
152  (
154  scalarField& G0,
155  scalarField& omega0
156  );
157 
158  //- Calculate the omega and G
159  virtual void calculate
160  (
162  const List<scalar>& cornerWeights,
163  const fvPatch& patch,
164  scalarField& G,
166  );
167 
168  //- Return non-const access to the master patch ID
169  virtual label& master()
170  {
171  return master_;
172  }
173 
175 
176 
177 public:
178 
179  //- Runtime type information
180  TypeName("omegaWallFunction");
181 
182 
183  // Constructors
184 
185  //- Construct from patch and internal field
187  (
188  const fvPatch&,
190  );
192  //- Construct from patch, internal field and dictionary
194  (
195  const fvPatch&,
197  const dictionary&
198  );
199 
200  //- Construct by mapping given
201  // omegaWallFunctionFvPatchScalarField
202  // onto a new patch
204  (
206  const fvPatch&,
208  const fvPatchFieldMapper&
209  );
210 
211  //- Copy constructor
213  (
215  );
216 
217  //- Construct and return a clone
218  virtual tmp<fvPatchScalarField> clone() const
219  {
221  (
223  );
224  }
225 
226  //- Copy constructor setting internal field reference
228  (
231  );
232 
233  //- Construct and return a clone setting internal field reference
235  (
237  ) const
238  {
240  (
242  );
243  }
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 omega field
254  scalarField& omega(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 // ************************************************************************* //
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:158
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:59
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.
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(dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar(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
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:325
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
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
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.
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.