yPlusTemplates.C
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) 2013-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
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 \*---------------------------------------------------------------------------*/
25 
26 #include "wallFvPatch.H"
28 #include "nearWallDist.H"
29 
30 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
31 
32 template<class TurbulenceModel>
34 (
36  const fvMesh& mesh,
38 )
39 {
41 
44 
47 
50 
51  const fvPatchList& patches = mesh.boundary();
52 
54  {
55  const fvPatch& patch = patches[patchi];
56 
57  if (isA<nutWallFunctionFvPatchScalarField>(nutBf[patchi]))
58  {
60  dynamic_cast<const nutWallFunctionFvPatchScalarField&>
61  (
62  nutBf[patchi]
63  );
64 
65  yPlus.boundaryField()[patchi] = nutPf.yPlus();
66  const scalarField& yPlusp = yPlus.boundaryField()[patchi];
67 
68  const scalar minYplus = gMin(yPlusp);
69  const scalar maxYplus = gMax(yPlusp);
70  const scalar avgYplus = gAverage(yPlusp);
71 
72  if (log_) Info<< " patch " << patch.name()
73  << " y+ : min = " << minYplus << ", max = " << maxYplus
74  << ", average = " << avgYplus << nl;
75 
76  writeTime(file());
77  file()
78  << token::TAB << patch.name()
79  << token::TAB << minYplus
80  << token::TAB << maxYplus
81  << token::TAB << avgYplus
82  << endl;
83  }
84  else if (isA<wallFvPatch>(patch))
85  {
86  yPlus.boundaryField()[patchi] =
87  d[patchi]
88  *sqrt
89  (
90  nuEffBf[patchi]
91  *mag(turbulenceModel.U().boundaryField()[patchi].snGrad())
92  )/nuBf[patchi];
93  const scalarField& yPlusp = yPlus.boundaryField()[patchi];
94 
95  const scalar minYplus = gMin(yPlusp);
96  const scalar maxYplus = gMax(yPlusp);
97  const scalar avgYplus = gAverage(yPlusp);
98 
99  if (log_) Info
100  << " patch " << patch.name()
101  << " y+ : min = " << minYplus << ", max = " << maxYplus
102  << ", average = " << avgYplus << nl;
103 
104  writeTime(file());
105  file()
106  << token::TAB << patch.name()
107  << token::TAB << minYplus
108  << token::TAB << maxYplus
109  << token::TAB << avgYplus
110  << endl;
111  }
112  }
113 }
114 
115 
116 // ************************************************************************* //
Foam::turbulenceModel::nu
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::turbulenceModel::nut
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:571
nutWallFunctionFvPatchScalarField.H
Foam::nutWallFunctionFvPatchScalarField::yPlus
virtual tmp< scalarField > yPlus() const =0
Calculate and return the yPlus at the boundary.
wallFvPatch.H
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::fvPatch::name
const word & name() const
Return name.
Definition: fvPatch.H:149
Foam::TurbulenceModel
Templated abstract base class for turbulence models.
Definition: TurbulenceModel.H:57
nearWallDist.H
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::nearWallDist::y
const volScalarField::GeometricBoundaryField & y() const
Definition: nearWallDist.H:87
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:60
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
boundaryField
cellIbMask *cellIbMaskExt *faceIbMask *cellIbMask boundaryField().evaluateCoupled()
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::yPlus
This function object evaluates and outputs turbulence y+ for turbulence models. The field is stored o...
Definition: yPlus.H:110
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::nearWallDist
Distance calculation for cells with face on a wall. Searches pointNeighbours to find closest.
Definition: nearWallDist.H:51
patches
patches[0]
Definition: createSingleCellMesh.H:36
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:563
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::turbulenceModel::nuEff
virtual tmp< volScalarField > nuEff() const =0
Return the effective viscosity.
Foam::nutWallFunctionFvPatchScalarField
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions,...
Definition: nutWallFunctionFvPatchScalarField.H:101
Foam::GeometricField::GeometricBoundaryField
Definition: GeometricField.H:105
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:562
Foam::turbulenceModel::U
const volVectorField & U() const
Access function to velocity field.
Definition: turbulenceModel.H:142
Foam::yPlus::calcYPlus
void calcYPlus(const TurbulenceModel &turbulenceModel, const fvMesh &mesh, volScalarField &yPlus)
Calculate y+.
Definition: yPlusTemplates.C:34