nutkAtmRoughWallFunctionFvPatchScalarField.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) 2011-2016 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 "turbulenceModel.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 turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
44  (
46  (
48  internalField().group()
49  )
50  );
51  const scalarField& y = turbModel.y()[patchi];
52  const tmp<volScalarField> tk = turbModel.k();
53  const volScalarField& k = tk();
54  const tmp<scalarField> tnuw = turbModel.nu(patchi);
55  const scalarField& nuw = tnuw();
56 
57  const scalar Cmu25 = pow025(Cmu_);
58 
59  tmp<scalarField> tnutw(new scalarField(*this));
60  scalarField& nutw = tnutw.ref();
61 
62  forAll(nutw, facei)
63  {
64  label faceCelli = patch().faceCells()[facei];
65 
66  scalar uStar = Cmu25*sqrt(k[faceCelli]);
67  scalar yPlus = uStar*y[facei]/nuw[facei];
68 
69  scalar Edash = (y[facei] + z0_[facei])/z0_[facei];
70 
71  nutw[facei] =
72  nuw[facei]*(yPlus*kappa_/log(max(Edash, 1+1e-4)) - 1);
73 
74  if (debug)
75  {
76  Info<< "yPlus = " << yPlus
77  << ", Edash = " << Edash
78  << ", nutw = " << nutw[facei]
79  << endl;
80  }
81  }
82 
83  return tnutw;
84 }
85 
86 
87 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
88 
91 (
92  const fvPatch& p,
94 )
95 :
97  z0_(p.size(), 0.0)
98 {}
99 
100 
103 (
105  const fvPatch& p,
107  const fvPatchFieldMapper& mapper
108 )
109 :
110  nutkWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
111  z0_(ptf.z0_, mapper)
112 {}
113 
114 
117 (
118  const fvPatch& p,
120  const dictionary& dict
121 )
122 :
124  z0_("z0", dict, p.size())
125 {}
126 
127 
130 (
132 )
133 :
135  z0_(rwfpsf.z0_)
136 {}
137 
138 
141 (
144 )
145 :
147  z0_(rwfpsf.z0_)
148 {}
149 
150 
151 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
152 
154 (
155  const fvPatchFieldMapper& m
156 )
157 {
158  nutkWallFunctionFvPatchScalarField::autoMap(m);
159  z0_.autoMap(m);
160 }
161 
162 
164 (
165  const fvPatchScalarField& ptf,
166  const labelList& addr
167 )
168 {
169  nutkWallFunctionFvPatchScalarField::rmap(ptf, addr);
170 
172  refCast<const nutkAtmRoughWallFunctionFvPatchScalarField>(ptf);
173 
174  z0_.rmap(nrwfpsf.z0_, addr);
175 }
176 
177 
179 {
181  writeLocalEntries(os);
182  z0_.writeEntry("z0", os);
183  writeEntry("value", os);
184 }
185 
186 
187 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
188 
190 (
193 );
194 
195 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
196 
197 } // End namespace Foam
198 
199 // ************************************************************************* //
const char *const group
Group name for atomic constants.
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:428
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
dimensionedScalar log(const dimensionedScalar &ds)
const double e
Elementary charge.
Definition: doubleFloat.H:78
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
dimensionedScalar pow025(const dimensionedScalar &ds)
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:719
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
label k
Boltzmann constant.
Abstract base class for turbulence models (RAS, LES and laminar).
Macros for easy insertion into run-time selection tables.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
scalar y
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
static const word propertiesName
Default name of the turbulence properties dictionary.
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:161
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
void autoMap(const FieldMapper &map, const bool applyFlip=true)
Map from self.
Definition: Field.C:514
virtual tmp< scalarField > calcNut() const
Calculate the turbulence viscosity.
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition: Field.C:577
label patchi
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
messageStream Info
A class for managing temporary objects.
Definition: PtrList.H:54
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
This boundary condition provides a turbulent kinematic viscosity for atmospheric velocity profiles...
const nearWallDist & y() const
Return the near wall distances.
nutkWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Namespace for OpenFOAM.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:363