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 auto& turbModel =
48  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  auto tnutw = tmp<scalarField>::New(*this);
65  auto& nutw = tnutw.ref();
66 
67  const scalar Cmu25 = pow025(Cmu_);
68 
69  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
70  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
71 
72  const scalar t = db().time().timeOutputValue();
73  const scalarField z0(z0_->value(t));
74 
75  #ifdef FULLDEBUG
76  for (const scalar z : z0)
77  {
78  if (z < VSMALL)
79  {
81  << "z0 field can only contain positive values. "
82  << "Please check input field z0."
83  << exit(FatalError);
84  }
85  }
86  #endif
87 
88  const labelList& faceCells = patch().faceCells();
89 
90  forAll(nutw, facei)
91  {
92  const label celli = faceCells[facei];
93 
94  // (RH:Eq. 6)
95  const scalar Edash = (y[facei] + z0[facei])/(z0[facei] + z0Min_);
96 
97  // (RH:Eq. 6)
98  const scalar uStarU = magUp[facei]*kappa_/log(max(Edash, 1 + SMALL));
99 
100  // (RH:Eq. 7)
101  const scalar uStarK = Cmu25*Foam::sqrt(k[celli]);
102 
103  // (SBJM:Eq. 7; SM:Eq. 25)
104  const scalar tauw = uStarU*uStarK;
105 
106  nutw[facei] =
107  max(tauw*y[facei]/(max(magUp[facei], SMALL)) - nuw[facei], 0.0);
108  }
109 
110  return tnutw;
111 }
112 
113 
114 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
115 
117 (
118  const fvPatch& p,
120 )
121 :
123  z0Min_(SMALL),
124  z0_(nullptr)
125 {}
126 
127 
129 (
131  const fvPatch& p,
133  const fvPatchFieldMapper& mapper
134 )
135 :
136  nutkWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
137  z0Min_(ptf.z0Min_),
138  z0_(ptf.z0_.clone(p.patch()))
139 {}
140 
141 
143 (
144  const fvPatch& p,
146  const dictionary& dict
147 )
148 :
150  z0Min_
151  (
152  dict.getCheckOrDefault<scalar>
153  (
154  "z0Min",
155  SMALL,
157  )
158  ),
159  z0_(PatchFunction1<scalar>::New(p.patch(), "z0", dict))
160 {}
161 
162 
164 (
166 )
167 :
169  z0Min_(rwfpsf.z0Min_),
170  z0_(rwfpsf.z0_.clone(this->patch().patch()))
171 {}
172 
173 
175 (
178 )
179 :
181  z0Min_(rwfpsf.z0Min_),
182  z0_(rwfpsf.z0_.clone(this->patch().patch()))
183 {}
184 
185 
186 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
187 
189 (
190  const fvPatchFieldMapper& m
191 )
192 {
193  nutkWallFunctionFvPatchScalarField::autoMap(m);
194  z0_->autoMap(m);
195 }
196 
197 
199 (
200  const fvPatchScalarField& ptf,
201  const labelList& addr
202 )
203 {
204  nutkWallFunctionFvPatchScalarField::rmap(ptf, addr);
205 
206  const atmNutWallFunctionFvPatchScalarField& nrwfpsf =
207  refCast<const atmNutWallFunctionFvPatchScalarField>(ptf);
208 
209  z0_->rmap(nrwfpsf.z0_(), addr);
210 }
211 
212 
214 {
217  os.writeEntry("z0Min", z0Min_);
218  z0_->writeData(os);
219  writeEntry("value", os);
220 }
221 
222 
223 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
224 
226 (
229 );
230 
231 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232 
233 } // End namespace Foam
234 
235 // ************************************************************************* //
Foam::fvPatchField< vector >
Foam::atmNutWallFunctionFvPatchScalarField::rmap
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: atmNutWallFunctionFvPatchScalarField.C:199
volFields.H
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
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:61
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
fvPatchFieldMapper.H
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::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::fvPatchField::patchInternalField
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:233
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:209
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:117
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
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::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::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::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:213
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:189
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