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-2022 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
27Class
28 Foam::DEShybrid
29
30Description
31 Hybrid convection scheme of Travin et al. for hybrid RAS/LES calculations.
32
33 The scheme provides a blend between two convection schemes, based on local
34 properties including the wall distance, velocity gradient and eddy
35 viscosity. The scheme was originally developed for DES calculations to
36 blend a low-dissipative scheme, e.g. linear, in the vorticity-dominated,
37 finely-resolved regions and a numerically more robust, e.g. upwind-biased,
38 convection scheme in irrotational or coarsely-resolved regions.
39
40 The routine calculates the blending factor denoted as "sigma" in the
41 literature reference, where 0 <= sigma <= sigmaMax, which is then employed
42 to set the weights:
43 \f[
44 weight = (1-sigma) w_1 + sigma w_2
45 \f]
46
47 where
48 \vartable
49 sigma | blending factor
50 w_1 | scheme 1 weights
51 w_2 | scheme 2 weights
52 \endvartable
53
54 First published in:
55 \verbatim
56 A. Travin, M. Shur, M. Strelets, P. Spalart (2000).
57 Physical and numerical upgrades in the detached-eddy simulation of
58 complex turbulent flows.
59 In Proceedings of the 412th Euromech Colloquium on LES and Complex
60 Transition and Turbulent Flows, Munich, Germany
61 \endverbatim
62
63 Original publication contained a typo for C_H3 constant. Corrected version
64 with minor changes for 2 lower limiters published in:
65 \verbatim
66 P. Spalart, M. Shur, M. Strelets, A. Travin (2012).
67 Sensitivity of Landing-Gear Noise Predictions by Large-Eddy
68 Simulation to Numerics and Resolution.
69 AIAA Paper 2012-1174, 50th AIAA Aerospace Sciences Meeting,
70 Nashville / TN, Jan. 2012
71 \endverbatim
72
73 Example of the DEShybrid scheme specification using linear within the LES
74 region and linearUpwind within the RAS region:
75 \verbatim
76 divSchemes
77 {
78 .
79 .
80 div(phi,U) Gauss DEShybrid
81 linear // scheme 1
82 linearUpwind grad(U) // scheme 2
83 hmax // LES delta name, e.g. 'delta', 'hmax'
84 0.65 // DES coefficient, typically = 0.65
85 30 // Reference velocity scale
86 2 // Reference length scale
87 0 // Minimum sigma limit (0-1)
88 1 // Maximum sigma limit (0-1)
89 1.0e-03; // Limiter of B function, typically 1e-03
90 .
91 .
92 }
93 \endverbatim
94
95Notes
96 - Scheme 1 should be linear (or other low-dissipative schemes) which will
97 be used in the detached/vortex shedding regions.
98 - Scheme 2 should be an upwind/deferred correction/TVD scheme which will
99 be used in the free-stream/Euler/boundary layer regions.
100 - the scheme is compiled into a separate library, and not available to
101 solvers by default. In order to use the scheme, add the library as a
102 run-time loaded library in the \$FOAM\_CASE/system/controlDict
103 dictionary, e.g.:
104 \verbatim
105 libs ("libturbulenceModelSchemes.so");
106 \endverbatim
107
108SourceFiles
109 DEShybrid.C
110
111\*---------------------------------------------------------------------------*/
112
113#ifndef DEShybrid_H
114#define DEShybrid_H
115
117#include "surfaceInterpolate.H"
118#include "fvcGrad.H"
119#include "blendedSchemeBase.H"
122
123// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
124
125namespace Foam
126{
127
128/*---------------------------------------------------------------------------*\
129 Class DEShybrid Declaration
130\*---------------------------------------------------------------------------*/
131
132template<class Type>
133class DEShybrid
134:
135 public surfaceInterpolationScheme<Type>,
136 public blendedSchemeBase<Type>
137{
138 // Private Data
139
140 //- Scheme 1
141 tmp<surfaceInterpolationScheme<Type>> tScheme1_;
142
143 //- Scheme 2
145
146 //- Name of the LES delta field
147 word deltaName_;
148
149 //- DES Coefficient
150 scalar CDES_;
151
152 //- Reference velocity scale [m/s]
154
155 //- Reference length scale [m]
157
158 //- Minimum bound for sigma (0 <= sigmaMin <= 1)
159 scalar sigmaMin_;
160
161 //- Maximum bound for sigma (0 <= sigmaMax <= 1)
162 scalar sigmaMax_;
163
164 //- Limiter of B function
165 scalar OmegaLim_;
166
167 //- Scheme constants
168 scalar CH1_;
169 scalar CH2_;
170 scalar CH3_;
171
172 //- No copy construct
173 DEShybrid(const DEShybrid&) = delete;
174
175 //- No copy assignment
176 void operator=(const DEShybrid&) = delete;
177
178
179 // Private Member Functions
180
181 //- Calculate the blending factor
182 tmp<surfaceScalarField> calcBlendingFactor
183 (
185 const volScalarField& nuEff,
186 const volVectorField& U,
187 const volScalarField& delta
188 ) const
189 {
191 const volScalarField S(sqrt(2.0)*mag(symm(gradU())));
192 const volScalarField Omega(sqrt(2.0)*mag(skew(gradU())));
193 const dimensionedScalar tau0_ = L0_/U0_;
194
195 const volScalarField B
196 (
197 CH3_*Omega*max(S, Omega)
198 /max(0.5*(sqr(S) + sqr(Omega)), sqr(OmegaLim_/tau0_))
199 );
200
201 const volScalarField K
202 (
203 max(Foam::sqrt(0.5*(sqr(S) + sqr(Omega))), 0.1/tau0_)
204 );
205
206 const volScalarField lTurb
207 (
209 (
210 max
211 (
212 nuEff/(pow(0.09, 1.5)*K),
214 )
215 )
216 );
217
218 const volScalarField g(tanh(pow4(B)));
219
220 const volScalarField A
221 (
222 CH2_*max
223 (
224 scalar(0),
225 CDES_*delta/max(lTurb*g, SMALL*L0_) - 0.5
226 )
227 );
228
229
230 const word factorName(IOobject::scopedName(typeName, "Factor"));
231 const fvMesh& mesh = this->mesh();
232
233 const IOobject factorIO
234 (
235 factorName,
236 mesh.time().timeName(),
237 mesh,
240 );
241
242 if (blendedSchemeBaseName::debug)
243 {
244 auto* factorPtr = mesh.getObjectPtr<volScalarField>(factorName);
245
246 if (!factorPtr)
247 {
248 factorPtr =
250 (
251 factorIO,
252 mesh,
254 );
255
256 const_cast<fvMesh&>(mesh).objectRegistry::store(factorPtr);
257 }
258
259 auto& factor = *factorPtr;
260
261 factor = max(sigmaMax_*tanh(pow(A, CH1_)), sigmaMin_);
262
264 (
265 vf.name() + "BlendingFactor",
266 fvc::interpolate(factor)
267 );
268 }
269 else
270 {
271 const volScalarField factor
272 (
273 factorIO,
274 max(sigmaMax_*tanh(pow(A, CH1_)), sigmaMin_)
275 );
276
278 (
279 vf.name() + "BlendingFactor",
280 fvc::interpolate(factor)
281 );
282 }
283 }
284
285
286public:
287
288 //- Runtime type information
289 TypeName("DEShybrid");
290
291
292 // Constructors
293
294 //- Construct from mesh and Istream.
295 // The name of the flux field is read from the Istream and looked-up
296 // from the mesh objectRegistry
297 DEShybrid(const fvMesh& mesh, Istream& is)
298 :
300 tScheme1_
301 (
303 ),
304 tScheme2_
305 (
307 ),
308 deltaName_(is),
309 CDES_(readScalar(is)),
310 U0_("U0", dimLength/dimTime, readScalar(is)),
311 L0_("L0", dimLength, readScalar(is)),
312 sigmaMin_(readScalar(is)),
313 sigmaMax_(readScalar(is)),
314 OmegaLim_(readScalar(is)),
315 CH1_(3.0),
316 CH2_(1.0),
317 CH3_(2.0)
318 {
319 if (U0_.value() <= 0)
320 {
322 << "U0 coefficient must be > 0. "
323 << "Current value: " << U0_ << exit(FatalError);
324 }
325 if (L0_.value() <= 0)
326 {
328 << "L0 coefficient must be > 0. "
329 << "Current value: " << L0_ << exit(FatalError);
330 }
331 if (sigmaMin_ < 0)
332 {
334 << "sigmaMin coefficient must be >= 0. "
335 << "Current value: " << sigmaMin_ << exit(FatalError);
336 }
337 if (sigmaMax_ < 0)
338 {
340 << "sigmaMax coefficient must be >= 0. "
341 << "Current value: " << sigmaMax_ << exit(FatalError);
342 }
343 if (sigmaMin_ > 1)
344 {
346 << "sigmaMin coefficient must be <= 1. "
347 << "Current value: " << sigmaMin_ << exit(FatalError);
348 }
349 if (sigmaMax_ > 1)
350 {
352 << "sigmaMax coefficient must be <= 1. "
353 << "Current value: " << sigmaMax_ << exit(FatalError);
354 }
355 }
356
357 //- Construct from mesh, faceFlux and Istream
359 (
360 const fvMesh& mesh,
361 const surfaceScalarField& faceFlux,
362 Istream& is
363 )
364 :
366 tScheme1_
367 (
368 surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
369 ),
370 tScheme2_
371 (
372 surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
373 ),
374 deltaName_(is),
375 CDES_(readScalar(is)),
376 U0_("U0", dimLength/dimTime, readScalar(is)),
377 L0_("L0", dimLength, readScalar(is)),
378 sigmaMin_(readScalar(is)),
379 sigmaMax_(readScalar(is)),
380 OmegaLim_(readScalar(is)),
381 CH1_(3.0),
382 CH2_(1.0),
383 CH3_(2.0)
384 {
385 if (U0_.value() <= 0)
386 {
388 << "U0 coefficient must be > 0. "
389 << "Current value: " << U0_ << exit(FatalError);
390 }
391 if (L0_.value() <= 0)
392 {
394 << "L0 coefficient must be > 0. "
395 << "Current value: " << U0_ << exit(FatalError);
396 }
397 if (sigmaMin_ < 0)
398 {
400 << "sigmaMin coefficient must be >= 0. "
401 << "Current value: " << sigmaMin_ << exit(FatalError);
402 }
403 if (sigmaMax_ < 0)
404 {
406 << "sigmaMax coefficient must be >= 0. "
407 << "Current value: " << sigmaMax_ << exit(FatalError);
408 }
409 if (sigmaMin_ > 1)
410 {
412 << "sigmaMin coefficient must be <= 1. "
413 << "Current value: " << sigmaMin_ << exit(FatalError);
414 }
415 if (sigmaMax_ > 1)
416 {
418 << "sigmaMax coefficient must be <= 1. "
419 << "Current value: " << sigmaMax_ << exit(FatalError);
420 }
421 }
422
423
424 // Member Functions
425
426 //- Return the face-based blending factor
428 (
430 ) const
431 {
432 const fvMesh& mesh = this->mesh();
433
434 typedef compressible::turbulenceModel cmpModel;
435 typedef incompressible::turbulenceModel icoModel;
436
437 // Lookup the LES delta from the mesh database
438 const volScalarField& delta = this->mesh().template
439 lookupObject<const volScalarField>(deltaName_);
440
441 // Could avoid the compressible/incompressible case by looking
442 // up all fields from the database - but retrieving from model
443 // ensures consistent fields are being employed e.g. for multiphase
444 // where group name is used
445
446 if (mesh.foundObject<icoModel>(icoModel::propertiesName))
447 {
448 const icoModel& model =
449 mesh.lookupObject<icoModel>(icoModel::propertiesName);
450
451 return calcBlendingFactor(vf, model.nuEff(), model.U(), delta);
452 }
453 else if (mesh.foundObject<cmpModel>(cmpModel::propertiesName))
454 {
455 const cmpModel& model =
456 mesh.lookupObject<cmpModel>(cmpModel::propertiesName);
457
458 return calcBlendingFactor
459 (
460 vf, model.muEff()/model.rho(), model.U(), delta
461 );
462 }
463
465 << "Scheme requires a turbulence model to be present. "
466 << "Unable to retrieve turbulence model from the mesh "
467 << "database" << exit(FatalError);
468
469 return nullptr;
470 }
471
472
473 //- Return the interpolation weighting factors
474 tmp<surfaceScalarField> weights
475 (
476 const GeometricField<Type, fvPatchField, volMesh>& vf
477 ) const
478 {
480
481 return
482 (scalar(1) - bf)*tScheme1_().weights(vf)
483 + bf*tScheme2_().weights(vf);
484 }
486
487 //- Return the face-interpolate of the given cell field
488 // with explicit correction
491 (
493 ) const
494 {
496
497 return
498 (scalar(1) - bf)*tScheme1_().interpolate(vf)
499 + bf*tScheme2_().interpolate(vf);
500 }
502
503 //- Return true if this scheme uses an explicit correction
504 virtual bool corrected() const
505 {
506 return tScheme1_().corrected() || tScheme2_().corrected();
507 }
508
509
510 //- Return the explicit correction to the face-interpolate
511 // for the given field
514 (
516 ) const
517 {
519
520 if (tScheme1_().corrected())
521 {
522 if (tScheme2_().corrected())
523 {
524 return
525 (
526 (scalar(1) - bf)
527 * tScheme1_().correction(vf)
528 + bf
529 * tScheme2_().correction(vf)
530 );
531 }
532 else
533 {
534 return
535 (
536 (scalar(1) - bf)
537 * tScheme1_().correction(vf)
538 );
539 }
540 }
541 else if (tScheme2_().corrected())
542 {
543 return (bf*tScheme2_().correction(vf));
544 }
545
546 return nullptr;
547 }
548};
549
550
551// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
552
553} // End namespace Foam
554
555// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
556
557#endif
558
559// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
CGAL::Exact_predicates_exact_constructions_kernel K
scalar delta
const uniformDimensionedVectorField & g
Hybrid convection scheme of Travin et al. for hybrid RAS/LES calculations.
Definition: DEShybrid.H:148
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
Definition: DEShybrid.H:515
tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the interpolation weighting factors.
Definition: DEShybrid.H:486
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
Definition: DEShybrid.H:525
DEShybrid(const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &is)
Construct from mesh, faceFlux and Istream.
Definition: DEShybrid.H:370
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the face-interpolate of the given cell field.
Definition: DEShybrid.H:502
virtual tmp< surfaceScalarField > blendingFactor(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the face-based blending factor.
Definition: DEShybrid.H:439
TypeName("DEShybrid")
Runtime type information.
Generic GeometricField class.
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
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:47
Templated abstract base class for single-phase incompressible turbulence models.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
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
Base class for blended schemes to provide access to the blending factor surface field.
const Type & value() const
Return const reference to value.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
const Type & lookupObject(const word &name, const bool recursive=false) const
Type * getObjectPtr(const word &name, const bool recursive=false) const
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:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Calculate the gradient of the given field.
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: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.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
dimensionedScalar tanh(const dimensionedScalar &ds)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensionedScalar pow4(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dimensionedTensor skew(const dimensionedTensor &dt)
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73