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 auto& turbModel =
47  db().lookupObject<turbulenceModel>
48  (
50  (
52  internalField().group()
53  )
54  );
55 
56  const scalarField& y = turbModel.y()[patchi];
57 
58  const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
59  const vectorField Up(Uw.patchInternalField() - Uw);
60  const scalarField magUpn(mag(Up - (Up & patch().nf())*patch().nf()));
61 
62  const tmp<scalarField> tnuw = turbModel.nu(patchi);
63  const scalarField& nuw = tnuw();
64 
65  const scalar t = db().time().timeOutputValue();
66  const scalarField z0(z0_->value(t));
67 
68  #ifdef FULLDEBUG
69  for (const scalar z : z0)
70  {
71  if (z < VSMALL)
72  {
74  << "z0 field can only contain positive values. "
75  << "Please check input field z0."
76  << exit(FatalError);
77  }
78  }
79  #endif
80 
81  auto tnutw = tmp<scalarField>::New(*this);
82  auto& nutw = tnutw.ref();
83 
84  forAll(nutw, facei)
85  {
86  const scalar Edash = (y[facei] + z0[facei])/(z0[facei] + 1e-4);
87  const scalar uStar = magUpn[facei]*kappa_/log(max(Edash, 1.0 + 1e-4));
88 
89  nutw[facei] = sqr(uStar)/max(magUpn[facei], 1e-6)*y[facei] - nuw[facei];
90  }
91 
92  if (boundNut_)
93  {
94  nutw = max(nutw, scalar(0));
95  }
96 
97  return tnutw;
98 }
99 
100 
101 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
102 
104 (
105  const fvPatch& p,
107 )
108 :
110  boundNut_(true),
111  z0_(nullptr)
112 {}
113 
114 
116 (
118  const fvPatch& p,
120  const fvPatchFieldMapper& mapper
121 )
122 :
123  nutUWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
124  boundNut_(ptf.boundNut_),
125  z0_(ptf.z0_.clone(p.patch()))
126 {}
127 
128 
130 (
131  const fvPatch& p,
133  const dictionary& dict
134 )
135 :
137  boundNut_(dict.getOrDefault<bool>("boundNut", true)),
138  z0_(PatchFunction1<scalar>::New(p.patch(), "z0", dict))
139 {}
140 
141 
143 (
145 )
146 :
148  boundNut_(rwfpsf.boundNut_),
149  z0_(rwfpsf.z0_.clone(this->patch().patch()))
150 {}
151 
152 
154 (
157 )
158 :
160  boundNut_(rwfpsf.boundNut_),
161  z0_(rwfpsf.z0_.clone(this->patch().patch()))
162 {}
163 
164 
165 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
166 
168 (
169  const fvPatchFieldMapper& m
170 )
171 {
172  nutUWallFunctionFvPatchScalarField::autoMap(m);
173  z0_->autoMap(m);
174 }
175 
176 
178 (
179  const fvPatchScalarField& ptf,
180  const labelList& addr
181 )
182 {
183  nutUWallFunctionFvPatchScalarField::rmap(ptf, addr);
184 
186  refCast<const atmNutUWallFunctionFvPatchScalarField>(ptf);
187 
188  z0_->rmap(nrwfpsf.z0_(), addr);
189 }
190 
191 
193 {
196  os.writeEntry("boundNut", boundNut_);
197  z0_->writeData(os);
198  writeEntry("value", os);
199 }
200 
201 
202 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
203 
205 (
208 );
209 
210 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
211 
212 } // End namespace Foam
213 
214 // ************************************************************************* //
Foam::fvPatchField< vector >
volFields.H
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
Foam::atmNutUWallFunctionFvPatchScalarField::calcNut
virtual tmp< scalarField > calcNut() const
Calculate the turbulent viscosity.
Definition: atmNutUWallFunctionFvPatchScalarField.C:42
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::atmNutUWallFunctionFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: atmNutUWallFunctionFvPatchScalarField.C:192
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::atmNutUWallFunctionFvPatchScalarField::atmNutUWallFunctionFvPatchScalarField
atmNutUWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: atmNutUWallFunctionFvPatchScalarField.C:104
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:178
fvPatchFieldMapper.H
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::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:233
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::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
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:453
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:168
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::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
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
turbulenceModel.H
Foam::nutWallFunctionFvPatchScalarField::writeLocalEntries
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
Definition: nutWallFunctionFvPatchScalarField.C:91
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