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-------------------------------------------------------------------------------
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/>.
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
39namespace Foam
40{
41namespace 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 (
69 (
72 mesh_,
75 true
76 ),
77 mesh_,
80 );
81
82 result1->store();
83
84 result2 = new volScalarField
85 (
87 (
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 (
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// ************************************************************************* //
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
#define Log
Definition: PDRblock.C:35
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
virtual bool read()
Re-read model coefficients if they have changed.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const dimensionSet & dimensions() const
Return const reference to dimensions.
Abstract base-class for Time/database function objects.
Computes the Obukhov length field and associated friction velocity field.
scalar rhoRef_
Reference density (to convert from kinematic to static pressure)
scalar kappa_
von Kármán constant [-]
word UName_
Name of velocity field.
void removeObukhovLength()
Remove (checkOut) the output fields from the object registry.
word resultName2_
Name of the output field for Ustar.
const dimensionedVector g_
Gravitational acceleration vector [m/s2].
word resultName1_
Name of the output field for ObukhovLength.
bool calcOL()
Hard-coded Obukhov length field and friction velocity.
Definition: ObukhovLength.C:51
virtual bool execute()
Calculate the output fields.
virtual bool write()
Write the output fields.
dimensionedScalar beta_
Thermal expansion coefficient [1/K].
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
Computes the magnitude of the square of an input field.
Definition: magSqr.H:153
Computes the magnitude of an input field.
Definition: mag.H:153
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:162
const polyMesh & mesh() const
Return polyMesh.
Definition: mapPolyMesh.H:363
void movePoints()
Update for new mesh geometry.
void updateMesh()
Update for new mesh topology.
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
const Type & lookupObject(const word &name, const bool recursive=false) const
Type * getObjectPtr(const word &name, const bool recursive=false) const
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:76
virtual bool write(const bool valid=true) const
Write using setting from DB.
A class for handling words, derived from Foam::string.
Definition: word.H:68
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
const volScalarField & T
engineTime & runTime
scalar nut
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
const dimensionSet dimless
Dimensionless.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
dimensionedScalar sign(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
const dimensionSet dimVelocity
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:54
dimensionedScalar sqrt(const dimensionedScalar &ds)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
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
dictionary dict
volScalarField & e
Definition: createFields.H:11