epsilonWallFunctionFvPatchScalarField.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::epsilonWallFunctionFvPatchScalarField
26 
27 Group
28  grpWallFunctions
29 
30 Description
31  This boundary condition provides a turbulence dissipation wall function
32  condition for high Reynolds number, turbulent flow cases.
33 
34  The condition can be applied to wall boundaries, whereby it
35  - calculates \c epsilon and \c G
36  - inserts near wall epsilon values directly into the epsilon equation
37  to act as a constraint
38 
39  where
40 
41  \vartable
42  epsilon | turblence dissipation field
43  G | turblence generation field
44  \endvartable
45 
46  \heading Patch usage
47 
48  \table
49  Property | Description | Required | Default value
50  Cmu | model coefficient | no | 0.09
51  kappa | Von Karman constant | no | 0.41
52  E | model coefficient | no | 9.8
53  \endtable
54 
55  Example of the boundary condition specification:
56  \verbatim
57  myPatch
58  {
59  type epsilonWallFunction;
60  }
61  \endverbatim
62 
63 SeeAlso
64  Foam::fixedInternalValueFvPatchField
65 
66 SourceFiles
67  epsilonWallFunctionFvPatchScalarField.C
68 
69 \*---------------------------------------------------------------------------*/
70 
71 #ifndef epsilonWallFunctionFvPatchScalarField_H
72 #define epsilonWallFunctionFvPatchScalarField_H
73 
74 #include "fixedValueFvPatchField.H"
75 
76 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
77 
78 namespace Foam
79 {
80 
81 class turbulenceModel;
82 
83 /*---------------------------------------------------------------------------*\
84  Class epsilonWallFunctionFvPatchScalarField Declaration
85 \*---------------------------------------------------------------------------*/
86 
87 class epsilonWallFunctionFvPatchScalarField
88 :
89  public fixedValueFvPatchField<scalar>
90 {
91 protected:
92 
93  // Protected data
94 
95  //- Tolerance used in weighted calculations
96  static scalar tolerance_;
97 
98  //- Cmu coefficient
99  scalar Cmu_;
100 
101  //- Von Karman constant
102  scalar kappa_;
103 
104  //- E coefficient
105  scalar E_;
106 
107  //- Local copy of turbulence G field
108  scalarField G_;
109 
110  //- Local copy of turbulence epsilon field
112 
113  //- Initialised flag
115 
116  //- Master patch ID
117  label master_;
118 
119  //- List of averaging corner weights
121 
122 
123  // Protected Member Functions
124 
125  //- Check the type of the patch
126  virtual void checkType();
127 
128  //- Write local wall function variables
129  virtual void writeLocalEntries(Ostream&) const;
130 
131  //- Set the master patch - master is responsible for updating all
132  // wall function patches
133  virtual void setMaster();
134 
135  //- Create the averaging weights for cells which are bounded by
136  // multiple wall function faces
137  virtual void createAveragingWeights();
139  //- Helper function to return non-const access to an epsilon patch
141  (
142  const label patchi
143  );
145  //- Main driver to calculate the turbulence fields
146  virtual void calculateTurbulenceFields
147  (
149  scalarField& G0,
151  );
152 
153  //- Calculate the epsilon and G
154  virtual void calculate
155  (
156  const turbulenceModel& turbulence,
157  const List<scalar>& cornerWeights,
158  const fvPatch& patch,
159  scalarField& G,
161  );
162 
163  //- Return non-const access to the master patch ID
164  virtual label& master()
165  {
166  return master_;
167  }
168 
169 
170 public:
171 
172  //- Runtime type information
173  TypeName("epsilonWallFunction");
174 
175 
176  // Constructors
177 
178  //- Construct from patch and internal field
180  (
181  const fvPatch&,
183  );
184 
185  //- Construct from patch, internal field and dictionary
187  (
188  const fvPatch&,
190  const dictionary&
191  );
192 
193  //- Construct by mapping given
194  // epsilonWallFunctionFvPatchScalarField
195  // onto a new patch
197  (
199  const fvPatch&,
201  const fvPatchFieldMapper&
202  );
203 
204  //- Construct as copy
206  (
208  );
209 
210  //- Construct and return a clone
211  virtual tmp<fvPatchScalarField> clone() const
212  {
214  (
216  );
217  }
218 
219  //- Construct as copy setting internal field reference
221  (
224  );
225 
226  //- Construct and return a clone setting internal field reference
228  (
230  ) const
231  {
233  (
235  );
236  }
237 
238  //- Destructor
240  {}
241 
242 
243  // Member functions
244 
245  // Access
246 
247  //- Return non-const access to the master's G field
248  scalarField& G(bool init = false);
249 
250  //- Return non-const access to the master's epsilon field
251  scalarField& epsilon(bool init = false);
252 
253 
254  // Evaluation functions
256  //- Update the coefficients associated with the patch field
257  virtual void updateCoeffs();
258 
259  //- Update the coefficients associated with the patch field
260  virtual void updateCoeffs(const scalarField& weights);
261 
262  //- Manipulate matrix
263  virtual void manipulateMatrix(fvMatrix<scalar>& matrix);
264 
265  //- Manipulate matrix with given weights
266  virtual void manipulateMatrix
267  (
268  fvMatrix<scalar>& matrix,
269  const scalarField& weights
270  );
271 
272 
273  // I-O
274 
275  //- Write
276  virtual void write(Ostream&) const;
277 };
278 
279 
280 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
281 
282 } // End namespace Foam
283 
284 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
285 
286 #endif
287 
288 // ************************************************************************* //
virtual tmp< fvPatchScalarField > clone() const
Construct and return a clone.
scalarField & G(bool init=false)
Return non-const access to the master&#39;s G field.
This boundary condition provides a turbulence dissipation wall function condition for high Reynolds n...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionedScalar G0
Conductance quantum: default SI units: [S].
virtual void calculateTurbulenceFields(const turbulenceModel &turbulence, scalarField &G0, scalarField &epsilon0)
Main driver to calculate the turbulence fields.
static scalar tolerance_
Tolerance used in weighted calculations.
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::fvPatchFieldMapper.
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:68
virtual void calculate(const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &epsilon)
Calculate the epsilon and G.
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
List< List< scalar > > cornerWeights_
List of averaging corner weights.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Namespace for OpenFOAM.
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
Manipulate matrix.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
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
label patchi
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Abstract base class for turbulence models (RAS, LES and laminar).
TypeName("epsilonWallFunction")
Runtime type information.
scalarField G_
Local copy of turbulence G field.
scalarField & epsilon(bool init=false)
Return non-const access to the master&#39;s epsilon field.
virtual void setMaster()
Set the master patch - master is responsible for updating all.
virtual epsilonWallFunctionFvPatchScalarField & epsilonPatch(const label patchi)
Helper function to return non-const access to an epsilon patch.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
virtual void checkType()
Check the type of the patch.
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 label & master()
Return non-const access to the master patch ID.
scalarField epsilon_
Local copy of turbulence epsilon field.
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:300
epsilonWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A class for managing temporary objects.
Definition: PtrList.H:118
const dimensionedScalar epsilon0
Electric constant: default SI units: [F/m].