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 | Copyright (C) 2013-2015 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 SourceFiles
37  LESModel.C
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #ifndef LESModel_H
42 #define LESModel_H
43 
44 #include "TurbulenceModel.H"
45 #include "LESdelta.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 /*---------------------------------------------------------------------------*\
53  Class LESModel Declaration
54 \*---------------------------------------------------------------------------*/
55 
56 template<class BasicTurbulenceModel>
57 class LESModel
58 :
59  public BasicTurbulenceModel
60 {
61 
62 protected:
63 
64  // Protected data
65 
66  //- LES coefficients dictionary
68 
69  //- Turbulence on/off flag
71 
72  //- Flag to print the model coeffs at run-time
74 
75  //- Model coefficients dictionary
77 
78  //- Lower limit of k
80 
81  //- Lower limit of omega
83 
84  //- Run-time selectable delta model
86 
87 
88  // Protected Member Functions
89 
90  //- Print model coefficients
91  virtual void printCoeffs(const word& type);
92 
93 
94 private:
95 
96  // Private Member Functions
97 
98  //- Disallow default bitwise copy construct
99  LESModel(const LESModel&);
100 
101  //- Disallow default bitwise assignment
102  void operator=(const LESModel&);
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()
169  {}
170 
171 
172  // Member Functions
173 
174  //- Read model coefficients if they have changed
175  virtual bool read();
176 
177 
178  // Access
179 
180  //- Const access to the coefficients dictionary
181  virtual const dictionary& coeffDict() const
182  {
183  return coeffDict_;
184  }
185 
186  //- Return the lower allowable limit for k (default: SMALL)
187  const dimensionedScalar& kMin() const
188  {
189  return kMin_;
190  }
191 
192  //- Allow kMin to be changed
194  {
195  return kMin_;
196  }
197 
198  //- Access function to filter width
199  inline const volScalarField& delta() const
200  {
201  return delta_();
202  }
203 
204 
205  //- Return the effective viscosity
206  virtual tmp<volScalarField> nuEff() const
207  {
208  return tmp<volScalarField>
209  (
210  new volScalarField
211  (
212  IOobject::groupName("nuEff", this->U_.group()),
213  this->nut() + this->nu()
214  )
215  );
216  }
217 
218  //- Return the effective viscosity on patch
219  virtual tmp<scalarField> nuEff(const label patchi) const
220  {
221  return this->nut(patchi) + this->nu(patchi);
222  }
223 
224  //- Solve the turbulence equations and correct the turbulence viscosity
225  virtual void correct();
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 // ************************************************************************* //
Foam::LESModel::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: LESModel.H:107
Foam::IOobject::groupName
static word groupName(Name name, const word &group)
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:60
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 string.
Definition: word.H:59
Foam::LESModel::coeffDict
virtual const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: LESModel.H:180
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::LESModel::turbulence_
Switch turbulence_
Turbulence on/off flag.
Definition: LESModel.H:69
LESModel.C
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
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:115
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:97
Foam::LESModel::operator=
void operator=(const LESModel &)
Disallow default bitwise assignment.
Foam::LESModel::kMin
dimensionedScalar & kMin()
Allow kMin to be changed.
Definition: LESModel.H:192
Foam::LESModel::omegaMin_
dimensionedScalar omegaMin_
Lower limit of omega.
Definition: LESModel.H:81
Foam::LESModel::delta
const volScalarField & delta() const
Access function to filter width.
Definition: LESModel.H:198
Foam::LESModel::delta_
autoPtr< Foam::LESdelta > delta_
Run-time selectable delta model.
Definition: LESModel.H:84
U
U
Definition: pEqn.H:46
nu
volScalarField & nu
Definition: readMechanicalProperties.H:179
Foam::LESModel::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: LESModel.C:194
Foam::LESModel::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: LESModel.H:108
Foam::LESModel::LESModel
LESModel(const LESModel &)
Disallow default bitwise copy construct.
Foam::LESModel::kMin_
dimensionedScalar kMin_
Lower limit of k.
Definition: LESModel.H:78
Foam::label
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
Foam::LESModel::read
virtual bool read()
Read model coefficients if they have changed.
Definition: LESModel.C:168
nut
nut
Definition: createTDFields.H:71
Foam::LESModel::printCoeffs
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: LESModel.C:31
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
rho
rho
Definition: pEqn.H:3
LESdelta.H
Foam::autoPtr< Foam::LESdelta >
Foam::LESModel::~LESModel
virtual ~LESModel()
Destructor.
Definition: LESModel.H:167
Foam::LESModel::nuEff
virtual tmp< scalarField > nuEff(const label patchi) const
Return the effective viscosity on patch.
Definition: LESModel.H:218
TurbulenceModel.H
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::LESModel::LESDict_
dictionary LESDict_
LES coefficients dictionary.
Definition: LESModel.H:66
Foam::LESModel
Templated abstract base class for LES SGS models.
Definition: LESModel.H:56
Foam::LESModel::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: LESModel.H:106
Foam::LESModel::TypeName
TypeName("LES")
Runtime type information.
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::LESModel::kMin
const dimensionedScalar & kMin() const
Return the lower allowable limit for k (default: SMALL)
Definition: LESModel.H:186
Foam::LESModel::printCoeffs_
Switch printCoeffs_
Flag to print the model coeffs at run-time.
Definition: LESModel.H:72
Foam::LESModel::coeffDict_
dictionary coeffDict_
Model coefficients dictionary.
Definition: LESModel.H:75
Foam::LESModel::nuEff
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
Definition: LESModel.H:205