Loading...
Searching...
No Matches
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,2023 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"
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
41namespace Foam
42{
43namespace functionObjects
44{
47}
48}
49
50
51// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
52
54{
55 const auto* rhoPtr = mesh_.findObject<volScalarField>("rho");
58 const volScalarField& alphat = mesh_.lookupObject<volScalarField>("alphat");
60
63
64 const bool isNew = !result1;
65
66 if (!result1)
67 {
68 result1 = new volScalarField
69 (
71 (
74 mesh_,
78 ),
79 mesh_,
82 );
83
84 result1->store();
85 }
86
87 if (!result2)
88 {
89 result2 = new volScalarField
90 (
91 IOobject
92 (
95 mesh_,
99 ),
100 mesh_,
103 );
104
105 result2->store();
106 }
107
108 volScalarField B(alphat*beta_*(fvc::grad(T) & g_));
109 if (rhoPtr)
110 {
111 const auto& rho = *rhoPtr;
112 B /= rho;
113 }
114 else
115 {
117 B /= rho;
118 }
119
120 *result2 = // Ustar
121 sqrt
122 (
123 max
124 (
126 dimensionedScalar(sqr(U.dimensions()), VSMALL)
127 )
128 );
129
130
131 // Special bit of coding here to handle cyclics
132 // - sign(B) can easily be positive in one cell but negative in its
133 // coupled cell
134 // - so cyclic value might evaluate to 0
135 // - which upsets the division
136 // - note that instead of sign() any other field might have the same
137 // positive and negative coupled values - it is just a lot less likely
138 // - problem is that overridden patch types do not propagate - use in-place
139 // operations only
140
141 volScalarField denom
142 (
143 IOobject
144 (
145 "denom",
146 mesh_.time().timeName(),
147 mesh_,
151 ),
152 sign(B)
153 );
154
155 // Override interpolated value on interpolated coupled patches
156 for (auto& pfld : denom.boundaryFieldRef())
157 {
158 if
159 (
160 isA<coupledFvPatchField<scalar>>(pfld)
161 && !isA<processorFvPatchField<scalar>>(pfld)
162 )
163 {
164 pfld = pfld.patchInternalField();
165 }
166 }
167
168 denom *= kappa_*max(mag(B), dimensionedScalar(B.dimensions(), VSMALL));
169
170 // (O:Eq. 26)
171 *result1 = // ObukhovLength
172 -min
173 (
174 dimensionedScalar(dimLength, ROOTVGREAT), // neutral stratification
175 pow3(*result2)/denom
176 );
178 return isNew;
179}
180
181
182// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
183
185(
186 const word& name,
187 const Time& runTime,
188 const dictionary& dict
189)
190:
192 UName_("U"),
193 resultName1_("ObukhovLength"),
194 resultName2_("Ustar"),
195 rhoRef_(1.0),
196 kappa_(0.40),
197 beta_
198 (
200 (
202 dict.getCheckOrDefault<scalar>
203 (
204 "beta",
205 3e-3,
206 [&](const scalar x){ return x > SMALL; }
207 )
208 )
209 ),
210 g_
211 (
212 "g",
214 meshObjects::gravity::New(mesh_.time()).value()
215 )
217 read(dict);
218}
219
220
221// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
222
224{
226
227 UName_ = dict.getOrDefault<word>("U", "U");
228 resultName1_ = dict.getOrDefault<word>("ObukhovLength", "ObukhovLength");
229 resultName2_ = dict.getOrDefault<word>("Ustar", "Ustar");
230
231 if (UName_ != "U" && resultName1_ == "ObukhovLength")
232 {
233 resultName1_ += '(' + UName_ + ')';
234 }
235
236 if (UName_ != "U" && resultName1_ == "Ustar")
237 {
238 resultName2_ += '(' + UName_ + ')';
239 }
240
241 rhoRef_ = dict.getOrDefault<scalar>("rhoRef", 1.0);
242 kappa_ = dict.getOrDefault<scalar>("kappa", 0.4);
243 beta_.value() = dict.getOrDefault<scalar>("beta", 3e-3);
244
245 return true;
246}
247
248
250{
251 Log << type() << " " << name() << " execute:" << endl;
252
253 bool isNew = false;
254
255 isNew = calcOL();
257 if (isNew) Log << " (new)" << nl << endl;
258
259 return true;
260}
261
262
264{
265 const auto* ioptr1 = mesh_.cfindObject<regIOobject>(resultName1_);
266 const auto* ioptr2 = mesh_.cfindObject<regIOobject>(resultName2_);
267
268 if (ioptr1)
269 {
270 Log << type() << " " << name() << " write:" << nl
271 << " writing field " << ioptr1->name() << nl
272 << " writing field " << ioptr2->name() << endl;
273
274 ioptr1->write();
275 ioptr2->write();
276 }
277
278 return true;
279}
280
281
283{
284 mesh_.thisDb().checkOut(resultName1_);
285 mesh_.thisDb().checkOut(resultName2_);
286}
287
288
290{
291 if (&mpm.mesh() == &mesh_)
292 {
294 }
295}
296
297
299{
300 if (&m == &mesh_)
301 {
302 removeObukhovLength();
303 }
304}
305
306
307// ************************************************************************* //
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
#define Log
Definition PDRblock.C:28
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
@ NO_REGISTER
Do not request registration (bool: false).
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
Abstract base class for coupled patches.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const dimensionSet & dimensions() const noexcept
Return const reference to dimensions.
Abstract base-class for Time/database function objects.
virtual bool read(const dictionary &dict)
Read and set the function object if its data have changed.
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.
virtual void movePoints(const polyMesh &m)
Update for mesh point-motion.
virtual bool read(const dictionary &dict)
Read the data.
ObukhovLength(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
virtual void updateMesh(const mapPolyMesh &mpm)
Update for changes of mesh.
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.
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.
fvMeshFunctionObject(const fvMeshFunctionObject &)=delete
No copy construct.
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const polyMesh & mesh() const noexcept
Return polyMesh.
static const gravity & New(const word &name, const Time &runTime)
Return named gravity field cached or construct on Time.
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
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
Type * getObjectPtr(const word &name, const bool recursive=false) const
Return non-const pointer to the object of the given Type, using a const-cast to have it behave like a...
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
This boundary condition enables processor communication across patches.
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition regIOobject.H:71
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
const volScalarField & T
engineTime & runTime
scalar nut
auto & name
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
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
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition int32.H:127
dimensionedScalar sign(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
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 POSIX.C:801
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
dimensionedScalar sqrt(const dimensionedScalar &ds)
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
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
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
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