phaseModel.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) 2015-2020 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"
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 bool referencePhase,
99  const label index
100  ),
101  (fluid, phaseName, referencePhase, index)
102  );
103 
104 
105  // Constructors
106 
107  phaseModel
108  (
109  const phaseSystem& fluid,
110  const word& phaseName,
111  const bool referencePhase,
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 bool referencePhase,
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  const word& referencePhaseName_;
135  mutable label indexCounter_;
136 
137  public:
138 
140  (
141  const phaseSystem& fluid,
142  const word& referencePhaseName
143  )
144  :
145  fluid_(fluid),
146  referencePhaseName_(referencePhaseName),
147  indexCounter_(-1)
148  {}
151  {
152  indexCounter_++;
153 
154  const word phaseName(is);
155 
156  return autoPtr<phaseModel>
157  (
159  (
160  fluid_,
161  phaseName,
162  phaseName == referencePhaseName_,
163  indexCounter_
164  )
165  );
166  }
167  };
168 
169 
170  //- Destructor
171  virtual ~phaseModel();
172 
173 
174  // Member Functions
175 
176  //- Return the name of this phase
177  const word& name() const;
178 
179  //- Return the name of the phase for use as the keyword in PtrDictionary
180  const word& keyword() const;
181 
182  //- Return the index of the phase
183  label index() const;
184 
185  //- Return the system to which this phase belongs
186  const phaseSystem& fluid() const;
187 
188  //- Return the residual phase-fraction for given phase
189  // Used to stabilize the phase momentum as the phase-fraction -> 0
190  const dimensionedScalar& residualAlpha() const;
191 
192  //- Return the maximum phase-fraction (e.g. packing limit)
193  scalar alphaMax() const;
194 
195  //- Return the Sauter-mean diameter
196  tmp<volScalarField> d() const;
197 
198  //- Return const-reference to diameterModel of the phase
199  const autoPtr<diameterModel>& dPtr() const;
200 
201  //- Correct the phase properties
202  virtual void correct();
203 
204  //- Correct the continuity error
205  virtual void correctContinuityError(const volScalarField& source);
206 
207  //- Correct the kinematics
208  virtual void correctKinematics();
209 
210  //- Correct the thermodynamics
211  virtual void correctThermo();
212 
213  //- Correct the reactions
214  virtual void correctReactions();
215 
216  //- Correct the species concentrations
217  virtual void correctSpecies();
218 
219  //- Correct the turbulence
220  virtual void correctTurbulence();
221 
222  //- Correct the energy transport
223  virtual void correctEnergyTransport();
224 
225  //- Ensure that the flux at inflow/outflow BCs is preserved
227 
228  //- Read phase properties dictionary
229  virtual bool read();
230 
231 
232  // Density variation and compressibility
233 
234  //- Return true if the phase is incompressible otherwise false
235  virtual bool incompressible() const = 0;
236 
237  //- Return true if the phase is constant density otherwise false
238  virtual bool isochoric() const = 0;
239 
240  //- Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
241  virtual tmp<volScalarField> divU() const = 0;
242 
243  //- Set the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
244  virtual void divU(tmp<volScalarField> divU) = 0;
245 
246 
247  // Thermo
248 
249  //- Return the thermophysical model
250  virtual const rhoThermo& thermo() const = 0;
251 
252  //- Access the thermophysical model
253  virtual rhoThermo& thermoRef() = 0;
254 
255  //- Return the density field
256  virtual tmp<volScalarField> rho() const = 0;
257 
258  //- Return whether the phase is isothermal
259  virtual bool isothermal() const = 0;
260 
261  //- Return the enthalpy equation
262  virtual tmp<fvScalarMatrix> heEqn() = 0;
263 
264 
265  // Species
266 
267  //- Return whether the phase is pure (i.e., not multi-component)
268  virtual bool pure() const = 0;
269 
270  //- Return the species fraction equation
271  virtual tmp<fvScalarMatrix> YiEqn(volScalarField& Yi) = 0;
272 
273  //- Return the species mass fractions
274  virtual const PtrList<volScalarField>& Y() const = 0;
275 
276  //- Return a species mass fraction by name
277  virtual const volScalarField& Y(const word& name) const = 0;
278 
279  //- Access the species mass fractions
280  virtual PtrList<volScalarField>& YRef() = 0;
281 
282  //- Return the active species mass fractions
283  virtual const UPtrList<volScalarField>& YActive() const = 0;
284 
285  //- Access the active species mass fractions
286  virtual UPtrList<volScalarField>& YActiveRef() = 0;
287 
288  //- Return the fuel consumption rate matrix
289  virtual tmp<fvScalarMatrix> R(volScalarField& Yi) const = 0;
290 
291 
292  // Momentum
293 
294  //- Return whether the phase is stationary
295  virtual bool stationary() const = 0;
296 
297  //- Return the momentum equation
298  virtual tmp<fvVectorMatrix> UEqn() = 0;
299 
300  //- Return the momentum equation for the face-based algorithm
301  virtual tmp<fvVectorMatrix> UfEqn() = 0;
302 
303  //- Return the velocity
304  virtual tmp<volVectorField> U() const = 0;
305 
306  //- Access the velocity
307  virtual volVectorField& URef() = 0;
308 
309  //- Return the volumetric flux
310  virtual tmp<surfaceScalarField> phi() const = 0;
311 
312  //- Access the volumetric flux
313  virtual surfaceScalarField& phiRef() = 0;
314 
315  //- Return the volumetric flux of the phase
316  virtual tmp<surfaceScalarField> alphaPhi() const = 0;
317 
318  //- Access the volumetric flux of the phase
319  virtual surfaceScalarField& alphaPhiRef() = 0;
320 
321  //- Return the mass flux of the phase
322  virtual tmp<surfaceScalarField> alphaRhoPhi() const = 0;
323 
324  //- Access the mass flux of the phase
325  virtual surfaceScalarField& alphaRhoPhiRef() = 0;
326 
327  //- Return the substantive acceleration
328  virtual tmp<volVectorField> DUDt() const = 0;
329 
330  //- Return the substantive acceleration on the faces
331  virtual tmp<surfaceScalarField> DUDtf() const = 0;
332 
333  //- Return the continuity error
334  virtual tmp<volScalarField> continuityError() const = 0;
335 
336  //- Return the phase kinetic energy
337  virtual tmp<volScalarField> K() const = 0;
338 
339 
340  // Transport
341 
342  //- Effective thermal turbulent diffusivity for temperature
343  // of mixture for patch [W/m/K]
344  virtual tmp<scalarField> kappaEff(const label patchi) const = 0;
345 
346  //- Return the turbulent kinetic energy
347  virtual tmp<volScalarField> k() const = 0;
348 
349  //- Return the phase-pressure'
350  // (derivative of phase-pressure w.r.t. phase-fraction)
351  virtual tmp<volScalarField> pPrime() const = 0;
352 };
353 
354 
355 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
356 
357 } // End namespace Foam
358 
359 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
360 
361 #endif
362 
363 // ************************************************************************* //
Foam::surfaceFields.
virtual surfaceScalarField & phiRef()=0
Access the volumetric flux.
virtual UPtrList< volScalarField > & YActiveRef()=0
Access the active species mass fractions.
void correctInflowOutflow(surfaceScalarField &alphaPhi) const
Ensure that the flux at inflow/outflow BCs is preserved.
label index() const
Return the index of the phase.
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
Return the continuity error.
virtual surfaceScalarField & alphaRhoPhiRef()=0
Access the mass flux of the phase.
virtual tmp< volScalarField > divU() const =0
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
virtual bool incompressible() const =0
Return true if the phase is incompressible otherwise false.
static autoPtr< phaseModel > New(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
virtual tmp< fvScalarMatrix > R(volScalarField &Yi) const =0
Return the fuel consumption rate matrix.
virtual tmp< volVectorField > U() const =0
Return the velocity.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
declareRunTimeSelectionTable(autoPtr, phaseModel, phaseSystem,(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index),(fluid, phaseName, referencePhase, index))
autoPtr< phaseModel > operator()(Istream &is) const
Definition: phaseModel.H:100
virtual void correctKinematics()
Correct the kinematics.
virtual tmp< volScalarField > rho() const =0
Return the density field.
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
Return the Sauter-mean diameter.
phaseModel(const word &phaseName, const volScalarField &p, const volScalarField &T)
Construct from components.
virtual bool isothermal() const =0
Return whether the phase is isothermal.
iNew(const volScalarField &p, const volScalarField &T)
Definition: phaseModel.H:91
virtual tmp< fvScalarMatrix > YiEqn(volScalarField &Yi)=0
Return the species fraction equation.
virtual tmp< surfaceScalarField > alphaPhi() const =0
Return the volumetric flux of the phase.
virtual void correctSpecies()
Correct the species concentrations.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
virtual tmp< fvScalarMatrix > heEqn()=0
Return the enthalpy equation.
virtual tmp< surfaceScalarField > phi() const =0
Return the volumetric flux.
virtual ~phaseModel()
Destructor.
virtual tmp< scalarField > kappaEff(const label patchi) const =0
Effective thermal turbulent diffusivity for temperature.
virtual PtrList< volScalarField > & YRef()=0
Access the species mass fractions.
virtual bool read()
Read phase properties dictionary.
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:67
virtual tmp< volScalarField > k() const =0
Return the turbulent kinetic energy.
A class for handling words, derived from string.
Definition: word.H:59
virtual const UPtrList< volScalarField > & YActive() const =0
Return the active species mass fractions.
const rhoThermo & thermo() const
Return const-access to phase rhoThermo.
Definition: phaseModel.H:120
virtual bool pure() const =0
Return whether the phase is pure (i.e., not multi-component)
const word & keyword() const
Definition: phaseModel.H:114
virtual tmp< volScalarField > K() const =0
Return the phase kinetic energy.
virtual tmp< volScalarField > pPrime() const =0
Return the phase-pressure&#39;.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:54
const dimensionedScalar & residualAlpha() const
Return the residual phase-fraction for given phase.
virtual const PtrList< volScalarField > & Y() const =0
Return the species mass fractions.
scalar alphaMax() const
Return the maximum phase-fraction (e.g. packing limit)
virtual tmp< fvVectorMatrix > UEqn()=0
Return the momentum equation.
virtual tmp< surfaceScalarField > alphaRhoPhi() const =0
Return the mass flux of the phase.
virtual tmp< surfaceScalarField > DUDtf() const =0
Return the substantive acceleration on the faces.
Forward declarations of fvMatrix specializations.
label patchi
virtual surfaceScalarField & alphaPhiRef()=0
Access the volumetric flux of the phase.
ClassName("phaseModel")
Runtime type information.
virtual bool isochoric() const =0
Return true if the phase is constant density otherwise false.
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:70
const phaseSystem & fluid() const
Return the system to which this phase belongs.
const autoPtr< diameterModel > & dPtr() const
Return const-reference to diameterModel of the phase.
virtual bool stationary() const =0
Return whether the phase is stationary.
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.
virtual void correctTurbulence()
Correct the turbulence.
virtual void correctReactions()
Correct the reactions.
virtual tmp< fvVectorMatrix > UfEqn()=0
Return the momentum equation for the face-based algorithm.
virtual volVectorField & URef()=0
Access the velocity.
Namespace for OpenFOAM.
virtual void correctContinuityError(const volScalarField &source)
Correct the continuity error.
virtual rhoThermo & thermoRef()=0
Access the thermophysical model.