phaseModel.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) 2015-2017 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::phaseModel
26 
27 SourceFiles
28  phaseModel.C
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #ifndef phaseModel_H
33 #define phaseModel_H
34 
35 #include "dictionary.H"
36 #include "dimensionedScalar.H"
37 #include "volFields.H"
38 #include "surfaceFields.H"
39 #include "fvMatricesFwd.H"
40 #include "rhoThermo.H"
41 #include "phaseCompressibleTurbulenceModelFwd.H"
42 #include "runTimeSelectionTables.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 class phaseSystem;
50 class diameterModel;
51 
52 /*---------------------------------------------------------------------------*\
53  Class phaseModel Declaration
54 \*---------------------------------------------------------------------------*/
55 
56 class phaseModel
57 :
58  public volScalarField
59 {
60  // Private data
61 
62  //- Reference to the phaseSystem to which this phase belongs
63  const phaseSystem& fluid_;
64 
65  //- Name of phase
66  word name_;
67 
68  //- Index of phase
69  label index_;
70 
71  //- Return the residual phase-fraction for given phase
72  // Used to stabilize the phase momentum as the phase-fraction -> 0
73  dimensionedScalar residualAlpha_;
74 
75  //- Optional maximum phase-fraction (e.g. packing limit)
76  scalar alphaMax_;
77 
78  //- Diameter model
79  autoPtr<diameterModel> diameterModel_;
80 
81 
82 public:
83 
84  //- Runtime type information
85  ClassName("phaseModel");
86 
87 
88  // Declare runtime construction
89 
91  (
92  autoPtr,
93  phaseModel,
94  phaseSystem,
95  (
96  const phaseSystem& fluid,
97  const word& phaseName,
98  const label index
99  ),
100  (fluid, phaseName, index)
101  );
102 
103 
104  // Constructors
105 
106  phaseModel
107  (
108  const phaseSystem& fluid,
109  const word& phaseName,
110  const label index
111  );
112 
113  //- Return clone
114  autoPtr<phaseModel> clone() const;
115 
116 
117  // Selectors
118 
119  static autoPtr<phaseModel> New
120  (
121  const phaseSystem& fluid,
122  const word& phaseName,
123  const label index
124  );
125 
126  //- Return a pointer to a new phase created on freestore
127  // from Istream
128  class iNew
129  {
130  const phaseSystem& fluid_;
131  mutable label indexCounter_;
132 
133  public:
134 
136  (
137  const phaseSystem& fluid
138  )
139  :
140  fluid_(fluid),
141  indexCounter_(-1)
142  {}
145  {
146  indexCounter_++;
147  return autoPtr<phaseModel>
148  (
149  phaseModel::New(fluid_, word(is), indexCounter_)
150  );
151  }
152  };
153 
154 
155  //- Destructor
156  virtual ~phaseModel();
157 
158 
159  // Member Functions
160 
161  //- Return the name of this phase
162  const word& name() const;
163 
164  //- Return the name of the phase for use as the keyword in PtrDictionary
165  const word& keyword() const;
166 
167  //- Return the index of the phase
168  label index() const;
169 
170  //- Return the system to which this phase belongs
171  const phaseSystem& fluid() const;
172 
173  //- Return the residual phase-fraction for given phase
174  // Used to stabilize the phase momentum as the phase-fraction -> 0
175  const dimensionedScalar& residualAlpha() const;
176 
177  //- Return the maximum phase-fraction (e.g. packing limit)
178  scalar alphaMax() const;
179 
180  //- Return the Sauter-mean diameter
181  tmp<volScalarField> d() const;
182 
183  //- Correct the phase properties
184  virtual void correct();
185 
186  //- Correct the kinematics
187  virtual void correctKinematics();
188 
189  //- Correct the thermodynamics
190  virtual void correctThermo();
191 
192  //- Correct the turbulence
193  virtual void correctTurbulence();
194 
195  //- Correct the energy transport e.g. alphat
196  virtual void correctEnergyTransport();
197 
198  //- Return the momentum equation
199  virtual tmp<fvVectorMatrix> UEqn() = 0;
200 
201  //- Return the enthalpy equation
202  virtual tmp<fvScalarMatrix> heEqn() = 0;
203 
204  //- Return the species fraction equation
205  virtual tmp<fvScalarMatrix> YiEqn(volScalarField& Yi) = 0;
206 
207  //- Read phase properties dictionary
208  virtual bool read();
209 
210 
211  // Compressibility (variable density)
212 
213  //- Return true if the phase is compressible otherwise false
214  virtual bool compressible() const;
215 
216  //- Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
217  virtual const tmp<volScalarField>& divU() const;
218 
219  //- Set the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
220  virtual void divU(const tmp<volScalarField>& divU);
221 
222  //- Return the phase kinetic energy
223  virtual const volScalarField& K() const;
224 
225 
226  // Implicit phase pressure and dispersion support
227 
228  //- Return the phase diffusivity divided by the momentum coefficient
229  virtual const surfaceScalarField& DbyA() const;
230 
231  //- Set the phase diffusivity divided by the momentum coefficient
232  virtual void DbyA(const tmp<surfaceScalarField>& DbyA);
233 
234 
235  // Thermo
236 
237  //- Return const access to the thermophysical model
238  virtual const rhoThermo& thermo() const = 0;
239 
240  //- Return non-const access to the thermophysical model
241  // for correction
242  virtual rhoThermo& thermo() = 0;
243 
244  //- Return the density field
245  virtual tmp<volScalarField> rho() const = 0;
246 
247  //- Constant access the species mass fractions
248  virtual const PtrList<volScalarField>& Y() const = 0;
249 
250  //- Access the species mass fractions
251  virtual PtrList<volScalarField>& Y() = 0;
252 
253 
254  // Momentum
255 
256  //- Constant access the velocity
257  virtual tmp<volVectorField> U() const = 0;
258 
259  //- Access the velocity
260  virtual volVectorField& U() = 0;
261 
262  //- Return the substantive acceleration
263  virtual tmp<volVectorField> DUDt() const = 0;
264 
265  //- Constant access the continuity error
266  virtual tmp<volScalarField> continuityError() const = 0;
267 
268  //- Constant access the volumetric flux
269  virtual tmp<surfaceScalarField> phi() const = 0;
270 
271  //- Access the volumetric flux
272  virtual surfaceScalarField& phi() = 0;
273 
274  //- Constant access the volumetric flux of the phase
275  virtual tmp<surfaceScalarField> alphaPhi() const = 0;
276 
277  //- Access the volumetric flux of the phase
278  virtual surfaceScalarField& alphaPhi() = 0;
279 
280  //- Constant access the mass flux of the phase
281  virtual tmp<surfaceScalarField> alphaRhoPhi() const = 0;
282 
283  //- Access the mass flux of the phase
284  virtual surfaceScalarField& alphaRhoPhi() = 0;
285 
286  //- Ensure that the flux at inflow/outflow BCs is preserved
288 
289 
290  // Transport
291 
292  //- Return the laminar dynamic viscosity
293  virtual tmp<volScalarField> mu() const = 0;
294 
295  //- Return the laminar dynamic viscosity on a patch
296  virtual tmp<scalarField> mu(const label patchi) const = 0;
297 
298  //- Return the laminar kinematic viscosity
299  virtual tmp<volScalarField> nu() const = 0;
300 
301  //- Return the laminar kinematic viscosity on a patch
302  virtual tmp<scalarField> nu(const label patchi) const = 0;
303 
304  //- Return the laminar thermal conductivity
305  virtual tmp<volScalarField> kappa() const = 0;
306 
307  //- Return the laminar thermal conductivity on a patch
308  virtual tmp<scalarField> kappa(const label patchi) const = 0;
309 
310  //- Return the effective thermal conductivity
312  (
313  const volScalarField& alphat
314  ) const = 0;
315 
316  //- Access the effective thermal conductivity
317  virtual tmp<scalarField> kappaEff
318  (
319  const scalarField& alphat,
320  const label patchi
321  ) const = 0;
322 
323  //- Return the laminar thermal diffusivity for enthalpy
324  virtual tmp<volScalarField> alpha() const = 0;
325 
326  //- Return the laminar thermal diffusivity for enthalpy on a patch
327  virtual tmp<scalarField> alpha(const label patchi) const = 0;
328 
329  //- Return the effective thermal diffusivity for enthalpy
331  (
332  const volScalarField& alphat
333  ) const = 0;
334 
335  //- Return the effective thermal diffusivity for enthalpy on a patch
336  virtual tmp<scalarField> alphaEff
337  (
338  const scalarField& alphat,
339  const label patchi
340  ) const = 0;
341 
342 
343  // Turbulence
344 
345  //- Return the turbulence model
346  virtual const phaseCompressibleTurbulenceModel&
347  turbulence() const = 0;
348 };
349 
350 
351 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
352 
353 } // End namespace Foam
354 
355 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
356 
357 #endif
358 
359 // ************************************************************************* //
Foam::surfaceFields.
void correctInflowOutflow(surfaceScalarField &alphaPhi) const
Ensure that the flux at inflow/outflow BCs is preserved.
label index() const
Return the index of the phase.
static autoPtr< phaseModel > New(const phaseSystem &fluid, const word &phaseName, const label index)
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 > continuityError() const =0
Constant access the continuity error.
virtual tmp< volScalarField > alpha() const =0
Return the laminar thermal diffusivity for enthalpy.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
autoPtr< phaseModel > operator()(Istream &is) const
Definition: phaseModel.H:100
virtual void correctKinematics()
Correct the kinematics.
virtual tmp< volScalarField > kappaEff(const volScalarField &alphat) const =0
Return the effective thermal conductivity.
autoPtr< phaseModel > clone() const
Return clone.
const word & name() const
Definition: phaseModel.H:109
virtual tmp< volVectorField > DUDt() const =0
Return the substantive acceleration.
tmp< volScalarField > d() const
phaseModel(const word &phaseName, const volScalarField &p, const volScalarField &T)
Construct from components.
const surfaceScalarField & phi() const
Definition: phaseModel.H:190
virtual const volScalarField & K() const
Return the phase kinetic energy.
iNew(const volScalarField &p, const volScalarField &T)
Definition: phaseModel.H:91
virtual tmp< fvScalarMatrix > YiEqn(volScalarField &Yi)=0
Return the species fraction equation.
const surfaceScalarField & alphaPhi() const
Definition: phaseModel.H:200
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
virtual tmp< fvScalarMatrix > heEqn()=0
Return the enthalpy equation.
virtual const phaseCompressibleTurbulenceModel & turbulence() const =0
Return the turbulence model.
const volVectorField & U() const
Definition: phaseModel.H:170
virtual const tmp< volScalarField > & divU() const
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
virtual tmp< volScalarField > alphaEff(const volScalarField &alphat) const =0
Return the effective thermal diffusivity for enthalpy.
virtual ~phaseModel()
Destructor.
virtual tmp< volScalarField > mu() const =0
Return the laminar dynamic viscosity.
virtual bool read()
Read phase properties dictionary.
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:67
A class for handling words, derived from string.
Definition: word.H:59
const rhoThermo & thermo() const
Return const-access to phase rhoThermo.
Definition: phaseModel.H:120
const word & keyword() const
Definition: phaseModel.H:114
virtual bool compressible() const
Return true if the phase is compressible otherwise false.
const dimensionedScalar & residualAlpha() const
Return the residual phase-fraction for given phase.
void correct()
Correct the laminar viscosity.
virtual const PtrList< volScalarField > & Y() const =0
Constant access the species mass fractions.
const dimensionedScalar & kappa() const
Definition: phaseModel.H:155
scalar alphaMax() const
Return the maximum phase-fraction (e.g. packing limit)
const dimensionedScalar & nu() const
Return the laminar viscosity.
Definition: phaseModel.H:150
virtual tmp< fvVectorMatrix > UEqn()=0
Return the momentum equation.
virtual tmp< surfaceScalarField > alphaRhoPhi() const =0
Constant access the mass flux of the phase.
virtual const surfaceScalarField & DbyA() const
Return the phase diffusivity divided by the momentum coefficient.
Forward declarations of fvMatrix specializations.
label patchi
ClassName("phaseModel")
Runtime type information.
const dimensionedScalar & rho() const
Definition: phaseModel.H:165
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:63
const phaseSystem & fluid() const
Return the system to which this phase belongs.
Basic thermodynamic properties based on density.
Definition: rhoThermo.H:49
virtual void correctThermo()
Correct the thermodynamics.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Macros to ease declaration of run-time selection tables.
A class for managing temporary objects.
Definition: PtrList.H:53
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat.
virtual void correctTurbulence()
Correct the turbulence.
declareRunTimeSelectionTable(autoPtr, phaseModel, phaseSystem,(const phaseSystem &fluid, const word &phaseName, const label index),(fluid, phaseName, index))
Namespace for OpenFOAM.