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 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 =
49  db().lookupObject<turbulenceModel>
50  (
52  (
54  internalField().group()
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(Cmu_);
69 
70  const scalar t = db().time().timeOutputValue();
71  const scalarField z0(z0_->value(t));
72 
73  #ifdef FULLDEBUG
74  for (const scalar z : z0)
75  {
76  if (z < VSMALL)
77  {
79  << "z0 field can only contain positive values. "
80  << "Please check input field z0."
81  << exit(FatalError);
82  }
83  }
84  #endif
85 
86  const labelList& faceCells = patch().faceCells();
87 
88  // (HW:Eq. 5)
89  forAll(nutw, facei)
90  {
91  const label celli = faceCells[facei];
92 
93  const scalar uStar = Cmu25*sqrt(k[celli]);
94  const scalar yPlus = uStar*y[facei]/nuw[facei];
95  const scalar Edash = (y[facei] + z0[facei])/z0[facei];
96 
97  nutw[facei] = nuw[facei]*(yPlus*kappa_/log(max(Edash, 1 + 1e-4)) - 1);
98  }
99 
100  if (boundNut_)
101  {
102  nutw = max(nutw, scalar(0));
103  }
104 
105  return tnutw;
106 }
107 
108 
109 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
110 
112 (
113  const fvPatch& p,
115 )
116 :
118  boundNut_(true),
119  z0_(nullptr)
120 {}
121 
122 
124 (
126  const fvPatch& p,
128  const fvPatchFieldMapper& mapper
129 )
130 :
131  nutkWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
132  boundNut_(ptf.boundNut_),
133  z0_(ptf.z0_.clone(p.patch()))
134 {}
135 
136 
138 (
139  const fvPatch& p,
141  const dictionary& dict
142 )
143 :
145  boundNut_(dict.getOrDefault<bool>("boundNut", false)),
146  z0_(PatchFunction1<scalar>::New(p.patch(), "z0", dict))
147 {}
148 
149 
151 (
153 )
154 :
156  boundNut_(rwfpsf.boundNut_),
157  z0_(rwfpsf.z0_.clone(this->patch().patch()))
158 {}
159 
160 
162 (
165 )
166 :
168  boundNut_(rwfpsf.boundNut_),
169  z0_(rwfpsf.z0_.clone(this->patch().patch()))
170 {}
171 
172 
173 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
174 
176 (
177  const fvPatchFieldMapper& m
178 )
179 {
180  nutkWallFunctionFvPatchScalarField::autoMap(m);
181  z0_->autoMap(m);
182 }
183 
184 
186 (
187  const fvPatchScalarField& ptf,
188  const labelList& addr
189 )
190 {
191  nutkWallFunctionFvPatchScalarField::rmap(ptf, addr);
192 
194  refCast<const atmNutkWallFunctionFvPatchScalarField>(ptf);
195 
196  z0_->rmap(nrwfpsf.z0_(), addr);
197 }
198 
199 
201 {
204  os.writeEntry("boundNut", boundNut_);
205  z0_->writeData(os);
206  writeEntry("value", os);
207 }
208 
209 
210 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
211 
213 (
216 );
217 
218 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
219 
220 } // End namespace Foam
221 
222 // ************************************************************************* //
Foam::fvPatchField< scalar >
volFields.H
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
Foam::atmNutkWallFunctionFvPatchScalarField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: atmNutkWallFunctionFvPatchScalarField.C:176
atmNutkWallFunctionFvPatchScalarField.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::constant::atomic::group
constexpr const char *const group
Group name for atomic constants.
Definition: atomicConstants.H:52
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::atmNutkWallFunctionFvPatchScalarField::calcNut
virtual tmp< scalarField > calcNut() const
Calculate the turbulent viscosity.
Definition: atmNutkWallFunctionFvPatchScalarField.C:44
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
fvPatchFieldMapper.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::nutWallFunctionFvPatchScalarField::kappa_
scalar kappa_
von Kármán constant
Definition: nutWallFunctionFvPatchScalarField.H:290
Foam::pow025
dimensionedScalar pow025(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:133
Foam::Field< scalar >
bound.H
Bound the given scalar field if it has gone unbounded.
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::atmNutkWallFunctionFvPatchScalarField::atmNutkWallFunctionFvPatchScalarField
atmNutkWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: atmNutkWallFunctionFvPatchScalarField.C:112
Foam::nutkWallFunctionFvPatchScalarField::yPlus
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
Definition: nutkWallFunctionFvPatchScalarField.C:140
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::atmNutkWallFunctionFvPatchScalarField
This boundary condition provides a wall constraint on the turbulent viscosity (i.e....
Definition: atmNutkWallFunctionFvPatchScalarField.H:173
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:63
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:262
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::PatchFunction1
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: PatchFunction1.H:60
Foam::nutWallFunctionFvPatchScalarField::Cmu_
scalar Cmu_
Empirical model coefficient.
Definition: nutWallFunctionFvPatchScalarField.H:287
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::atmNutkWallFunctionFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: atmNutkWallFunctionFvPatchScalarField.C:200
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::nutWallFunctionFvPatchScalarField::nutw
static const nutWallFunctionFvPatchScalarField & nutw(const turbulenceModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
Definition: nutWallFunctionFvPatchScalarField.C:235
Foam::List< label >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::atmNutkWallFunctionFvPatchScalarField::rmap
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: atmNutkWallFunctionFvPatchScalarField.C:186
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::faceCells
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:56
turbulenceModel.H
Foam::nutWallFunctionFvPatchScalarField::writeLocalEntries
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
Definition: nutWallFunctionFvPatchScalarField.C:91
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::nutkWallFunctionFvPatchScalarField
This boundary condition provides a wall constraint on the turbulent viscosity, i.e....
Definition: nutkWallFunctionFvPatchScalarField.H:94