atmNutWallFunctionFvPatchScalarField.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 CENER
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"
33 #include "bound.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 
41 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
42 
44 {
45  const label patchi = patch().index();
46 
47  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
48  (
50  (
52  internalField().group()
53  )
54  );
55  const scalarField& y = turbModel.y()[patchi];
56 
57  const tmp<volScalarField> tk = turbModel.k();
58  const volScalarField& k = tk();
59 
60  const tmp<scalarField> tnuw = turbModel.nu(patchi);
61  const scalarField& nuw = tnuw();
62 
63  tmp<scalarField> tnutw(new scalarField(*this));
64  scalarField& nutw = tnutw.ref();
65 
66  const scalar Cmu25 = pow025(Cmu_);
67 
68  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
69  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
70 
71  const scalar t = db().time().timeOutputValue();
72  const scalarField z0(z0_->value(t));
73 
74  #ifdef FULLDEBUG
75  for (const auto& z : z0)
76  {
77  if (z < VSMALL)
78  {
80  << "z0 field can only contain positive values. "
81  << "Please check input field z0."
82  << exit(FatalIOError);
83  }
84  }
85  #endif
86 
87  const labelList& faceCells = patch().faceCells();
88 
89  forAll(nutw, facei)
90  {
91  const label celli = faceCells[facei];
92 
93  // (RH:Eq. 6)
94  const scalar Edash = (y[facei] + z0[facei])/(z0[facei] + z0Min_);
95 
96  // (RH:Eq. 6)
97  const scalar uStarU = magUp[facei]*kappa_/log(max(Edash, 1 + SMALL));
98 
99  // (RH:Eq. 7)
100  const scalar uStarK = Cmu25*Foam::sqrt(k[celli]);
101 
102  // (SBJM:Eq. 7; SM:Eq. 25)
103  const scalar tauw = uStarU*uStarK;
104 
105  nutw[facei] =
106  max(tauw*y[facei]/(max(magUp[facei], SMALL)) - nuw[facei], 0.0);
107  }
108 
109  return tnutw;
110 }
111 
112 
113 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
114 
116 (
117  const fvPatch& p,
119 )
120 :
122  z0Min_(SMALL),
123  z0_(nullptr)
124 {}
125 
126 
128 (
130  const fvPatch& p,
132  const fvPatchFieldMapper& mapper
133 )
134 :
135  nutkWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
136  z0Min_(ptf.z0Min_),
137  z0_(ptf.z0_.clone(p.patch()))
138 {}
139 
140 
142 (
143  const fvPatch& p,
145  const dictionary& dict
146 )
147 :
149  z0Min_
150  (
151  dict.getCheckOrDefault<scalar>
152  (
153  "z0Min",
154  SMALL,
156  )
157  ),
158  z0_(PatchFunction1<scalar>::New(p.patch(), "z0", dict))
159 {}
160 
161 
163 (
165 )
166 :
168  z0Min_(rwfpsf.z0Min_),
169  z0_(rwfpsf.z0_.clone(this->patch().patch()))
170 {}
171 
172 
174 (
177 )
178 :
180  z0Min_(rwfpsf.z0Min_),
181  z0_(rwfpsf.z0_.clone(this->patch().patch()))
182 {}
183 
184 
185 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
186 
188 (
189  const fvPatchFieldMapper& m
190 )
191 {
192  nutkWallFunctionFvPatchScalarField::autoMap(m);
193  z0_->autoMap(m);
194 }
195 
196 
198 (
199  const fvPatchScalarField& ptf,
200  const labelList& addr
201 )
202 {
203  nutkWallFunctionFvPatchScalarField::rmap(ptf, addr);
204 
205  const atmNutWallFunctionFvPatchScalarField& nrwfpsf =
206  refCast<const atmNutWallFunctionFvPatchScalarField>(ptf);
207 
208  z0_->rmap(nrwfpsf.z0_(), addr);
209 }
210 
211 
213 {
215  os.writeEntry("z0Min", z0Min_);
216  z0_->writeData(os);
217  writeEntry("value", os);
218 }
219 
220 
221 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
222 
224 (
227 );
228 
229 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
230 
231 } // End namespace Foam
232 
233 // ************************************************************************* //
Foam::fvPatchField< vector >
Foam::atmNutWallFunctionFvPatchScalarField::rmap
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: atmNutWallFunctionFvPatchScalarField.C:198
volFields.H
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::atmNutWallFunctionFvPatchScalarField::calcNut
virtual tmp< scalarField > calcNut() const
Calculate the turbulent viscosity.
Definition: atmNutWallFunctionFvPatchScalarField.C:43
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::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
magUp
scalar magUp
Definition: evaluateNearWall.H:10
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
atmNutWallFunctionFvPatchScalarField.H
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::fvPatchField::patchInternalField
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:225
tauw
scalar tauw
Definition: evaluateNearWall.H:12
Foam::dictionary::getCheckOrDefault
T getCheckOrDefault(const word &keyword, const T &deflt, const Predicate &pred, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:203
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::atmNutWallFunctionFvPatchScalarField::atmNutWallFunctionFvPatchScalarField
atmNutWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: atmNutWallFunctionFvPatchScalarField.C:116
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
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::MinMax::ge
static MinMax< T > ge(const T &minVal)
A semi-infinite range from minVal to the type max.
Definition: MinMaxI.H:31
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
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
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::atmNutWallFunctionFvPatchScalarField
This boundary condition provides a wall constraint on the turbulent viscosity, i.e....
Definition: atmNutWallFunctionFvPatchScalarField.H:211
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::atmNutWallFunctionFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: atmNutWallFunctionFvPatchScalarField.C:212
Foam::faceCells
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:56
turbulenceModel.H
Foam::atmNutWallFunctionFvPatchScalarField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: atmNutWallFunctionFvPatchScalarField.C:188
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
Foam::turbulenceModel::U
const volVectorField & U() const
Access function to velocity field.
Definition: turbulenceModel.H:144
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