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 constraint
32  for low- and high-Reynolds number turbulence models.
33 
34  The condition can be applied to wall boundaries for which it
35  - calculates \c epsilon and \c G
36  - specifies the near-wall epsilon value
37 
38  where
39 
40  \vartable
41  epsilon | turblence dissipation field
42  G | turblence generation field
43  \endvartable
44 
45  The model switches between laminar and turbulent functions based on the
46  laminar-to-turbulent y+ value derived from kappa and E.
47 
48 Usage
49  \table
50  Property | Description | Required | Default value
51  Cmu | model coefficient | no | 0.09
52  kappa | Von Karman constant | no | 0.41
53  E | model coefficient | no | 9.8
54  \endtable
55 
56  Example of the boundary condition specification:
57  \verbatim
58  <patchName>
59  {
60  type epsilonWallFunction;
61  }
62  \endverbatim
63 
64 See also
65  Foam::fixedInternalValueFvPatchField
66  Foam::omegaWallFunctionFvPatchScalarField
67 
68 SourceFiles
69  epsilonWallFunctionFvPatchScalarField.C
70 
71 \*---------------------------------------------------------------------------*/
72 
73 #ifndef epsilonWallFunctionFvPatchScalarField_H
74 #define epsilonWallFunctionFvPatchScalarField_H
75 
76 #include "fixedValueFvPatchField.H"
77 
78 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
79 
80 namespace Foam
81 {
82 
83 class turbulenceModel;
84 
85 /*---------------------------------------------------------------------------*\
86  Class epsilonWallFunctionFvPatchScalarField Declaration
87 \*---------------------------------------------------------------------------*/
88 
89 class epsilonWallFunctionFvPatchScalarField
90 :
91  public fixedValueFvPatchField<scalar>
92 {
93 protected:
94 
95  // Protected data
96 
97  //- Tolerance used in weighted calculations
98  static scalar tolerance_;
99 
100  //- Cmu coefficient
101  scalar Cmu_;
102 
103  //- Von Karman constant
104  scalar kappa_;
105 
106  //- E coefficient
107  scalar E_;
108 
109  //- y+ at the edge of the laminar sublayer
110  scalar yPlusLam_;
111 
112  //- Local copy of turbulence G field
113  scalarField G_;
114 
115  //- Local copy of turbulence epsilon field
117 
118  //- Initialised flag
119  bool initialised_;
120 
121  //- Master patch ID
122  label master_;
123 
124  //- List of averaging corner weights
126 
127 
128  // Protected Member Functions
129 
130  //- Check the type of the patch
131  virtual void checkType();
132 
133  //- Write local wall function variables
134  virtual void writeLocalEntries(Ostream&) const;
135 
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();
144  //- Helper function to return non-const access to an epsilon patch
146  (
147  const label patchi
148  );
150  //- Main driver to calculate the turbulence fields
151  virtual void calculateTurbulenceFields
152  (
154  scalarField& G0,
156  );
157 
158  //- Calculate the epsilon and G
159  virtual void calculate
160  (
161  const turbulenceModel& turbulence,
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 
174 
175 public:
176 
177  //- Runtime type information
178  TypeName("epsilonWallFunction");
179 
180 
181  // Constructors
182 
183  //- Construct from patch and internal field
185  (
186  const fvPatch&,
188  );
189 
190  //- Construct from patch, internal field and dictionary
192  (
193  const fvPatch&,
195  const dictionary&
196  );
197 
198  //- Construct by mapping given
199  // epsilonWallFunctionFvPatchScalarField
200  // onto a new patch
202  (
204  const fvPatch&,
206  const fvPatchFieldMapper&
207  );
208 
209  //- Construct as copy
211  (
213  );
214 
215  //- Construct and return a clone
216  virtual tmp<fvPatchScalarField> clone() const
217  {
219  (
221  );
222  }
223 
224  //- Construct as copy setting internal field reference
226  (
229  );
230 
231  //- Construct and return a clone setting internal field reference
233  (
235  ) const
236  {
238  (
240  );
241  }
242 
243  //- Destructor
245  {}
246 
247 
248  // Member functions
249 
250  // Access
251 
252  //- Return non-const access to the master's G field
253  scalarField& G(bool init = false);
254 
255  //- Return non-const access to the master's epsilon field
256  scalarField& epsilon(bool init = false);
257 
258 
259  // Evaluation functions
261  //- Update the coefficients associated with the patch field
262  virtual void updateCoeffs();
263 
264  //- Update the coefficients associated with the patch field
265  virtual void updateWeightedCoeffs(const scalarField& weights);
266 
267  //- Manipulate matrix
268  virtual void manipulateMatrix(fvMatrix<scalar>& matrix);
269 
270  //- Manipulate matrix with given weights
271  virtual void manipulateMatrix
272  (
273  fvMatrix<scalar>& matrix,
274  const scalarField& weights
275  );
276 
277 
278  // I-O
279 
280  //- Write
281  virtual void write(Ostream&) const;
282 };
283 
284 
285 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
286 
287 } // End namespace Foam
288 
289 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
290 
291 #endif
292 
293 // ************************************************************************* //
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 writeLocalEntries(Ostream &) const
Write local wall function variables.
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.
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:60
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].
Abstract base class for turbulence models (RAS, LES and laminar).
virtual label & master()
Return non-const access to the master patch ID.
virtual tmp< fvPatchScalarField > clone() const
Construct and return a clone.
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.
scalar yPlusLam_
y+ at the edge of the laminar sublayer
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:72
This boundary condition provides a turbulence dissipation wall constraint for low- and high-Reynolds ...
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:340
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
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("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", 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
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:53
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.