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-2024 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  //- Turbulence on/off flag
72 
73  //- Lower limit of k
75 
76  //- Upper limit coefficient for nut
78 
79  //- Run-time selectable generalised Newtonian viscosity model
82 
83  //- Run-time selectable delta model
85 
86 
87  // Protected member functions
88 
89  //- Const access to the LES dictionary
90  const dictionary& LESDict() const;
91 
92  //- Const access to the coefficients dictionary
93  const dictionary& coeffDict() const;
94 
95 
96 public:
97 
98  typedef typename BasicMomentumTransportModel::alphaField alphaField;
99  typedef typename BasicMomentumTransportModel::rhoField rhoField;
100 
101 
102  //- Runtime type information
103  TypeName("LES");
104 
105 
106  // Declare run-time constructor selection table
107 
109  (
110  autoPtr,
111  LESModel,
112  dictionary,
113  (
114  const alphaField& alpha,
115  const rhoField& rho,
116  const volVectorField& U,
117  const surfaceScalarField& alphaRhoPhi,
118  const surfaceScalarField& phi,
119  const viscosity& viscosity
120  ),
121  (alpha, rho, U, alphaRhoPhi, phi, viscosity)
122  );
123 
124 
125  // Constructors
126 
127  //- Construct from components
128  LESModel
129  (
130  const word& type,
131  const alphaField& alpha,
132  const rhoField& rho,
133  const volVectorField& U,
134  const surfaceScalarField& alphaRhoPhi,
135  const surfaceScalarField& phi,
136  const viscosity& viscosity
137  );
138 
139  //- Disallow default bitwise copy construction
140  LESModel(const LESModel&) = delete;
141 
142 
143  // Selectors
144 
145  //- Return a reference to the selected LES model
146  static autoPtr<LESModel> New
147  (
148  const alphaField& alpha,
149  const rhoField& rho,
150  const volVectorField& U,
151  const surfaceScalarField& alphaRhoPhi,
152  const surfaceScalarField& phi,
153  const viscosity& viscosity
154  );
155 
156 
157  //- Destructor
158  virtual ~LESModel()
159  {}
160 
161 
162  // Member Functions
163 
164  //- Read model coefficients if they have changed
165  virtual bool read();
166 
167  //- Access function to filter width
168  inline const volScalarField& delta() const
169  {
170  return delta_();
171  }
172 
173  //- Return the laminar viscosity
174  virtual tmp<volScalarField> nu() const
175  {
176  return viscosityModel_->nu();
177  }
178 
179  //- Return the laminar viscosity on patchi
180  virtual tmp<scalarField> nu(const label patchi) const
181  {
182  return viscosityModel_->nu(patchi);
183  }
184 
185  //- Return the effective viscosity
186  virtual tmp<volScalarField> nuEff() const
187  {
188  return volScalarField::New
189  (
190  this->groupName("nuEff"),
191  this->nut() + this->nu()
192  );
193  }
194 
195  //- Return the effective viscosity on patch
196  virtual tmp<scalarField> nuEff(const label patchi) const
197  {
198  return this->nut(patchi) + this->nu(patchi);
199  }
200 
201  //- Predict the turbulence transport coefficients if possible
202  // without solving turbulence transport model equations
203  virtual void predict()
204  {}
205 
206  //- Solve the turbulence equations and correct the turbulence viscosity
207  virtual void correct();
208 
209 
210  // Member Operators
211 
212  //- Disallow default bitwise assignment
213  void operator=(const LESModel&) = delete;
214 };
215 
216 
217 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
218 
219 } // End namespace Foam
220 
221 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
222 
223 #ifdef NoRepository
224  #include "LESModel.C"
225 #endif
226 
227 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
228 
229 #endif
230 
231 // ************************************************************************* //
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
Templated abstract base class for LES SGS models.
Definition: LESModel.H:63
BasicMomentumTransportModel::alphaField alphaField
Definition: LESModel.H:97
const dictionary & LESDict() const
Const access to the LES dictionary.
Definition: LESModel.C: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:33
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
Definition: LESModel.H:185
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: LESModel.C:191
dimensionedScalar nutMaxCoeff_
Upper limit coefficient for nut.
Definition: LESModel.H:76
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))
autoPtr< laminarModels::generalisedNewtonianViscosityModel > viscosityModel_
Run-time selectable generalised Newtonian viscosity model.
Definition: LESModel.H:80
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
Definition: LESModel.H:173
virtual void predict()
Predict the turbulence transport coefficients if possible.
Definition: LESModel.H:202
autoPtr< Foam::LESdelta > delta_
Run-time selectable delta model.
Definition: LESModel.H:83
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: LESModel.C:164
Switch turbulence_
Turbulence on/off flag.
Definition: LESModel.H:70
TypeName("LES")
Runtime type information.
virtual ~LESModel()
Destructor.
Definition: LESModel.H:157
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:99
dimensionedScalar kMin_
Lower limit of k.
Definition: LESModel.H:73
const volScalarField & delta() const
Access function to filter width.
Definition: LESModel.H:167
virtual bool read()
Read model coefficients if they have changed.
Definition: LESModel.C:171
void operator=(const LESModel &)=delete
Disallow default bitwise assignment.
BasicMomentumTransportModel::rhoField rhoField
Definition: LESModel.H:98
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
A class for managing temporary objects.
Definition: tmp.H:55
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
const scalar nut
label patchi
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Namespace for OpenFOAM.
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
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488