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-2015 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  \heading Patch usage
54 
55  \table
56  Property | Description | Required | Default value
57  Cmu | model coefficient | no | 0.09
58  kappa | Von Karman constant | no | 0.41
59  E | model coefficient | no | 9.8
60  beta1 | model coefficient | no | 0.075
61  \endtable
62 
63  Example of the boundary condition specification:
64  \verbatim
65  myPatch
66  {
67  type omegaWallFunction;
68  }
69  \endverbatim
70 
71 SourceFiles
72  omegaWallFunctionFvPatchScalarField.C
73 
74 \*---------------------------------------------------------------------------*/
75 
76 #ifndef omegaWallFunctionFvPatchScalarField_H
77 #define omegaWallFunctionFvPatchScalarField_H
78 
79 #include "fixedValueFvPatchField.H"
80 
81 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
82 
83 namespace Foam
84 {
85 
86 class turbulenceModel;
87 
88 /*---------------------------------------------------------------------------*\
89  Class omegaWallFunctionFvPatchScalarField Declaration
90 \*---------------------------------------------------------------------------*/
91 
92 class omegaWallFunctionFvPatchScalarField
93 :
94  public fixedValueFvPatchField<scalar>
95 {
96 protected:
97 
98  // Protected data
99 
100  //- Tolerance used in weighted calculations
101  static scalar tolerance_;
102 
103  //- Cmu coefficient
104  scalar Cmu_;
105 
106  //- Von Karman constant
107  scalar kappa_;
108 
109  //- E coefficient
110  scalar E_;
111 
112  //- beta1 coefficient
113  scalar beta1_;
114 
115  //- Y+ at the edge of the laminar sublayer
116  scalar yPlusLam_;
117 
118  //- Local copy of turbulence G field
119  scalarField G_;
120 
121  //- Local copy of turbulence omega field
123 
124  //- Initialised flag
125  bool initialised_;
126 
127  //- Master patch ID
128  label master_;
129 
130  //- List of averaging corner weights
132 
133 
134  // Protected Member Functions
135 
136  //- Check the type of the patch
137  virtual void checkType();
138 
139  //- Write local wall function variables
140  virtual void writeLocalEntries(Ostream&) const;
141 
142  //- Set the master patch - master is responsible for updating all
143  // wall function patches
144  virtual void setMaster();
145 
146  //- Create the averaging weights for cells which are bounded by
147  // multiple wall function faces
148  virtual void createAveragingWeights();
149 
150  //- Helper function to return non-const access to an omega patch
152  (
153  const label patchi
154  );
155 
156  //- Main driver to calculate the turbulence fields
158  (
161  scalarField& omega0
162  );
163 
164  //- Calculate the omega and G
165  virtual void calculate
166  (
168  const List<scalar>& cornerWeights,
169  const fvPatch& patch,
170  scalarField& G,
172  );
173 
174  //- Return non-const access to the master patch ID
175  virtual label& master()
176  {
177  return master_;
178  }
179 
180 
181 public:
182 
183  //- Runtime type information
184  TypeName("omegaWallFunction");
185 
186 
187  // Constructors
188 
189  //- Construct from patch and internal field
191  (
192  const fvPatch&,
194  );
195 
196  //- Construct from patch, internal field and dictionary
198  (
199  const fvPatch&,
201  const dictionary&
202  );
203 
204  //- Construct by mapping given
205  // omegaWallFunctionFvPatchScalarField
206  // onto a new patch
208  (
210  const fvPatch&,
212  const fvPatchFieldMapper&
213  );
214 
215  //- Construct as copy
217  (
219  );
220 
221  //- Construct and return a clone
222  virtual tmp<fvPatchScalarField> clone() const
223  {
225  (
227  );
228  }
229 
230  //- Construct as copy setting internal field reference
232  (
235  );
236 
237  //- Construct and return a clone setting internal field reference
239  (
241  ) const
242  {
244  (
246  );
247  }
248 
249 
250  // Member functions
251 
252  // Access
253 
254  //- Return non-const access to the master's G field
255  scalarField& G(bool init = false);
256 
257  //- Return non-const access to the master's omega field
258  scalarField& omega(bool init = false);
259 
260 
261  // Evaluation functions
262 
263  //- Update the coefficients associated with the patch field
264  virtual void updateCoeffs();
265 
266  //- Update the coefficients associated with the patch field
267  virtual void updateCoeffs(const scalarField& weights);
268 
269  //- Manipulate matrix
270  virtual void manipulateMatrix(fvMatrix<scalar>& matrix);
271 
272  //- Manipulate matrix with given weights
273  virtual void manipulateMatrix
274  (
275  fvMatrix<scalar>& matrix,
276  const scalarField& weights
277  );
278 
279 
280  // I-O
281 
282  //- Write
283  virtual void write(Ostream&) const;
284 };
285 
286 
287 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
288 
289 } // End namespace Foam
290 
291 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
292 
293 #endif
294 
295 // ************************************************************************* //
Foam::omegaWallFunctionFvPatchScalarField::createAveragingWeights
virtual void createAveragingWeights()
Create the averaging weights for cells which are bounded by.
Definition: omegaWallFunctionFvPatchScalarField.C:99
Foam::omegaWallFunctionFvPatchScalarField::omegaWallFunctionFvPatchScalarField
omegaWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: omegaWallFunctionFvPatchScalarField.C:261
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:48
Foam::omegaWallFunctionFvPatchScalarField::manipulateMatrix
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
Manipulate matrix.
Definition: omegaWallFunctionFvPatchScalarField.C:520
Foam::omegaWallFunctionFvPatchScalarField::calculate
virtual void calculate(const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &omega)
Calculate the omega and G.
Definition: omegaWallFunctionFvPatchScalarField.C:208
Foam::omegaWallFunctionFvPatchScalarField::yPlusLam_
scalar yPlusLam_
Y+ at the edge of the laminar sublayer.
Definition: omegaWallFunctionFvPatchScalarField.H:148
Foam::omegaWallFunctionFvPatchScalarField
This boundary condition provides a wall function constraint on turbulnce specific dissipation,...
Definition: omegaWallFunctionFvPatchScalarField.H:124
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::omegaWallFunctionFvPatchScalarField::kappa_
scalar kappa_
Von Karman constant.
Definition: omegaWallFunctionFvPatchScalarField.H:139
Foam::omegaWallFunctionFvPatchScalarField::G_
scalarField G_
Local copy of turbulence G field.
Definition: omegaWallFunctionFvPatchScalarField.H:151
Foam::omegaWallFunctionFvPatchScalarField::beta1_
scalar beta1_
beta1 coefficient
Definition: omegaWallFunctionFvPatchScalarField.H:145
Foam::omegaWallFunctionFvPatchScalarField::writeLocalEntries
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
Definition: omegaWallFunctionFvPatchScalarField.C:60
Foam::omegaWallFunctionFvPatchScalarField::cornerWeights_
List< List< scalar > > cornerWeights_
List of averaging corner weights.
Definition: omegaWallFunctionFvPatchScalarField.H:163
Foam::omegaWallFunctionFvPatchScalarField::master_
label master_
Master patch ID.
Definition: omegaWallFunctionFvPatchScalarField.H:160
Foam::omegaWallFunctionFvPatchScalarField::omega
scalarField & omega(bool init=false)
Return non-const access to the master's omega field.
Definition: omegaWallFunctionFvPatchScalarField.C:393
Foam::fixedValueFvPatchField
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
Definition: fixedValueFvPatchField.H:79
Foam::omegaWallFunctionFvPatchScalarField::omega_
scalarField omega_
Local copy of turbulence omega field.
Definition: omegaWallFunctionFvPatchScalarField.H:154
Foam::omegaWallFunctionFvPatchScalarField::TypeName
TypeName("omegaWallFunction")
Runtime type information.
Foam::label
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
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::compressible::turbulenceModel
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Definition: turbulentFluidThermoModel.H:60
turbulence
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
Foam::constant::electromagnetic::G0
const dimensionedScalar G0
Conductance quantum: default SI units: [S].
Foam::omegaWallFunctionFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: omegaWallFunctionFvPatchScalarField.C:588
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::omegaWallFunctionFvPatchScalarField::Cmu_
scalar Cmu_
Cmu coefficient.
Definition: omegaWallFunctionFvPatchScalarField.H:136
Foam::omegaWallFunctionFvPatchScalarField::initialised_
bool initialised_
Initialised flag.
Definition: omegaWallFunctionFvPatchScalarField.H:157
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:60
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::omegaWallFunctionFvPatchScalarField::calculateTurbulenceFields
virtual void calculateTurbulenceFields(const turbulenceModel &turbulence, scalarField &G0, scalarField &omega0)
Main driver to calculate the turbulence fields.
Definition: omegaWallFunctionFvPatchScalarField.C:175
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::omegaWallFunctionFvPatchScalarField::checkType
virtual void checkType()
Check the type of the patch.
Definition: omegaWallFunctionFvPatchScalarField.C:46
Foam::List
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
Foam::omegaWallFunctionFvPatchScalarField::clone
virtual tmp< fvPatchScalarField > clone() const
Construct and return a clone.
Definition: omegaWallFunctionFvPatchScalarField.H:254
Foam::omegaWallFunctionFvPatchScalarField::G
scalarField & G(bool init=false)
Return non-const access to the master's G field.
Definition: omegaWallFunctionFvPatchScalarField.C:377
Foam::omegaWallFunctionFvPatchScalarField::E_
scalar E_
E coefficient.
Definition: omegaWallFunctionFvPatchScalarField.H:142
Foam::omegaWallFunctionFvPatchScalarField::omegaPatch
virtual omegaWallFunctionFvPatchScalarField & omegaPatch(const label patchi)
Helper function to return non-const access to an omega patch.
Definition: omegaWallFunctionFvPatchScalarField.C:160
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::omegaWallFunctionFvPatchScalarField::setMaster
virtual void setMaster()
Set the master patch - master is responsible for updating all.
Definition: omegaWallFunctionFvPatchScalarField.C:69
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::omegaWallFunctionFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: omegaWallFunctionFvPatchScalarField.C:409
Foam::omegaWallFunctionFvPatchScalarField::master
virtual label & master()
Return non-const access to the master patch ID.
Definition: omegaWallFunctionFvPatchScalarField.H:207
fixedValueFvPatchField.H
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::fvPatchField< scalar >::patch
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:300
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::omegaWallFunctionFvPatchScalarField::tolerance_
static scalar tolerance_
Tolerance used in weighted calculations.
Definition: omegaWallFunctionFvPatchScalarField.H:133
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51