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-2024 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 "rhoFluidThermo.H"
41 #include "runTimeSelectionTables.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47 
48 class phaseSystem;
49 class diameterModel;
50 
51 /*---------------------------------------------------------------------------*\
52  Class phaseModel Declaration
53 \*---------------------------------------------------------------------------*/
54 
55 class phaseModel
56 :
57  public volScalarField
58 {
59  // Private Data
60 
61  //- Reference to the phaseSystem to which this phase belongs
62  const phaseSystem& fluid_;
63 
64  //- Name of phase
65  word name_;
66 
67  //- Index of phase
68  label index_;
69 
70  //- Return the residual phase-fraction for given phase
71  // Used to stabilise the phase momentum as the phase-fraction -> 0
72  dimensionedScalar residualAlpha_;
73 
74  //- Optional maximum phase-fraction (e.g. packing limit)
75  scalar alphaMax_;
76 
77  //- Diameter model
78  autoPtr<diameterModel> diameterModel_;
79 
80 
81 public:
82 
83  //- Runtime type information
84  ClassName("phaseModel");
85 
86 
87  // Declare runtime construction
88 
90  (
91  autoPtr,
92  phaseModel,
94  (
95  const phaseSystem& fluid,
96  const word& phaseName,
97  const bool referencePhase,
98  const label index
99  ),
100  (fluid, phaseName, referencePhase, index)
101  );
102 
103 
104  // Constructors
105 
106  phaseModel
107  (
108  const phaseSystem& fluid,
109  const word& phaseName,
110  const bool referencePhase,
111  const label index
112  );
113 
114  //- Return clone
115  autoPtr<phaseModel> clone() const;
116 
117 
118  // Selectors
119 
120  static autoPtr<phaseModel> New
121  (
122  const phaseSystem& fluid,
123  const word& phaseName,
124  const bool referencePhase,
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  const word& referencePhaseName_;
134  mutable label indexCounter_;
135 
136  public:
137 
138  iNew
139  (
140  const phaseSystem& fluid,
141  const word& referencePhaseName
142  )
143  :
144  fluid_(fluid),
145  referencePhaseName_(referencePhaseName),
146  indexCounter_(-1)
147  {}
148 
150  {
151  indexCounter_++;
152 
153  const word phaseName(is);
154 
155  return autoPtr<phaseModel>
156  (
158  (
159  fluid_,
160  phaseName,
161  phaseName == referencePhaseName_,
162  indexCounter_
163  )
164  );
165  }
166  };
167 
168 
169  //- Destructor
170  virtual ~phaseModel();
171 
172 
173  // Member Functions
174 
175  //- Return the name of this phase
176  const word& name() const;
177 
178  //- Return the name of the phase for use as the keyword in PtrDictionary
179  const word& keyword() const;
180 
181  //- Return the index of the phase
182  label index() const;
183 
184  //- Return the system to which this phase belongs
185  const phaseSystem& fluid() const;
186 
187  //- Return the residual phase-fraction for given phase
188  // Used to stabilise the phase momentum as the phase-fraction -> 0
189  const dimensionedScalar& residualAlpha() const;
190 
191  //- Return the maximum phase-fraction (e.g. packing limit)
192  scalar alphaMax() const;
193 
194  //- Return the Sauter-mean diameter
195  tmp<volScalarField> d() const;
196 
197  //- Return a reference to the diameterModel of the phase
198  const diameterModel& diameter() const;
199 
200  //- Correct the phase properties
201  virtual void correct();
202 
203  //- Correct the continuity error
204  virtual void correctContinuityError(const volScalarField& source);
205 
206  //- Correct the kinematics
207  virtual void correctKinematics();
208 
209  //- Correct the thermodynamics
210  virtual void correctThermo();
211 
212  //- Correct the reactions
213  virtual void correctReactions();
214 
215  //- Correct the species concentrations
216  virtual void correctSpecies();
217 
218  //- Predict the momentumTransport
219  virtual void predictMomentumTransport();
220 
221  //- Predict the energy transport
222  virtual void predictThermophysicalTransport();
223 
224  //- Correct the momentumTransport
225  virtual void correctMomentumTransport();
226 
227  //- Correct the energy transport
228  virtual void correctThermophysicalTransport();
229 
230  //- Correct the face velocity for moving meshes
231  virtual void correctUf();
232 
233  //- Ensure that the flux at inflow/outflow BCs is preserved
235 
236  //- Read phase properties dictionary
237  virtual bool read();
238 
239 
240  // Density variation and compressibility
241 
242  //- Return true if the phase is incompressible otherwise false
243  virtual bool incompressible() const = 0;
244 
245  //- Return true if the phase is constant density otherwise false
246  virtual bool isochoric() const = 0;
247 
248  //- Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
249  virtual const autoPtr<volScalarField>& divU() const = 0;
250 
251  //- Set the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
252  virtual void divU(tmp<volScalarField> divU) = 0;
253 
254 
255  // Thermo
256 
257  //- Return the thermophysical model
258  virtual const rhoThermo& thermo() const = 0;
259 
260  //- Access the thermophysical model
261  virtual rhoThermo& thermo() = 0;
262 
263  //- Return the thermophysical model
264  virtual const rhoFluidThermo& fluidThermo() const = 0;
265 
266  //- Access the thermophysical model
267  virtual rhoFluidThermo& fluidThermo() = 0;
268 
269  //- Return the density field
270  virtual const volScalarField& rho() const = 0;
271 
272  //- Access the density field
273  virtual volScalarField& rho() = 0;
274 
275  //- Return whether the phase is isothermal
276  virtual bool isothermal() const = 0;
277 
278  //- Return the enthalpy equation
279  virtual tmp<fvScalarMatrix> heEqn() = 0;
280 
281 
282  // Species
283 
284  //- Return whether the phase is pure (i.e., not multi-component)
285  virtual bool pure() const = 0;
286 
287  //- Return the species fraction equation
288  virtual tmp<fvScalarMatrix> YiEqn(volScalarField& Yi) = 0;
289 
290  //- Return the species mass fractions
291  virtual const PtrList<volScalarField>& Y() const = 0;
292 
293  //- Return a species mass fraction by name
294  virtual const volScalarField& Y(const word& name) const = 0;
295 
296  //- Access the species mass fractions
297  virtual PtrList<volScalarField>& YRef() = 0;
298 
299  //- Should the given specie be solved for? I.e., is it active and
300  // not the default specie?
301  virtual bool solveSpecie(const label speciei) const = 0;
302 
303  //- Return the fuel consumption rate matrix
304  virtual tmp<fvScalarMatrix> R(volScalarField& Yi) const = 0;
305 
306 
307  // Momentum
308 
309  //- Return whether the phase is stationary
310  virtual bool stationary() const = 0;
311 
312  //- Return the momentum equation
313  virtual tmp<fvVectorMatrix> UEqn() = 0;
314 
315  //- Return the momentum equation for the face-based algorithm
316  virtual tmp<fvVectorMatrix> UfEqn() = 0;
317 
318  //- Return the velocity
319  virtual tmp<volVectorField> U() const = 0;
320 
321  //- Access the velocity
322  virtual volVectorField& URef() = 0;
323 
324  //- Access the velocity
325  virtual const volVectorField& URef() const = 0;
326 
327  //- Return the volumetric flux
328  virtual tmp<surfaceScalarField> phi() const = 0;
329 
330  //- Access the volumetric flux
331  virtual surfaceScalarField& phiRef() = 0;
332 
333  //- Access the volumetric flux
334  virtual const surfaceScalarField& phiRef() const = 0;
335 
336  //- Return the face velocity
337  // Required for moving mesh cases
338  virtual const autoPtr<surfaceVectorField>& Uf() const = 0;
339 
340  //- Access the face velocity
341  // Required for moving mesh cases
342  virtual surfaceVectorField& UfRef() = 0;
343 
344  //- Access the face velocity
345  // Required for moving mesh cases
346  virtual const surfaceVectorField& UfRef() const = 0;
347 
348  //- Return the volumetric flux of the phase
349  virtual tmp<surfaceScalarField> alphaPhi() const = 0;
350 
351  //- Access the volumetric flux of the phase
352  virtual surfaceScalarField& alphaPhiRef() = 0;
353 
354  //- Access the volumetric flux of the phase
355  virtual const surfaceScalarField& alphaPhiRef() const = 0;
356 
357  //- Return the mass flux of the phase
358  virtual tmp<surfaceScalarField> alphaRhoPhi() const = 0;
359 
360  //- Access the mass flux of the phase
361  virtual surfaceScalarField& alphaRhoPhiRef() = 0;
362 
363  //- Access the mass flux of the phase
364  virtual const surfaceScalarField& alphaRhoPhiRef() const = 0;
365 
366  //- Return the velocity transport matrix
367  virtual tmp<fvVectorMatrix> UgradU() const = 0;
368 
369  //- Return the substantive acceleration matrix
370  virtual tmp<fvVectorMatrix> DUDt() const = 0;
371 
372  //- Return the continuity error
373  virtual tmp<volScalarField> continuityError() const = 0;
374 
375  //- Return the phase kinetic energy
376  virtual tmp<volScalarField> K() const = 0;
377 
378 
379  // Transport
380 
381  //- Effective thermal turbulent conductivity
382  // of mixture for patch [W/m/K]
383  virtual tmp<scalarField> kappaEff(const label patchi) const = 0;
384 
385  //- Return the turbulent kinetic energy
386  virtual tmp<volScalarField> k() const = 0;
387 
388  //- Return the face-phase-pressure'
389  // (derivative of phase-pressure w.r.t. phase-fraction)
390  virtual tmp<surfaceScalarField> pPrimef() const = 0;
391 };
392 
393 
394 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
395 
396 } // End namespace Foam
397 
398 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
399 
400 #endif
401 
402 // ************************************************************************* //
Generic GeometricField class.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Abstract base-class for dispersed-phase particle diameter models.
Definition: diameterModel.H:52
Return a pointer to a new phase created on freestore.
Definition: phaseModel.H:130
iNew(const phaseSystem &fluid, const word &referencePhaseName)
Definition: phaseModel.H:138
autoPtr< phaseModel > operator()(Istream &is) const
Definition: phaseModel.H:148
virtual void correctUf()
Correct the face velocity for moving meshes.
Definition: phaseModel.C:235
virtual void correctKinematics()
Correct the kinematics.
Definition: phaseModel.C:203
virtual const autoPtr< surfaceVectorField > & Uf() const =0
Return the face velocity.
ClassName("phaseModel")
Runtime type information.
const phaseSystem & fluid() const
Return the system to which this phase belongs.
Definition: phaseModel.C:163
virtual void correctReactions()
Correct the reactions.
Definition: phaseModel.C:211
virtual bool pure() const =0
Return whether the phase is pure (i.e., not multi-component)
const dimensionedScalar & residualAlpha() const
Return the residual phase-fraction for given phase.
Definition: phaseModel.C:169
virtual void correctSpecies()
Correct the species concentrations.
Definition: phaseModel.C:215
virtual tmp< fvVectorMatrix > DUDt() const =0
Return the substantive acceleration matrix.
virtual void correctThermo()
Correct the thermodynamics.
Definition: phaseModel.C:207
virtual void correct()
Correct the phase properties.
Definition: phaseModel.C:193
virtual ~phaseModel()
Destructor.
Definition: phaseModel.C:139
virtual bool solveSpecie(const label speciei) const =0
Should the given specie be solved for? I.e., is it active and.
virtual tmp< scalarField > kappaEff(const label patchi) const =0
Effective thermal turbulent conductivity.
virtual surfaceScalarField & phiRef()=0
Access the volumetric flux.
virtual tmp< fvVectorMatrix > UfEqn()=0
Return the momentum equation for the face-based algorithm.
virtual void predictThermophysicalTransport()
Predict the energy transport.
Definition: phaseModel.C:223
virtual const rhoThermo & thermo() const =0
Return the thermophysical model.
virtual tmp< volScalarField > continuityError() const =0
Return the continuity error.
virtual tmp< volScalarField > k() const =0
Return the turbulent kinetic energy.
declareRunTimeSelectionTable(autoPtr, phaseModel, phaseSystem,(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index),(fluid, phaseName, referencePhase, index))
virtual void correctMomentumTransport()
Correct the momentumTransport.
Definition: phaseModel.C:227
virtual tmp< volScalarField > K() const =0
Return the phase kinetic energy.
scalar alphaMax() const
Return the maximum phase-fraction (e.g. packing limit)
Definition: phaseModel.C:175
virtual const PtrList< volScalarField > & Y() const =0
Return the species mass fractions.
virtual tmp< surfaceScalarField > alphaRhoPhi() const =0
Return the mass flux of the phase.
virtual bool stationary() const =0
Return whether the phase is stationary.
const diameterModel & diameter() const
Return a reference to the diameterModel of the phase.
Definition: phaseModel.C:187
virtual const autoPtr< volScalarField > & divU() const =0
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
tmp< volScalarField > d() const
Return the Sauter-mean diameter.
Definition: phaseModel.C:181
virtual const rhoFluidThermo & fluidThermo() const =0
Return the thermophysical model.
virtual tmp< surfaceScalarField > phi() const =0
Return the volumetric flux.
virtual void correctThermophysicalTransport()
Correct the energy transport.
Definition: phaseModel.C:231
virtual bool incompressible() const =0
Return true if the phase is incompressible otherwise false.
virtual tmp< fvVectorMatrix > UEqn()=0
Return the momentum equation.
label index() const
Return the index of the phase.
Definition: phaseModel.C:157
virtual tmp< fvScalarMatrix > heEqn()=0
Return the enthalpy equation.
virtual surfaceScalarField & alphaRhoPhiRef()=0
Access the mass flux of the phase.
const word & keyword() const
Return the name of the phase for use as the keyword in PtrDictionary.
Definition: phaseModel.C:151
phaseModel(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
Definition: phaseModel.C:42
virtual tmp< volVectorField > U() const =0
Return the velocity.
virtual bool isochoric() const =0
Return true if the phase is constant density otherwise false.
virtual tmp< fvScalarMatrix > R(volScalarField &Yi) const =0
Return the fuel consumption rate matrix.
autoPtr< phaseModel > clone() const
Return clone.
Definition: phaseModel.C:130
virtual bool isothermal() const =0
Return whether the phase is isothermal.
virtual surfaceScalarField & alphaPhiRef()=0
Access the volumetric flux of the phase.
virtual volVectorField & URef()=0
Access the velocity.
static autoPtr< phaseModel > New(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
Definition: phaseModelNew.C:32
void correctInflowOutflow(surfaceScalarField &alphaPhi) const
Ensure that the flux at inflow/outflow BCs is preserved.
Definition: phaseModel.C:245
virtual PtrList< volScalarField > & YRef()=0
Access the species mass fractions.
virtual tmp< fvScalarMatrix > YiEqn(volScalarField &Yi)=0
Return the species fraction equation.
virtual tmp< fvVectorMatrix > UgradU() const =0
Return the velocity transport matrix.
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:145
virtual tmp< surfaceScalarField > pPrimef() const =0
Return the face-phase-pressure'.
virtual const volScalarField & rho() const =0
Return the density field.
virtual tmp< surfaceScalarField > alphaPhi() const =0
Return the volumetric flux of the phase.
virtual void predictMomentumTransport()
Predict the momentumTransport.
Definition: phaseModel.C:219
virtual void correctContinuityError(const volScalarField &source)
Correct the continuity error.
Definition: phaseModel.C:199
virtual bool read()
Read phase properties dictionary.
Definition: phaseModel.C:239
virtual surfaceVectorField & UfRef()=0
Access the face velocity.
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:73
Base-class for fluid thermodynamic properties based on density.
Base-class for thermodynamic properties based on density.
Definition: rhoThermo.H:54
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Forward declarations of fvMatrix specialisations.
label patchi
Namespace for OpenFOAM.
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
Macros to ease declaration of run-time selection tables.
Foam::surfaceFields.