omegaWallFunctionFvPatchScalarField.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::omegaWallFunctionFvPatchScalarField
29 
30 Group
31  grpWallFunctions
32 
33 Description
34  This boundary condition provides a wall function for the specific
35  dissipation rate (i.e. \c omega) and the turbulent kinetic energy
36  production contribution (i.e. \c G) for low- and high-Reynolds number
37  applications.
38 
39 Usage
40  Example of the boundary condition specification:
41  \verbatim
42  <patchName>
43  {
44  // Mandatory entries
45  type omegaWallFunction;
46 
47  // Optional entries
48  beta1 0.075;
49 
50  // Inherited entries
51  ...
52  }
53  \endverbatim
54 
55  \table
56  Property | Description | Type | Reqd | Deflt
57  type | Type name: omegaWallFunction | word | yes | -
58  beta1 | Model coefficient | scalar | no | 0.075
59  \endtable
60 
61  The inherited entries are elaborated in:
62  - \link wallFunctionCoefficients.H \endlink
63  - \link wallFunctionBlenders.H \endlink
64 
65 SourceFiles
66  omegaWallFunctionFvPatchScalarField.C
67 
68 \*---------------------------------------------------------------------------*/
69 
70 #ifndef omegaWallFunctionFvPatchScalarField_H
71 #define omegaWallFunctionFvPatchScalarField_H
72 
73 #include "fixedValueFvPatchField.H"
75 #include "wallFunctionBlenders.H"
76 
77 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
78 
79 namespace Foam
80 {
81 
82 class turbulenceModel;
83 
84 /*---------------------------------------------------------------------------*\
85  Class omegaWallFunctionFvPatchScalarField Declaration
86 \*---------------------------------------------------------------------------*/
87 
88 class omegaWallFunctionFvPatchScalarField
89 :
90  public fixedValueFvPatchField<scalar>,
91  private wallFunctionBlenders
92 {
93 protected:
94 
95  // Protected Data
96 
97  //- Tolerance used in weighted calculations
98  static scalar tolerance_;
99 
100  //- Initialised flag
102 
103  //- Master patch ID
104  label master_;
105 
106  //- beta1 coefficient
107  scalar beta1_;
108 
109  //- Wall-function coefficients
111 
112  //- Local copy of turbulence G field
114 
115  //- Local copy of turbulence omega field
117 
118  //- List of averaging corner weights
120 
121 
122  // Protected Member Functions
124  //- Set the master patch - master is responsible for updating all
125  //- wall function patches
126  virtual void setMaster();
127 
128  //- Create the averaging weights for cells which are bounded by
129  //- multiple wall function faces
130  virtual void createAveragingWeights();
131 
132  //- Helper function to return non-const access to an omega patch
134  (
135  const label patchi
136  );
137 
138  //- Main driver to calculate the turbulence fields
139  virtual void calculateTurbulenceFields
140  (
142  scalarField& G0,
143  scalarField& omega0
144  );
145 
146  //- Calculate the omega and G
147  virtual void calculate
148  (
150  const List<scalar>& cornerWeights,
151  const fvPatch& patch,
152  scalarField& G,
154  );
155 
156  //- Return non-const access to the master patch ID
157  virtual label& master()
158  {
159  return master_;
160  }
161 
162  //- Write local wall function variables
163  void writeLocalEntries(Ostream&) const;
164 
165 
166 public:
167 
168  //- Runtime type information
169  TypeName("omegaWallFunction");
170 
171 
172  // Constructors
173 
174  //- Construct from patch and internal field
176  (
177  const fvPatch&,
179  );
180 
181  //- Construct from patch, internal field and dictionary
183  (
184  const fvPatch&,
186  const dictionary&
187  );
188 
189  //- Construct by mapping given
190  //- omegaWallFunctionFvPatchScalarField
191  //- onto a new patch
193  (
195  const fvPatch&,
197  const fvPatchFieldMapper&
198  );
199 
200  //- Construct as copy
202  (
204  );
205 
206  //- Construct as copy setting internal field reference
208  (
211  );
212 
213  //- Return a clone
214  virtual tmp<fvPatchField<scalar>> clone() const
215  {
216  return fvPatchField<scalar>::Clone(*this);
217  }
218 
219  //- Clone with an internal field reference
221  (
223  ) const
224  {
225  return fvPatchField<scalar>::Clone(*this, iF);
226  }
227 
228 
229  //- Destructor
230  virtual ~omegaWallFunctionFvPatchScalarField() = default;
231 
232 
233  // Member Functions
234 
235  // Access
236 
237  //- Return non-const access to the master's G field
238  scalarField& G(bool init = false);
239 
240  //- Return non-const access to the master's omega field
241  scalarField& omega(bool init = false);
242 
243 
244  // Evaluation
245 
246  //- Update the coefficients associated with the patch field
247  virtual void updateCoeffs();
248 
249  //- Update the coefficients associated with the patch field
250  virtual void updateWeightedCoeffs(const scalarField& weights);
251 
252  //- Manipulate matrix
253  virtual void manipulateMatrix(fvMatrix<scalar>& matrix);
254 
255  //- Manipulate matrix with given weights
256  virtual void manipulateMatrix
257  (
258  fvMatrix<scalar>& matrix,
259  const scalarField& weights
260  );
261 
262 
263  // I-O
264 
265  //- Write
266  virtual void write(Ostream&) const;
267 };
268 
269 
270 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
272 } // End namespace Foam
273 
274 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
275 
276 #endif
277 
278 // ************************************************************************* //
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
const dimensionedScalar G0
Conductance quantum: default SI units: [S].
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
omegaWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
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
virtual void updateWeightedCoeffs(const scalarField &weights)
Update the coefficients associated with the patch field.
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
virtual void createAveragingWeights()
Create the averaging weights for cells which are bounded by multiple wall function faces...
scalarField & omega(bool init=false)
Return non-const access to the master&#39;s omega field.
Abstract base class for turbulence models (RAS, LES and laminar).
virtual omegaWallFunctionFvPatchScalarField & omegaPatch(const label patchi)
Helper function to return non-const access to an omega patch.
static tmp< fvPatchField< Type > > Clone(const DerivedPatchField &pf, Args &&... args)
Clone a patch field, optionally with internal field reference etc.
Definition: fvPatchField.H:597
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
static scalar tolerance_
Tolerance used in weighted calculations.
virtual void calculateTurbulenceFields(const turbulenceModel &turbulence, scalarField &G0, scalarField &omega0)
Main driver to calculate the turbulence fields.
scalarField & G(bool init=false)
Return non-const access to the master&#39;s G field.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
virtual tmp< fvPatchField< scalar > > clone() const
Return a clone.
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
wallFunctionCoefficients wallCoeffs_
Wall-function coefficients.
List< List< scalar > > cornerWeights_
List of averaging corner weights.
scalarField G_
Local copy of turbulence G field.
void writeLocalEntries(Ostream &) const
Write local wall function variables.
This boundary condition provides a wall function for the specific dissipation rate (i...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual ~omegaWallFunctionFvPatchScalarField()=default
Destructor.
virtual void setMaster()
Set the master patch - master is responsible for updating all wall function patches.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
TypeName("omegaWallFunction")
Runtime type information.
Class to host the wall-function coefficients being used in the wall function boundary conditions...
virtual label & master()
Return non-const access to the master patch ID.
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
Manipulate matrix.
scalarField omega_
Local copy of turbulence omega field.
Namespace for OpenFOAM.
virtual void calculate(const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &omega)
Calculate the omega and G.