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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2015-2018 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 Class
27  Foam::phaseModel
28 
29 SourceFiles
30  phaseModel.C
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #ifndef phaseModel_H
35 #define phaseModel_H
36 
37 #include "dictionary.H"
38 #include "dimensionedScalar.H"
39 #include "volFields.H"
40 #include "surfaceFields.H"
41 #include "fvMatricesFwd.H"
42 #include "rhoThermo.H"
43 #include "phaseCompressibleTurbulenceModelFwd.H"
44 #include "runTimeSelectionTables.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 class phaseSystem;
52 class diameterModel;
53 
54 /*---------------------------------------------------------------------------*\
55  Class phaseModel Declaration
56 \*---------------------------------------------------------------------------*/
57 
58 class phaseModel
59 :
60  public volScalarField
61 {
62  // Private data
63 
64  //- Reference to the phaseSystem to which this phase belongs
65  const phaseSystem& fluid_;
66 
67  //- Name of phase
68  word name_;
69 
70  //- Index of phase
71  label index_;
72 
73  //- The residual phase-fraction for given phase
74  // Used to stabilize the phase momentum as the phase-fraction -> 0
75  dimensionedScalar residualAlpha_;
76 
77  //- Optional maximum phase-fraction (e.g. packing limit)
78  scalar alphaMax_;
79 
80  //- Diameter model
81  autoPtr<diameterModel> diameterModel_;
82 
83 
84 public:
85 
86  //- Runtime type information
87  ClassName("phaseModel");
88 
89 
90  // Declare runtime construction
91 
93  (
94  autoPtr,
95  phaseModel,
97  (
98  const phaseSystem& fluid,
99  const word& phaseName,
100  const label index
101  ),
102  (fluid, phaseName, index)
103  );
104 
105 
106  // Constructors
107 
108  phaseModel
109  (
110  const phaseSystem& fluid,
111  const word& phaseName,
112  const label index
113  );
114 
115  //- Return clone
116  autoPtr<phaseModel> clone() const;
117 
118 
119  // Selectors
120 
121  static autoPtr<phaseModel> New
122  (
123  const phaseSystem& fluid,
124  const word& phaseName,
125  const label index
126  );
127 
128  //- Return a pointer to a new phase created on freestore
129  // from Istream
130  class iNew
131  {
132  const phaseSystem& fluid_;
133  mutable label indexCounter_;
134 
135  public:
136 
138  (
139  const phaseSystem& fluid
140  )
141  :
142  fluid_(fluid),
143  indexCounter_(-1)
144  {}
145 
147  {
148  indexCounter_++;
149  return autoPtr<phaseModel>
150  (
151  phaseModel::New(fluid_, word(is), indexCounter_)
152  );
153  }
154  };
155 
156 
157  //- Destructor
158  virtual ~phaseModel();
159 
160 
161  // Member Functions
162 
163  //- Return the name of this phase
164  const word& name() const;
165 
166  //- Return the name of the phase for use as the keyword in PtrDictionary
167  const word& keyword() const;
168 
169  //- Return the index of the phase
170  label index() const;
171 
172  //- Return the system to which this phase belongs
173  const phaseSystem& fluid() const;
174 
175  //- Return the residual phase-fraction for given phase
176  // Used to stabilize the phase momentum as the phase-fraction -> 0
177  const dimensionedScalar& residualAlpha() const;
178 
179  //- Return the maximum phase-fraction (e.g. packing limit)
180  scalar alphaMax() const;
181 
182  //- Return the Sauter-mean diameter
183  tmp<volScalarField> d() const;
184 
185  //- Return const-reference to diameterModel of the phase
186  const autoPtr<diameterModel>& dPtr() const;
187 
188  //- Correct the phase properties
189  virtual void correct();
190 
191  //- Correct the kinematics
192  virtual void correctKinematics();
193 
194  //- Correct the thermodynamics
195  virtual void correctThermo();
196 
197  //- Correct the turbulence
198  virtual void correctTurbulence();
199 
200  //- Correct the energy transport
201  virtual void correctEnergyTransport();
202 
203  //- Ensure that the flux at inflow/outflow BCs is preserved
205 
206  //- Read phase properties dictionary
207  virtual bool read();
208 
209 
210  // Compressibility (variable density)
211 
212  //- Return true if the phase is compressible otherwise false
213  virtual bool compressible() const = 0;
214 
215  //- Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
216  virtual tmp<volScalarField> divU() const = 0;
217 
218  //- Set the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
219  virtual void divU(tmp<volScalarField> divU) = 0;
220 
221 
222  // Thermo
223 
224  //- Return whether the phase is isothermal
225  virtual bool isothermal() const = 0;
226 
227  //- Return the enthalpy equation
228  virtual tmp<fvScalarMatrix> heEqn() = 0;
229 
230  //- Return the thermophysical model
231  virtual const rhoThermo& thermo() const = 0;
232 
233  //- Access the thermophysical model
234  virtual rhoThermo& thermoRef() = 0;
235 
236  //- Return the density field
237  virtual tmp<volScalarField> rho() const = 0;
238 
239 
240  // Species
241 
242  //- Return whether the phase is pure (i.e., not multi-component)
243  virtual bool pure() const = 0;
244 
245  //- Return the species fraction equation
246  virtual tmp<fvScalarMatrix> YiEqn(volScalarField& Yi) = 0;
247 
248  //- Return the species mass fractions
249  virtual const PtrList<volScalarField>& Y() const = 0;
250 
251  //- Return a species mass fraction by name
252  virtual const volScalarField& Y(const word& name) const = 0;
253 
254  //- Access the species mass fractions
255  virtual PtrList<volScalarField>& YRef() = 0;
256 
257  //- Return the active species mass fractions
258  virtual const UPtrList<volScalarField>& YActive() const = 0;
259 
260  //- Access the active species mass fractions
261  virtual UPtrList<volScalarField>& YActiveRef() = 0;
262 
263 
264  // Momentum
265 
266  //- Return whether the phase is stationary
267  virtual bool stationary() const = 0;
268 
269  //- Return the momentum equation
270  virtual tmp<fvVectorMatrix> UEqn() = 0;
271 
272  //- Return the momentum equation for the face-based algorithm
273  virtual tmp<fvVectorMatrix> UfEqn() = 0;
274 
275  //- Return the velocity
276  virtual tmp<volVectorField> U() const = 0;
277 
278  //- Access the velocity
279  virtual volVectorField& URef() = 0;
280 
281  //- Return the volumetric flux
282  virtual tmp<surfaceScalarField> phi() const = 0;
283 
284  //- Access the volumetric flux
285  virtual surfaceScalarField& phiRef() = 0;
286 
287  //- Return the volumetric flux of the phase
288  virtual tmp<surfaceScalarField> alphaPhi() const = 0;
289 
290  //- Access the volumetric flux of the phase
291  virtual surfaceScalarField& alphaPhiRef() = 0;
292 
293  //- Return the mass flux of the phase
294  virtual tmp<surfaceScalarField> alphaRhoPhi() const = 0;
295 
296  //- Access the mass flux of the phase
297  virtual surfaceScalarField& alphaRhoPhiRef() = 0;
298 
299  //- Return the substantive acceleration
300  virtual tmp<volVectorField> DUDt() const = 0;
301 
302  //- Return the substantive acceleration on the faces
303  virtual tmp<surfaceScalarField> DUDtf() const = 0;
304 
305  //- Return the continuity error
306  virtual tmp<volScalarField> continuityError() const = 0;
307 
308  //- Return the continuity error due to the flow field
309  virtual tmp<volScalarField> continuityErrorFlow() const = 0;
310 
311  //- Return the continuity error due to any sources
312  virtual tmp<volScalarField> continuityErrorSources() const = 0;
313 
314  //- Return the phase kinetic energy
315  virtual tmp<volScalarField> K() const = 0;
316 
317 
318  // Transport
319 
320  //- Return the laminar dynamic viscosity
321  virtual tmp<volScalarField> mu() const = 0;
322 
323  //- Return the laminar dynamic viscosity on a patch
324  virtual tmp<scalarField> mu(const label patchi) const = 0;
325 
326  //- Return the laminar kinematic viscosity
327  virtual tmp<volScalarField> nu() const = 0;
328 
329  //- Return the laminar kinematic viscosity on a patch
330  virtual tmp<scalarField> nu(const label patchi) const = 0;
331 
332  //- Thermal diffusivity for enthalpy of mixture [kg/m/s]
333  virtual tmp<volScalarField> alpha() const = 0;
334 
335  //- Thermal diffusivity for enthalpy of mixture for patch [kg/m/s]
336  virtual tmp<scalarField> alpha(const label patchi) const = 0;
337 
338  //- Thermal diffusivity for temperature of mixture [J/m/s/K]
339  virtual tmp<volScalarField> kappa() const = 0;
340 
341  //- Thermal diffusivity for temperature of mixture
342  // for patch [J/m/s/K]
343  virtual tmp<scalarField> kappa(const label patchi) const = 0;
344 
345  //- Thermal diffusivity for energy of mixture [kg/m/s]
346  virtual tmp<volScalarField> alphahe() const = 0;
347 
348  //- Thermal diffusivity for energy of mixture for patch [kg/m/s]
349  virtual tmp<scalarField> alphahe(const label patchi) const = 0;
350 
351  //- Effective thermal turbulent diffusivity for temperature
352  // of mixture [J/m/s/K]
354  (
355  const volScalarField& alphat
356  ) const = 0;
357 
358  //- Effective thermal turbulent diffusivity for temperature
359  // of mixture for patch [J/m/s/K]
360  virtual tmp<scalarField> kappaEff
361  (
362  const scalarField& alphat,
363  const label patchi
364  ) const = 0;
365 
366  //- Effective thermal turbulent diffusivity of mixture [kg/m/s]
368  (
369  const volScalarField& alphat
370  ) const = 0;
371 
372  //- Effective thermal turbulent diffusivity of mixture
373  // for patch [kg/m/s]
374  virtual tmp<scalarField> alphaEff
375  (
376  const scalarField& alphat,
377  const label patchi
378  ) const = 0;
379 
380 
381  // Turbulence
382 
383  //- Return the turbulent dynamic viscosity
384  virtual tmp<volScalarField> mut() const = 0;
385 
386  //- Return the effective dynamic viscosity
387  virtual tmp<volScalarField> muEff() const = 0;
388 
389  //- Return the turbulent kinematic viscosity
390  virtual tmp<volScalarField> nut() const = 0;
391 
392  //- Return the effective kinematic viscosity
393  virtual tmp<volScalarField> nuEff() const = 0;
394 
395  //- Effective thermal turbulent diffusivity for temperature
396  // of mixture [J/m/s/K]
397  virtual tmp<volScalarField> kappaEff() const = 0;
398 
399  //- Effective thermal turbulent diffusivity for temperature
400  // of mixture for patch [J/m/s/K]
401  virtual tmp<scalarField> kappaEff(const label patchi) const = 0;
402 
403  //- Effective thermal turbulent diffusivity of mixture [kg/m/s]
404  virtual tmp<volScalarField> alphaEff() const = 0;
405 
406  //- Effective thermal turbulent diffusivity of mixture
407  // for patch [kg/m/s]
408  virtual tmp<scalarField> alphaEff(const label patchi) const = 0;
409 
410  //- Return the turbulent kinetic energy
411  virtual tmp<volScalarField> k() const = 0;
412 
413  //- Return the phase-pressure'
414  // (derivative of phase-pressure w.r.t. phase-fraction)
415  virtual tmp<volScalarField> pPrime() const = 0;
416 };
417 
418 
419 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
420 
421 } // End namespace Foam
422 
423 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
424 
425 #endif
426 
427 // ************************************************************************* //
volFields.H
Foam::phaseModel::DUDtf
virtual tmp< surfaceScalarField > DUDtf() const =0
Return the substantive acceleration on the faces.
Foam::phaseModel
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:57
Foam::phaseModel::correctTurbulence
virtual void correctTurbulence()
Correct the turbulence.
Definition: phaseModel.C:156
Foam::phaseModel::compressible
virtual bool compressible() const =0
Return true if the phase is compressible otherwise false.
Foam::phaseModel::index
label index() const
Return the index of the phase.
Definition: phaseModel.C:106
Foam::phaseModel::heEqn
virtual tmp< fvScalarMatrix > heEqn()=0
Return the enthalpy equation.
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::phaseModel::k
virtual tmp< volScalarField > k() const =0
Return the turbulent kinetic energy.
Foam::phaseModel::YiEqn
virtual tmp< fvScalarMatrix > YiEqn(volScalarField &Yi)=0
Return the species fraction equation.
Foam::tmp< volScalarField >
Foam::phaseModel::kappaEff
virtual tmp< volScalarField > kappaEff() const =0
Effective thermal turbulent diffusivity for temperature.
Foam::phaseModel::rho
virtual tmp< volScalarField > rho() const =0
Return the density field.
Foam::phaseModel::d
tmp< volScalarField > d() const
Return the Sauter-mean diameter.
Definition: phaseModel.C:130
Foam::phaseModel::alphaMax
scalar alphaMax() const
Return the maximum phase-fraction (e.g. packing limit)
Definition: phaseModel.C:124
Foam::phaseModel::correctEnergyTransport
virtual void correctEnergyTransport()
Correct the energy transport.
Definition: phaseModel.C:160
Foam::phaseModel::keyword
const word & keyword() const
Return the name of the phase for use as the keyword in PtrDictionary.
Definition: phaseModel.C:100
Foam::phaseModel::continuityErrorFlow
virtual tmp< volScalarField > continuityErrorFlow() const =0
Return the continuity error due to the flow field.
Foam::phaseModel::correct
virtual void correct()
Correct the phase properties.
Definition: phaseModel.C:142
fvMatricesFwd.H
Forward declarations of fvMatrix specializations.
surfaceFields.H
Foam::surfaceFields.
Foam::phaseModel::thermo
virtual const rhoThermo & thermo() const =0
Return the thermophysical model.
Foam::phaseModel::UEqn
virtual tmp< fvVectorMatrix > UEqn()=0
Return the momentum equation.
Foam::phaseModel::URef
virtual volVectorField & URef()=0
Access the velocity.
Foam::phaseModel::New
static autoPtr< phaseModel > New(const phaseSystem &fluid, const word &phaseName, const label index)
Definition: newPhaseModel.C:35
Foam::phaseModel::Y
virtual const PtrList< volScalarField > & Y() const =0
Return the species mass fractions.
Foam::phaseModel::declareRunTimeSelectionTable
declareRunTimeSelectionTable(autoPtr, phaseModel, phaseSystem,(const phaseSystem &fluid, const word &phaseName, const label index),(fluid, phaseName, index))
Foam::phaseModel::divU
virtual tmp< volScalarField > divU() const =0
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
Foam::phaseModel::~phaseModel
virtual ~phaseModel()
Destructor.
Definition: phaseModel.C:88
Foam::phaseModel::thermoRef
virtual rhoThermo & thermoRef()=0
Access the thermophysical model.
Foam::phaseModel::kappa
virtual tmp< volScalarField > kappa() const =0
Thermal diffusivity for temperature of mixture [J/m/s/K].
rhoThermo.H
Foam::phaseModel::correctThermo
virtual void correctThermo()
Correct the thermodynamics.
Definition: phaseModel.C:152
Foam::rhoThermo
Basic thermodynamic properties based on density.
Definition: rhoThermo.H:55
Foam::phaseModel::nu
virtual tmp< volScalarField > nu() const =0
Return the laminar kinematic viscosity.
Foam::phaseModel::alpha
virtual tmp< volScalarField > alpha() const =0
Thermal diffusivity for enthalpy of mixture [kg/m/s].
Foam::phaseModel::alphaRhoPhi
virtual tmp< surfaceScalarField > alphaRhoPhi() const =0
Return the mass flux of the phase.
Foam::phaseModel::continuityErrorSources
virtual tmp< volScalarField > continuityErrorSources() const =0
Return the continuity error due to any sources.
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:62
Foam::Field< scalar >
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::phaseModel::nuEff
virtual tmp< volScalarField > nuEff() const =0
Return the effective kinematic viscosity.
Foam::phaseModel::iNew::operator()
autoPtr< phaseModel > operator()(Istream &is) const
Definition: phaseModel.H:145
Foam::phaseModel::alphahe
virtual tmp< volScalarField > alphahe() const =0
Thermal diffusivity for energy of mixture [kg/m/s].
Foam::phaseModel::ClassName
ClassName("phaseModel")
Runtime type information.
Foam::UPtrList
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: UPtrList.H:63
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:65
Foam::phaseModel::iNew
Return a pointer to a new phase created on freestore.
Definition: phaseModel.H:129
Foam::phaseModel::isothermal
virtual bool isothermal() const =0
Return whether the phase is isothermal.
Foam::phaseModel::mu
virtual tmp< volScalarField > mu() const =0
Return the laminar dynamic viscosity.
Foam::phaseModel::YRef
virtual PtrList< volScalarField > & YRef()=0
Access the species mass fractions.
Foam::dimensioned< scalar >
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::phaseModel::YActive
virtual const UPtrList< volScalarField > & YActive() const =0
Return the active species mass fractions.
Foam::phaseModel::dPtr
const autoPtr< diameterModel > & dPtr() const
Return const-reference to diameterModel of the phase.
Definition: phaseModel.C:136
Foam::phaseModel::phi
virtual tmp< surfaceScalarField > phi() const =0
Return the volumetric flux.
Foam::phaseModel::alphaPhi
virtual tmp< surfaceScalarField > alphaPhi() const =0
Return the volumetric flux of the phase.
Foam::phaseModel::name
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:94
dimensionedScalar.H
Foam::phaseModel::mut
virtual tmp< volScalarField > mut() const =0
Return the turbulent dynamic viscosity.
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::phaseModel::phiRef
virtual surfaceScalarField & phiRef()=0
Access the volumetric flux.
Foam::phaseModel::UfEqn
virtual tmp< fvVectorMatrix > UfEqn()=0
Return the momentum equation for the face-based algorithm.
Foam::phaseModel::correctInflowOutflow
void correctInflowOutflow(surfaceScalarField &alphaPhi) const
Ensure that the flux at inflow/outflow BCs is preserved.
Definition: phaseModel.C:170
runTimeSelectionTables.H
Macros to ease declaration of run-time selection tables.
Foam::phaseModel::alphaRhoPhiRef
virtual surfaceScalarField & alphaRhoPhiRef()=0
Access the mass flux of the phase.
Foam::phaseModel::read
virtual bool read()
Read phase properties dictionary.
Definition: phaseModel.C:164
dictionary.H
Foam::phaseModel::continuityError
virtual tmp< volScalarField > continuityError() const =0
Return the continuity error.
Foam::phaseModel::fluid
const phaseSystem & fluid() const
Return the system to which this phase belongs.
Definition: phaseModel.C:112
Foam::phaseModel::U
virtual tmp< volVectorField > U() const =0
Return the velocity.
Foam::phaseModel::iNew::iNew
iNew(const phaseSystem &fluid)
Definition: phaseModel.H:137
Foam::phaseModel::alphaEff
virtual tmp< volScalarField > alphaEff() const =0
Effective thermal turbulent diffusivity of mixture [kg/m/s].
Foam::phaseModel::pure
virtual bool pure() const =0
Return whether the phase is pure (i.e., not multi-component)
Foam::phaseSystem
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:69
Foam::phaseModel::alphaPhiRef
virtual surfaceScalarField & alphaPhiRef()=0
Access the volumetric flux of the phase.
Foam::phaseModel::K
virtual tmp< volScalarField > K() const =0
Return the phase kinetic energy.
Foam::phaseModel::pPrime
virtual tmp< volScalarField > pPrime() const =0
Return the phase-pressure'.
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::phaseModel::clone
autoPtr< phaseModel > clone() const
Return clone.
Definition: phaseModel.C:79
Foam::phaseModel::residualAlpha
const dimensionedScalar & residualAlpha() const
Return the residual phase-fraction for given phase.
Definition: phaseModel.C:118
Foam::phaseModel::muEff
virtual tmp< volScalarField > muEff() const =0
Return the effective dynamic viscosity.
Foam::phaseModel::YActiveRef
virtual UPtrList< volScalarField > & YActiveRef()=0
Access the active species mass fractions.
Foam::phaseModel::stationary
virtual bool stationary() const =0
Return whether the phase is stationary.
Foam::phaseModel::correctKinematics
virtual void correctKinematics()
Correct the kinematics.
Definition: phaseModel.C:148
Foam::phaseModel::phaseModel
phaseModel(const phaseSystem &fluid, const word &phaseName, const label index)
Definition: phaseModel.C:44
Foam::phaseModel::DUDt
virtual tmp< volVectorField > DUDt() const =0
Return the substantive acceleration.
Foam::phaseModel::nut
virtual tmp< volScalarField > nut() const =0
Return the turbulent kinematic viscosity.