Loading...
Searching...
No Matches
DeltaOmegaTildeDelta.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) 2022 Upstream CFD GmbH
9 Copyright (C) 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::LESModels::DeltaOmegaTildeDelta
29
30Description
31 Delta formulation that accounts for the orientation of the vorticity
32 vector. In "2D-regions" (i.e. xy-plane), delta is of the order
33 max(delta_x,delta_y), so that the influence of delta_z is discarded.
34 This can help to accelerate the transition from RANS to LES in hybrid
35 RANS/LES simulations.
36
37 Reference:
38 \verbatim
39 Shur, M. L., Spalart, P. R., Strelets, M. K., & Travin, A. K. (2015).
40 An enhanced version of DES with rapid transition
41 from RANS to LES in separated flows.
42 Flow, turbulence and combustion, 95(4), 709-737.
43 DOI:10.1007/s10494-015-9618-0
44 \endverbatim
45
46SourceFiles
47 DeltaOmegaTildeDelta.C
48
49\*---------------------------------------------------------------------------*/
50
51#ifndef LESModels_DeltaOmegaTildeDelta_H
52#define LESModels_DeltaOmegaTildeDelta_H
53
54#include "LESdelta.H"
55
56// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57
58namespace Foam
59{
60namespace LESModels
61{
63/*---------------------------------------------------------------------------*\
64 Class DeltaOmegaTildeDelta Declaration
65\*---------------------------------------------------------------------------*/
66
67class DeltaOmegaTildeDelta
68:
69 public LESdelta
70{
71 // Private Data
72
73 //- Run-time selectable delta for hmax
74 // Defaults to the maxDeltaxyz model if not supplied
75 autoPtr<LESdelta> hmaxPtr_;
76
77 //- Model coefficient
78 scalar deltaCoeff_;
79
80 //- Flag to indicate whether hmax requires updating
81 bool requireUpdate_;
82
83
84 // Private Member Functions
85
86 //- Calculate the delta values
87 void calcDelta();
88
89 //- No copy construct
90 DeltaOmegaTildeDelta(const DeltaOmegaTildeDelta&) = delete;
91
92 //- No copy assignment
93 void operator=(const DeltaOmegaTildeDelta&) = delete;
94
95
96public:
97
98 //- Runtime type information
99 TypeName("DeltaOmegaTilde");
100
101
102 // Constructors
103
104 //- Construct from name, turbulenceModel and dictionary
105 DeltaOmegaTildeDelta
106 (
107 const word& name,
110 );
111
112
113 //- Destructor
114 virtual ~DeltaOmegaTildeDelta() = default;
115
116
117 // Member Functions
118
119 //- Read the LESdelta dictionary
120 virtual void read(const dictionary&);
121
122 //- Return the hmax delta field
123 const volScalarField& hmax() const
124 {
125 return hmaxPtr_();
126 }
127
128 // Correct values
129 void correct();
130};
131
132
133// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
134
135} // End namespace LESModels
136} // End namespace Foam
137
138// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
139
140#endif
142// ************************************************************************* //
virtual ~DeltaOmegaTildeDelta()=default
Destructor.
TypeName("DeltaOmegaTilde")
Runtime type information.
const volScalarField & hmax() const
Return the hmax delta field.
virtual void read(const dictionary &)
Read the LESdelta dictionary.
LESdelta(const LESdelta &)=delete
No copy construct.
const turbulenceModel & turbulence() const
Return turbulenceModel reference.
Definition LESdelta.H:146
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Abstract base class for turbulence models (RAS, LES and laminar).
A class for handling words, derived from Foam::string.
Definition word.H:66
Namespace for LES SGS models.
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68