atmBoundaryLayer.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2014-2021 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "atmBoundaryLayer.H"
27 
28 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
29 
30 const Foam::scalar Foam::atmBoundaryLayer::kappaDefault_ = 0.41;
31 
32 const Foam::scalar Foam::atmBoundaryLayer::CmuDefault_ = 0.09;
33 
34 
35 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36 
37 void Foam::atmBoundaryLayer::init()
38 {
39  if (mag(flowDir_) < small || mag(zDir_) < small)
40  {
42  << "magnitude of n or z must be greater than zero"
43  << abort(FatalError);
44  }
45 
46  // Ensure direction vectors are normalised
47  flowDir_ /= mag(flowDir_);
48  zDir_ /= mag(zDir_);
49 
50  Ustar_ = kappa_*Uref_/(log((Zref_ + z0_)/z0_));
51 }
52 
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
57 :
58  flowDir_(Zero),
59  zDir_(Zero),
60  kappa_(0.41),
61  Cmu_(0.09),
62  Uref_(0),
63  Zref_(0),
64  z0_(0),
65  zGround_(0),
66  Ustar_(0),
67  offset_(false),
68  Ulower_(0),
69  kLower_(0),
70  epsilonLower_(0)
71 {}
72 
73 
75 (
76  const vector& flowDir,
77  const vector& zDir,
78  const scalar Uref,
79  const scalar Zref,
80  const scalarField& z0,
81  const scalarField& zGround,
82  const scalar kappa,
83  const scalar Cmu,
84  const scalar Ulower,
85  const scalar kLower,
86  const scalar epsilonLower
87 )
88 :
89  flowDir_(flowDir),
90  zDir_(zDir),
91  kappa_(kappa),
92  Cmu_(Cmu),
93  Uref_(Uref),
94  Zref_(Zref),
95  z0_(z0),
96  zGround_(zGround),
97  Ustar_(z0.size()),
98  offset_(Ulower != 0),
99  Ulower_(Ulower),
100  kLower_(kLower),
101  epsilonLower_(epsilonLower)
102 {
103  init();
104 }
105 
106 
108 (
109  const vectorField& p,
110  const dictionary& dict
111 )
112 :
113  flowDir_(dict.lookup("flowDir")),
114  zDir_(dict.lookup("zDir")),
115  kappa_(dict.lookupOrDefault<scalar>("kappa", kappaDefault_)),
116  Cmu_(dict.lookupOrDefault<scalar>("Cmu", CmuDefault_)),
117  Uref_(dict.lookup<scalar>("Uref")),
118  Zref_(dict.lookup<scalar>("Zref")),
119  z0_("z0", dict, p.size()),
120  zGround_("zGround", dict, p.size()),
121  Ustar_(p.size()),
122  offset_(dict.found("Ulower")),
123  Ulower_(dict.lookupOrDefault<scalar>("Ulower", 0)),
124  kLower_(dict.lookupOrDefault<scalar>("kLower", 0)),
125  epsilonLower_(dict.lookupOrDefault<scalar>("epsilonLower", 0))
126 {
127  init();
128 }
129 
130 
132 (
133  const atmBoundaryLayer& abl,
134  const fvPatchFieldMapper& mapper
135 )
136 :
137  flowDir_(abl.flowDir_),
138  zDir_(abl.zDir_),
139  kappa_(abl.kappa_),
140  Cmu_(abl.Cmu_),
141  Uref_(abl.Uref_),
142  Zref_(abl.Zref_),
143  z0_(mapper(abl.z0_)),
144  zGround_(mapper(abl.zGround_)),
145  Ustar_(mapper(abl.Ustar_)),
146  offset_(abl.offset_),
147  Ulower_(abl.Ulower_),
148  kLower_(abl.kLower_),
149  epsilonLower_(abl.epsilonLower_)
150 {}
151 
152 
154 :
155  flowDir_(abl.flowDir_),
156  zDir_(abl.zDir_),
157  kappa_(abl.kappa_),
158  Cmu_(abl.Cmu_),
159  Uref_(abl.Uref_),
160  Zref_(abl.Zref_),
161  z0_(abl.z0_),
162  zGround_(abl.zGround_),
163  Ustar_(abl.Ustar_),
164  offset_(abl.offset_),
165  Ulower_(abl.Ulower_),
166  kLower_(abl.kLower_),
167  epsilonLower_(abl.epsilonLower_)
168 {}
169 
170 
171 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
172 
174 {
175  m(z0_, z0_);
176  m(zGround_, zGround_);
177  m(Ustar_, Ustar_);
178 }
179 
180 
182 (
183  const atmBoundaryLayer& blptf,
184  const labelList& addr
185 )
186 {
187  z0_.rmap(blptf.z0_, addr);
188  zGround_.rmap(blptf.zGround_, addr);
189  Ustar_.rmap(blptf.Ustar_, addr);
190 }
191 
192 
194 (
195  const vectorField& p
196 ) const
197 {
198  const scalarField Un
199  (
200  (Ustar_/kappa_)
201  *log(max((zDir_ & p) - zGround_ + z0_, z0_)/z0_)
202  );
203 
204  if (offset_)
205  {
206  return flowDir_*Un + flowDir_*Ulower_;
207  }
208  else
209  {
210  return flowDir_*Un;
211  }
212 }
213 
214 
216 (
217  const vectorField& p
218 ) const
219 {
221  (
222  sqr(Ustar_)/sqrt(Cmu_)
223  );
224 
225  if (offset_)
226  {
227  const scalarField z((zDir_ & p) - zGround_);
228  tk.ref() = pos0(z)*tk() + neg(z)*kLower_;
229  }
230 
231  return tk;
232 }
233 
234 
236 (
237  const vectorField& p
238 ) const
239 {
240  tmp<scalarField> tepsilon
241  (
242  pow3(Ustar_)/(kappa_*((zDir_ & p) - zGround_ + z0_))
243  );
244 
245  if (offset_)
246  {
247  const scalarField z((zDir_ & p) - zGround_);
248  tepsilon.ref() = pos0(z)*tepsilon() + neg(z)*epsilonLower_;
249  }
250 
251  return tepsilon;
252 }
253 
254 
256 {
257  writeEntry(os, "z0", z0_) ;
258  writeEntry(os, "flowDir", flowDir_);
259  writeEntry(os, "zDir", zDir_);
260  writeEntry(os, "kappa", kappa_);
261  writeEntry(os, "Cmu", Cmu_);
262  writeEntry(os, "Uref", Uref_);
263  writeEntry(os, "Zref", Zref_);
264 
265  if (offset_)
266  {
267  writeEntry(os, "Ulower", Ulower_);
268  writeEntry(os, "kLower", kLower_);
269  writeEntry(os, "epsilonLower", epsilonLower_);
270  }
271 
272  writeEntry(os, "zGround", zGround_);
273 }
274 
275 
276 // ************************************************************************* //
This class provides functions to evaluate the velocity and turbulence distributions appropriate for a...
tmp< scalarField > epsilon(const vectorField &p) const
Return the turbulent dissipation rate distribution for the ATM.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:643
dimensionedScalar log(const dimensionedScalar &ds)
void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
dimensionedScalar sqrt(const dimensionedScalar &ds)
tmp< vectorField > U(const vectorField &p) const
Return the velocity distribution for the ATM.
tmp< scalarField > k(const vectorField &p) const
Return the turbulent kinetic energy distribution for the ATM.
dimensionedScalar neg(const dimensionedScalar &ds)
Foam::fvPatchFieldMapper.
void write(Ostream &) const
Write.
static const zero Zero
Definition: zero.H:97
errorManip< error > abort(error &err)
Definition: errorManip.H:131
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
dimensionedScalar pos0(const dimensionedScalar &ds)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar pow3(const dimensionedScalar &ds)
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition: Field.C:362
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
void rmap(const atmBoundaryLayer &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
atmBoundaryLayer()
Construct null.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844