kineticTheoryModel.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) 2011-2019 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::RASModels::kineticTheoryModel
26 
27 Description
28  Kinetic theory particle phase RAS model
29 
30  Reference:
31  \verbatim
32  van Wachem, B. G. M. (2000).
33  Derivation, implementation, and validation of computer simulation models
34  for gas-solid fluidized beds.
35  PhD Thesis, TU Delft.
36  \endverbatim
37 
38  There are no default model coefficients.
39 
40 SourceFiles
41  kineticTheoryModel.C
42 
43 \*---------------------------------------------------------------------------*/
44 
45 #ifndef kineticTheoryModel_H
46 #define kineticTheoryModel_H
47 
48 #include "RASModel.H"
49 #include "eddyViscosity.H"
50 #include "phaseCompressibleTurbulenceModel.H"
51 #include "EddyDiffusivity.H"
52 #include "phaseModel.H"
53 #include "dragModel.H"
54 #include "viscosityModel.H"
55 #include "conductivityModel.H"
56 #include "radialModel.H"
57 #include "granularPressureModel.H"
58 #include "frictionalStressModel.H"
59 
60 
61 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62 
63 namespace Foam
64 {
65 namespace RASModels
66 {
67 
68 /*---------------------------------------------------------------------------*\
69  Class kineticTheoryModel Declaration
70 \*---------------------------------------------------------------------------*/
71 
73 :
74  public eddyViscosity
75  <
76  RASModel<EddyDiffusivity<phaseCompressibleTurbulenceModel>>
77  >
78 {
79  // Private Data
80 
81  // Input Fields
82 
83  const phaseModel& phase_;
84 
85 
86  // Sub-models
87 
88  //- Run-time selected viscosity model
90 
91  //- Run-time selected conductivity model
93 
94  //- Run-time selected radial distribution model
96 
97  //- Run-time selected granular pressure model
99  granularPressureModel_;
100 
101  //- Run-time selected frictional stress model
103  frictionalStressModel_;
104 
105 
106  // Kinetic Theory Model coefficients
107 
108  //- Use equilibrium approximation: generation == dissipation
109  Switch equilibrium_;
110 
111  //- Coefficient of restitution
113 
114  //- Maximum packing phase-fraction
115  dimensionedScalar alphaMax_;
116 
117  //- Min value for which the frictional stresses are zero
118  dimensionedScalar alphaMinFriction_;
119 
120  //- Residual phase fraction
121  dimensionedScalar residualAlpha_;
122 
123  //- Maximum turbulent viscosity
124  dimensionedScalar maxNut_;
125 
126 
127  // Kinetic Theory Model Fields
128 
129  //- The granular energy/temperature
130  volScalarField Theta_;
131 
132  //- The granular bulk viscosity
133  volScalarField lambda_;
134 
135  //- The granular radial distribution
136  volScalarField gs0_;
137 
138  //- The granular "thermal" conductivity
139  volScalarField kappa_;
140 
141  //- The frictional viscosity
142  volScalarField nuFric_;
143 
144 
145  // Private Member Functions
146 
147  void correctNut()
148  {}
149 
150 
151 public:
152 
153  //- Runtime type information
154  TypeName("kineticTheory");
155 
156 
157  // Constructors
158 
159  //- Construct from components
161  (
162  const volScalarField& alpha,
163  const volScalarField& rho,
164  const volVectorField& U,
166  const surfaceScalarField& phi,
167  const phaseModel& transport,
169  const word& type = typeName
170  );
171 
172  //- Disallow default bitwise copy construction
173  kineticTheoryModel(const kineticTheoryModel&) = delete;
174 
175 
176  //- Destructor
177  virtual ~kineticTheoryModel();
178 
179 
180  // Member Functions
181 
182  //- Re-read model coefficients if they have changed
183  virtual bool read();
184 
185  //- Return the effective viscosity
186  virtual tmp<volScalarField> nuEff() const
187  {
188  return this->nut();
189  }
190 
191  //- Return the effective viscosity on patch
192  virtual tmp<scalarField> nuEff(const label patchi) const
193  {
194  return this->nut(patchi);
195  }
196 
197  //- Return the turbulence kinetic energy
198  virtual tmp<volScalarField> k() const;
199 
200  //- Return the turbulence kinetic energy dissipation rate
201  virtual tmp<volScalarField> epsilon() const;
202 
203  //- Return the Reynolds stress tensor
204  virtual tmp<volSymmTensorField> R() const;
205 
206  //- Return the phase-pressure'
207  // (derivative of phase-pressure w.r.t. phase-fraction)
208  virtual tmp<volScalarField> pPrime() const;
209 
210  //- Return the face-phase-pressure'
211  // (derivative of phase-pressure w.r.t. phase-fraction)
212  virtual tmp<surfaceScalarField> pPrimef() const;
213 
214  //- Return the effective stress tensor
215  virtual tmp<volSymmTensorField> devRhoReff() const;
216 
217  //- Return the source term for the momentum equation
219 
220  //- Solve the kinetic theory equations and correct the viscosity
221  virtual void correct();
222 
223 
224  // Member Operators
225 
226  //- Disallow default bitwise assignment
227  void operator=(const kineticTheoryModel&) = delete;
228 };
229 
230 
231 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232 
233 } // End namespace RASModels
234 } // End namespace Foam
235 
236 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
237 
238 #endif
239 
240 // ************************************************************************* //
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
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
virtual tmp< volScalarField > nut() const
Return the turbulence viscosity.
const volVectorField & U() const
Access function to velocity field.
virtual tmp< surfaceScalarField > pPrimef() const
Return the face-phase-pressure&#39;.
kineticTheoryModel(const volScalarField &alpha, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const phaseModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
TypeName("kineticTheory")
Runtime type information.
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/any.
Definition: Switch.H:60
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:49
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
const transportModel & transport() const
Access function to incompressible transport model.
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
const surfaceScalarField & alphaRhoPhi() const
Access function to phase flux field.
virtual void correct()
Solve the kinetic theory equations and correct the viscosity.
void operator=(const kineticTheoryModel &)=delete
Disallow default bitwise assignment.
virtual ~kineticTheoryModel()
Destructor.
const volScalarField & rho() const
Return the density field.
Kinetic theory particle phase RAS model.
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
label patchi
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure&#39;.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
virtual bool read()
Re-read model coefficients if they have changed.
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
virtual tmp< surfaceScalarField > phi() const
Return the volumetric flux field.
A class for managing temporary objects.
Definition: PtrList.H:53
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
virtual tmp< volScalarField > alpha() const
Thermal diffusivity for enthalpy of mixture [kg/m/s].
Namespace for OpenFOAM.