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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Namespace
27  Foam::LESModels
28 
29 Description
30  Namespace for LES SGS models.
31 
32 Class
33  Foam::LESModel
34 
35 Description
36  Templated abstract base class for LES SGS models
37 
38 SourceFiles
39  LESModel.C
40 
41 \*---------------------------------------------------------------------------*/
42 
43 #ifndef LESModel_H
44 #define LESModel_H
45 
46 #include "TurbulenceModel.H"
47 #include "LESdelta.H"
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 namespace Foam
52 {
53 
54 /*---------------------------------------------------------------------------*\
55  Class LESModel Declaration
56 \*---------------------------------------------------------------------------*/
57 
58 template<class BasicTurbulenceModel>
59 class LESModel
60 :
61  public BasicTurbulenceModel
62 {
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  //- Lower limit of k
82 
83  //- Lower limit of epsilon
85 
86  //- Lower limit for omega
88 
89  //- Run-time selectable delta model
91 
92 
93  // Protected Member Functions
94 
95  //- Print model coefficients
96  virtual void printCoeffs(const word& type);
97 
98  //- No copy construct
99  LESModel(const LESModel&) = delete;
100 
101  //- No copy assignment
102  void operator=(const LESModel&) = delete;
103 
104 
105 public:
106 
107  typedef typename BasicTurbulenceModel::alphaField alphaField;
108  typedef typename BasicTurbulenceModel::rhoField rhoField;
109  typedef typename BasicTurbulenceModel::transportModel transportModel;
110 
111 
112  //- Runtime type information
113  TypeName("LES");
114 
115 
116  // Declare run-time constructor selection table
117 
119  (
120  autoPtr,
121  LESModel,
122  dictionary,
123  (
124  const alphaField& alpha,
125  const rhoField& rho,
126  const volVectorField& U,
127  const surfaceScalarField& alphaRhoPhi,
128  const surfaceScalarField& phi,
129  const transportModel& transport,
130  const word& propertiesName
131  ),
132  (alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName)
133  );
134 
135 
136  // Constructors
137 
138  //- Construct from components
139  LESModel
140  (
141  const word& type,
142  const alphaField& alpha,
143  const rhoField& rho,
144  const volVectorField& U,
145  const surfaceScalarField& alphaRhoPhi,
146  const surfaceScalarField& phi,
147  const transportModel& transport,
148  const word& propertiesName
149  );
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 transportModel& transport,
163  const word& propertiesName = turbulenceModel::propertiesName
164  );
165 
166 
167  //- Destructor
168  virtual ~LESModel() = default;
169 
170 
171  // Member Functions
172 
173  //- Read model coefficients if they have changed
174  virtual bool read();
175 
176 
177  // Access
178 
179  //- Const access to the coefficients dictionary
180  virtual const dictionary& coeffDict() const
181  {
182  return coeffDict_;
183  }
184 
185  //- Return the lower allowable limit for k (default: SMALL)
186  const dimensionedScalar& kMin() const
187  {
188  return kMin_;
189  }
190 
191  //- Allow kMin to be changed
193  {
194  return kMin_;
195  }
196 
197  //- Access function to filter width
198  inline const volScalarField& delta() const
199  {
200  return *delta_;
201  }
202 
203 
204  //- Return the effective viscosity
205  virtual tmp<volScalarField> nuEff() const
206  {
207  return tmp<volScalarField>
208  (
209  new volScalarField
210  (
211  IOobject::groupName("nuEff", this->alphaRhoPhi_.group()),
212  this->nut() + this->nu()
213  )
214  );
215  }
216 
217  //- Return the effective viscosity on patch
218  virtual tmp<scalarField> nuEff(const label patchi) const
219  {
220  return this->nut(patchi) + this->nu(patchi);
221  }
222 
223  //- Solve the turbulence equations and correct the turbulence viscosity
224  virtual void correct();
225 };
226 
227 
228 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
229 
230 } // End namespace Foam
231 
232 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
233 
234 #ifdef NoRepository
235  #include "LESModel.C"
236 #endif
237 
238 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
239 
240 #endif
241 
242 // ************************************************************************* //
Foam::LESModel::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: LESModel.H:107
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:62
Foam::LESModel::coeffDict
virtual const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: LESModel.H:179
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:129
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:191
Foam::LESModel::omegaMin_
dimensionedScalar omegaMin_
Lower limit for omega.
Definition: LESModel.H:86
Foam::LESModel::delta
const volScalarField & delta() const
Access function to filter width.
Definition: LESModel.H:197
Foam::LESModel::delta_
autoPtr< Foam::LESdelta > delta_
Run-time selectable delta model.
Definition: LESModel.H:89
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:206
Foam::LESModel::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: LESModel.H:108
Foam::LESModel::kMin_
dimensionedScalar kMin_
Lower limit of k.
Definition: LESModel.H:80
Foam::LESModel::read
virtual bool read()
Read model coefficients if they have changed.
Definition: LESModel.C:185
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:121
Foam::LESModel::epsilonMin_
dimensionedScalar epsilonMin_
Lower limit of epsilon.
Definition: LESModel.H:83
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:217
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:58
Foam::LESModel::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: LESModel.H:106
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:185
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:204