epsilonWallFunctionFvPatchScalarField.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 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 Description
28  This boundary condition provides a turbulence dissipation wall constraint
29  for low- and high-Reynolds number turbulence models.
30 
31  The condition can be applied to wall boundaries for which it
32  - calculates \c epsilon and \c G
33  - specifies the near-wall epsilon value
34 
35  where
36 
37  \vartable
38  epsilon | turbulence dissipation field
39  G | turbulence generation field
40  \endvartable
41 
42  The model switches between laminar and turbulent functions based on the
43  laminar-to-turbulent y+ value derived from the kappa and E specified in the
44  corresponding nutWallFunction.
45 
46 Usage
47  Example of the boundary condition specification:
48  \verbatim
49  <patchName>
50  {
51  type epsilonWallFunction;
52  }
53  \endverbatim
54 
55 See also
56  Foam::fixedInternalValueFvPatchField
57  Foam::omegaWallFunctionFvPatchScalarField
58 
59 SourceFiles
60  epsilonWallFunctionFvPatchScalarField.C
61 
62 \*---------------------------------------------------------------------------*/
63 
64 #ifndef epsilonWallFunctionFvPatchScalarField_H
65 #define epsilonWallFunctionFvPatchScalarField_H
66 
67 #include "fixedValueFvPatchField.H"
68 
69 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
70 
71 namespace Foam
72 {
73 
75 
76 /*---------------------------------------------------------------------------*\
77  Class epsilonWallFunctionFvPatchScalarField Declaration
78 \*---------------------------------------------------------------------------*/
79 
80 class epsilonWallFunctionFvPatchScalarField
81 :
82  public fixedValueFvPatchField<scalar>
83 {
84 protected:
85 
86  // Protected data
87 
88  //- Tolerance used in weighted calculations
89  static scalar tolerance_;
90 
91  //- Local copy of turbulence G field
93 
94  //- Local copy of turbulence epsilon field
96 
97  //- Initialised flag
98  bool initialised_;
99 
100  //- Master patch ID
101  label master_;
103  //- List of averaging corner weights
106 
107  // Protected Member Functions
109  //- Set the master patch - master is responsible for updating all
110  // wall function patches
111  virtual void setMaster();
112 
113  //- Create the averaging weights for cells which are bounded by
114  // multiple wall function faces
115  virtual void createAveragingWeights();
116 
117  //- Helper function to return non-const access to an epsilon patch
119  (
120  const label patchi
121  );
122 
123  //- Main driver to calculate the turbulence fields
124  virtual void calculateTurbulenceFields
125  (
126  const momentumTransportModel& turbModel,
127  scalarField& G0,
129  );
130 
131  //- Calculate the epsilon and G
132  virtual void calculate
133  (
134  const momentumTransportModel& turbModel,
135  const List<scalar>& cornerWeights,
136  const fvPatch& patch,
137  scalarField& G,
139  );
140 
141  //- Return non-const access to the master patch ID
142  virtual label& master()
143  {
144  return master_;
145  }
146 
147 
148 public:
150  //- Runtime type information
151  TypeName("epsilonWallFunction");
152 
153 
154  // Constructors
155 
156  //- Construct from patch and internal field
158  (
159  const fvPatch&,
161  );
162 
163  //- Construct from patch, internal field and dictionary
165  (
166  const fvPatch&,
168  const dictionary&
169  );
170 
171  //- Construct by mapping given
172  // epsilonWallFunctionFvPatchScalarField
173  // onto a new patch
175  (
177  const fvPatch&,
179  const fvPatchFieldMapper&
180  );
181 
182  //- Disallow copy without setting internal field reference
184  (
186  ) = delete;
187 
188  //- Copy constructor setting internal field reference
190  (
193  );
194 
195  //- Construct and return a clone setting internal field reference
197  (
199  ) const
200  {
202  (
204  );
205  }
206 
207 
208  //- Destructor
210  {}
211 
212 
213  // Member Functions
214 
215  // Access
217  //- Return non-const access to the master's G field
218  scalarField& G(bool init = false);
219 
220  //- Return non-const access to the master's epsilon field
221  scalarField& epsilon(bool init = false);
222 
223 
224  // Evaluation functions
225 
226  //- Update the coefficients associated with the patch field
227  virtual void updateCoeffs();
228 
229  //- Manipulate matrix
230  virtual void manipulateMatrix(fvMatrix<scalar>& matrix);
231 };
232 
233 
234 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
235 
236 } // End namespace Foam
237 
238 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
239 
240 #endif
241 
242 // ************************************************************************* //
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
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.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
Manipulate matrix.
static scalar tolerance_
Tolerance used in weighted calculations.
virtual void calculate(const momentumTransportModel &turbModel, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &epsilon)
Calculate the epsilon and G.
virtual void calculateTurbulenceFields(const momentumTransportModel &turbModel, scalarField &G0, scalarField &epsilon0)
Main driver to calculate the turbulence fields.
compressibleMomentumTransportModel momentumTransportModel
virtual label & master()
Return non-const access to the master patch ID.
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:351
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.
Abstract base class for turbulence models (RAS, LES and laminar).
label patchi
List< List< scalar > > cornerWeights_
List of averaging corner weights.
const dimensionedScalar G0
Conductance quantum: default SI units: [S].
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionedScalar epsilon0
Electric constant: default SI units: [F/m].
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.
tmp< fvPatchField< Type > > clone() const
Disallow clone without setting internal field reference.
Definition: fvPatchField.H:199
Namespace for OpenFOAM.
scalarField & epsilon(bool init=false)
Return non-const access to the master&#39;s epsilon field.