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-2023 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  //- Upper limit coefficient for nut
87 
88  //- Run-time selectable generalised Newtonian viscosity model
91 
92  //- Run-time selectable delta model
94 
95 
96  // Protected Member Functions
97 
98  //- Print model coefficients
99  virtual void printCoeffs(const word& type);
100 
101 
102 public:
103 
104  typedef typename BasicMomentumTransportModel::alphaField alphaField;
105  typedef typename BasicMomentumTransportModel::rhoField rhoField;
106 
107 
108  //- Runtime type information
109  TypeName("LES");
110 
111 
112  // Declare run-time constructor selection table
113 
115  (
116  autoPtr,
117  LESModel,
118  dictionary,
119  (
120  const alphaField& alpha,
121  const rhoField& rho,
122  const volVectorField& U,
123  const surfaceScalarField& alphaRhoPhi,
124  const surfaceScalarField& phi,
125  const viscosity& viscosity
126  ),
127  (alpha, rho, U, alphaRhoPhi, phi, viscosity)
128  );
129 
130 
131  // Constructors
132 
133  //- Construct from components
134  LESModel
135  (
136  const word& type,
137  const alphaField& alpha,
138  const rhoField& rho,
139  const volVectorField& U,
140  const surfaceScalarField& alphaRhoPhi,
141  const surfaceScalarField& phi,
142  const viscosity& viscosity
143  );
144 
145  //- Disallow default bitwise copy construction
146  LESModel(const LESModel&) = delete;
147 
148 
149  // Selectors
150 
151  //- Return a reference to the selected LES model
152  static autoPtr<LESModel> New
153  (
154  const alphaField& alpha,
155  const rhoField& rho,
156  const volVectorField& U,
157  const surfaceScalarField& alphaRhoPhi,
158  const surfaceScalarField& phi,
159  const viscosity& viscosity
160  );
161 
162 
163  //- Destructor
164  virtual ~LESModel()
165  {}
166 
167 
168  // Member Functions
169 
170  //- Read model coefficients if they have changed
171  virtual bool read();
172 
173  //- Const access to the coefficients dictionary
174  virtual const dictionary& coeffDict() const
175  {
176  return coeffDict_;
177  }
178 
179  //- Access function to filter width
180  inline const volScalarField& delta() const
181  {
182  return delta_();
183  }
184 
185  //- Return the laminar viscosity
186  virtual tmp<volScalarField> nu() const
187  {
188  return viscosityModel_->nu();
189  }
190 
191  //- Return the laminar viscosity on patchi
192  virtual tmp<scalarField> nu(const label patchi) const
193  {
194  return viscosityModel_->nu(patchi);
195  }
196 
197  //- Return the effective viscosity
198  virtual tmp<volScalarField> nuEff() const
199  {
200  return volScalarField::New
201  (
202  this->groupName("nuEff"),
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  //- Predict the turbulence transport coefficients if possible
214  // without solving turbulence transport model equations
215  virtual void predict()
216  {}
217 
218  //- Solve the turbulence equations and correct the turbulence viscosity
219  virtual void correct();
220 
221 
222  // Member Operators
223 
224  //- Disallow default bitwise assignment
225  void operator=(const LESModel&) = delete;
226 };
227 
228 
229 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
230 
231 } // End namespace Foam
232 
233 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
234 
235 #ifdef NoRepository
236  #include "LESModel.C"
237 #endif
238 
239 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
240 
241 #endif
242 
243 // ************************************************************************* //
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, 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:103
dictionary coeffDict_
Model coefficients dictionary.
Definition: LESModel.H:79
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: LESModel.C:32
virtual const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: LESModel.H:173
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
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
Definition: LESModel.H:197
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: LESModel.C:207
dimensionedScalar nutMaxCoeff_
Upper limit coefficient for nut.
Definition: LESModel.H:85
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:89
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
Definition: LESModel.H:185
virtual void predict()
Predict the turbulence transport coefficients if possible.
Definition: LESModel.H:214
autoPtr< Foam::LESdelta > delta_
Run-time selectable delta model.
Definition: LESModel.H:92
dictionary LESDict_
LES coefficients dictionary.
Definition: LESModel.H:70
Switch turbulence_
Turbulence on/off flag.
Definition: LESModel.H:73
TypeName("LES")
Runtime type information.
virtual ~LESModel()
Destructor.
Definition: LESModel.H:163
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:134
dimensionedScalar kMin_
Lower limit of k.
Definition: LESModel.H:82
const volScalarField & delta() const
Access function to filter width.
Definition: LESModel.H:179
Switch printCoeffs_
Flag to print the model coeffs at run-time.
Definition: LESModel.H:76
virtual bool read()
Read model coefficients if they have changed.
Definition: LESModel.C:183
void operator=(const LESModel &)=delete
Disallow default bitwise assignment.
BasicMomentumTransportModel::rhoField rhoField
Definition: LESModel.H:104
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 keyword definitions, which are a keyword followed by any number of values (e....
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