atmNutUWallFunctionFvPatchScalarField.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) 2020 ENERCON GmbH
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 #include "turbulenceModel.H"
31 #include "fvPatchFieldMapper.H"
32 #include "volFields.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 
40 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
41 
43 {
44  const label patchi = patch().index();
45 
46  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
47  (
49  (
51  internalField().group()
52  )
53  );
54 
55  const scalarField& y = turbModel.y()[patchi];
56 
57  const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
58  const vectorField Up(Uw.patchInternalField() - Uw);
59  const scalarField magUpn(mag(Up - (Up & patch().nf())*patch().nf()));
60 
61  const tmp<scalarField> tnuw = turbModel.nu(patchi);
62  const scalarField& nuw = tnuw();
63 
64  const scalar t = db().time().timeOutputValue();
65  const scalarField z0(z0_->value(t));
66 
67  #ifdef FULLDEBUG
68  for (const auto& z : z0)
69  {
70  if (z < VSMALL)
71  {
73  << "z0 field can only contain positive values. "
74  << "Please check input field z0."
75  << exit(FatalIOError);
76  }
77  }
78  #endif
79 
80  tmp<scalarField> tnutw(new scalarField(*this));
81  scalarField& nutw = tnutw.ref();
82 
83  forAll(nutw, facei)
84  {
85  const scalar Edash = (y[facei] + z0[facei])/(z0[facei] + 1e-4);
86  const scalar uStar = magUpn[facei]*kappa_/log(max(Edash, 1.0 + 1e-4));
87 
88  nutw[facei] = sqr(uStar)/max(magUpn[facei], 1e-6)*y[facei] - nuw[facei];
89  }
90 
91  if (boundNut_)
92  {
93  nutw = max(nutw, scalar(0.0));
94  }
95 
96  return tnutw;
97 }
98 
99 
100 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
101 
103 (
104  const fvPatch& p,
106 )
107 :
109  boundNut_(true),
110  z0_(nullptr)
111 {}
112 
113 
115 (
117  const fvPatch& p,
119  const fvPatchFieldMapper& mapper
120 )
121 :
122  nutUWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
123  boundNut_(ptf.boundNut_),
124  z0_(ptf.z0_.clone(p.patch()))
125 {}
126 
127 
129 (
130  const fvPatch& p,
132  const dictionary& dict
133 )
134 :
136  boundNut_(dict.getOrDefault<Switch>("boundNut", true)),
137  z0_(PatchFunction1<scalar>::New(p.patch(), "z0", dict))
138 {}
139 
140 
142 (
144 )
145 :
147  boundNut_(rwfpsf.boundNut_),
148  z0_(rwfpsf.z0_.clone(this->patch().patch()))
149 {}
150 
151 
153 (
156 )
157 :
159  boundNut_(rwfpsf.boundNut_),
160  z0_(rwfpsf.z0_.clone(this->patch().patch()))
161 {}
162 
163 
164 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
165 
167 (
168  const fvPatchFieldMapper& m
169 )
170 {
171  nutUWallFunctionFvPatchScalarField::autoMap(m);
172  z0_->autoMap(m);
173 }
174 
175 
177 (
178  const fvPatchScalarField& ptf,
179  const labelList& addr
180 )
181 {
182  nutUWallFunctionFvPatchScalarField::rmap(ptf, addr);
183 
185  refCast<const atmNutUWallFunctionFvPatchScalarField>(ptf);
186 
187  z0_->rmap(nrwfpsf.z0_(), addr);
188 }
189 
190 
192 {
194  os.writeEntry("boundNut", boundNut_);
195  z0_->writeData(os);
196  writeEntry("value", os);
197 }
198 
199 
200 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
201 
203 (
206 );
207 
208 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
209 
210 } // End namespace Foam
211 
212 // ************************************************************************* //
Foam::fvPatchField< vector >
volFields.H
Foam::atmNutUWallFunctionFvPatchScalarField::calcNut
virtual tmp< scalarField > calcNut() const
Calculate the turbulent viscosity.
Definition: atmNutUWallFunctionFvPatchScalarField.C:42
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.
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::atmNutUWallFunctionFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: atmNutUWallFunctionFvPatchScalarField.C:191
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::atmNutUWallFunctionFvPatchScalarField::atmNutUWallFunctionFvPatchScalarField
atmNutUWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: atmNutUWallFunctionFvPatchScalarField.C:103
atmNutUWallFunctionFvPatchScalarField.H
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
Foam::nutWallFunctionFvPatchScalarField::U
virtual const volVectorField & U(const turbulenceModel &turb) const
Definition: nutWallFunctionFvPatchScalarField.C:75
Foam::atmNutUWallFunctionFvPatchScalarField::rmap
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: atmNutUWallFunctionFvPatchScalarField.C:177
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::atmNutUWallFunctionFvPatchScalarField
This boundary condition provides a wall constraint on the turbulent viscosity, i.e....
Definition: atmNutUWallFunctionFvPatchScalarField.H:153
Foam::nutWallFunctionFvPatchScalarField::kappa_
scalar kappa_
von Kármán constant
Definition: nutWallFunctionFvPatchScalarField.H:290
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
Foam::Field< scalar >
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::fvPatchField::patchInternalField
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:225
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
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
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:372
Foam::nutUWallFunctionFvPatchScalarField
This boundary condition provides a wall constraint on the turbulent viscosity, i.e....
Definition: nutUWallFunctionFvPatchScalarField.H:98
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::atmNutUWallFunctionFvPatchScalarField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: atmNutUWallFunctionFvPatchScalarField.C:167
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::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::nutUWallFunctionFvPatchScalarField::write
virtual void write(Ostream &os) const
Write.
Definition: nutUWallFunctionFvPatchScalarField.C:202
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
turbulenceModel.H
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
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