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-2019 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 public:
98 
99  typedef typename BasicTurbulenceModel::alphaField alphaField;
100  typedef typename BasicTurbulenceModel::rhoField rhoField;
101  typedef typename BasicTurbulenceModel::transportModel transportModel;
102 
103 
104  //- Runtime type information
105  TypeName("LES");
106 
107 
108  // Declare run-time constructor selection table
109 
111  (
112  autoPtr,
113  LESModel,
114  dictionary,
115  (
116  const alphaField& alpha,
117  const rhoField& rho,
118  const volVectorField& U,
119  const surfaceScalarField& alphaRhoPhi,
120  const surfaceScalarField& phi,
121  const transportModel& transport,
122  const word& propertiesName
123  ),
124  (alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName)
125  );
126 
127 
128  // Constructors
129 
130  //- Construct from components
131  LESModel
132  (
133  const word& type,
134  const alphaField& alpha,
135  const rhoField& rho,
136  const volVectorField& U,
137  const surfaceScalarField& alphaRhoPhi,
138  const surfaceScalarField& phi,
139  const transportModel& transport,
140  const word& propertiesName
141  );
142 
143  //- Disallow default bitwise copy construction
144  LESModel(const LESModel&) = delete;
145 
146 
147  // Selectors
148 
149  //- Return a reference to the selected LES model
150  static autoPtr<LESModel> New
151  (
152  const alphaField& alpha,
153  const rhoField& rho,
154  const volVectorField& U,
155  const surfaceScalarField& alphaRhoPhi,
156  const surfaceScalarField& phi,
157  const transportModel& transport,
158  const word& propertiesName = turbulenceModel::propertiesName
159  );
160 
161 
162  //- Destructor
163  virtual ~LESModel()
164  {}
165 
166 
167  // Member Functions
168 
169  //- Read model coefficients if they have changed
170  virtual bool read();
171 
172 
173  // Access
174 
175  //- Const access to the coefficients dictionary
176  virtual const dictionary& coeffDict() const
177  {
178  return coeffDict_;
179  }
180 
181  //- Return the lower allowable limit for k (default: small)
182  const dimensionedScalar& kMin() const
183  {
184  return kMin_;
185  }
186 
187  //- Allow kMin to be changed
189  {
190  return kMin_;
191  }
192 
193  //- Access function to filter width
194  inline const volScalarField& delta() const
195  {
196  return delta_();
197  }
198 
199 
200  //- Return the effective viscosity
201  virtual tmp<volScalarField> nuEff() const
202  {
203  return volScalarField::New
204  (
205  IOobject::groupName("nuEff", this->alphaRhoPhi_.group()),
206  this->nut() + this->nu()
207  );
208  }
209 
210  //- Return the effective viscosity on patch
211  virtual tmp<scalarField> nuEff(const label patchi) const
212  {
213  return this->nut(patchi) + this->nu(patchi);
214  }
215 
216  //- Solve the turbulence equations and correct the turbulence viscosity
217  virtual void correct();
218 
219 
220  // Member Operators
221 
222  //- Disallow default bitwise assignment
223  void operator=(const LESModel&) = delete;
224 };
225 
226 
227 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
228 
229 } // End namespace Foam
230 
231 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232 
233 #ifdef NoRepository
234  #include "LESModel.C"
235 #endif
236 
237 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
238 
239 #endif
240 
241 // ************************************************************************* //
BasicTurbulenceModel::transportModel transportModel
Definition: LESModel.H:100
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:158
BasicTurbulenceModel::alphaField alphaField
Definition: LESModel.H:98
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:200
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/any.
Definition: Switch.H:60
dimensionedScalar epsilonMin_
Lower limit of epsilon.
Definition: LESModel.H:81
void operator=(const LESModel &)=delete
Disallow default bitwise assignment.
dictionary coeffDict_
Model coefficients dictionary.
Definition: LESModel.H:75
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
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:162
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:193
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:99
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:181
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
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
virtual const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: LESModel.H:175
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
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.
Definition: PtrList.H:53
volScalarField & nu
LESModel(const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
Construct from components.
Definition: LESModel.C:44
scalar nut
Namespace for OpenFOAM.