atmAlphatkWallFunctionFvPatchScalarField.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-2022 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 "wallFvPatch.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
45 
46 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
47 
49 {
50  if (!isA<wallFvPatch>(patch()))
51  {
53  << "Invalid wall function specification" << nl
54  << " Patch type for patch " << patch().name()
55  << " must be wall" << nl
56  << " Current patch type is " << patch().type() << nl << endl
57  << abort(FatalError);
58  }
59 }
60 
61 
63 (
64  Ostream& os
65 ) const
66 {
67  os.writeEntryIfDifferent<scalar>("Cmu", 0.09, Cmu_);
68  os.writeEntryIfDifferent<scalar>("kappa", 0.41, kappa_);
69 
70  if (Pr_)
71  {
72  Pr_->writeData(os);
73  }
74  if (Prt_)
75  {
76  Prt_->writeData(os);
77  }
78  if (z0_)
79  {
80  z0_->writeData(os);
81  }
82 }
83 
84 
85 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
86 
89 (
90  const fvPatch& p,
92 )
93 :
94  fixedValueFvPatchScalarField(p, iF),
95  Cmu_(0.09),
96  kappa_(0.41),
97  Pr_(nullptr),
98  Prt_(nullptr),
99  z0_(nullptr)
100 {
101  checkType();
102 }
103 
104 
107 (
109  const fvPatch& p,
111  const fvPatchFieldMapper& mapper
112 )
113 :
114  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
115  Cmu_(ptf.Cmu_),
116  kappa_(ptf.kappa_),
117  Pr_(ptf.Pr_.clone()),
118  Prt_(ptf.Prt_.clone(p.patch())),
119  z0_(ptf.z0_.clone(p.patch()))
120 {
121  checkType();
122 }
123 
124 
127 (
128  const fvPatch& p,
130  const dictionary& dict
131 )
132 :
133  fixedValueFvPatchScalarField(p, iF, dict),
134  Cmu_
135  (
136  dict.getCheckOrDefault<scalar>
137  (
138  "Cmu",
139  0.09,
140  scalarMinMax::ge(SMALL)
141  )
142  ),
143  kappa_
144  (
145  dict.getCheckOrDefault<scalar>
146  (
147  "kappa",
148  0.41,
149  scalarMinMax::ge(SMALL)
150  )
151  ),
152  Pr_(Function1<scalar>::New("Pr", dict, &db())),
153  Prt_(PatchFunction1<scalar>::New(p.patch(), "Prt", dict)),
154  z0_(PatchFunction1<scalar>::New(p.patch(), "z0", dict))
155 {
156  checkType();
157 }
158 
159 
162 (
164 )
165 :
166  fixedValueFvPatchScalarField(wfpsf),
167  Cmu_(wfpsf.Cmu_),
168  kappa_(wfpsf.kappa_),
169  Pr_(wfpsf.Pr_.clone()),
170  Prt_(wfpsf.Prt_.clone(this->patch().patch())),
171  z0_(wfpsf.z0_.clone(this->patch().patch()))
172 {
173  checkType();
174 }
175 
176 
179 (
182 )
183 :
184  fixedValueFvPatchScalarField(wfpsf, iF),
185  Cmu_(wfpsf.Cmu_),
186  kappa_(wfpsf.kappa_),
187  Pr_(wfpsf.Pr_.clone()),
188  Prt_(wfpsf.Prt_.clone(this->patch().patch())),
189  z0_(wfpsf.z0_.clone(this->patch().patch()))
190 {
191  checkType();
192 }
193 
194 
195 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
196 
198 {
199  if (updated())
200  {
201  return;
202  }
203 
204  const label patchi = patch().index();
205 
206  // Retrieve turbulence properties from model
207  const auto& turbModel = db().lookupObject<turbulenceModel>
208  (
210  (
212  internalField().group()
213  )
214  );
215 
216  const scalarField& y = turbModel.y()[patchi];
217 
218  const tmp<scalarField> tnuw = turbModel.nu(patchi);
219  const scalarField& nuw = tnuw();
220 
221  const tmp<volScalarField> tk = turbModel.k();
222  const volScalarField& k = tk();
223 
224  const scalar Cmu25 = pow025(Cmu_);
225 
226  const scalar t = db().time().timeOutputValue();
227  const scalar Pr = Pr_->value(t);
228 
229  #ifdef FULLDEBUG
230  if (Pr < VSMALL)
231  {
233  << "Pr cannot be negative or zero. "
234  << "Please check input Pr = " << Pr
235  << exit(FatalError);
236  }
237  #endif
238 
239  const scalarField Prt(Prt_->value(t));
240  const scalarField z0(z0_->value(t));
241 
242  #ifdef FULLDEBUG
243  forAll(Prt, i)
244  {
245  if (Prt[i] < VSMALL || z0[i] < VSMALL)
246  {
248  << "Elements of input surface fields can only be positive. "
249  << "Please check input fields z0 and Prt."
250  << exit(FatalError);
251  }
252  }
253  #endif
254 
255  const labelUList& faceCells = patch().faceCells();
256 
257  scalarField& alphatw = *this;
258 
259  forAll(alphatw, facei)
260  {
261  const label celli = faceCells[facei];
262 
263  const scalar uStar = Cmu25*Foam::sqrt(k[celli]);
264  const scalar Edash = (y[facei] + z0[facei])/(z0[facei] + 1e-4);
265 
266  // Update turbulent thermal conductivity
267  alphatw[facei] =
268  uStar*kappa_*y[facei]/(Prt[facei]*log(max(Edash, 1 + 1e-4)))
269  + nuw[facei]/Pr;
270  }
271 
272  // lower bound values to avoid unrealistic
273  // negative temperatures on the ground
274  alphatw = max(alphatw, scalar(0.01));
275 
277 }
278 
279 
281 (
282  const fvPatchFieldMapper& m
283 )
284 {
285  fixedValueFvPatchScalarField::autoMap(m);
286 
287  if (Prt_)
288  {
289  Prt_->autoMap(m);
290  }
291  if (z0_)
292  {
293  z0_->autoMap(m);
294  }
295 }
296 
297 
299 (
300  const fvPatchScalarField& ptf,
301  const labelList& addr
302 )
303 {
304  fixedValueFvPatchScalarField::rmap(ptf, addr);
305 
306  const auto& nrwfpsf =
307  refCast<const atmAlphatkWallFunctionFvPatchScalarField>(ptf);
308 
309  if (Prt_)
310  {
311  Prt_->rmap(nrwfpsf.Prt_(), addr);
312  }
313  if (z0_)
314  {
315  z0_->rmap(nrwfpsf.z0_(), addr);
316  }
317 }
318 
319 
321 {
325 }
326 
327 
328 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
329 
331 (
333  atmAlphatkWallFunctionFvPatchScalarField
334 );
335 
336 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
337 
338 } // End namespace Foam
339 
340 // ************************************************************************* //
dimensionedScalar Pr("Pr", dimless, laminarTransport)
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
dictionary dict
dimensionedScalar log(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
autoPtr< PatchFunction1< scalar > > Prt_
Turbulent Prandtl number field.
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
dimensionedScalar pow025(const dimensionedScalar &ds)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
This boundary condition provides a wall constraint on the kinematic turbulent thermal conductivity (i...
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
label k
Boltzmann constant.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
Macros for easy insertion into run-time selection tables.
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
autoPtr< PatchFunction1< scalar > > z0_
Surface roughness length [m].
atmAlphatkWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
constexpr const char *const group
Group name for atomic constants.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
scalar y
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
fvPatchField< scalar > fvPatchScalarField
static const word propertiesName
Default name of the turbulence properties dictionary.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
errorManip< error > abort(error &err)
Definition: errorManip.H:139
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
OBJstream os(runTime.globalPath()/outputName)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
dimensionedScalar Prt("Prt", dimless, laminarTransport)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:309
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
autoPtr< Function1< scalar > > Pr_
Molecular Prandtl number.
Namespace for OpenFOAM.
void writeLocalEntries(Ostream &) const
Write local wall function variables.