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-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::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 Usage
47  \table
48  Property | Description | Required | Default value
49  Cmu | model coefficient | no | 0.09
50  kappa | Von Karman constant | no | 0.41
51  E | model coefficient | no | 9.8
52  \endtable
53 
54  Example of the boundary condition specification:
55  \verbatim
56  <patchName>
57  {
58  type epsilonWallFunction;
59  }
60  \endverbatim
61 
62 See also
63  Foam::fixedInternalValueFvPatchField
64 
65 SourceFiles
66  epsilonWallFunctionFvPatchScalarField.C
67 
68 \*---------------------------------------------------------------------------*/
69 
70 #ifndef epsilonWallFunctionFvPatchScalarField_H
71 #define epsilonWallFunctionFvPatchScalarField_H
72 
73 #include "fixedValueFvPatchField.H"
74 
75 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
76 
77 namespace Foam
78 {
79 
80 class turbulenceModel;
81 
82 /*---------------------------------------------------------------------------*\
83  Class epsilonWallFunctionFvPatchScalarField Declaration
84 \*---------------------------------------------------------------------------*/
85 
86 class epsilonWallFunctionFvPatchScalarField
87 :
88  public fixedValueFvPatchField<scalar>
89 {
90 protected:
91 
92  // Protected data
93 
94  //- Tolerance used in weighted calculations
95  static scalar tolerance_;
96 
97  //- Cmu coefficient
98  scalar Cmu_;
99 
100  //- Von Karman constant
101  scalar kappa_;
102 
103  //- E coefficient
104  scalar E_;
105 
106  //- Local copy of turbulence G field
107  scalarField G_;
108 
109  //- Local copy of turbulence epsilon field
111 
112  //- Initialised flag
114 
115  //- Master patch ID
116  label master_;
117 
118  //- List of averaging corner weights
120 
121 
122  // Protected Member Functions
123 
124  //- Check the type of the patch
125  virtual void checkType();
126 
127  //- Write local wall function variables
128  virtual void writeLocalEntries(Ostream&) const;
129 
130  //- Set the master patch - master is responsible for updating all
131  // wall function patches
132  virtual void setMaster();
133 
134  //- Create the averaging weights for cells which are bounded by
135  // multiple wall function faces
136  virtual void createAveragingWeights();
138  //- Helper function to return non-const access to an epsilon patch
140  (
141  const label patchi
142  );
144  //- Main driver to calculate the turbulence fields
145  virtual void calculateTurbulenceFields
146  (
148  scalarField& G0,
150  );
151 
152  //- Calculate the epsilon and G
153  virtual void calculate
154  (
155  const turbulenceModel& turbulence,
156  const List<scalar>& cornerWeights,
157  const fvPatch& patch,
158  scalarField& G,
160  );
161 
162  //- Return non-const access to the master patch ID
163  virtual label& master()
164  {
165  return master_;
166  }
167 
168 
169 public:
170 
171  //- Runtime type information
172  TypeName("epsilonWallFunction");
173 
174 
175  // Constructors
176 
177  //- Construct from patch and internal field
179  (
180  const fvPatch&,
182  );
183 
184  //- Construct from patch, internal field and dictionary
186  (
187  const fvPatch&,
189  const dictionary&
190  );
191 
192  //- Construct by mapping given
193  // epsilonWallFunctionFvPatchScalarField
194  // onto a new patch
196  (
198  const fvPatch&,
200  const fvPatchFieldMapper&
201  );
202 
203  //- Construct as copy
205  (
207  );
208 
209  //- Construct and return a clone
210  virtual tmp<fvPatchScalarField> clone() const
211  {
213  (
215  );
216  }
217 
218  //- Construct as copy setting internal field reference
220  (
223  );
224 
225  //- Construct and return a clone setting internal field reference
227  (
229  ) const
230  {
232  (
234  );
235  }
236 
237  //- Destructor
239  {}
240 
241 
242  // Member functions
243 
244  // Access
245 
246  //- Return non-const access to the master's G field
247  scalarField& G(bool init = false);
248 
249  //- Return non-const access to the master's epsilon field
250  scalarField& epsilon(bool init = false);
251 
252 
253  // Evaluation functions
255  //- Update the coefficients associated with the patch field
256  virtual void updateCoeffs();
257 
258  //- Update the coefficients associated with the patch field
259  virtual void updateWeightedCoeffs(const scalarField& weights);
260 
261  //- Manipulate matrix
262  virtual void manipulateMatrix(fvMatrix<scalar>& matrix);
263 
264  //- Manipulate matrix with given weights
265  virtual void manipulateMatrix
266  (
267  fvMatrix<scalar>& matrix,
268  const scalarField& weights
269  );
270 
271 
272  // I-O
273 
274  //- Write
275  virtual void write(Ostream&) const;
276 };
277 
278 
279 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
280 
281 } // End namespace Foam
282 
283 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
284 
285 #endif
286 
287 // ************************************************************************* //
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
const dimensionedScalar G0
Conductance quantum: default SI units: [S].
virtual void calculate(const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &epsilon)
Calculate the epsilon and G.
virtual tmp< fvPatchScalarField > clone() const
Construct and return a clone.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
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
scalarField G_
Local copy of turbulence G field.
TypeName("epsilonWallFunction")
Runtime type information.
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
Manipulate matrix.
static scalar tolerance_
Tolerance used in weighted calculations.
const dimensionedScalar epsilon0
Electric constant: default SI units: [F/m].
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
Abstract base class for turbulence models (RAS, LES and laminar).
virtual label & master()
Return non-const access to the master patch ID.
virtual void calculateTurbulenceFields(const turbulenceModel &turbulence, scalarField &G0, scalarField &epsilon0)
Main driver to calculate the turbulence fields.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Foam::fvPatchFieldMapper.
scalarField epsilon_
Local copy of turbulence epsilon field.
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
This boundary condition provides a turbulence dissipation wall function condition for high Reynolds n...
scalarField & G(bool init=false)
Return non-const access to the master&#39;s G field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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
label patchi
List< List< scalar > > cornerWeights_
List of averaging corner weights.
virtual void updateWeightedCoeffs(const scalarField &weights)
Update the coefficients associated with the patch field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual void checkType()
Check the type of the patch.
virtual void setMaster()
Set the master patch - master is responsible for updating all.
A class for managing temporary objects.
Definition: PtrList.H:54
virtual void createAveragingWeights()
Create the averaging weights for cells which are bounded by.
epsilonWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Namespace for OpenFOAM.
scalarField & epsilon(bool init=false)
Return non-const access to the master&#39;s epsilon field.