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 | Copyright (C) 2011-2016 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 Group
28  grpWallFunctions
29 
30 Description
31  This boundary condition provides a wall function constraint on turbulnce
32  specific dissipation, omega. The values are computed 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  Model described by 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 Usage
54  \table
55  Property | Description | Required | Default value
56  Cmu | model coefficient | no | 0.09
57  kappa | Von Karman constant | no | 0.41
58  E | model coefficient | no | 9.8
59  beta1 | model coefficient | no | 0.075
60  \endtable
61 
62  Example of the boundary condition specification:
63  \verbatim
64  <patchName>
65  {
66  type omegaWallFunction;
67  }
68  \endverbatim
69 
70 SourceFiles
71  omegaWallFunctionFvPatchScalarField.C
72 
73 \*---------------------------------------------------------------------------*/
74 
75 #ifndef omegaWallFunctionFvPatchScalarField_H
76 #define omegaWallFunctionFvPatchScalarField_H
77 
78 #include "fixedValueFvPatchField.H"
79 
80 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
81 
82 namespace Foam
83 {
84 
85 class turbulenceModel;
86 
87 /*---------------------------------------------------------------------------*\
88  Class omegaWallFunctionFvPatchScalarField Declaration
89 \*---------------------------------------------------------------------------*/
90 
91 class omegaWallFunctionFvPatchScalarField
92 :
93  public fixedValueFvPatchField<scalar>
94 {
95 protected:
96 
97  // Protected data
98 
99  //- Tolerance used in weighted calculations
100  static scalar tolerance_;
101 
102  //- Cmu coefficient
103  scalar Cmu_;
104 
105  //- Von Karman constant
106  scalar kappa_;
107 
108  //- E coefficient
109  scalar E_;
110 
111  //- beta1 coefficient
112  scalar beta1_;
113 
114  //- Y+ at the edge of the laminar sublayer
115  scalar yPlusLam_;
116 
117  //- Local copy of turbulence G field
118  scalarField G_;
119 
120  //- Local copy of turbulence omega field
122 
123  //- Initialised flag
124  bool initialised_;
125 
126  //- Master patch ID
127  label master_;
128 
129  //- List of averaging corner weights
131 
133  // Protected Member Functions
134 
135  //- Check the type of the patch
136  virtual void checkType();
137 
138  //- Write local wall function variables
139  virtual void writeLocalEntries(Ostream&) const;
140 
141  //- Set the master patch - master is responsible for updating all
142  // wall function patches
143  virtual void setMaster();
145  //- Create the averaging weights for cells which are bounded by
146  // multiple wall function faces
147  virtual void createAveragingWeights();
148 
149  //- Helper function to return non-const access to an omega patch
151  (
152  const label patchi
153  );
154 
155  //- Main driver to calculate the turbulence fields
157  (
160  scalarField& omega0
161  );
163  //- Calculate the omega and G
164  virtual void calculate
165  (
167  const List<scalar>& cornerWeights,
168  const fvPatch& patch,
169  scalarField& G,
171  );
172 
173  //- Return non-const access to the master patch ID
174  virtual label& master()
175  {
176  return master_;
177  }
178 
179 
180 public:
181 
182  //- Runtime type information
183  TypeName("omegaWallFunction");
184 
185 
186  // Constructors
187 
188  //- Construct from patch and internal field
190  (
191  const fvPatch&,
193  );
194 
195  //- Construct from patch, internal field and dictionary
197  (
198  const fvPatch&,
200  const dictionary&
201  );
202 
203  //- Construct by mapping given
204  // omegaWallFunctionFvPatchScalarField
205  // onto a new patch
207  (
209  const fvPatch&,
211  const fvPatchFieldMapper&
212  );
213 
214  //- Construct as copy
216  (
218  );
219 
220  //- Construct and return a clone
221  virtual tmp<fvPatchScalarField> clone() const
222  {
224  (
226  );
227  }
228 
229  //- Construct as copy setting internal field reference
231  (
234  );
235 
236  //- Construct and return a clone setting internal field reference
238  (
240  ) const
241  {
243  (
245  );
246  }
247 
248 
249  // Member functions
250 
251  // Access
252 
253  //- Return non-const access to the master's G field
254  scalarField& G(bool init = false);
255 
256  //- Return non-const access to the master's omega field
257  scalarField& omega(bool init = false);
258 
259 
260  // Evaluation functions
261 
262  //- Update the coefficients associated with the patch field
263  virtual void updateCoeffs();
264 
265  //- Update the coefficients associated with the patch field
266  virtual void updateWeightedCoeffs(const scalarField& weights);
267 
268  //- Manipulate matrix
269  virtual void manipulateMatrix(fvMatrix<scalar>& matrix);
271  //- Manipulate matrix with given weights
272  virtual void manipulateMatrix
273  (
274  fvMatrix<scalar>& matrix,
275  const scalarField& weights
276  );
277 
278 
279  // I-O
280 
281  //- Write
282  virtual void write(Ostream&) const;
283 };
284 
285 
286 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
287 
288 } // End namespace Foam
289 
290 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
291 
292 #endif
293 
294 // ************************************************************************* //
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:339
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:59
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 tmp< fvPatchScalarField > clone() const
Construct and return a clone.
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.
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
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:71
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 function constraint on turbulnce specific dissipation...
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:54
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.