SLADelta.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 -------------------------------------------------------------------------------
11 License
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 Class
28  Foam::LESModels::SLADelta
29 
30 Description
31  Delta formulation that accounts for the orientation of the vorticity vector
32  and a flow-sensitised function to detect 2D flow regions. Delta value is
33  strongly reduced in these regions to accelerate the transition from RANS
34  to LES in hybrid RANS/LES simulations.
35 
36  Reference:
37  \verbatim
38  Shur, M. L., Spalart, P. R., Strelets, M. K., & Travin, A. K. (2015).
39  An enhanced version of DES with rapid transition
40  from RANS to LES in separated flows.
41  Flow, turbulence and combustion, 95(4), 709-737.
42  DOI:10.1007/s10494-015-9618-0
43  \endverbatim
44 
45 SourceFiles
46  SLADelta.C
47 
48 \*---------------------------------------------------------------------------*/
49 
50 #ifndef LESModels_SLADelta_H
51 #define LESModels_SLADelta_H
52 
53 #include "LESdelta.H"
54 
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 
57 namespace Foam
58 {
59 namespace LESModels
60 {
61 
62 /*---------------------------------------------------------------------------*\
63  Class SLADelta Declaration
64 \*---------------------------------------------------------------------------*/
65 
66 class SLADelta
67 :
68  public LESdelta
69 {
70  // Private Data
71 
72  //- Run-time selectable delta for hmax
73  // Defaults to the maxDeltaxyz model if not supplied
74  autoPtr<LESdelta> hmaxPtr_;
75 
76  // Model coefficients
77 
78  scalar deltaCoeff_;
79  bool requireUpdate_;
80  bool shielding_;
81  scalar FKHmin_;
82  scalar FKHmax_;
83  scalar a1_;
84  scalar a2_;
85  scalar epsilon_;
86  scalar kappa_;
87  scalar Cd1_;
88  scalar Cd2_;
89 
90 
91  // Private Member Functions
92 
93  // Calculate the delta values
94  void calcDelta();
95 
96  //- No copy construct
97  SLADelta(const SLADelta&) = delete;
98 
99  //- No copy assignment
100  void operator=(const SLADelta&) = delete;
101 
102 
103 public:
104 
105  //- Runtime type information
106  TypeName("SLADelta");
107 
108 
109  // Constructors
110 
111  //- Construct from name, turbulenceModel and IOdictionary
112  SLADelta
113  (
114  const word& name,
116  const dictionary&
117  );
118 
119 
120  //- Destructor
121  virtual ~SLADelta() = default;
122 
123 
124  // Member Functions
125 
126  //- Read the LESdelta dictionary
127  void read(const dictionary&);
128 
129  //- Return the hmax delta field
130  const volScalarField& hmax() const
131  {
132  return hmaxPtr_();
133  }
134 
135  // Correct values
136  void correct();
137 };
138 
139 
140 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
141 
142 } // End namespace LESModels
143 } // End namespace Foam
144 
145 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
146 
147 #endif
148 
149 // ************************************************************************* //
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
const turbulenceModel & turbulence() const
Return turbulenceModel reference.
Definition: LESdelta.H:146
Abstract base class for LES deltas.
Definition: LESdelta.H:49
virtual ~SLADelta()=default
Destructor.
Abstract base class for turbulence models (RAS, LES and laminar).
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
A class for handling words, derived from Foam::string.
Definition: word.H:63
Delta formulation that accounts for the orientation of the vorticity vector and a flow-sensitised fun...
Definition: SLADelta.H:61
TypeName("SLADelta")
Runtime type information.
void read(const dictionary &)
Read the LESdelta dictionary.
Definition: SLADelta.C:375
const volScalarField & hmax() const
Return the hmax delta field.
Definition: SLADelta.H:142
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Namespace for OpenFOAM.