All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2020 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 "MomentumTransportModel.H"
45 #include "LESdelta.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 /*---------------------------------------------------------------------------*\
53  Class LESModel Declaration
54 \*---------------------------------------------------------------------------*/
55 
56 template<class BasicMomentumTransportModel>
57 class LESModel
58 :
59  public BasicMomentumTransportModel
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 BasicMomentumTransportModel::alphaField alphaField;
100  typedef typename BasicMomentumTransportModel::rhoField rhoField;
101  typedef typename BasicMomentumTransportModel::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  ),
123  (alpha, rho, U, alphaRhoPhi, phi, transport)
124  );
125 
126 
127  // Constructors
128 
129  //- Construct from components
130  LESModel
131  (
132  const word& type,
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  );
140 
141  //- Disallow default bitwise copy construction
142  LESModel(const LESModel&) = delete;
143 
144 
145  // Selectors
146 
147  //- Return a reference to the selected LES model
148  static autoPtr<LESModel> New
149  (
150  const alphaField& alpha,
151  const rhoField& rho,
152  const volVectorField& U,
153  const surfaceScalarField& alphaRhoPhi,
154  const surfaceScalarField& phi,
155  const transportModel& transport
156  );
157 
158 
159  //- Destructor
160  virtual ~LESModel()
161  {}
162 
163 
164  // Member Functions
165 
166  //- Read model coefficients if they have changed
167  virtual bool read();
168 
169 
170  // Access
171 
172  //- Const access to the coefficients dictionary
173  virtual const dictionary& coeffDict() const
174  {
175  return coeffDict_;
176  }
177 
178  //- Return the lower allowable limit for k (default: small)
179  const dimensionedScalar& kMin() const
180  {
181  return kMin_;
182  }
183 
184  //- Allow kMin to be changed
186  {
187  return kMin_;
188  }
189 
190  //- Access function to filter width
191  inline const volScalarField& delta() const
192  {
193  return delta_();
194  }
195 
196 
197  //- Return the effective viscosity
198  virtual tmp<volScalarField> nuEff() const
199  {
200  return volScalarField::New
201  (
202  IOobject::groupName("nuEff", this->alphaRhoPhi_.group()),
203  this->nut() + this->nu()
204  );
205  }
206 
207  //- Return the effective viscosity on patch
208  virtual tmp<scalarField> nuEff(const label patchi) const
209  {
210  return this->nut(patchi) + this->nu(patchi);
211  }
212 
213  //- Solve the turbulence equations and correct the turbulence viscosity
214  virtual void correct();
215 
216 
217  // Member Operators
218 
219  //- Disallow default bitwise assignment
220  void operator=(const LESModel&) = delete;
221 };
222 
223 
224 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
225 
226 } // End namespace Foam
227 
228 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
229 
230 #ifdef NoRepository
231  #include "LESModel.C"
232 #endif
233 
234 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
235 
236 #endif
237 
238 // ************************************************************************* //
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
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
declareRunTimeSelectionTable(autoPtr, LESModel, dictionary,(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport),(alpha, rho, U, alphaRhoPhi, phi, transport))
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:197
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
phi
Definition: pEqn.H:104
Switch turbulence_
Turbulence on/off flag.
Definition: LESModel.H:69
virtual ~LESModel()
Destructor.
Definition: LESModel.H:159
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:190
autoPtr< Foam::LESdelta > delta_
Run-time selectable delta model.
Definition: LESModel.H:87
LESModel(const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
Construct from components.
Definition: LESModel.C:44
Switch printCoeffs_
Flag to print the model coeffs at run-time.
Definition: LESModel.H:72
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: LESModel.C:197
TypeName("LES")
Runtime type information.
BasicMomentumTransportModel::rhoField rhoField
Definition: LESModel.H:99
const dimensionedScalar & kMin() const
Return the lower allowable limit for k (default: small)
Definition: LESModel.H:178
dimensionedScalar omegaMin_
Lower limit for omega.
Definition: LESModel.H:84
BasicMomentumTransportModel::transportModel transportModel
Definition: LESModel.H:100
label patchi
U
Definition: pEqn.H:72
virtual bool read()
Read model coefficients if they have changed.
Definition: LESModel.C:174
dimensionedScalar kMin_
Lower limit of k.
Definition: LESModel.H:78
static autoPtr< LESModel > New(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
Return a reference to the selected LES model.
Definition: LESModel.C:124
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:172
BasicMomentumTransportModel::alphaField alphaField
Definition: LESModel.H:98
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
scalar nut
Namespace for OpenFOAM.