nutLowReWallFunctionFvPatchScalarField.C
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-2021 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 \*---------------------------------------------------------------------------*/
25 
27 #include "momentumTransportModel.H"
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
38 
40 {
41  return tmp<scalarField>(new scalarField(patch().size(), 0.0));
42 }
43 
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
48 (
49  const fvPatch& p,
51 )
52 :
54 {}
55 
56 
58 (
59  const fvPatch& p,
61  const dictionary& dict
62 )
63 :
65 {}
66 
67 
69 (
71  const fvPatch& p,
73  const fvPatchFieldMapper& mapper
74 )
75 :
76  nutWallFunctionFvPatchScalarField(ptf, p, iF, mapper)
77 {}
78 
79 
81 (
84 )
85 :
87 {}
88 
89 
90 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
91 
93 {
94  const label patchi = patch().index();
95  const momentumTransportModel& turbModel =
96  db().lookupObject<momentumTransportModel>
97  (
99  (
100  momentumTransportModel::typeName,
101  internalField().group()
102  )
103  );
104  const scalarField& y = turbModel.y()[patchi];
105  const tmp<scalarField> tnuw = turbModel.nu(patchi);
106  const scalarField& nuw = tnuw();
107  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
108 
109  return y*sqrt(nuw*mag(Uw.snGrad()))/nuw;
110 }
111 
112 
113 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
114 
116 (
119 );
120 
121 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
122 
123 } // End namespace Foam
124 
125 // ************************************************************************* //
This boundary condition provides a turbulent kinematic viscosity condition for use with low Reynolds ...
const char *const group
Group name for atomic constants.
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:166
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
const Boundary & boundaryField() const
Return const-reference to the boundary field.
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
const volScalarField::Boundary & y() const
Return the near wall distance.
Macros for easy insertion into run-time selection tables.
scalar y
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
Foam::fvPatchFieldMapper.
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
Abstract base class for turbulence models (RAS, LES and laminar).
const volVectorField & U() const
Access function to velocity field.
label patchi
virtual tmp< scalarField > nut() const
Calculate the turbulence viscosity.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensioned< scalar > mag(const dimensioned< Type > &)
nutWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A class for managing temporary objects.
Definition: PtrList.H:53
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
nutLowReWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Namespace for OpenFOAM.