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