LESModel.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) 2013-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 Namespace
25  Foam::LESModels
26 
27 Description
28  Namespace for LES SGS models.
29 
30 Class
31  Foam::LESModel
32 
33 Description
34  Templated abstract base class for LES SGS models
35 
36 SourceFiles
37  LESModel.C
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #ifndef LESModel_H
42 #define LESModel_H
43 
44 #include "TurbulenceModel.H"
45 #include "LESdelta.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 /*---------------------------------------------------------------------------*\
53  Class LESModel Declaration
54 \*---------------------------------------------------------------------------*/
55 
56 template<class BasicTurbulenceModel>
57 class LESModel
58 :
59  public BasicTurbulenceModel
60 {
61 
62 protected:
63 
64  // Protected data
65 
66  //- LES coefficients dictionary
68 
69  //- Turbulence on/off flag
71 
72  //- Flag to print the model coeffs at run-time
74 
75  //- Model coefficients dictionary
77 
78  //- Lower limit of k
80 
81  //- Lower limit of epsilon
83 
84  //- Lower limit for omega
86 
87  //- Run-time selectable delta model
89 
90 
91  // Protected Member Functions
92 
93  //- Print model coefficients
94  virtual void printCoeffs(const word& type);
95 
96 
97 private:
98 
99  // Private Member Functions
100 
101  //- Disallow default bitwise copy construct
102  LESModel(const LESModel&);
103 
104  //- Disallow default bitwise assignment
105  void operator=(const LESModel&);
106 
107 
108 public:
110  typedef typename BasicTurbulenceModel::alphaField alphaField;
111  typedef typename BasicTurbulenceModel::rhoField rhoField;
112  typedef typename BasicTurbulenceModel::transportModel transportModel;
113 
114 
115  //- Runtime type information
116  TypeName("LES");
117 
118 
119  // Declare run-time constructor selection table
120 
122  (
123  autoPtr,
124  LESModel,
125  dictionary,
126  (
127  const alphaField& alpha,
128  const rhoField& rho,
129  const volVectorField& U,
130  const surfaceScalarField& alphaRhoPhi,
131  const surfaceScalarField& phi,
132  const transportModel& transport,
133  const word& propertiesName
134  ),
135  (alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName)
136  );
137 
138 
139  // Constructors
140 
141  //- Construct from components
142  LESModel
143  (
144  const word& type,
145  const alphaField& alpha,
146  const rhoField& rho,
147  const volVectorField& U,
148  const surfaceScalarField& alphaRhoPhi,
149  const surfaceScalarField& phi,
150  const transportModel& transport,
151  const word& propertiesName
152  );
153 
154 
155  // Selectors
156 
157  //- Return a reference to the selected LES model
158  static autoPtr<LESModel> New
159  (
160  const alphaField& alpha,
161  const rhoField& rho,
162  const volVectorField& U,
163  const surfaceScalarField& alphaRhoPhi,
164  const surfaceScalarField& phi,
165  const transportModel& transport,
166  const word& propertiesName = turbulenceModel::propertiesName
167  );
168 
169 
170  //- Destructor
171  virtual ~LESModel()
172  {}
173 
174 
175  // Member Functions
176 
177  //- Read model coefficients if they have changed
178  virtual bool read();
179 
180 
181  // Access
182 
183  //- Const access to the coefficients dictionary
184  virtual const dictionary& coeffDict() const
185  {
186  return coeffDict_;
187  }
188 
189  //- Return the lower allowable limit for k (default: small)
190  const dimensionedScalar& kMin() const
191  {
192  return kMin_;
193  }
194 
195  //- Allow kMin to be changed
197  {
198  return kMin_;
199  }
200 
201  //- Access function to filter width
202  inline const volScalarField& delta() const
203  {
204  return delta_();
205  }
206 
207 
208  //- Return the effective viscosity
209  virtual tmp<volScalarField> nuEff() const
210  {
211  return tmp<volScalarField>
212  (
213  new volScalarField
214  (
215  IOobject::groupName("nuEff", this->alphaRhoPhi_.group()),
216  this->nut() + this->nu()
217  )
218  );
219  }
220 
221  //- Return the effective viscosity on patch
222  virtual tmp<scalarField> nuEff(const label patchi) const
223  {
224  return this->nut(patchi) + this->nu(patchi);
225  }
226 
227  //- Solve the turbulence equations and correct the turbulence viscosity
228  virtual void correct();
229 };
230 
231 
232 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
233 
234 } // End namespace Foam
235 
236 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
237 
238 #ifdef NoRepository
239  #include "LESModel.C"
240 #endif
241 
242 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
243 
244 #endif
245 
246 // ************************************************************************* //
BasicTurbulenceModel::transportModel transportModel
Definition: LESModel.H:111
Templated abstract base class for LES SGS models.
Definition: LESModel.H:56
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
surfaceScalarField & phi
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
BasicTurbulenceModel::alphaField alphaField
Definition: LESModel.H:109
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: LESModel.C:31
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
Definition: LESModel.H:208
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none.
Definition: Switch.H:60
dimensionedScalar epsilonMin_
Lower limit of epsilon.
Definition: LESModel.H:81
dictionary coeffDict_
Model coefficients dictionary.
Definition: LESModel.H:75
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
declareRunTimeSelectionTable(autoPtr, LESModel, dictionary,(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName),(alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName))
Switch turbulence_
Turbulence on/off flag.
Definition: LESModel.H:69
virtual ~LESModel()
Destructor.
Definition: LESModel.H:170
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
dictionary LESDict_
LES coefficients dictionary.
Definition: LESModel.H:66
const volScalarField & delta() const
Access function to filter width.
Definition: LESModel.H:201
autoPtr< Foam::LESdelta > delta_
Run-time selectable delta model.
Definition: LESModel.H:87
Switch printCoeffs_
Flag to print the model coeffs at run-time.
Definition: LESModel.H:72
BasicTurbulenceModel::rhoField rhoField
Definition: LESModel.H:110
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: LESModel.C:202
TypeName("LES")
Runtime type information.
const dimensionedScalar & kMin() const
Return the lower allowable limit for k (default: small)
Definition: LESModel.H:189
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
dimensionedScalar omegaMin_
Lower limit for omega.
Definition: LESModel.H:84
label patchi
U
Definition: pEqn.H:72
virtual bool read()
Read model coefficients if they have changed.
Definition: LESModel.C:179
dimensionedScalar kMin_
Lower limit of k.
Definition: LESModel.H:78
virtual const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: LESModel.H:183
static autoPtr< LESModel > New(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName)
Return a reference to the selected LES model.
Definition: LESModel.C:126
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
volScalarField & nu
scalar nut
Namespace for OpenFOAM.