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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016, 2019 OpenFOAM Foundation
9  Copyright (C) 2019-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::epsilonWallFunctionFvPatchScalarField
29 
30 Group
31  grpWallFunctions
32 
33 Description
34  This boundary condition provides wall functions for the turbulent kinetic
35  energy dissipation rate (i.e. \c epsilon) and the turbulent kinetic
36  energy production contribution (i.e. \c G) for low- and high-Reynolds
37  number simulations.
38 
39 Usage
40  Example of the boundary condition specification:
41  \verbatim
42  <patchName>
43  {
44  // Mandatory entries
45  type epsilonWallFunction;
46 
47  // Optional entries
48  lowReCorrection false;
49 
50  // Inherited entries
51  ...
52  }
53  \endverbatim
54 
55  where the entries mean:
56  \table
57  Property | Description | Type | Reqd | Deflt
58  type | Type name: epsilonWallFunction | word | yes | -
59  lowReCorrection | Flag: apply low-Re correction | bool | no | false
60  \endtable
61 
62  The inherited entries are elaborated in:
63  - \link wallFunctionCoefficients.H \endlink
64  - \link wallFunctionBlenders.H \endlink
65 
66 Note
67  - \c lowReCorrection operates with only \c stepwise blending treatment to
68  ensure the backward compatibility.
69  - If \c lowReCorrection is \c on, \c stepwise blending treatment is fully
70  active.
71  - If \c lowReCorrection is \c off, only the inertial sublayer prediction
72  is used in the wall function, hence high-Re mode operation.
73 
74 SourceFiles
75  epsilonWallFunctionFvPatchScalarField.C
76 
77 \*---------------------------------------------------------------------------*/
78 
79 #ifndef epsilonWallFunctionFvPatchScalarField_H
80 #define epsilonWallFunctionFvPatchScalarField_H
81 
82 #include "fixedValueFvPatchField.H"
84 #include "wallFunctionBlenders.H"
85 
86 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
87 
88 namespace Foam
89 {
90 
91 class turbulenceModel;
92 
93 /*---------------------------------------------------------------------------*\
94  Class epsilonWallFunctionFvPatchScalarField Declaration
95 \*---------------------------------------------------------------------------*/
96 
97 class epsilonWallFunctionFvPatchScalarField
98 :
99  public fixedValueFvPatchField<scalar>,
100  private wallFunctionBlenders
101 {
102 protected:
103 
104  // Protected Data
105 
106  //- Tolerance used in weighted calculations
107  static scalar tolerance_;
108 
109  //- Apply low-Re correction term (default = no)
110  const bool lowReCorrection_;
111 
112  //- Initialised flag
113  bool initialised_;
114 
115  //- Master patch ID
116  label master_;
117 
118  //- Wall-function coefficients
120 
121  //- Local copy of turbulence G field
123 
124  //- Local copy of turbulence epsilon field
126 
127  //- List of averaging corner weights
129 
130 
131  // Protected Member Functions
133  //- Set the master patch - master is responsible for updating all
134  //- wall function patches
135  virtual void setMaster();
136 
137  //- Create the averaging weights for cells which are bounded by
138  //- multiple wall function faces
139  virtual void createAveragingWeights();
140 
141  //- Helper function to return non-const access to an epsilon patch
143  (
144  const label patchi
145  );
146 
147  //- Main driver to calculate the turbulence fields
148  virtual void calculateTurbulenceFields
149  (
151  scalarField& G0,
153  );
154 
155  //- Calculate the epsilon and G
156  virtual void calculate
157  (
159  const List<scalar>& cornerWeights,
160  const fvPatch& patch,
161  scalarField& G,
163  );
164 
165  //- Return non-const access to the master patch ID
166  virtual label& master()
167  {
168  return master_;
169  }
170 
171  //- Write local wall function variables
172  void writeLocalEntries(Ostream&) const;
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 
244  //- Destructor
245  virtual ~epsilonWallFunctionFvPatchScalarField() = default;
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
260 
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
293 // ************************************************************************* //
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 a number of values (eg...
Definition: dictionary.H:129
const fvPatch & patch() const noexcept
Return the patch.
Definition: fvPatchField.H:269
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
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:70
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 ~epsilonWallFunctionFvPatchScalarField()=default
Destructor.
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.
A FieldMapper for finite-volume patch fields.
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:64
This boundary condition provides wall functions for the turbulent kinetic energy dissipation rate (i...
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:56
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
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...
Definition: areaFieldsFwd.H:42
virtual void setMaster()
Set the master patch - master is responsible for updating all wall function patches.
wallFunctionCoefficients wallCoeffs_
Wall-function coefficients.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Class to host the wall-function coefficients being used in the wall function boundary conditions...
virtual void createAveragingWeights()
Create the averaging weights for cells which are bounded by multiple wall function faces...
epsilonWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
const bool lowReCorrection_
Apply low-Re correction term (default = no)
Namespace for OpenFOAM.
scalarField & epsilon(bool init=false)
Return non-const access to the master&#39;s epsilon field.