multiphaseMixtureThermo.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) 2013-2017 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::multiphaseMixtureThermo
29 
30 Description
31 
32 SourceFiles
33  multiphaseMixtureThermo.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef multiphaseMixtureThermo_H
38 #define multiphaseMixtureThermo_H
39 
40 #include "phaseModel.H"
41 #include "PtrDictionary.H"
42 #include "volFields.H"
43 #include "surfaceFields.H"
44 #include "rhoThermo.H"
45 #include "psiThermo.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 /*---------------------------------------------------------------------------*\
53  Class multiphaseMixtureThermo Declaration
54 \*---------------------------------------------------------------------------*/
55 
57 :
58  public psiThermo
59 {
60 public:
61 
62  class interfacePair
63  :
64  public Pair<word>
65  {
66  public:
67 
68  struct hash
69  {
70  label operator()(const interfacePair& key) const
71  {
72  return word::hash()(key.first()) + word::hash()(key.second());
73  }
74  };
75 
76 
77  // Constructors
78 
79  interfacePair() {} // = default
80 
81  interfacePair(const word& alpha1Name, const word& alpha2Name)
82  :
83  Pair<word>(alpha1Name, alpha2Name)
84  {}
85 
87  :
89  {}
90 
91 
92  // Friend Operators
93 
94  friend bool operator==
95  (
96  const interfacePair& a,
97  const interfacePair& b
98  )
99  {
100  return
101  (
102  ((a.first() == b.first()) && (a.second() == b.second()))
103  || ((a.first() == b.second()) && (a.second() == b.first()))
104  );
105  }
106 
107  friend bool operator!=
108  (
109  const interfacePair& a,
110  const interfacePair& b
111  )
112  {
113  return (!(a == b));
114  }
115  };
116 
117 
118 private:
119 
120  // Private data
121 
122  //- Dictionary of phases
124 
125  const fvMesh& mesh_;
126  const volVectorField& U_;
127  const surfaceScalarField& phi_;
128 
129  surfaceScalarField rhoPhi_;
130 
131  volScalarField alphas_;
132 
134  sigmaTable;
135 
136  sigmaTable sigmas_;
137  dimensionSet dimSigma_;
138 
139  //- Stabilisation for normalisation of the interface normal
140  const dimensionedScalar deltaN_;
141 
142 
143  // Private member functions
144 
145  void calcAlphas();
146 
147  void solveAlphas(const scalar cAlpha);
148 
150  (
151  const volScalarField& alpha1,
152  const volScalarField& alpha2
153  ) const;
154 
156  (
157  const volScalarField& alpha1,
158  const volScalarField& alpha2
159  ) const;
160 
161  void correctContactAngle
162  (
163  const phaseModel& alpha1,
164  const phaseModel& alpha2,
165  surfaceVectorField::Boundary& nHatb
166  ) const;
167 
169  (
170  const phaseModel& alpha1,
171  const phaseModel& alpha2
172  ) const;
173 
174 
175 public:
176 
177  //- Runtime type information
178  TypeName("multiphaseMixtureThermo");
179 
180 
181  // Constructors
182 
183  //- Construct from components
185  (
186  const volVectorField& U,
187  const surfaceScalarField& phi
188  );
189 
190 
191  //- Destructor
192  virtual ~multiphaseMixtureThermo() = default;
193 
194 
195  // Member Functions
196 
197  //- Return the phases
198  const PtrDictionary<phaseModel>& phases() const
199  {
200  return phases_;
201  }
202 
203  //- Return non-const access to the phases
205  {
206  return phases_;
207  }
208 
209  //- Return the velocity
210  const volVectorField& U() const
211  {
212  return U_;
213  }
214 
215  //- Return the volumetric flux
216  const surfaceScalarField& phi() const
217  {
218  return phi_;
219  }
220 
221  const surfaceScalarField& rhoPhi() const
222  {
223  return rhoPhi_;
224  }
225 
226  //- Update properties
227  virtual void correct();
228 
229  //- Update densities for given pressure change
230  void correctRho(const volScalarField& dp);
231 
232  //- Return the name of the thermo physics
233  virtual word thermoName() const;
234 
235  //- Return true if the equation of state is incompressible
236  // i.e. rho != f(p)
237  virtual bool incompressible() const;
238 
239  //- Return true if the equation of state is isochoric
240  // i.e. rho = const
241  virtual bool isochoric() const;
242 
243 
244  // Access to thermodynamic state variables
245 
246  //- Enthalpy/Internal energy [J/kg]
247  // Non-const access allowed for transport equations
248  virtual volScalarField& he()
249  {
251  return phases_[0].thermo().he();
252  }
253 
254  //- Enthalpy/Internal energy [J/kg]
255  virtual const volScalarField& he() const
256  {
258  return phases_[0].thermo().he();
259  }
260 
261  //- Enthalpy/Internal energy
262  // for given pressure and temperature [J/kg]
263  virtual tmp<volScalarField> he
264  (
265  const volScalarField& p,
266  const volScalarField& T
267  ) const;
268 
269  //- Enthalpy/Internal energy for cell-set [J/kg]
270  virtual tmp<scalarField> he
271  (
272  const scalarField& p,
273  const scalarField& T,
274  const labelList& cells
275  ) const;
276 
277  //- Enthalpy/Internal energy for patch [J/kg]
278  virtual tmp<scalarField> he
279  (
280  const scalarField& p,
281  const scalarField& T,
282  const label patchi
283  ) const;
284 
285  //- Chemical enthalpy [J/kg]
286  virtual tmp<volScalarField> hc() const;
287 
288  //- Temperature from enthalpy/internal energy for cell-set
289  virtual tmp<scalarField> THE
290  (
291  const scalarField& h,
292  const scalarField& p,
293  const scalarField& T0, // starting temperature
294  const labelList& cells
295  ) const;
296 
297  //- Temperature from enthalpy/internal energy for patch
298  virtual tmp<scalarField> THE
299  (
300  const scalarField& h,
301  const scalarField& p,
302  const scalarField& T0, // starting temperature
303  const label patchi
304  ) const;
305 
306 
307  // Fields derived from thermodynamic state variables
308 
309  //- Density [kg/m^3]
310  virtual tmp<volScalarField> rho() const;
311 
312  //- Density for patch [kg/m^3]
313  virtual tmp<scalarField> rho(const label patchi) const;
314 
315  //- Heat capacity at constant pressure [J/kg/K]
316  virtual tmp<volScalarField> Cp() const;
317 
318  //- Heat capacity at constant pressure for patch [J/kg/K]
319  virtual tmp<scalarField> Cp
320  (
321  const scalarField& p,
322  const scalarField& T,
323  const label patchi
324  ) const;
325 
326  //- Heat capacity at constant volume [J/kg/K]
327  virtual tmp<volScalarField> Cv() const;
328 
329  //- Heat capacity at constant volume for patch [J/kg/K]
330  virtual tmp<scalarField> Cv
331  (
332  const scalarField& p,
333  const scalarField& T,
334  const label patchi
335  ) const;
336 
337  //- Gamma = Cp/Cv []
338  virtual tmp<volScalarField> gamma() const;
339 
340  //- Gamma = Cp/Cv for patch []
341  virtual tmp<scalarField> gamma
342  (
343  const scalarField& p,
344  const scalarField& T,
345  const label patchi
346  ) const;
347 
348  //- Heat capacity at constant pressure/volume [J/kg/K]
349  virtual tmp<volScalarField> Cpv() const;
350 
351  //- Heat capacity at constant pressure/volume for patch [J/kg/K]
352  virtual tmp<scalarField> Cpv
353  (
354  const scalarField& p,
355  const scalarField& T,
356  const label patchi
357  ) const;
358 
359  //- Heat capacity ratio []
360  virtual tmp<volScalarField> CpByCpv() const;
361 
362  //- Heat capacity ratio for patch []
363  virtual tmp<scalarField> CpByCpv
364  (
365  const scalarField& p,
366  const scalarField& T,
367  const label patchi
368  ) const;
369 
370  //- Molecular weight [kg/kmol]
371  virtual tmp<volScalarField> W() const;
372 
373 
374  // Fields derived from transport state variables
375 
376  //- Kinematic viscosity of mixture [m^2/s]
377  virtual tmp<volScalarField> nu() const;
378 
379  //- Kinematic viscosity of mixture for patch [m^2/s]
380  virtual tmp<scalarField> nu(const label patchi) const;
381 
382  //- Thermal diffusivity for temperature of mixture [J/m/s/K]
383  virtual tmp<volScalarField> kappa() const;
384 
385  //- Thermal diffusivity of mixture for patch [J/m/s/K]
386  virtual tmp<scalarField> kappa
387  (
388  const label patchi
389  ) const;
390 
391 
392  //- Thermal diffusivity for energy of mixture [kg/m/s]
393  virtual tmp<volScalarField> alphahe() const;
394 
395  //- Thermal diffusivity for energy of mixture for patch [kg/m/s]
396  virtual tmp<scalarField> alphahe(const label patchi) const;
397 
398  //- Effective thermal diffusivity of mixture [J/m/s/K]
400  (
401  const volScalarField& alphat
402  ) const;
403 
404  //- Effective thermal diffusivity of mixture for patch [J/m/s/K]
405  virtual tmp<scalarField> kappaEff
406  (
407  const scalarField& alphat,
408  const label patchi
409  ) const;
410 
411  //- Effective thermal diffusivity of mixture [J/m/s/K]
413  (
414  const volScalarField& alphat
415  ) const;
416 
417  //- Effective thermal diffusivity of mixture for patch [J/m/s/K]
418  virtual tmp<scalarField> alphaEff
419  (
420  const scalarField& alphat,
421  const label patchi
422  ) const;
423 
424 
425  //- Return the phase-averaged reciprocal Cv
426  tmp<volScalarField> rCv() const;
427 
429 
430  //- Indicator of the proximity of the interface
431  // Field values are 1 near and 0 away for the interface.
433 
434  //- Solve for the mixture phase-fractions
435  void solve();
436 };
437 
438 
439 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
440 
441 } // End namespace Foam
442 
443 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
444 
445 #endif
446 
447 // ************************************************************************* //
Foam::Pair::second
const T & second() const noexcept
Return second element, which is also the last element.
Definition: PairI.H:122
Foam::multiphaseMixtureThermo::nu
virtual tmp< volScalarField > nu() const
Kinematic viscosity of mixture [m^2/s].
volFields.H
Foam::phaseModel
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:54
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::multiphaseMixtureThermo::alphahe
virtual tmp< volScalarField > alphahe() const
Thermal diffusivity for energy of mixture [kg/m/s].
Foam::multiphaseMixtureThermo::CpByCpv
virtual tmp< volScalarField > CpByCpv() const
Heat capacity ratio [].
Foam::multiphaseMixtureThermo::U
const volVectorField & U() const
Return the velocity.
Definition: multiphaseMixtureThermo.H:209
Foam::basicThermo::p
virtual volScalarField & p()
Pressure [Pa].
Definition: basicThermo.C:512
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::multiphaseMixtureThermo::interfacePair
Definition: multiphaseMixtureThermo.H:61
alpha2
const volScalarField & alpha2
Definition: setRegionFluidFields.H:9
Foam::PtrDictionary
Template dictionary class which manages the storage associated with it.
Definition: PtrDictionary.H:55
PtrDictionary.H
Foam::multiphaseMixtureThermo::~multiphaseMixtureThermo
virtual ~multiphaseMixtureThermo()=default
Destructor.
Foam::multiphaseMixtureThermo::he
virtual const volScalarField & he() const
Enthalpy/Internal energy [J/kg].
Definition: multiphaseMixtureThermo.H:254
Foam::multiphaseMixtureThermo::interfacePair::hash::operator()
label operator()(const interfacePair &key) const
Definition: multiphaseMixtureThermo.H:69
surfaceFields.H
Foam::surfaceFields.
Foam::multiphaseMixtureThermo::Cpv
virtual tmp< volScalarField > Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
alpha1
const volScalarField & alpha1
Definition: setRegionFluidFields.H:8
Foam::multiphaseMixtureThermo::interfacePair::interfacePair
interfacePair()
Definition: multiphaseMixtureThermo.H:78
Foam::multiphaseMixtureThermo::interfacePair::interfacePair
interfacePair(const phaseModel &alpha1, const phaseModel &alpha2)
Definition: multiphaseMixtureThermo.H:85
Foam::multiphaseMixtureThermo::isochoric
virtual bool isochoric() const
Return true if the equation of state is isochoric.
Foam::dimensionSet
Dimension set for the base types.
Definition: dimensionSet.H:65
Foam::multiphaseMixtureThermo
Definition: multiphaseMixtureThermo.H:55
Foam::baseIOdictionary::name
const word & name() const
Definition: baseIOdictionary.C:82
Foam::multiphaseMixtureThermo::nearInterface
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
Foam::multiphaseMixtureThermo::surfaceTensionForce
tmp< surfaceScalarField > surfaceTensionForce() const
Foam::multiphaseMixtureThermo::phi
const surfaceScalarField & phi() const
Return the volumetric flux.
Definition: multiphaseMixtureThermo.H:215
K
CGAL::Exact_predicates_exact_constructions_kernel K
Definition: CGALTriangulation3DKernel.H:58
rhoThermo.H
Foam::multiphaseMixtureThermo::kappaEff
virtual tmp< volScalarField > kappaEff(const volScalarField &alphat) const
Effective thermal diffusivity of mixture [J/m/s/K].
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:445
Foam::multiphaseMixtureThermo::interfacePair::interfacePair
interfacePair(const word &alpha1Name, const word &alpha2Name)
Definition: multiphaseMixtureThermo.H:80
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::Field< scalar >
Foam::constant::universal::h
const dimensionedScalar h
Planck constant.
Definition: setRegionSolidFields.H:33
Foam::string::hash
Hashing function for string and derived string classes.
Definition: string.H:151
Foam::psiThermo
Basic thermodynamic properties based on compressibility.
Definition: psiThermo.H:55
Foam::multiphaseMixtureThermo::Cv
virtual tmp< volScalarField > Cv() const
Heat capacity at constant volume [J/kg/K].
Foam::multiphaseMixtureThermo::rho
virtual tmp< volScalarField > rho() const
Density [kg/m^3].
Foam::multiphaseMixtureThermo::he
virtual volScalarField & he()
Enthalpy/Internal energy [J/kg].
Definition: multiphaseMixtureThermo.H:247
Foam::multiphaseMixtureThermo::thermoName
virtual word thermoName() const
Return the name of the thermo physics.
Foam::dimensioned< scalar >
Foam::multiphaseMixtureThermo::TypeName
TypeName("multiphaseMixtureThermo")
Runtime type information.
Foam::multiphaseMixtureThermo::W
virtual tmp< volScalarField > W() const
Molecular weight [kg/kmol].
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:83
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::multiphaseMixtureThermo::phases
const PtrDictionary< phaseModel > & phases() const
Return the phases.
Definition: multiphaseMixtureThermo.H:197
Foam::multiphaseMixtureThermo::interfacePair::hash
Definition: multiphaseMixtureThermo.H:67
Foam::multiphaseMixtureThermo::hc
virtual tmp< volScalarField > hc() const
Chemical enthalpy [J/kg].
Foam::HashTable< scalar, interfacePair, interfacePair::hash >
psiThermo.H
Foam::multiphaseMixtureThermo::rCv
tmp< volScalarField > rCv() const
Return the phase-averaged reciprocal Cv.
Foam::Pair
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:54
Foam::basicThermo::T
virtual const volScalarField & T() const
Temperature [K].
Definition: basicThermo.C:524
Foam::multiphaseMixtureThermo::THE
virtual tmp< scalarField > THE(const scalarField &h, const scalarField &p, const scalarField &T0, const labelList &cells) const
Temperature from enthalpy/internal energy for cell-set.
Foam::multiphaseMixtureThermo::phases
PtrDictionary< phaseModel > & phases()
Return non-const access to the phases.
Definition: multiphaseMixtureThermo.H:203
Foam::List< label >
Foam::multiphaseMixtureThermo::Cp
virtual tmp< volScalarField > Cp() const
Heat capacity at constant pressure [J/kg/K].
Foam::multiphaseMixtureThermo::alphaEff
virtual tmp< volScalarField > alphaEff(const volScalarField &alphat) const
Effective thermal diffusivity of mixture [J/m/s/K].
Foam::multiphaseMixtureThermo::multiphaseMixtureThermo
multiphaseMixtureThermo(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
Foam::multiphaseMixtureThermo::kappa
virtual tmp< volScalarField > kappa() const
Thermal diffusivity for temperature of mixture [J/m/s/K].
Foam::multiphaseMixtureThermo::gamma
virtual tmp< volScalarField > gamma() const
Gamma = Cp/Cv [].
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::multiphaseMixtureThermo::solve
void solve()
Solve for the mixture phase-fractions.
Foam::multiphaseMixtureThermo::incompressible
virtual bool incompressible() const
Return true if the equation of state is incompressible.
Foam::multiphaseMixtureThermo::correct
virtual void correct()
Update properties.
Foam::GeometricField< vector, fvPatchField, volMesh >
T0
scalar T0
Definition: createFields.H:22
Foam::multiphaseMixtureThermo::rhoPhi
const surfaceScalarField & rhoPhi() const
Definition: multiphaseMixtureThermo.H:220
Foam::multiphaseMixtureThermo::correctRho
void correctRho(const volScalarField &dp)
Update densities for given pressure change.