DeardorffDiffStress.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) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2019-2021 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::DeardorffDiffStress
29
30Group
31 grpLESTurbulence
32
33Description
34 Differential SGS Stress Equation Model for incompressible and
35 compressible flows
36
37 Reference:
38 \verbatim
39 Deardorff, J. W. (1973).
40 The use of subgrid transport equations in a three-dimensional model
41 of atmospheric turbulence.
42 Journal of Fluids Engineering, 95(3), 429-438.
43 \endverbatim
44
45 This SGS model uses a full balance equation for the SGS stress tensor to
46 simulate the behaviour of B.
47
48 This implementation is as described in the above paper except that the
49 triple correlation model of Donaldson is replaced with the generalized
50 gradient diffusion model of Daly and Harlow:
51 \verbatim
52 Daly, B. J., & Harlow, F. H. (1970).
53 Transport equations in turbulence.
54 Physics of Fluids (1958-1988), 13(11), 2634-2649.
55 \endverbatim
56 with the default value for the coefficient Cs of 0.25 from
57 \verbatim
58 Launder, B. E., Reece, G. J., & Rodi, W. (1975).
59 Progress in the development of a Reynolds-stress turbulence closure.
60 Journal of fluid mechanics, 68(03), 537-566.
61 \endverbatim
62
63SourceFiles
64 DeardorffDiffStress.C
65
66\*---------------------------------------------------------------------------*/
67
68#ifndef DeardorffDiffStress_H
69#define DeardorffDiffStress_H
70
71#include "LESModel.H"
72#include "ReynoldsStress.H"
73
74// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
75
76namespace Foam
77{
78namespace LESModels
79{
80
81/*---------------------------------------------------------------------------*\
82 Class DeardorffDiffStress Declaration
83\*---------------------------------------------------------------------------*/
84
85template<class BasicTurbulenceModel>
87:
88 public ReynoldsStress<LESModel<BasicTurbulenceModel>>
89{
90 // Private Member Functions
91
92 //- No copy construct
94
95 //- No copy assignment
96 void operator=(const DeardorffDiffStress&) = delete;
97
98
99protected:
100
101 // Protected data
102
103 // Model constants
109
110
111 // Protected Member Functions
112
113 //- Update the eddy-viscosity
114 virtual void correctNut();
115
116
117public:
119 typedef typename BasicTurbulenceModel::alphaField alphaField;
120 typedef typename BasicTurbulenceModel::rhoField rhoField;
121 typedef typename BasicTurbulenceModel::transportModel transportModel;
122
123
124 //- Runtime type information
125 TypeName("DeardorffDiffStress");
126
127
128 // Constructors
129
130 //- Constructor from components
132 (
133 const alphaField& alpha,
134 const rhoField& rho,
135 const volVectorField& U,
136 const surfaceScalarField& alphaRhoPhi,
137 const surfaceScalarField& phi,
138 const transportModel& transport,
139 const word& propertiesName = turbulenceModel::propertiesName,
140 const word& type = typeName
141 );
142
143
144 //- Destructor
145 virtual ~DeardorffDiffStress() = default;
146
147
148 // Member Functions
149
150 //- Read model coefficients if they have changed
151 virtual bool read();
152
153 //- Correct sub-grid stress, eddy-Viscosity and related properties
154 virtual void correct();
155};
156
157
158// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
159
160} // End namespace LESModels
161} // End namespace Foam
162
163// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
164
165#ifdef NoRepository
166 #include "DeardorffDiffStress.C"
167#endif
168
169// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
170
171#endif
172
173// ************************************************************************* //
surfaceScalarField & phi
Differential SGS Stress Equation Model for incompressible and compressible flows.
BasicTurbulenceModel::alphaField alphaField
BasicTurbulenceModel::rhoField rhoField
virtual void correct()
Correct sub-grid stress, eddy-Viscosity and related properties.
virtual ~DeardorffDiffStress()=default
Destructor.
TypeName("DeardorffDiffStress")
Runtime type information.
virtual void correctNut()
Update the eddy-viscosity.
BasicTurbulenceModel::transportModel transportModel
virtual bool read()
Read model coefficients if they have changed.
Reynolds-stress turbulence model base class.
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
Namespace for OpenFOAM.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
volScalarField & alpha
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73