nutkAtmRoughWallFunctionFvPatchScalarField.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-2024 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 "fieldMapper.H"
29 #include "volFields.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
38 
40 {
41  const label patchi = patch().index();
42 
43  const momentumTransportModel& turbModel =
44  db().lookupType<momentumTransportModel>(internalField().group());
45 
46  const scalarField& y = turbModel.yb()[patchi];
47  const tmp<volScalarField> tk = turbModel.k();
48  const volScalarField& k = tk();
49  const tmp<scalarField> tnuw = turbModel.nu(patchi);
50  const scalarField& nuw = tnuw();
51 
52  const scalar Cmu25 = pow025(Cmu_);
53 
54  tmp<scalarField> tnutw(new scalarField(*this));
55  scalarField& nutw = tnutw.ref();
56 
57  forAll(nutw, facei)
58  {
59  label celli = patch().faceCells()[facei];
60 
61  scalar uStar = Cmu25*sqrt(k[celli]);
62  scalar yPlus = uStar*y[facei]/nuw[facei];
63 
64  scalar Edash = (y[facei] + z0_[facei])/z0_[facei];
65 
66  nutw[facei] =
67  nuw[facei]*(yPlus*kappa_/log(max(Edash, 1+1e-4)) - 1);
68 
69  if (debug)
70  {
71  Info<< "yPlus = " << yPlus
72  << ", Edash = " << Edash
73  << ", nutw = " << nutw[facei]
74  << endl;
75  }
76  }
77 
78  return tnutw;
79 }
80 
81 
82 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
83 
86 (
87  const fvPatch& p,
89  const dictionary& dict
90 )
91 :
93  z0_("z0", dimLength, dict, p.size())
94 {}
95 
96 
99 (
101  const fvPatch& p,
103  const fieldMapper& mapper
104 )
105 :
106  nutkWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
107  z0_(mapper(ptf.z0_))
108 {}
109 
110 
113 (
116 )
117 :
119  z0_(rwfpsf.z0_)
120 {}
121 
122 
123 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
124 
126 (
127  const fvPatchScalarField& ptf,
128  const fieldMapper& mapper
129 )
130 {
132 
134  refCast<const nutkAtmRoughWallFunctionFvPatchScalarField>(ptf);
135 
136  mapper(z0_, nrwfpsf.z0_);
137 }
138 
139 
141 (
142  const fvPatchScalarField& ptf
143 )
144 {
145  nutkWallFunctionFvPatchScalarField::reset(ptf);
146 
148  refCast<const nutkAtmRoughWallFunctionFvPatchScalarField>(ptf);
149 
150  z0_.reset(nrwfpsf.z0_);
151 }
152 
153 
155 {
157  writeLocalEntries(os);
158  writeEntry(os, "z0", z0_);
159  writeEntry(os, "value", *this);
160 }
161 
162 
163 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
164 
166 (
169 );
170 
171 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172 
173 } // End namespace Foam
174 
175 // ************************************************************************* //
scalar y
label k
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
void reset(const Field< Type > &)
Reset the field values to the given field.
Definition: Field.C:456
Generic GeometricField class.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:229
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
Abstract base class for turbulence models (RAS, LES and laminar).
const volScalarField::Boundary & yb() const
Return the near wall distance.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
virtual void map(const fvPatchScalarField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
static const nutWallFunctionFvPatchScalarField & nutw(const momentumTransportModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
This boundary condition provides a turbulent kinematic viscosity for atmospheric velocity profiles....
virtual tmp< scalarField > nut() const
Calculate the turbulence viscosity.
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
nutkAtmRoughWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void map(const fvPatchScalarField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
Turbulent viscosity wall-function boundary condition for high Reynolds number flows based on near-wal...
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
label patchi
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:106
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
messageStream Info
const dimensionSet dimLength
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar sqrt(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dimensionedScalar pow025(const dimensionedScalar &ds)
dictionary dict
volScalarField & p