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-2021 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  with support for generalised Newtonian viscosity models including
37  strain-rate dependency.
38 
39 SourceFiles
40  LESModel.C
41 
42 \*---------------------------------------------------------------------------*/
43 
44 #ifndef LESModel_H
45 #define LESModel_H
46 
47 #include "momentumTransportModel.H"
49 #include "LESdelta.H"
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 namespace Foam
54 {
55 
56 /*---------------------------------------------------------------------------*\
57  Class LESModel Declaration
58 \*---------------------------------------------------------------------------*/
59 
60 template<class BasicMomentumTransportModel>
61 class LESModel
62 :
63  public BasicMomentumTransportModel
64 {
65 
66 protected:
67 
68  // Protected data
69 
70  //- LES coefficients dictionary
72 
73  //- Turbulence on/off flag
75 
76  //- Flag to print the model coeffs at run-time
78 
79  //- Model coefficients dictionary
81 
82  //- Lower limit of k
84 
85  //- Lower limit of epsilon
87 
88  //- Lower limit for omega
90 
91  //- Run-time selectable generalised Newtonian viscosity model
94 
95  //- Run-time selectable delta model
97 
98 
99  // Protected Member Functions
100 
101  //- Print model coefficients
102  virtual void printCoeffs(const word& type);
103 
104 
105 public:
107  typedef typename BasicMomentumTransportModel::alphaField alphaField;
108  typedef typename BasicMomentumTransportModel::rhoField rhoField;
109 
110 
111  //- Runtime type information
112  TypeName("LES");
113 
114 
115  // Declare run-time constructor selection table
116 
118  (
119  autoPtr,
120  LESModel,
121  dictionary,
122  (
123  const alphaField& alpha,
124  const rhoField& rho,
125  const volVectorField& U,
126  const surfaceScalarField& alphaRhoPhi,
127  const surfaceScalarField& phi,
128  const viscosity& viscosity
129  ),
130  (alpha, rho, U, alphaRhoPhi, phi, viscosity)
131  );
132 
133 
134  // Constructors
135 
136  //- Construct from components
137  LESModel
138  (
139  const word& type,
140  const alphaField& alpha,
141  const rhoField& rho,
142  const volVectorField& U,
143  const surfaceScalarField& alphaRhoPhi,
144  const surfaceScalarField& phi,
145  const viscosity& viscosity
146  );
147 
148  //- Disallow default bitwise copy construction
149  LESModel(const LESModel&) = delete;
150 
151 
152  // Selectors
153 
154  //- Return a reference to the selected LES model
155  static autoPtr<LESModel> New
156  (
157  const alphaField& alpha,
158  const rhoField& rho,
159  const volVectorField& U,
160  const surfaceScalarField& alphaRhoPhi,
161  const surfaceScalarField& phi,
162  const viscosity& viscosity
163  );
164 
165 
166  //- Destructor
167  virtual ~LESModel()
168  {}
169 
170 
171  // Member Functions
172 
173  //- Read model coefficients if they have changed
174  virtual bool read();
175 
176  //- Const access to the coefficients dictionary
177  virtual const dictionary& coeffDict() const
178  {
179  return coeffDict_;
180  }
181 
182  //- Return the lower allowable limit for k (default: small)
183  const dimensionedScalar& kMin() const
184  {
185  return kMin_;
186  }
187 
188  //- Allow kMin to be changed
190  {
191  return kMin_;
192  }
193 
194  //- Access function to filter width
195  inline const volScalarField& delta() const
196  {
197  return delta_();
198  }
199 
200  //- Return the laminar viscosity
201  virtual tmp<volScalarField> nu() const
202  {
203  return viscosityModel_->nu();
204  }
205 
206  //- Return the laminar viscosity on patchi
207  virtual tmp<scalarField> nu(const label patchi) const
208  {
209  return viscosityModel_->nu(patchi);
210  }
211 
212  //- Return the effective viscosity
213  virtual tmp<volScalarField> nuEff() const
214  {
215  return volScalarField::New
216  (
217  IOobject::groupName("nuEff", this->alphaRhoPhi_.group()),
218  this->nut() + this->nu()
219  );
220  }
221 
222  //- Return the effective viscosity on patch
223  virtual tmp<scalarField> nuEff(const label patchi) const
224  {
225  return this->nut(patchi) + this->nu(patchi);
226  }
227 
228  //- Solve the turbulence equations and correct the turbulence viscosity
229  virtual void correct();
230 
231 
232  // Member Operators
233 
234  //- Disallow default bitwise assignment
235  void operator=(const LESModel&) = delete;
236 };
237 
238 
239 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
240 
241 } // End namespace Foam
242 
243 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
244 
245 #ifdef NoRepository
246  #include "LESModel.C"
247 #endif
248 
249 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
250 
251 #endif
252 
253 // ************************************************************************* //
Templated abstract base class for LES SGS models.
Definition: LESModel.H:60
U
Definition: pEqn.H:72
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
LESModel(const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
Construct from components.
Definition: LESModel.C:45
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: LESModel.C:32
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
Definition: LESModel.H:212
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:85
void operator=(const LESModel &)=delete
Disallow default bitwise assignment.
dictionary coeffDict_
Model coefficients dictionary.
Definition: LESModel.H:79
Switch turbulence_
Turbulence on/off flag.
Definition: LESModel.H:73
virtual ~LESModel()
Destructor.
Definition: LESModel.H:166
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
Definition: LESModel.H:200
phi
Definition: correctPhi.H:3
static autoPtr< LESModel > New(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
Return a reference to the selected LES model.
Definition: LESModel.C:145
dictionary LESDict_
LES coefficients dictionary.
Definition: LESModel.H:70
autoPtr< laminarModels::generalisedNewtonianViscosityModel > viscosityModel_
Run-time selectable generalised Newtonian viscosity model.
Definition: LESModel.H:92
const volScalarField & delta() const
Access function to filter width.
Definition: LESModel.H:194
autoPtr< Foam::LESdelta > delta_
Run-time selectable delta model.
Definition: LESModel.H:95
Switch printCoeffs_
Flag to print the model coeffs at run-time.
Definition: LESModel.H:76
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: LESModel.C:217
TypeName("LES")
Runtime type information.
BasicMomentumTransportModel::rhoField rhoField
Definition: LESModel.H:107
const dimensionedScalar & kMin() const
Return the lower allowable limit for k (default: small)
Definition: LESModel.H:182
Abstract base class for all fluid physical properties.
Definition: viscosity.H:49
dimensionedScalar omegaMin_
Lower limit for omega.
Definition: LESModel.H:88
label patchi
virtual bool read()
Read model coefficients if they have changed.
Definition: LESModel.C:194
dimensionedScalar kMin_
Lower limit of k.
Definition: LESModel.H:82
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:176
BasicMomentumTransportModel::alphaField alphaField
Definition: LESModel.H:106
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
declareRunTimeSelectionTable(autoPtr, LESModel, dictionary,(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity),(alpha, rho, U, alphaRhoPhi, phi, viscosity))
A class for managing temporary objects.
Definition: PtrList.H:53
scalar nut
Namespace for OpenFOAM.