Loading...
Searching...
No Matches
DEShybrid.H
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) 2015 OpenFOAM Foundation
9 Copyright (C) 2015-2024 OpenCFD Ltd.
10 Copyright (C) 2022 Upstream CFD GmbH
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28Class
29 Foam::DEShybrid
30
31Description
32 Improved hybrid convection scheme of Travin et al. for hybrid RAS/LES
33 calculations with enhanced Grey Area Mitigation (GAM) behaviour.
34
35 The scheme provides a blend between two convection schemes, based on local
36 properties including the wall distance, velocity gradient and eddy
37 viscosity. The scheme was originally developed for DES calculations to
38 blend a low-dissipative scheme, e.g. linear, in the vorticity-dominated,
39 finely-resolved regions and a numerically more robust, e.g. upwind-biased,
40 convection scheme in irrotational or coarsely-resolved regions.
41
42 The routine calculates the blending factor denoted as "sigma" in the
43 literature reference, where 0 <= sigma <= sigmaMax, which is then employed
44 to set the weights:
45 \f[
46 weight = (1-sigma) w_1 + sigma w_2
47 \f]
48
49 where
50 \vartable
51 sigma | blending factor
52 w_1 | scheme 1 weights
53 w_2 | scheme 2 weights
54 \endvartable
55
56 First published in:
57 \verbatim
58 Travin, A., Shur, M., Strelets, M., & Spalart, P. R. (2000).
59 Physical and numerical upgrades in the detached-eddy
60 simulation of complex turbulent flows.
61 In LES of Complex Transitional and Turbulent Flows.
62 Proceedings of the Euromech Colloquium 412. Munich, Germany
63 \endverbatim
64
65 Original publication contained a typo for \c C_H3 constant.
66 Corrected version with minor changes for 2 lower limiters published in:
67 \verbatim
68 Spalart, P., Shur, M., Strelets, M., & Travin, A. (2012).
69 Sensitivity of landing-gear noise predictions by large-eddy
70 simulation to numerics and resolution.
71 In 50th AIAA Aerospace Sciences Meeting Including the
72 New Horizons Forum and Aerospace Exposition. Nashville, US.
73 DOI:10.2514/6.2012-1174
74 \endverbatim
75
76 Example of the \c DEShybrid scheme specification using \c linear
77 within the LES region and \c linearUpwind within the RAS region:
78 \verbatim
79 divSchemes
80 {
81 .
82 .
83 div(phi,U) Gauss DEShybrid
84 linear // scheme 1
85 linearUpwind grad(U) // scheme 2
86 delta // LES delta name, e.g. 'delta', 'hmax'
87 0.65 // CDES coefficient
88 30 // Reference velocity scale
89 2 // Reference length scale
90 0 // Minimum sigma limit (0-1)
91 1 // Maximum sigma limit (0-1)
92 1.0e-03 // Limiter of B function, typically 1e-03
93 1.0; // nut limiter (if > 1, GAM extension is active)
94 .
95 .
96 }
97 \endverbatim
98
99Notes
100 - Scheme 1 should be linear (or other low-dissipative schemes) which will
101 be used in the detached/vortex shedding regions.
102 - Scheme 2 should be an upwind/deferred correction/TVD scheme which will
103 be used in the free-stream/Euler/boundary layer regions.
104 - The scheme is compiled into a separate library, and not available to
105 solvers by default. In order to use the scheme, add the library as a
106 run-time loaded library in the \$FOAM\_CASE/system/controlDict
107 dictionary, e.g.:
108 \verbatim
109 libs (turbulenceModelSchemes);
110 \endverbatim
111
112SourceFiles
113 DEShybrid.C
114
115\*---------------------------------------------------------------------------*/
116
117#ifndef Foam_DEShybrid_H
118#define Foam_DEShybrid_H
119
121#include "surfaceInterpolate.H"
122#include "fvcGrad.H"
123#include "blendedSchemeBase.H"
124#include "turbulenceModel.H"
125
126// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
127
128namespace Foam
129{
130
131/*---------------------------------------------------------------------------*\
132 Class DEShybrid Declaration
133\*---------------------------------------------------------------------------*/
134
135template<class Type>
136class DEShybrid
137:
138 public surfaceInterpolationScheme<Type>,
139 public blendedSchemeBase<Type>
140{
141 typedef GeometricField<Type, fvPatchField, volMesh> VolFieldType;
142 typedef GeometricField<Type, fvsPatchField, surfaceMesh> SurfaceFieldType;
144 // Private Data
145
146 //- Scheme 1
148
149 //- Scheme 2
151
152 //- Name of the LES delta field
153 word deltaName_;
154
155 //- DES coefficient
156 scalar CDES_;
157
158 //- Reference velocity scale [m/s]
160
161 //- Reference length scale [m]
163
164 //- Minimum bound for sigma (0 <= sigmaMin <= 1)
165 scalar sigmaMin_;
166
167 //- Maximum bound for sigma (0 <= sigmaMax <= 1)
168 scalar sigmaMax_;
169
170 //- Limiter of B function
171 scalar OmegaLim_;
172
173 //- Limiter for modified GAM behaviour
174 scalar nutLim_;
175
176 //- Scheme constants
177 scalar CH1_;
178 scalar CH2_;
179 scalar CH3_;
180 scalar Cs_;
181
182 //- No copy construct
183 DEShybrid(const DEShybrid&) = delete;
184
185 //- No copy assignment
186 void operator=(const DEShybrid&) = delete;
187
188
189 // Private Member Functions
190
191 //- Check the scheme coefficients
192 void checkValues()
193 {
194 if (U0_.value() <= 0)
195 {
197 << "U0 coefficient must be > 0. "
198 << "Current value: " << U0_ << exit(FatalError);
199 }
200 if (L0_.value() <= 0)
201 {
203 << "L0 coefficient must be > 0. "
204 << "Current value: " << L0_ << exit(FatalError);
205 }
206 if (sigmaMin_ < 0)
207 {
209 << "sigmaMin coefficient must be >= 0. "
210 << "Current value: " << sigmaMin_ << exit(FatalError);
211 }
212 if (sigmaMax_ < 0)
213 {
215 << "sigmaMax coefficient must be >= 0. "
216 << "Current value: " << sigmaMax_ << exit(FatalError);
217 }
218 if (sigmaMin_ > 1)
219 {
221 << "sigmaMin coefficient must be <= 1. "
222 << "Current value: " << sigmaMin_ << exit(FatalError);
223 }
224 if (sigmaMax_ > 1)
225 {
227 << "sigmaMax coefficient must be <= 1. "
228 << "Current value: " << sigmaMax_ << exit(FatalError);
229 }
230
231 if (debug)
232 {
233 Info<< type() << "coefficients:" << nl
234 << " delta : " << deltaName_ << nl
235 << " CDES : " << CDES_ << nl
236 << " U0 : " << U0_.value() << nl
237 << " L0 : " << L0_.value() << nl
238 << " sigmaMin : " << sigmaMin_ << nl
239 << " sigmaMax : " << sigmaMax_ << nl
240 << " OmegaLim : " << OmegaLim_ << nl
241 << " nutLim : " << nutLim_ << nl
242 << " CH1 : " << CH1_ << nl
243 << " CH2 : " << CH2_ << nl
244 << " CH3 : " << CH3_ << nl
245 << " Cs : " << Cs_ << nl
246 << endl;
247 }
248 }
249
250
251 //- Calculate the blending factor
252 tmp<surfaceScalarField> calcBlendingFactor
253 (
254 const VolFieldType& vf,
255 const volScalarField& nut,
256 const volScalarField& nu,
257 const volVectorField& U,
258 const volScalarField& delta
259 ) const
260 {
262 const volTensorField& gradU = tgradU.cref();
263 const volScalarField S(sqrt(2.0)*mag(symm(gradU)));
264 const volScalarField Omega(sqrt(2.0)*mag(skew(tgradU)));
265 const dimensionedScalar tau0_ = L0_/U0_;
266
268 CH3_*Omega*max(S, Omega)
269 /max(0.5*(sqr(S) + sqr(Omega)), sqr(OmegaLim_/tau0_));
270
271 tmp<volScalarField> tg = tanh(pow4(tB));
272
274 max(Foam::sqrt(0.5*(sqr(S) + sqr(Omega))), 0.1/tau0_);
275
276 tmp<volScalarField> tlTurb =
278 (
279 max
280 (
281 (max(nut, min(sqr(Cs_*delta)*S, nutLim_*nut)) + nu)
282 /(pow(0.09, 1.5)*tK),
284 )
285 );
286
287 const volScalarField A
288 (
289 CH2_*max
290 (
291 scalar(0),
292 CDES_*delta/max(tlTurb*tg, SMALL*L0_) - 0.5
293 )
294 );
295
296
297 const word factorName(IOobject::scopedName(typeName, "Factor"));
298 const fvMesh& mesh = this->mesh();
299
300 IOobject factorIO
301 (
302 factorName,
303 mesh.time().timeName(),
304 mesh.thisDb(),
308 );
309
310 if (blendedSchemeBaseName::debug)
311 {
312 auto* factorPtr = mesh.getObjectPtr<volScalarField>(factorName);
313
314 if (!factorPtr)
315 {
316 factorPtr =
318 (
319 factorIO,
320 mesh,
322 );
323
324 regIOobject::store(factorPtr);
325 }
326
327 auto& factor = *factorPtr;
328
329 factor = max(sigmaMax_*tanh(pow(A, CH1_)), sigmaMin_);
330
332 (
333 vf.name() + "BlendingFactor",
334 fvc::interpolate(factor)
335 );
336 }
337 else
338 {
339 factorIO.registerObject(IOobjectOption::NO_REGISTER);
340
341 volScalarField factor
342 (
343 factorIO,
344 max(sigmaMax_*tanh(pow(A, CH1_)), sigmaMin_)
345 );
346
348 (
349 vf.name() + "BlendingFactor",
350 fvc::interpolate(factor)
351 );
352 }
353 }
354
355
356public:
357
358 //- Runtime type information
359 TypeName("DEShybrid");
360
361
362 // Constructors
363
364 //- Construct from mesh and Istream.
365 // The name of the flux field is read from the Istream and looked-up
366 // from the mesh objectRegistry
367 DEShybrid(const fvMesh& mesh, Istream& is)
368 :
369 surfaceInterpolationScheme<Type>(mesh),
370 tScheme1_(surfaceInterpolationScheme<Type>::New(mesh, is)),
371 tScheme2_(surfaceInterpolationScheme<Type>::New(mesh, is)),
372 deltaName_(is),
373 CDES_(readScalar(is)),
374 U0_("U0", dimLength/dimTime, readScalar(is)),
375 L0_("L0", dimLength, readScalar(is)),
376 sigmaMin_(readScalar(is)),
377 sigmaMax_(readScalar(is)),
378 OmegaLim_(readScalar(is)),
379 nutLim_(readScalarOrDefault(is, scalar(1))),
380 CH1_(3.0),
381 CH2_(1.0),
382 CH3_(2.0),
383 Cs_(0.18)
384 {
385 checkValues();
386 }
387
388 //- Construct from mesh, faceFlux and Istream
389 DEShybrid
390 (
391 const fvMesh& mesh,
392 const surfaceScalarField& faceFlux,
393 Istream& is
394 )
395 :
396 surfaceInterpolationScheme<Type>(mesh),
397 tScheme1_
399 surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
400 ),
401 tScheme2_
402 (
403 surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
404 ),
405 deltaName_(is),
406 CDES_(readScalar(is)),
407 U0_("U0", dimLength/dimTime, readScalar(is)),
408 L0_("L0", dimLength, readScalar(is)),
409 sigmaMin_(readScalar(is)),
410 sigmaMax_(readScalar(is)),
411 OmegaLim_(readScalar(is)),
412 nutLim_(readScalarOrDefault(is, scalar(1))),
413 CH1_(3.0),
414 CH2_(1.0),
415 CH3_(2.0),
416 Cs_(0.18)
417 {
418 checkValues();
419 }
420
421
422 // Member Functions
423
424 //- Return the face-based blending factor
426 (
428 ) const
429 {
430 const fvMesh& mesh = this->mesh();
431
432 // Retrieve LES delta from the mesh database
433 const auto& delta =
434 mesh.lookupObject<const volScalarField>(deltaName_);
435
436 // Retrieve turbulence model from the mesh database
437 const auto* modelPtr =
438 mesh.cfindObject<turbulenceModel>
439 (
441 );
442
443 if (modelPtr)
444 {
445 const auto& model = *modelPtr;
446
447 return calcBlendingFactor
448 (
449 vf, model.nut(), model.nu(), model.U(), delta
450 );
451 }
452
454 << "Scheme requires a turbulence model to be present. "
455 << "Unable to retrieve turbulence model from the mesh "
456 << "database" << exit(FatalError);
457
458 return nullptr;
459 }
460
461
462 //- Return the interpolation weighting factors
463 tmp<surfaceScalarField> weights(const VolFieldType& vf) const
464 {
466
467 return
468 (scalar(1) - bf)*tScheme1_().weights(vf)
469 + bf*tScheme2_().weights(vf);
470 }
472
473 //- Return the face-interpolate of the given cell field
474 //- with explicit correction
475 tmp<SurfaceFieldType> interpolate(const VolFieldType& vf) const
476 {
478
479 return
480 (scalar(1) - bf)*tScheme1_().interpolate(vf)
481 + bf*tScheme2_().interpolate(vf);
482 }
483
484
485 //- Return true if this scheme uses an explicit correction
486 virtual bool corrected() const
487 {
488 return tScheme1_().corrected() || tScheme2_().corrected();
489 }
490
491
492 //- Return the explicit correction to the face-interpolate
493 //- for the given field
494 virtual tmp<SurfaceFieldType> correction(const VolFieldType& vf) const
495 {
497
498 if (tScheme1_().corrected())
499 {
500 if (tScheme2_().corrected())
501 {
502 return
503 (
504 (scalar(1) - bf)
505 * tScheme1_().correction(vf)
506 + bf
507 * tScheme2_().correction(vf)
508 );
509 }
510 else
512 return
513 (
514 (scalar(1) - bf)
515 * tScheme1_().correction(vf)
516 );
517 }
518 }
519 else if (tScheme2_().corrected())
520 {
521 return (bf*tScheme2_().correction(vf));
522 }
523
524 return nullptr;
526};
527
528
529// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
530
531} // End namespace Foam
532
533// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
534
535#endif
536
537// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
scalar delta
Improved hybrid convection scheme of Travin et al. for hybrid RAS/LES calculations with enhanced Grey...
Definition DEShybrid.H:147
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
Definition DEShybrid.H:538
DEShybrid(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
Definition DEShybrid.H:409
tmp< surfaceScalarField > weights(const VolFieldType &vf) const
Return the interpolation weighting factors.
Definition DEShybrid.H:511
virtual tmp< SurfaceFieldType > correction(const VolFieldType &vf) const
Return the explicit correction to the face-interpolate for the given field.
Definition DEShybrid.H:548
DEShybrid(const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &is)
Construct from mesh, faceFlux and Istream.
Definition DEShybrid.H:434
tmp< SurfaceFieldType > interpolate(const VolFieldType &vf) const
Return the face-interpolate of the given cell field with explicit correction.
Definition DEShybrid.H:525
virtual tmp< surfaceScalarField > blendingFactor(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the face-based blending factor.
Definition DEShybrid.H:472
TypeName("DEShybrid")
Runtime type information.
Generic GeometricField class.
@ NO_REGISTER
Do not request registration (bool: false).
@ REGISTER
Request registration (bool: true).
bool registerObject() const noexcept
Should objects created with this IOobject be registered?
@ 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
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
Base class for blended schemes to provide access to the blending factor surface field.
blendedSchemeBase()=default
Constructor.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
bool store()
Register object with its registry and transfer ownership to the registry.
Abstract base class for surface interpolation schemes.
const fvMesh & mesh() const
Return mesh reference.
static tmp< surfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
A class for managing temporary objects.
Definition tmp.H:75
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
Definition tmpI.H:221
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
Abstract base class for turbulence models (RAS, LES and laminar).
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
scalar nut
Calculate the gradient of the given field.
Namespace for handling debugging switches.
Definition debug.C:45
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
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.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar tanh(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
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
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
GeometricField< tensor, fvPatchField, volMesh > volTensorField
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
dimensionedScalar pow4(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensionedTensor skew(const dimensionedTensor &dt)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
volScalarField & nu
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68