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) 2011-2018 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 "transportModel.H"
40 #include "rhoThermo.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 // Forward declarations
48 class twoPhaseSystem;
49 class diameterModel;
50 
51 template<class Phase>
52 class PhaseCompressibleTurbulenceModel;
53 
54 
55 /*---------------------------------------------------------------------------*\
56  Class phaseModel Declaration
57 \*---------------------------------------------------------------------------*/
58 
59 class phaseModel
60 :
61  public volScalarField,
62  public transportModel
63 {
64  // Private data
65 
66  //- Reference to the twoPhaseSystem to which this phase belongs
67  const twoPhaseSystem& fluid_;
68 
69  //- Name of phase
70  word name_;
71 
72  dictionary phaseDict_;
73 
74  //- Return 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  //- Thermophysical properties
82  autoPtr<rhoThermo> thermo_;
83 
84  //- Velocity
85  volVectorField U_;
86 
87  //- Volumetric flux of the phase
88  surfaceScalarField alphaPhi_;
89 
90  //- Mass flux of the phase
91  surfaceScalarField alphaRhoPhi_;
92 
93  //- Volumetric flux of the phase
94  autoPtr<surfaceScalarField> phiPtr_;
95 
96  //- Diameter model
97  autoPtr<diameterModel> dPtr_;
98 
99  //- Turbulence model
100  autoPtr<PhaseCompressibleTurbulenceModel<phaseModel>> turbulence_;
101 
102 
103 public:
104 
105  // Constructors
106 
107  phaseModel
108  (
109  const twoPhaseSystem& fluid,
110  const dictionary& phaseProperties,
111  const word& phaseName
112  );
113 
114 
115  //- Destructor
116  virtual ~phaseModel();
117 
118 
119  // Member Functions
120 
121  //- Return the name of this phase
122  const word& name() const
123  {
124  return name_;
125  }
126 
127  //- Return the twoPhaseSystem to which this phase belongs
128  const twoPhaseSystem& fluid() const
129  {
130  return fluid_;
131  }
132 
133  //- Return the other phase in this two-phase system
134  const phaseModel& otherPhase() const;
135 
136  //- Return the residual phase-fraction for given phase
137  // Used to stabilize the phase momentum as the phase-fraction -> 0
138  const dimensionedScalar& residualAlpha() const
139  {
140  return residualAlpha_;
141  }
142 
143  //- Optional maximum phase-fraction (e.g. packing limit)
144  // Defaults to 1
145  scalar alphaMax() const
146  {
147  return alphaMax_;
148  }
149 
150  //- Return the Sauter-mean diameter
151  tmp<volScalarField> d() const;
152 
153  //- Return the turbulence model
155  turbulence() const;
156 
157  //- Return non-const access to the turbulence model
158  // for correction
160  turbulence();
161 
162  //- Return the thermophysical model
163  const rhoThermo& thermo() const
164  {
165  return thermo_();
166  }
167 
168  //- Return non-const access to the thermophysical model
169  // for correction
170  rhoThermo& thermo()
171  {
172  return thermo_();
173  }
174 
175  //- Return the laminar viscosity
176  tmp<volScalarField> nu() const
177  {
178  return thermo_->nu();
179  }
180 
181  //- Return the laminar viscosity for patch
182  tmp<scalarField> nu(const label patchi) const
183  {
184  return thermo_->nu(patchi);
185  }
186 
187  //- Return the laminar dynamic viscosity
188  tmp<volScalarField> mu() const
189  {
190  return thermo_->mu();
191  }
192 
193  //- Return the laminar dynamic viscosity for patch
194  tmp<scalarField> mu(const label patchi) const
195  {
196  return thermo_->mu(patchi);
197  }
198 
199  //- Thermal diffusivity for enthalpy of mixture [kg/m/s]
200  tmp<volScalarField> alpha() const
201  {
202  return thermo_->alpha();
203  }
204 
205  //- Thermal diffusivity for enthalpy of mixture for patch [kg/m/s]
206  tmp<scalarField> alpha(const label patchi) const
207  {
208  return thermo_->alpha(patchi);
209  }
210 
211  //- Thermal diffusivity for temperature of mixture [J/m/s/K]
212  tmp<scalarField> kappa(const label patchi) const
213  {
214  return thermo_->kappa(patchi);
215  }
216 
217  //- Thermal diffusivity for temperature of mixture
218  // for patch [J/m/s/K]
219  tmp<volScalarField> kappa() const
220  {
221  return thermo_->kappa();
222  }
223 
224  //- Thermal diffusivity for energy of mixture [kg/m/s]
226  {
227  return thermo_->alphahe();
228  }
229 
230  //- Thermal diffusivity for energy of mixture for patch [kg/m/s]
231  tmp<scalarField> alphahe(const label patchi) const
232  {
233  return thermo_->alphahe(patchi);
234  }
235 
236  //- Effective thermal turbulent diffusivity for temperature
237  // of mixture [J/m/s/K]
239  (
240  const volScalarField& alphat
241  ) const
242  {
243  return thermo_->kappaEff(alphat);
244  }
245 
246  //- Effective thermal turbulent diffusivity for temperature
247  // of mixture for patch [J/m/s/K]
249  (
250  const scalarField& alphat,
251  const label patchi
252  ) const
253  {
254  return thermo_->kappaEff(alphat, patchi);
255  }
256 
257  //- Effective thermal turbulent diffusivity of mixture [kg/m/s]
259  (
260  const volScalarField& alphat
261  ) const
262  {
263  return thermo_->alphaEff(alphat);
264  }
265 
266  //- Effective thermal turbulent diffusivity of mixture
267  // for patch [kg/m/s]
269  (
270  const scalarField& alphat,
271  const label patchi
272  ) const
273  {
274  return thermo_->alphaEff(alphat, patchi);
275  }
276 
277  //- Return the specific heat capacity
278  tmp<volScalarField> Cp() const
279  {
280  return thermo_->Cp();
281  }
282 
283  //- Return the density
284  const volScalarField& rho() const
285  {
286  return thermo_->rho();
287  }
288 
289  //- Return the velocity
290  const volVectorField& U() const
291  {
292  return U_;
293  }
294 
295  //- Return non-const access to the velocity
296  // Used in the momentum equation
297  volVectorField& U()
298  {
299  return U_;
300  }
301 
302  //- Return the volumetric flux
303  const surfaceScalarField& phi() const
304  {
305  return phiPtr_();
306  }
307 
308  //- Return non-const access to the volumetric flux
310  {
311  return phiPtr_();
312  }
313 
314  //- Return the volumetric flux of the phase
315  const surfaceScalarField& alphaPhi() const
316  {
317  return alphaPhi_;
318  }
319 
320  //- Return non-const access to the volumetric flux of the phase
322  {
323  return alphaPhi_;
324  }
325 
326  //- Return the mass flux of the phase
327  const surfaceScalarField& alphaRhoPhi() const
328  {
329  return alphaRhoPhi_;
330  }
331 
332  //- Return non-const access to the mass flux of the phase
334  {
335  return alphaRhoPhi_;
336  }
337 
338  //- Ensure that the flux at inflow/outflow BCs is preserved
340 
341  //- Correct the phase properties
342  // other than the thermodynamics and turbulence
343  // which have special treatment
344  void correct();
345 
346  //- Read phaseProperties dictionary
347  virtual bool read(const dictionary& phaseProperties);
348 
349  //- Dummy Read for transportModel
350  virtual bool read()
351  {
352  return true;
353  }
354 };
355 
356 
357 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
358 
359 } // End namespace Foam
360 
361 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
362 
363 #endif
364 
365 // ************************************************************************* //
Foam::surfaceFields.
const dimensionedScalar & Cp() const
Definition: phaseModel.H:160
void correctInflowOutflow(surfaceScalarField &alphaPhi) const
Ensure that the flux at inflow/outflow BCs is preserved.
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
const PhaseCompressibleTurbulenceModel< phaseModel > & turbulence() const
Return the turbulence model.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Helper class to manage multi-specie phase properties.
virtual tmp< volScalarField > alpha() const =0
Thermal diffusivity for enthalpy of mixture [kg/m/s].
virtual tmp< volScalarField > alphaEff() const =0
Effective thermal turbulent diffusivity of mixture [kg/m/s].
const word & name() const
Definition: phaseModel.H:109
tmp< volScalarField > d() const
phaseModel(const word &phaseName, const volScalarField &p, const volScalarField &T)
Construct from components.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
const surfaceScalarField & phi() const
Definition: phaseModel.H:190
Class which solves the volume fraction equations for two phases.
const surfaceScalarField & alphaPhi() const
Definition: phaseModel.H:200
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
const volVectorField & U() const
Definition: phaseModel.H:170
virtual ~phaseModel()
Destructor.
virtual tmp< volScalarField > mu() const =0
Return the laminar dynamic viscosity.
virtual bool read()
Read phase properties dictionary.
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 phaseModel & otherPhase() const
Return the other phase in this two-phase system.
const dimensionedScalar & residualAlpha() const
Return the residual phase-fraction for given phase.
void correct()
Correct the laminar viscosity.
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
Templated abstract base class for multiphase compressible turbulence models.
virtual tmp< surfaceScalarField > alphaRhoPhi() const =0
Return the mass flux of the phase.
const dimensionedScalar & rho() const
Definition: phaseModel.H:165
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const phaseSystem & fluid() const
Return the system to which this phase belongs.
virtual tmp< volScalarField > alphahe() const =0
Thermal diffusivity for energy of mixture [kg/m/s].
virtual tmp< volScalarField > kappaEff() const =0
Effective thermal turbulent diffusivity for temperature.
Basic thermodynamic properties based on density.
Definition: rhoThermo.H:49
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
Namespace for OpenFOAM.