atmNutkWallFunctionFvPatchScalarField.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2018 OpenFOAM Foundation
9  Copyright (C) 2020 ENERCON GmbH
10  Copyright (C) 2020-2022 OpenCFD Ltd.
11 -------------------------------------------------------------------------------
12 License
13  This file is part of OpenFOAM.
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 
28 \*---------------------------------------------------------------------------*/
29 
31 #include "turbulenceModel.H"
32 #include "fvPatchFieldMapper.H"
33 #include "volFields.H"
34 #include "bound.H"
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41 
42 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
43 
45 {
46  const label patchi = patch().index();
47 
48  const auto& turbModel = db().lookupObject<turbulenceModel>
49  (
51  (
53  internalField().group()
54  )
55  );
56 
57  const scalarField& y = turbModel.y()[patchi];
58 
59  const tmp<volScalarField> tk = turbModel.k();
60  const volScalarField& k = tk();
61 
62  const tmp<scalarField> tnuw = turbModel.nu(patchi);
63  const scalarField& nuw = tnuw();
64 
65  auto tnutw = tmp<scalarField>::New(*this);
66  auto& nutw = tnutw.ref();
67 
68  const scalar Cmu25 = pow025(wallCoeffs_.Cmu());
69  const scalar kappa = wallCoeffs_.kappa();
70 
71  const scalar t = db().time().timeOutputValue();
72  const scalarField z0(z0_->value(t));
73 
74  #ifdef FULLDEBUG
75  for (const scalar z : z0)
76  {
77  if (z < VSMALL)
78  {
80  << "z0 field can only contain positive values. "
81  << "Please check input field z0."
82  << exit(FatalError);
83  }
84  }
85  #endif
86 
87  const labelList& faceCells = patch().faceCells();
88 
89  // (HW:Eq. 5)
90  forAll(nutw, facei)
91  {
92  const label celli = faceCells[facei];
93 
94  const scalar uStar = Cmu25*sqrt(k[celli]);
95  const scalar yPlus = uStar*y[facei]/nuw[facei];
96  const scalar Edash = (y[facei] + z0[facei])/z0[facei];
97 
98  nutw[facei] = nuw[facei]*(yPlus*kappa/log(max(Edash, 1 + 1e-4)) - 1);
99  }
100 
101  if (boundNut_)
102  {
103  nutw = max(nutw, scalar(0));
104  }
105 
106  return tnutw;
107 }
108 
109 
111 (
112  Ostream& os
113 ) const
114 {
115  os.writeEntryIfDifferent<bool>("boundNut", false, boundNut_);
116 
117  if (z0_)
118  {
119  z0_->writeData(os);
120  }
121 }
122 
123 
124 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
125 
127 (
128  const fvPatch& p,
129  const DimensionedField<scalar, volMesh>& iF
130 )
131 :
133  boundNut_(true),
134  z0_(nullptr)
135 {}
136 
137 
139 (
141  const fvPatch& p,
143  const fvPatchFieldMapper& mapper
144 )
145 :
147  boundNut_(ptf.boundNut_),
148  z0_(ptf.z0_.clone(p.patch()))
149 {}
150 
151 
153 (
154  const fvPatch& p,
156  const dictionary& dict
157 )
158 :
160  boundNut_(dict.getOrDefault<bool>("boundNut", false)),
161  z0_(PatchFunction1<scalar>::New(p.patch(), "z0", dict))
162 {}
163 
164 
166 (
168 )
169 :
171  boundNut_(rwfpsf.boundNut_),
172  z0_(rwfpsf.z0_.clone(this->patch().patch()))
173 {}
174 
175 
177 (
180 )
181 :
183  boundNut_(rwfpsf.boundNut_),
184  z0_(rwfpsf.z0_.clone(this->patch().patch()))
185 {}
186 
187 
188 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
189 
191 (
192  const fvPatchFieldMapper& m
193 )
194 {
195  nutkWallFunctionFvPatchScalarField::autoMap(m);
196 
197  if (z0_)
198  {
199  z0_->autoMap(m);
200  }
201 }
202 
203 
205 (
206  const fvPatchScalarField& ptf,
207  const labelList& addr
208 )
209 {
210  nutkWallFunctionFvPatchScalarField::rmap(ptf, addr);
211 
212  const auto& nrwfpsf =
213  refCast<const atmNutkWallFunctionFvPatchScalarField>(ptf);
214 
215  if (z0_)
216  {
217  z0_->rmap(nrwfpsf.z0_(), addr);
218  }
219 }
220 
221 
223 {
227 }
228 
229 
230 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
231 
233 (
235  atmNutkWallFunctionFvPatchScalarField
236 );
237 
238 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
239 
240 } // End namespace Foam
241 
242 // ************************************************************************* //
dictionary dict
dimensionedScalar log(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
This boundary condition provides a wall constraint on the turbulent viscosity (i.e. nut) based on the turbulent kinetic energy (i.e. k) for atmospheric boundary layer modelling. It is designed to be used in conjunction with the atmBoundaryLayerInletVelocity boundary condition.
virtual tmp< scalarField > calcNut() const
Calculate the turbulent viscosity.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
dimensionedScalar pow025(const dimensionedScalar &ds)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
label k
Boltzmann constant.
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:52
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Abstract base class for turbulence models (RAS, LES and laminar).
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
constexpr const char *const group
Group name for atomic constants.
scalar y
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
This boundary condition provides a wall function for the turbulent viscosity (i.e. nut) based on the turbulent kinetic energy, (i.e. k) for for low- and high-Reynolds number applications.
fvPatchField< scalar > fvPatchScalarField
static const word propertiesName
Default name of the turbulence properties dictionary.
A FieldMapper for finite-volume patch fields.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:336
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Bound the given scalar field if it has gone unbounded.
OBJstream os(runTime.globalPath()/outputName)
scalar Cmu() const noexcept
Return the object: Cmu.
void writeLocalEntries(Ostream &) const
Write local wall function variables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
wallFunctionCoefficients wallCoeffs_
Wall-function coefficients.
List< label > labelList
A List of labels.
Definition: List.H:62
scalar kappa() const noexcept
Return the object: kappa.
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
static const nutWallFunctionFvPatchScalarField & nutw(const turbulenceModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
atmNutkWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Namespace for OpenFOAM.