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 -------------------------------------------------------------------------------
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::DeardorffDiffStress
29 
30 Group
31  grpLESTurbulence
32 
33 Description
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 
63 SourceFiles
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 
76 namespace Foam
77 {
78 namespace LESModels
79 {
80 
81 /*---------------------------------------------------------------------------*\
82  Class DeardorffDiffStress Declaration
83 \*---------------------------------------------------------------------------*/
84 
85 template<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 
99 protected:
100 
101  // Protected data
102 
103  // Model constants
104 
109 
110 
111  // Protected Member Functions
112 
113  //- Update the eddy-viscosity
114  virtual void correctNut();
115 
116 
117 public:
118 
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 // ************************************************************************* //
Foam::LESModels::DeardorffDiffStress
Differential SGS Stress Equation Model for incompressible and compressible flows.
Definition: DeardorffDiffStress.H:81
Foam::LESModels::DeardorffDiffStress::Cm_
dimensionedScalar Cm_
Definition: DeardorffDiffStress.H:101
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:63
Foam::LESModels::DeardorffDiffStress::correctNut
virtual void correctNut()
Definition: DeardorffDiffStress.C:36
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Definition: readThermalProperties.H:212
Foam::turbulenceModel::propertiesName
static const word propertiesName
Definition: turbulenceModel.H:96
Foam::LESModels::DeardorffDiffStress::read
virtual bool read()
Definition: DeardorffDiffStress.C:121
ReynoldsStress.H
Foam::ReynoldsStress
Reynolds-stress turbulence model base class.
Definition: ReynoldsStress.H:49
DeardorffDiffStress.C
rho
rho
Definition: readInitialConditions.H:88
Foam::LESModels::DeardorffDiffStress::Cs_
dimensionedScalar Cs_
Definition: DeardorffDiffStress.H:103
LESModel.H
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:36
Foam
Definition: atmBoundaryLayer.C:26
Foam::LESModels::DeardorffDiffStress::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: DeardorffDiffStress.H:114
U
U
Definition: pEqn.H:72
Foam::LESModels::DeardorffDiffStress::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: DeardorffDiffStress.H:116
Foam::LESModels::DeardorffDiffStress::Ck_
dimensionedScalar Ck_
Definition: DeardorffDiffStress.H:100
Foam::LESModels::DeardorffDiffStress::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: DeardorffDiffStress.H:115
Foam::LESModels::DeardorffDiffStress::TypeName
TypeName("DeardorffDiffStress")
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Definition: POSIX.C:717
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:49
Foam::LESModels::DeardorffDiffStress::correct
virtual void correct()
Definition: DeardorffDiffStress.C:138
Foam::LESModels::DeardorffDiffStress::Ce_
dimensionedScalar Ce_
Definition: DeardorffDiffStress.H:102
Foam::LESModels::DeardorffDiffStress::~DeardorffDiffStress
virtual ~DeardorffDiffStress()=default