DeardorffDiffStress.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::LESModels::DeardorffDiffStress
26 
27 Description
28  Differential SGS Stress Equation Model for incompressible and
29  compressible flows
30 
31  Reference:
32  \verbatim
33  Deardorff, J. W. (1973).
34  The use of subgrid transport equations in a three-dimensional model
35  of atmospheric turbulence.
36  Journal of Fluids Engineering, 95(3), 429-438.
37  \endverbatim
38 
39  This SGS model uses a full balance equation for the SGS stress tensor to
40  simulate the behaviour of B.
41 
42  This implementation is as described in the above paper except that the
43  triple correlation model of Donaldson is replaced with the generalized
44  gradient diffusion model of Daly and Harlow:
45  \verbatim
46  Daly, B. J., & Harlow, F. H. (1970).
47  Transport equations in turbulence.
48  Physics of Fluids (1958-1988), 13(11), 2634-2649.
49  \endverbatim
50  with the default value for the coefficient Cs of 0.25 from
51  \verbatim
52  Launder, B. E., Reece, G. J., & Rodi, W. (1975).
53  Progress in the development of a Reynolds-stress turbulence closure.
54  Journal of fluid mechanics, 68(03), 537-566.
55  \endverbatim
56 
57 SourceFiles
58  DeardorffDiffStress.C
59 
60 \*---------------------------------------------------------------------------*/
61 
62 #ifndef DeardorffDiffStress_H
63 #define DeardorffDiffStress_H
64 
65 #include "LESModel.H"
66 #include "ReynoldsStress.H"
67 
68 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
69 
70 namespace Foam
71 {
72 namespace LESModels
73 {
74 
75 /*---------------------------------------------------------------------------*\
76  Class DeardorffDiffStress Declaration
77 \*---------------------------------------------------------------------------*/
78 
79 template<class BasicTurbulenceModel>
81 :
82  public ReynoldsStress<LESModel<BasicTurbulenceModel>>
83 {
84  // Private Member Functions
85 
86  // Disallow default bitwise copy construct and assignment
88  void operator=(const DeardorffDiffStress&);
89 
90 
91 protected:
92 
93  // Protected data
94 
95  // Model constants
96 
101 
102 
103  // Protected Member Functions
104 
105  //- Update the eddy-viscosity
106  virtual void correctNut();
107 
108 
109 public:
111  typedef typename BasicTurbulenceModel::alphaField alphaField;
112  typedef typename BasicTurbulenceModel::rhoField rhoField;
113  typedef typename BasicTurbulenceModel::transportModel transportModel;
114 
115 
116  //- Runtime type information
117  TypeName("DeardorffDiffStress");
118 
119 
120  // Constructors
121 
122  //- Constructor from components
124  (
125  const alphaField& alpha,
126  const rhoField& rho,
127  const volVectorField& U,
128  const surfaceScalarField& alphaRhoPhi,
129  const surfaceScalarField& phi,
130  const transportModel& transport,
131  const word& propertiesName = turbulenceModel::propertiesName,
132  const word& type = typeName
133  );
134 
135 
136  //- Destructor
137  virtual ~DeardorffDiffStress()
138  {}
139 
140 
141  // Member Functions
142 
143  //- Read model coefficients if they have changed
144  virtual bool read();
145 
146  //- Return the turbulence kinetic energy dissipation rate
147  virtual tmp<volScalarField> epsilon() const;
148 
149  //- Correct sub-grid stress, eddy-Viscosity and related properties
150  virtual void correct();
151 };
152 
153 
154 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
155 
156 } // End namespace LESModels
157 } // End namespace Foam
158 
159 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
160 
161 #ifdef NoRepository
162  #include "DeardorffDiffStress.C"
163 #endif
164 
165 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
166 
167 #endif
168 
169 // ************************************************************************* //
surfaceScalarField & phi
Differential SGS Stress Equation Model for incompressible and compressible flows. ...
BasicTurbulenceModel::rhoField rhoField
BasicTurbulenceModel::alphaField alphaField
BasicTurbulenceModel::transportModel transportModel
TypeName("DeardorffDiffStress")
Runtime type information.
virtual void correct()
Correct sub-grid stress, eddy-Viscosity and related properties.
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
U
Definition: pEqn.H:72
virtual void correctNut()
Update the eddy-viscosity.
virtual bool read()
Read model coefficients if they have changed.
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Namespace for OpenFOAM.
Reynolds-stress turbulence model base class.