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-2022 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  const label patchi = patch().index();
42 
43  const momentumTransportModel& turbModel =
44  db().lookupObject<momentumTransportModel>
45  (
47  (
48  momentumTransportModel::typeName,
49  internalField().group()
50  )
51  );
52  const scalarField& y = turbModel.y()[patchi];
53  const tmp<volScalarField> tk = turbModel.k();
54  const volScalarField& k = tk();
55  const tmp<scalarField> tnuw = turbModel.nu(patchi);
56  const scalarField& nuw = tnuw();
57 
58  const scalar Cmu25 = pow025(Cmu_);
59 
60  tmp<scalarField> tnutw(new scalarField(*this));
61  scalarField& nutw = tnutw.ref();
62 
63  forAll(nutw, facei)
64  {
65  label celli = patch().faceCells()[facei];
66 
67  scalar uStar = Cmu25*sqrt(k[celli]);
68  scalar yPlus = uStar*y[facei]/nuw[facei];
69 
70  scalar Edash = (y[facei] + z0_[facei])/z0_[facei];
71 
72  nutw[facei] =
73  nuw[facei]*(yPlus*kappa_/log(max(Edash, 1+1e-4)) - 1);
74 
75  if (debug)
76  {
77  Info<< "yPlus = " << yPlus
78  << ", Edash = " << Edash
79  << ", nutw = " << nutw[facei]
80  << endl;
81  }
82  }
83 
84  return tnutw;
85 }
86 
87 
88 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
89 
92 (
93  const fvPatch& p,
95 )
96 :
98  z0_(p.size(), 0.0)
99 {}
100 
101 
104 (
105  const fvPatch& p,
107  const dictionary& dict
108 )
109 :
111  z0_("z0", dict, p.size())
112 {}
113 
114 
117 (
119  const fvPatch& p,
121  const fvPatchFieldMapper& mapper
122 )
123 :
124  nutkWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
125  z0_(mapper(ptf.z0_))
126 {}
127 
128 
131 (
134 )
135 :
137  z0_(rwfpsf.z0_)
138 {}
139 
140 
141 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
142 
144 (
145  const fvPatchFieldMapper& m
146 )
147 {
148  nutkWallFunctionFvPatchScalarField::autoMap(m);
149  m(z0_, z0_);
150 }
151 
152 
154 (
155  const fvPatchScalarField& ptf,
156  const labelList& addr
157 )
158 {
159  nutkWallFunctionFvPatchScalarField::rmap(ptf, addr);
160 
162  refCast<const nutkAtmRoughWallFunctionFvPatchScalarField>(ptf);
163 
164  z0_.rmap(nrwfpsf.z0_, addr);
165 }
166 
167 
169 (
170  const fvPatchScalarField& ptf
171 )
172 {
173  nutkWallFunctionFvPatchScalarField::reset(ptf);
174 
176  refCast<const nutkAtmRoughWallFunctionFvPatchScalarField>(ptf);
177 
178  z0_.reset(nrwfpsf.z0_);
179 }
180 
181 
183 {
185  writeLocalEntries(os);
186  writeEntry(os, "z0", z0_);
187  writeEntry(os, "value", *this);
188 }
189 
190 
191 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
192 
194 (
197 );
198 
199 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
200 
201 } // End namespace Foam
202 
203 // ************************************************************************* //
const char *const group
Group name for atomic constants.
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
nutkAtmRoughWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar log(const dimensionedScalar &ds)
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
dimensionedScalar pow025(const dimensionedScalar &ds)
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
label k
Boltzmann constant.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:243
const volScalarField::Boundary & y() const
Return the near wall distance.
Macros for easy insertion into run-time selection tables.
static const nutWallFunctionFvPatchScalarField & nutw(const momentumTransportModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
scalar y
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
Foam::fvPatchFieldMapper.
virtual label size() const
Return size.
Definition: fvPatch.H:157
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Abstract base class for turbulence models (RAS, LES and laminar).
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition: Field.C:362
label patchi
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
messageStream Info
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
A class for managing temporary objects.
Definition: PtrList.H:53
This boundary condition provides a turbulent kinematic viscosity for atmospheric velocity profiles...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
void reset(const Field< Type > &)
Reset the field values to the given field.
Definition: Field.C:432
virtual tmp< scalarField > nut() const
Calculate the turbulence viscosity.
nutkWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Namespace for OpenFOAM.