RASModel.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-2016 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 Class
25  Foam::RASModel
26 
27 Description
28  Templated abstract base class for RAS turbulence models
29 
30 SourceFiles
31  RASModel.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef RASModel_H
36 #define RASModel_H
37 
38 #include "TurbulenceModel.H"
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 
45 /*---------------------------------------------------------------------------*\
46  Class RASModel Declaration
47 \*---------------------------------------------------------------------------*/
48 
49 template<class BasicTurbulenceModel>
50 class RASModel
51 :
52  public BasicTurbulenceModel
53 {
54 
55 protected:
56 
57  // Protected data
58 
59  //- RAS coefficients dictionary
61 
62  //- Turbulence on/off flag
64 
65  //- Flag to print the model coeffs at run-time
67 
68  //- Model coefficients dictionary
70 
71  //- Lower limit of k
73 
74  //- Lower limit of epsilon
76 
77  //- Lower limit for omega
79 
80 
81  // Protected Member Functions
82 
83  //- Print model coefficients
84  virtual void printCoeffs(const word& type);
85 
86 
87 private:
88 
89  // Private Member Functions
90 
91  //- Disallow default bitwise copy construct
92  RASModel(const RASModel&);
93 
94  //- Disallow default bitwise assignment
95  void operator=(const RASModel&);
96 
97 
98 public:
99 
100  typedef typename BasicTurbulenceModel::alphaField alphaField;
101  typedef typename BasicTurbulenceModel::rhoField rhoField;
102  typedef typename BasicTurbulenceModel::transportModel transportModel;
103 
104 
105  //- Runtime type information
106  TypeName("RAS");
107 
108 
109  // Declare run-time constructor selection table
110 
112  (
113  autoPtr,
114  RASModel,
115  dictionary,
116  (
117  const alphaField& alpha,
118  const rhoField& rho,
119  const volVectorField& U,
120  const surfaceScalarField& alphaRhoPhi,
121  const surfaceScalarField& phi,
122  const transportModel& transport,
123  const word& propertiesName
124  ),
125  (alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName)
126  );
127 
128 
129  // Constructors
130 
131  //- Construct from components
132  RASModel
133  (
134  const word& type,
135  const alphaField& alpha,
136  const rhoField& rho,
137  const volVectorField& U,
138  const surfaceScalarField& alphaRhoPhi,
139  const surfaceScalarField& phi,
140  const transportModel& transport,
141  const word& propertiesName
142  );
143 
144 
145  // Selectors
146 
147  //- Return a reference to the selected RAS model
148  static autoPtr<RASModel> New
149  (
150  const alphaField& alpha,
151  const rhoField& rho,
152  const volVectorField& U,
153  const surfaceScalarField& alphaRhoPhi,
154  const surfaceScalarField& phi,
155  const transportModel& transport,
156  const word& propertiesName = turbulenceModel::propertiesName
157  );
158 
159 
160  //- Destructor
161  virtual ~RASModel()
162  {}
163 
164 
165  // Member Functions
166 
167  //- Read model coefficients if they have changed
168  virtual bool read();
169 
170 
171  // Access
172 
173  //- Return the lower allowable limit for k (default: SMALL)
174  const dimensionedScalar& kMin() const
175  {
176  return kMin_;
177  }
178 
179  //- Return the lower allowable limit for epsilon (default: SMALL)
180  const dimensionedScalar& epsilonMin() const
181  {
182  return epsilonMin_;
183  }
184 
185  //- Return the lower allowable limit for omega (default: SMALL)
186  const dimensionedScalar& omegaMin() const
187  {
188  return omegaMin_;
189  }
190 
191  //- Allow kMin to be changed
193  {
194  return kMin_;
195  }
196 
197  //- Allow epsilonMin to be changed
199  {
200  return epsilonMin_;
201  }
202 
203  //- Allow omegaMin to be changed
205  {
206  return omegaMin_;
207  }
208 
209  //- Const access to the coefficients dictionary
210  virtual const dictionary& coeffDict() const
211  {
212  return coeffDict_;
213  }
214 
215 
216  //- Return the effective viscosity
217  virtual tmp<volScalarField> nuEff() const
218  {
219  return tmp<volScalarField>
220  (
221  new volScalarField
222  (
223  IOobject::groupName("nuEff", this->U_.group()),
224  this->nut() + this->nu()
225  )
226  );
227  }
228 
229  //- Return the effective viscosity on patch
230  virtual tmp<scalarField> nuEff(const label patchi) const
231  {
232  return this->nut(patchi) + this->nu(patchi);
233  }
234 
235  //- Solve the turbulence equations and correct the turbulence viscosity
236  virtual void correct();
237 };
238 
239 
240 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
241 
242 } // End namespace Foam
243 
244 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
245 
246 #ifdef NoRepository
247  #include "RASModel.C"
248 #endif
249 
250 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
251 
252 #endif
253 
254 // ************************************************************************* //
surfaceScalarField & phi
virtual const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: RASModel.H:209
U
Definition: pEqn.H:83
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
dimensionedScalar kMin_
Lower limit of k.
Definition: RASModel.H:71
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dictionary RASDict_
RAS coefficients dictionary.
Definition: RASModel.H:59
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:100
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none.
Definition: Switch.H:60
const dimensionedScalar & epsilonMin() const
Return the lower allowable limit for epsilon (default: SMALL)
Definition: RASModel.H:179
declareRunTimeSelectionTable(autoPtr, RASModel, 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))
virtual bool read()
Read model coefficients if they have changed.
Definition: RASModel.C:169
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:101
virtual ~RASModel()
Destructor.
Definition: RASModel.H:160
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:49
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dictionary coeffDict_
Model coefficients dictionary.
Definition: RASModel.H:68
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: RASModel.C:31
const dimensionedScalar & omegaMin() const
Return the lower allowable limit for omega (default: SMALL)
Definition: RASModel.H:185
const dimensionedScalar & kMin() const
Return the lower allowable limit for k (default: SMALL)
Definition: RASModel.H:173
dimensionedScalar omegaMin_
Lower limit for omega.
Definition: RASModel.H:77
label patchi
TypeName("RAS")
Runtime type information.
static autoPtr< RASModel > 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 RAS model.
Definition: RASModel.C:116
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
Switch printCoeffs_
Flag to print the model coeffs at run-time.
Definition: RASModel.H:65
A class for managing temporary objects.
Definition: PtrList.H:54
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
Definition: RASModel.H:216
volScalarField & nu
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: RASModel.C:195
dimensionedScalar epsilonMin_
Lower limit of epsilon.
Definition: RASModel.H:74
scalar nut
Namespace for OpenFOAM.
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:99
Switch turbulence_
Turbulence on/off flag.
Definition: RASModel.H:62