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 turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
49  (
51  (
53  internalField().group()
54  )
55  );
56  const scalarField& y = turbModel.y()[patchi];
57 
58  const tmp<volScalarField> tk = turbModel.k();
59  const volScalarField& k = tk();
60 
61  const tmp<scalarField> tnuw = turbModel.nu(patchi);
62  const scalarField& nuw = tnuw();
63 
64  tmp<scalarField> tnutw(new scalarField(*this));
65  scalarField& nutw = tnutw.ref();
66 
67  const scalar Cmu25 = pow025(Cmu_);
68 
69  const scalar t = db().time().timeOutputValue();
70  const scalarField z0(z0_->value(t));
71 
72  #ifdef FULLDEBUG
73  for (const auto& z : z0)
74  {
75  if (z < VSMALL)
76  {
78  << "z0 field can only contain positive values. "
79  << "Please check input field z0."
80  << exit(FatalIOError);
81  }
82  }
83  #endif
84 
85  const labelList& faceCells = patch().faceCells();
86 
87  // (HW:Eq. 5)
88  forAll(nutw, facei)
89  {
90  const label celli = faceCells[facei];
91 
92  const scalar uStar = Cmu25*sqrt(k[celli]);
93  const scalar yPlus = uStar*y[facei]/nuw[facei];
94  const scalar Edash = (y[facei] + z0[facei])/z0[facei];
95 
96  nutw[facei] = nuw[facei]*(yPlus*kappa_/log(max(Edash, 1 + 1e-4)) - 1);
97  }
98 
99  if (boundNut_)
100  {
101  nutw = max(nutw, scalar(0.0));
102  }
103 
104  return tnutw;
105 }
106 
107 
108 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
109 
111 (
112  const fvPatch& p,
114 )
115 :
117  boundNut_(true),
118  z0_(nullptr)
119 {}
120 
121 
123 (
125  const fvPatch& p,
127  const fvPatchFieldMapper& mapper
128 )
129 :
130  nutkWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
131  boundNut_(ptf.boundNut_),
132  z0_(ptf.z0_.clone(p.patch()))
133 {}
134 
135 
137 (
138  const fvPatch& p,
140  const dictionary& dict
141 )
142 :
144  boundNut_(dict.getOrDefault<Switch>("boundNut", false)),
145  z0_(PatchFunction1<scalar>::New(p.patch(), "z0", dict))
146 {}
147 
148 
150 (
152 )
153 :
155  boundNut_(rwfpsf.boundNut_),
156  z0_(rwfpsf.z0_.clone(this->patch().patch()))
157 {}
158 
159 
161 (
164 )
165 :
167  boundNut_(rwfpsf.boundNut_),
168  z0_(rwfpsf.z0_.clone(this->patch().patch()))
169 {}
170 
171 
172 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
173 
175 (
176  const fvPatchFieldMapper& m
177 )
178 {
179  nutkWallFunctionFvPatchScalarField::autoMap(m);
180  z0_->autoMap(m);
181 }
182 
183 
185 (
186  const fvPatchScalarField& ptf,
187  const labelList& addr
188 )
189 {
190  nutkWallFunctionFvPatchScalarField::rmap(ptf, addr);
191 
193  refCast<const atmNutkWallFunctionFvPatchScalarField>(ptf);
194 
195  z0_->rmap(nrwfpsf.z0_(), addr);
196 }
197 
198 
200 {
202  os.writeEntry("boundNut", boundNut_);
203  z0_->writeData(os);
204  writeEntry("value", os);
205 }
206 
207 
208 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
209 
211 (
214 );
215 
216 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
217 
218 } // End namespace Foam
219 
220 // ************************************************************************* //
Foam::fvPatchField< scalar >
volFields.H
Foam::atmNutkWallFunctionFvPatchScalarField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: atmNutkWallFunctionFvPatchScalarField.C:175
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::turbulenceModel::nu
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
atmNutkWallFunctionFvPatchScalarField.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:76
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:59
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
Foam::FatalIOError
IOerror FatalIOError
fvPatchFieldMapper.H
Foam::turbulenceModel::y
const nearWallDist & y() const
Return the near wall distances.
Definition: turbulenceModel.H:159
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::tmp::ref
T & ref() const
Definition: tmpI.H:258
Foam::Field< scalar >
Foam::nutWallFunctionFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: nutWallFunctionFvPatchScalarField.C:344
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:111
Foam::nutkWallFunctionFvPatchScalarField::yPlus
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
Definition: nutkWallFunctionFvPatchScalarField.C:141
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::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:262
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:63
Foam::nutWallFunctionFvPatchScalarField::Cmu_
scalar Cmu_
Empirical model coefficient.
Definition: nutWallFunctionFvPatchScalarField.H:287
Foam::turbulenceModel::k
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
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:199
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:372
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::atmNutkWallFunctionFvPatchScalarField::rmap
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: atmNutkWallFunctionFvPatchScalarField.C:185
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:232
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:122
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::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