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