LESModel.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2013-2017 OpenFOAM Foundation
9  Copyright (C) 2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Namespace
28  Foam::LESModels
29 
30 Description
31  Namespace for LES SGS models.
32 
33 Class
34  Foam::LESModel
35 
36 Description
37  Templated abstract base class for LES SGS models
38 
39 SourceFiles
40  LESModel.C
41 
42 \*---------------------------------------------------------------------------*/
43 
44 #ifndef LESModel_H
45 #define LESModel_H
46 
47 #include "TurbulenceModel.H"
48 #include "LESdelta.H"
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 namespace Foam
53 {
54 
55 /*---------------------------------------------------------------------------*\
56  Class LESModel Declaration
57 \*---------------------------------------------------------------------------*/
58 
59 template<class BasicTurbulenceModel>
60 class LESModel
61 :
62  public BasicTurbulenceModel
63 {
64 protected:
65 
66  // Protected Data
67 
68  //- LES coefficients dictionary
70 
71  //- Turbulence on/off flag
73 
74  //- Flag to print the model coeffs at run-time
76 
77  //- Model coefficients dictionary
79 
80  //- Empirical model constant
82 
83  //- Lower limit of k
85 
86  //- Lower limit of epsilon
88 
89  //- Lower limit for omega
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  //- No copy construct
102  LESModel(const LESModel&) = delete;
103 
104  //- No copy assignment
105  void operator=(const LESModel&) = delete;
106 
107 
108 public:
109 
110  typedef typename BasicTurbulenceModel::alphaField alphaField;
111  typedef typename BasicTurbulenceModel::rhoField rhoField;
112  typedef typename BasicTurbulenceModel::transportModel transportModel;
113 
114 
115  //- Runtime type information
116  TypeName("LES");
117 
118 
119  // Declare run-time constructor selection table
120 
122  (
123  autoPtr,
124  LESModel,
125  dictionary,
126  (
127  const alphaField& alpha,
128  const rhoField& rho,
129  const volVectorField& U,
130  const surfaceScalarField& alphaRhoPhi,
131  const surfaceScalarField& phi,
132  const transportModel& transport,
133  const word& propertiesName
134  ),
135  (alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName)
136  );
137 
138 
139  // Constructors
140 
141  //- Construct from components
142  LESModel
143  (
144  const word& type,
145  const alphaField& alpha,
146  const rhoField& rho,
147  const volVectorField& U,
148  const surfaceScalarField& alphaRhoPhi,
149  const surfaceScalarField& phi,
150  const transportModel& transport,
151  const word& propertiesName
152  );
153 
154 
155  // Selectors
156 
157  //- Return a reference to the selected LES model
158  static autoPtr<LESModel> New
159  (
160  const alphaField& alpha,
161  const rhoField& rho,
162  const volVectorField& U,
163  const surfaceScalarField& alphaRhoPhi,
164  const surfaceScalarField& phi,
165  const transportModel& transport,
166  const word& propertiesName = turbulenceModel::propertiesName
167  );
168 
169 
170  //- Destructor
171  virtual ~LESModel() = default;
172 
173 
174  // Member Functions
175 
176  //- Read model coefficients if they have changed
177  virtual bool read();
178 
179 
180  // Access
181 
182  //- Const access to the coefficients dictionary
183  virtual const dictionary& coeffDict() const
184  {
185  return coeffDict_;
186  }
187 
188  // Return the empirical model constant
189  const dimensionedScalar& Ce() const noexcept
190  {
191  return Ce_;
192  }
193 
194  //- Return the lower allowable limit for k (default: SMALL)
195  const dimensionedScalar& kMin() const
196  {
197  return kMin_;
198  }
199 
200  //- Allow kMin to be changed
202  {
203  return kMin_;
204  }
205 
206  //- Access function to filter width
207  inline const volScalarField& delta() const
208  {
209  return *delta_;
210  }
211 
212 
213  //- Return the effective viscosity
214  virtual tmp<volScalarField> nuEff() const
215  {
216  return tmp<volScalarField>
217  (
218  new volScalarField
219  (
220  IOobject::groupName("nuEff", this->alphaRhoPhi_.group()),
221  this->nut() + this->nu()
222  )
223  );
224  }
225 
226  //- Return the effective viscosity on patch
227  virtual tmp<scalarField> nuEff(const label patchi) const
228  {
229  return this->nut(patchi) + this->nu(patchi);
230  }
231 
232  //- Return the turbulence kinetic energy dissipation rate
233  virtual tmp<volScalarField> epsilon() const;
234 
235  //- Return the specific dissipation rate
236  virtual tmp<volScalarField> omega() const;
237 
238  //- Solve the turbulence equations and correct the turbulence viscosity
239  virtual void correct();
240 };
241 
242 
243 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
244 
245 } // End namespace Foam
246 
247 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
248 
249 #ifdef NoRepository
250  #include "LESModel.C"
251 #endif
252 
253 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
254 
255 #endif
256 
257 // ************************************************************************* //
Foam::LESModel::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: LESModel.H:110
Foam::LESModel::epsilon
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
Definition: LESModel.C:213
Foam::LESModel::omega
virtual tmp< volScalarField > omega() const
Return the specific dissipation rate.
Definition: LESModel.C:230
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:77
Foam::LESModel::declareRunTimeSelectionTable
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))
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::LESModel::coeffDict
virtual const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: LESModel.H:182
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::LESModel::turbulence_
Switch turbulence_
Turbulence on/off flag.
Definition: LESModel.H:71
LESModel.C
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::LESModel::New
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:133
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
rho
rho
Definition: readInitialConditions.H:88
Foam::LESModel::kMin
dimensionedScalar & kMin()
Allow kMin to be changed.
Definition: LESModel.H:200
Foam::LESModel::omegaMin_
dimensionedScalar omegaMin_
Lower limit for omega.
Definition: LESModel.H:89
Foam::LESModel::delta
const volScalarField & delta() const
Access function to filter width.
Definition: LESModel.H:206
Foam::LESModel::delta_
autoPtr< Foam::LESdelta > delta_
Run-time selectable delta model.
Definition: LESModel.H:92
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::LESModel::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: LESModel.C:249
Foam::LESModel::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: LESModel.H:111
Foam::LESModel::kMin_
dimensionedScalar kMin_
Lower limit of k.
Definition: LESModel.H:83
Foam::LESModel::read
virtual bool read()
Read model coefficients if they have changed.
Definition: LESModel.C:189
Foam::LESModel::printCoeffs
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: LESModel.C:34
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::LESModel::LESModel
LESModel(const LESModel &)=delete
No copy construct.
Foam::LESModel::~LESModel
virtual ~LESModel()=default
Destructor.
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::LESModel::epsilonMin_
dimensionedScalar epsilonMin_
Lower limit of epsilon.
Definition: LESModel.H:86
Foam::dimensioned< scalar >
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::LESModel::operator=
void operator=(const LESModel &)=delete
No copy assignment.
LESdelta.H
Foam::autoPtr< Foam::LESdelta >
U
U
Definition: pEqn.H:72
Foam::LESModel::nuEff
virtual tmp< scalarField > nuEff(const label patchi) const
Return the effective viscosity on patch.
Definition: LESModel.H:226
Foam::LESModel::Ce
const dimensionedScalar & Ce() const noexcept
Definition: LESModel.H:188
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
TurbulenceModel.H
Foam::LESModel::LESDict_
dictionary LESDict_
LES coefficients dictionary.
Definition: LESModel.H:68
Foam::LESModel
Templated abstract base class for LES SGS models.
Definition: LESModel.H:59
Foam::LESModel::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: LESModel.H:109
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Foam::LESModel::TypeName
TypeName("LES")
Runtime type information.
nut
scalar nut
Definition: evaluateNearWall.H:5
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::LESModel::kMin
const dimensionedScalar & kMin() const
Return the lower allowable limit for k (default: SMALL)
Definition: LESModel.H:194
Foam::LESModel::Ce_
dimensionedScalar Ce_
Empirical model constant.
Definition: LESModel.H:80
Foam::LESModel::printCoeffs_
Switch printCoeffs_
Flag to print the model coeffs at run-time.
Definition: LESModel.H:74
Foam::LESModel::coeffDict_
dictionary coeffDict_
Model coefficients dictionary.
Definition: LESModel.H:77
Foam::LESModel::nuEff
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
Definition: LESModel.H:213