ObukhovLength.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 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 
29 #include "ObukhovLength.H"
30 #include "volFields.H"
31 #include "dictionary.H"
32 #include "Time.H"
33 #include "mapPolyMesh.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41 namespace functionObjects
42 {
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
50 
52 {
53  const auto* rhoPtr = mesh_.findObject<volScalarField>("rho");
56  const volScalarField& alphat = mesh_.lookupObject<volScalarField>("alphat");
58 
61 
62  const bool isNew = !result1;
63 
64  if (!result1)
65  {
66  result1 = new volScalarField
67  (
68  IOobject
69  (
71  mesh_.time().timeName(),
72  mesh_,
75  true
76  ),
77  mesh_,
78  dimLength,
80  );
81 
82  result1->store();
83 
84  result2 = new volScalarField
85  (
86  IOobject
87  (
89  mesh_.time().timeName(),
90  mesh_,
93  true
94  ),
95  mesh_,
98  );
99 
100  result2->store();
101  }
102 
103  volScalarField B(alphat*beta_*(fvc::grad(T) & g_));
104  if (rhoPtr)
105  {
106  const auto& rho = *rhoPtr;
107  B /= rho;
108  }
109  else
110  {
112  B /= rho;
113  }
114 
115  *result2 = // Ustar
116  sqrt
117  (
118  max
119  (
120  nut*sqrt(2*magSqr(symm(fvc::grad(U)))),
121  dimensionedScalar(sqr(U.dimensions()), VSMALL)
122  )
123  );
124 
125  // (O:Eq. 26)
126  *result1 = // ObukhovLength
127  -min
128  (
129  dimensionedScalar(dimLength, ROOTVGREAT), // neutral stratification
130  pow3(*result2)/
131  (
132  sign(B)*kappa_
133  *max(mag(B), dimensionedScalar(B.dimensions(), VSMALL))
134  )
135  );
136 
137  return isNew;
138 }
139 
140 
141 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
142 
144 (
145  const word& name,
146  const Time& runTime,
147  const dictionary& dict
148 )
149 :
151  UName_("U"),
152  resultName1_("ObukhovLength"),
153  resultName2_("Ustar"),
154  rhoRef_(1.0),
155  kappa_(0.40),
156  beta_
157  (
159  (
161  dict.getCheckOrDefault<scalar>
162  (
163  "beta",
164  3e-3,
165  [&](const scalar x){ return x > SMALL; }
166  )
167  )
168  ),
169  g_
170  (
171  "g",
173  meshObjects::gravity::New(mesh_.time()).value()
174  )
175 {
176  read(dict);
177 }
178 
179 
180 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
181 
183 {
185 
186  UName_ = dict.getOrDefault<word>("U", "U");
187  resultName1_ = dict.getOrDefault<word>("ObukhovLength", "ObukhovLength");
188  resultName2_ = dict.getOrDefault<word>("Ustar", "Ustar");
189 
190  if (UName_ != "U" && resultName1_ == "ObukhovLength")
191  {
192  resultName1_ += '(' + UName_ + ')';
193  }
194 
195  if (UName_ != "U" && resultName1_ == "Ustar")
196  {
197  resultName2_ += '(' + UName_ + ')';
198  }
199 
200  rhoRef_ = dict.getOrDefault<scalar>("rhoRef", 1.0);
201  kappa_ = dict.getOrDefault<scalar>("kappa", 0.4);
202  beta_.value() = dict.getOrDefault<scalar>("beta", 3e-3);
203 
204  return true;
205 }
206 
207 
209 {
210  Log << type() << " " << name() << " execute:" << endl;
211 
212  bool isNew = false;
213 
214  isNew = calcOL();
215 
216  if (isNew) Log << " (new)" << nl << endl;
217 
218  return true;
219 }
220 
221 
223 {
224  const auto* ioptr1 = mesh_.cfindObject<regIOobject>(resultName1_);
225  const auto* ioptr2 = mesh_.cfindObject<regIOobject>(resultName2_);
226 
227  if (ioptr1)
228  {
229  Log << type() << " " << name() << " write:" << nl
230  << " writing field " << ioptr1->name() << nl
231  << " writing field " << ioptr2->name() << endl;
232 
233  ioptr1->write();
234  ioptr2->write();
235  }
236 
237  return true;
238 }
239 
240 
242 {
243  mesh_.thisDb().checkOut(resultName1_);
244  mesh_.thisDb().checkOut(resultName2_);
245 }
246 
247 
249 {
250  if (&mpm.mesh() == &mesh_)
251  {
252  removeObukhovLength();
253  }
254 }
255 
256 
258 {
259  if (&m == &mesh_)
260  {
261  removeObukhovLength();
262  }
263 }
264 
265 
266 // ************************************************************************* //
Foam::functionObjects::ObukhovLength::resultName1_
word resultName1_
Name of the output field for ObukhovLength.
Definition: ObukhovLength.H:309
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
volFields.H
Foam::objectRegistry::getObjectPtr
Type * getObjectPtr(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:423
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::functionObjects::ObukhovLength::calcOL
bool calcOL()
Hard-coded Obukhov length field and friction velocity.
Definition: ObukhovLength.C:51
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::symm
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:84
Log
#define Log
Definition: PDRblock.C:35
Foam::functionObjects::ObukhovLength::resultName2_
word resultName2_
Name of the output field for Ustar.
Definition: ObukhovLength.H:312
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
ObukhovLength.H
Foam::dimVelocity
const dimensionSet dimVelocity
mapPolyMesh.H
Foam::functionObjects::ObukhovLength::read
virtual bool read(const dictionary &dict)
Read the data.
Definition: ObukhovLength.C:182
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
B
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
Foam::functionObjects::ObukhovLength::rhoRef_
scalar rhoRef_
Reference density (to convert from kinematic to static pressure)
Definition: ObukhovLength.H:315
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::sign
dimensionedScalar sign(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:166
zeroGradientFvPatchField.H
Foam::functionObjects::fvMeshFunctionObject
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Definition: fvMeshFunctionObject.H:64
rho
rho
Definition: readInitialConditions.H:88
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::functionObject
Abstract base-class for Time/database function objects.
Definition: functionObject.H:332
Foam::zeroGradientFvPatchField
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
Definition: zeroGradientFvPatchField.H:64
Foam::functionObjects::ObukhovLength::UName_
word UName_
Name of velocity field.
Definition: ObukhovLength.H:306
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
rhoPtr
Info<< "Reading mechanical properties\n"<< endl;IOdictionary mechanicalProperties(IOobject("mechanicalProperties", runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE));const dictionary &rhoDict(mechanicalProperties.subDict("rho"));word rhoType(rhoDict.get< word >"type"));autoPtr< volScalarField > rhoPtr
Definition: readMechanicalProperties.H:18
Foam::functionObjects::ObukhovLength::removeObukhovLength
void removeObukhovLength()
Remove (checkOut) the output fields from the object registry.
Definition: ObukhovLength.C:241
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::functionObjects::ObukhovLength::ObukhovLength
ObukhovLength(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: ObukhovLength.C:144
Foam::regIOobject::write
virtual bool write(const bool valid=true) const
Write using setting from DB.
Definition: regIOobjectWrite.C:132
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::functionObjects::ObukhovLength::write
virtual bool write()
Write the output fields.
Definition: ObukhovLength.C:222
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::functionObjects::ObukhovLength::kappa_
scalar kappa_
von Kármán constant [-]
Definition: ObukhovLength.H:318
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:434
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
Foam::functionObjects::ObukhovLength::movePoints
virtual void movePoints(const polyMesh &m)
Update for mesh point-motion.
Definition: ObukhovLength.C:257
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::functionObjects::regionFunctionObject::read
virtual bool read(const dictionary &dict)
Read optional controls.
Definition: regionFunctionObject.C:173
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned< scalar >
Foam::Ostream::write
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::IOobject::name
const word & name() const noexcept
Return name.
Definition: IOobjectI.H:65
Foam::functionObjects::ObukhovLength::beta_
dimensionedScalar beta_
Thermal expansion coefficient [1/K].
Definition: ObukhovLength.H:321
Time.H
U
U
Definition: pEqn.H:72
Foam::regIOobject
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:73
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::objectRegistry::findObject
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
Definition: objectRegistryTemplates.C:401
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::functionObjects::ObukhovLength
Computes the Obukhov length field and associated friction velocity field.
Definition: ObukhovLength.H:297
Foam::functionObjects::addToRunTimeSelectionTable
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
Foam::functionObjects::ObukhovLength::execute
virtual bool execute()
Calculate the output fields.
Definition: ObukhovLength.C:208
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::functionObjects::fvMeshFunctionObject::mesh_
const fvMesh & mesh_
Reference to the fvMesh.
Definition: fvMeshFunctionObject.H:73
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dictionary.H
Foam::functionObjects::defineTypeNameAndDebug
defineTypeNameAndDebug(ObukhovLength, 0)
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::dimTemperature
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:161
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
nut
scalar nut
Definition: evaluateNearWall.H:5
Foam::dimensioned::dimensions
const dimensionSet & dimensions() const
Return const reference to dimensions.
Definition: dimensionedType.C:420
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::mapPolyMesh::mesh
const polyMesh & mesh() const
Return polyMesh.
Definition: mapPolyMesh.H:363
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::meshObjects::gravity::New
static const gravity & New(const Time &runTime)
Construct on Time.
Definition: gravityMeshObject.H:93
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
Foam::functionObjects::ObukhovLength::updateMesh
virtual void updateMesh(const mapPolyMesh &mpm)
Update for changes of mesh.
Definition: ObukhovLength.C:248
Foam::functionObjects::ObukhovLength::g_
const dimensionedVector g_
Gravitational acceleration vector [m/s2].
Definition: ObukhovLength.H:324