Loading...
Searching...
No Matches
atmBoundaryLayer.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) 2014-2016 OpenFOAM Foundation
9 Copyright (C) 2018-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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/>.
27\*---------------------------------------------------------------------------*/
28
29#include "atmBoundaryLayer.H"
30
31// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32
33namespace Foam
34{
35
36// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37
39:
40 initABL_(false),
41 kappa_(0.41),
42 Cmu_(0.09),
43 C1_(0.0),
44 C2_(1.0),
45 ppMin_((boundBox(pp.localPoints())).min()),
46 time_(time),
47 patch_(pp),
48 flowDir_(nullptr),
49 zDir_(nullptr),
50 Uref_(nullptr),
51 Zref_(nullptr),
52 z0_(nullptr),
53 d_(nullptr)
54{}
55
56
58(
59 const Time& time,
60 const polyPatch& pp,
61 const dictionary& dict
62)
63:
64 initABL_(dict.getOrDefault<bool>("initABL", true)),
65 kappa_
66 (
67 dict.getCheckOrDefault<scalar>("kappa", 0.41, scalarMinMax::ge(SMALL))
68 ),
69 Cmu_(dict.getCheckOrDefault<scalar>("Cmu", 0.09, scalarMinMax::ge(SMALL))),
70 C1_(dict.getOrDefault("C1", 0.0)),
71 C2_(dict.getOrDefault("C2", 1.0)),
72 ppMin_((boundBox(pp.localPoints())).min()),
73 time_(time),
74 patch_(pp),
75 flowDir_(Function1<vector>::New("flowDir", dict, &time)),
76 zDir_(Function1<vector>::New("zDir", dict, &time)),
77 Uref_(Function1<scalar>::New("Uref", dict, &time)),
78 Zref_(Function1<scalar>::New("Zref", dict, &time)),
79 z0_(PatchFunction1<scalar>::New(pp, "z0", dict)),
80 d_(PatchFunction1<scalar>::New(pp, "d", dict))
81{}
82
83
85(
86 const atmBoundaryLayer& abl,
87 const fvPatch& patch,
88 const fvPatchFieldMapper& mapper
89)
90:
91 initABL_(abl.initABL_),
92 kappa_(abl.kappa_),
93 Cmu_(abl.Cmu_),
94 C1_(abl.C1_),
95 C2_(abl.C2_),
96 ppMin_(abl.ppMin_),
97 time_(abl.time_),
98 patch_(patch.patch()),
99 flowDir_(abl.flowDir_.clone()),
100 zDir_(abl.zDir_.clone()),
101 Uref_(abl.Uref_.clone()),
102 Zref_(abl.Zref_.clone()),
103 z0_(abl.z0_.clone(patch_)),
104 d_(abl.d_.clone(patch_))
105{}
106
107
109:
110 initABL_(abl.initABL_),
111 kappa_(abl.kappa_),
112 Cmu_(abl.Cmu_),
113 C1_(abl.C1_),
114 C2_(abl.C2_),
115 ppMin_(abl.ppMin_),
116 time_(abl.time_),
117 patch_(abl.patch_),
118 flowDir_(abl.flowDir_.clone()),
119 zDir_(abl.zDir_.clone()),
120 Uref_(abl.Uref_.clone()),
121 Zref_(abl.Zref_.clone()),
122 z0_(abl.z0_.clone(patch_)),
123 d_(abl.d_.clone(patch_))
124{}
125
126
127// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
128
130{
131 const scalar t = time_.timeOutputValue();
132 const vector dir(flowDir_->value(t));
133 const scalar magDir = mag(dir);
134
135 if (magDir < SMALL)
136 {
138 << "magnitude of " << flowDir_->name() << " = " << magDir
139 << " vector must be greater than zero"
141 }
142
143 return dir/magDir;
144}
145
146
148{
149 const scalar t = time_.timeOutputValue();
150 const vector dir(zDir_->value(t));
151 const scalar magDir = mag(dir);
152
153 if (magDir < SMALL)
154 {
156 << "magnitude of " << zDir_->name() << " = " << magDir
157 << " vector must be greater than zero"
159 }
160
161 return dir/magDir;
162}
163
164
165tmp<scalarField> atmBoundaryLayer::Ustar(const scalarField& z0) const
166{
167 const scalar t = time_.timeOutputValue();
168 const scalar Uref = Uref_->value(t);
169 const scalar Zref = Zref_->value(t);
170
171 if (Zref < 0)
172 {
174 << "Negative entry in " << Zref_->name() << " = " << Zref
175 << abort(FatalError);
177
178 // (derived from RH:Eq. 6, HW:Eq. 5)
179 return kappa_*Uref/log((Zref + z0)/z0);
180}
181
182
183void atmBoundaryLayer::autoMap(const fvPatchFieldMapper& mapper)
184{
185 if (z0_)
186 {
187 z0_->autoMap(mapper);
188 }
189 if (d_)
190 {
191 d_->autoMap(mapper);
192 }
193}
194
195
197(
198 const atmBoundaryLayer& abl,
199 const labelList& addr
200)
201{
202 if (z0_)
203 {
204 z0_->rmap(abl.z0_(), addr);
205 }
206 if (d_)
207 {
208 d_->rmap(abl.d_(), addr);
209 }
210}
211
212
213tmp<vectorField> atmBoundaryLayer::U(const vectorField& pCf) const
214{
215 const scalar t = time_.timeOutputValue();
216 const scalarField d(d_->value(t));
217 const scalarField z0(max(z0_->value(t), ROOTVSMALL));
218 const scalar groundMin = zDir() & ppMin_;
219
220 // (YGCJ:Table 1, RH:Eq. 6, HW:Eq. 5)
221 scalarField zEff(max((zDir() & pCf) - groundMin - d + z0, z0));
222
223 scalarField Un
224 (
225 (Ustar(z0)/kappa_)*log(zEff/z0)
226 );
227
228 return flowDir()*Un;
229}
230
231
233{
234 const scalar t = time_.timeOutputValue();
235 const scalarField d(d_->value(t));
236 const scalarField z0(max(z0_->value(t), ROOTVSMALL));
237 const scalar groundMin = zDir() & ppMin_;
238
239 // (YGCJ:Eq. 21; RH:Eq. 7, HW:Eq. 6 when C1=0 and C2=1)
240 scalarField zEff(max((zDir() & pCf) - groundMin - d + z0, z0));
241
242 return sqr(Ustar(z0))/sqrt(Cmu_)*sqrt(C1_*log(zEff/z0) + C2_);
243}
244
245
247{
248 const scalar t = time_.timeOutputValue();
249 const scalarField d(d_->value(t));
250 const scalarField z0(max(z0_->value(t), ROOTVSMALL));
251 const scalar groundMin = zDir() & ppMin_;
252
253 // (YGCJ:Eq. 22; RH:Eq. 8, HW:Eq. 7 when C1=0 and C2=1)
254 scalarField zEff(max((zDir() & pCf) - groundMin - d + z0, z0));
255
256 return pow3(Ustar(z0))/(kappa_*zEff)*sqrt(C1_*log(zEff/z0) + C2_);
257}
258
259
261{
262 const scalar t = time_.timeOutputValue();
263 const scalarField d(d_->value(t));
264 const scalarField z0(max(z0_->value(t), ROOTVSMALL));
265 const scalar groundMin = zDir() & ppMin_;
266
267 // (YGJ:Eq. 13)
268 scalarField zEff(max((zDir() & pCf) - groundMin - d + z0, z0));
269
270 return Ustar(z0)/(kappa_*sqrt(Cmu_)*zEff);
271}
272
273
275{
276 os.writeEntry("initABL", initABL_);
277 os.writeEntry("kappa", kappa_);
278 os.writeEntry("Cmu", Cmu_);
279 os.writeEntry("C1", C1_);
280 os.writeEntry("C2", C2_);
281 if (flowDir_)
282 {
283 flowDir_->writeData(os);
284 }
285 if (zDir_)
286 {
287 zDir_->writeData(os);
288 }
289 if (Uref_)
290 {
291 Uref_->writeData(os);
292 }
293 if (Zref_)
294 {
295 Zref_->writeData(os);
296 }
297 if (z0_)
298 {
299 z0_->writeData(os) ;
300 }
301 if (d_)
302 {
303 d_->writeData(os);
304 }
305}
306
307
308// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
309
310} // End namespace Foam
311
312// ************************************************************************* //
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition Function1.H:92
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
scalar timeOutputValue() const
Return the current user-time value. (ie, after applying any timeToUserTime() conversion).
Definition TimeStateI.H:24
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
Base class to set log-law type ground-normal inlet boundary conditions for wind velocity and turbulen...
void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
bool initABL_
Flag to initialise profiles with the theoretical ABL expressions, otherwise initialises by using "val...
void write(Ostream &) const
Write.
void rmap(const atmBoundaryLayer &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
tmp< scalarField > Ustar(const scalarField &z0) const
Return friction velocity.
tmp< scalarField > epsilon(const vectorField &pCf) const
Return the turbulent dissipation rate distribution for the ATM.
atmBoundaryLayer(const Time &time, const polyPatch &pp)
Construct null from time.
tmp< scalarField > omega(const vectorField &pCf) const
Return the specific dissipation rate distribution for the ATM.
tmp< vectorField > U(const vectorField &pCf) const
Return the velocity distribution for the ATM.
vector flowDir() const
Return flow direction.
tmp< scalarField > k(const vectorField &pCf) const
Return the turbulent kinetic energy distribution for the ATM.
vector zDir() const
Return the ground-normal direction.
A bounding box defined in terms of min/max extrema points.
Definition boundBox.H:71
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A FieldMapper for finite-volume patch fields.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
A class for managing temporary objects.
Definition tmp.H:75
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
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.
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition MinMax.H:97
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
Vector< scalar > vector
Definition vector.H:57
dictionary dict