multiphaseMixtureThermo.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) 2013-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::multiphaseMixtureThermo
26 
27 Description
28 
29 SourceFiles
30  multiphaseMixtureThermo.C
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #ifndef multiphaseMixtureThermo_H
35 #define multiphaseMixtureThermo_H
36 
37 #include "phaseModel.H"
38 #include "PtrDictionary.H"
39 #include "volFields.H"
40 #include "surfaceFields.H"
41 #include "rhoThermo.H"
42 #include "psiThermo.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 /*---------------------------------------------------------------------------*\
50  Class multiphaseMixtureThermo Declaration
51 \*---------------------------------------------------------------------------*/
52 
54 :
55  public psiThermo
56 {
57 public:
58 
59  class interfacePair
60  :
61  public Pair<word>
62  {
63  public:
64 
65  class hash
66  :
67  public Hash<interfacePair>
68  {
69  public:
70 
71  hash()
72  {}
73 
74  label operator()(const interfacePair& key) const
75  {
76  return word::hash()(key.first()) + word::hash()(key.second());
77  }
78  };
79 
80 
81  // Constructors
82 
84  {}
85 
86  interfacePair(const word& alpha1Name, const word& alpha2Name)
87  :
88  Pair<word>(alpha1Name, alpha2Name)
89  {}
90 
92  :
93  Pair<word>(alpha1.name(), alpha2.name())
94  {}
95 
96 
97  // Friend Operators
98 
99  friend bool operator==
100  (
101  const interfacePair& a,
102  const interfacePair& b
103  )
104  {
105  return
106  (
107  ((a.first() == b.first()) && (a.second() == b.second()))
108  || ((a.first() == b.second()) && (a.second() == b.first()))
109  );
110  }
111 
112  friend bool operator!=
113  (
114  const interfacePair& a,
115  const interfacePair& b
116  )
117  {
118  return (!(a == b));
119  }
120  };
121 
122 
123 private:
124 
125  // Private Data
126 
127  //- Dictionary of phases
129 
130  const fvMesh& mesh_;
131  const volVectorField& U_;
132  const surfaceScalarField& phi_;
133 
134  surfaceScalarField rhoPhi_;
135 
136  volScalarField alphas_;
137 
139  sigmaTable;
140 
141  sigmaTable sigmas_;
142  dimensionSet dimSigma_;
143 
144  //- Stabilisation for normalisation of the interface normal
145  const dimensionedScalar deltaN_;
146 
147 
148  // Private Member Functions
149 
150  void calcAlphas();
151 
152  void solveAlphas(const scalar cAlpha);
153 
155  (
156  const volScalarField& alpha1,
157  const volScalarField& alpha2
158  ) const;
159 
161  (
162  const volScalarField& alpha1,
163  const volScalarField& alpha2
164  ) const;
165 
166  void correctContactAngle
167  (
168  const phaseModel& alpha1,
169  const phaseModel& alpha2,
170  surfaceVectorField::Boundary& nHatb
171  ) const;
172 
174  (
175  const phaseModel& alpha1,
176  const phaseModel& alpha2
177  ) const;
178 
179 
180 public:
181 
182  //- Runtime type information
183  TypeName("multiphaseMixtureThermo");
184 
185 
186  // Constructors
187 
188  //- Construct from components
190  (
191  const volVectorField& U,
192  const surfaceScalarField& phi
193  );
194 
195 
196  //- Destructor
197  virtual ~multiphaseMixtureThermo()
198  {}
199 
200 
201  // Member Functions
202 
203  //- Return the phases
204  const PtrDictionary<phaseModel>& phases() const
205  {
206  return phases_;
207  }
208 
209  //- Return non-const access to the phases
211  {
212  return phases_;
213  }
214 
215  //- Return the velocity
216  const volVectorField& U() const
217  {
218  return U_;
219  }
220 
221  //- Return the volumetric flux
222  const surfaceScalarField& phi() const
223  {
224  return phi_;
225  }
227  const surfaceScalarField& rhoPhi() const
228  {
229  return rhoPhi_;
230  }
231 
232  //- Update properties
233  virtual void correct();
234 
235  //- Update densities for given pressure change
236  void correctRho(const volScalarField& dp);
237 
238  //- Return the name of the thermo physics
239  virtual word thermoName() const;
240 
241  //- Return true if the equation of state is incompressible
242  // i.e. rho != f(p)
243  virtual bool incompressible() const;
244 
245  //- Return true if the equation of state is isochoric
246  // i.e. rho = const
247  virtual bool isochoric() const;
248 
249 
250  // Access to thermodynamic state variables
251 
252  //- Enthalpy/Internal energy [J/kg]
253  // Non-const access allowed for transport equations
254  virtual volScalarField& he()
255  {
257  return phases_[0].thermo().he();
258  }
259 
260  //- Enthalpy/Internal energy [J/kg]
261  virtual const volScalarField& he() const
262  {
264  return phases_[0].thermo().he();
265  }
266 
267  //- Enthalpy/Internal energy
268  // for given pressure and temperature [J/kg]
269  virtual tmp<volScalarField> he
270  (
271  const volScalarField& p,
272  const volScalarField& T
273  ) const;
274 
275  //- Enthalpy/Internal energy for cell-set [J/kg]
276  virtual tmp<scalarField> he
277  (
278  const scalarField& T,
279  const labelList& cells
280  ) const;
281 
282  //- Enthalpy/Internal energy for patch [J/kg]
283  virtual tmp<scalarField> he
284  (
285  const scalarField& T,
286  const label patchi
287  ) const;
288 
289  //- Sensible enthalpy [J/kg]
290  virtual tmp<volScalarField> hs() const;
291 
292  //- Sensible enthalpy
293  // for given pressure and temperature [J/kg]
294  virtual tmp<volScalarField> hs
295  (
296  const volScalarField& p,
297  const volScalarField& T
298  ) const;
299 
300  //- Sensible enthalpy for cell-set [J/kg]
301  virtual tmp<scalarField> hs
302  (
303  const scalarField& T,
304  const labelList& cells
305  ) const;
306 
307  //- Sensible enthalpy for patch [J/kg]
308  virtual tmp<scalarField> hs
309  (
310  const scalarField& T,
311  const label patchi
312  ) const;
313 
314  //- Absolute enthalpy [J/kg]
315  virtual tmp<volScalarField> ha() const;
316 
317  //- Absolute enthalpy
318  // for given pressure and temperature [J/kg]
319  virtual tmp<volScalarField> ha
320  (
321  const volScalarField& p,
322  const volScalarField& T
323  ) const;
324 
325  //- Absolute enthalpy for cell-set [J/kg]
326  virtual tmp<scalarField> ha
327  (
328  const scalarField& T,
329  const labelList& cells
330  ) const;
331 
332  //- Absolute enthalpy for patch [J/kg]
333  virtual tmp<scalarField> ha
334  (
335  const scalarField& T,
336  const label patchi
337  ) const;
338 
339  //- Enthalpy of formation [J/kg]
340  virtual tmp<volScalarField> hc() const;
341 
342  //- Temperature from enthalpy/internal energy for cell-set
343  virtual tmp<scalarField> THE
344  (
345  const scalarField& h,
346  const scalarField& T0, // starting temperature
347  const labelList& cells
348  ) const;
349 
350  //- Temperature from enthalpy/internal energy for patch
351  virtual tmp<scalarField> THE
352  (
353  const scalarField& h,
354  const scalarField& T0, // starting temperature
355  const label patchi
356  ) const;
357 
358 
359  // Fields derived from thermodynamic state variables
360 
361  //- Density [kg/m^3]
362  virtual tmp<volScalarField> rho() const;
363 
364  //- Density for patch [kg/m^3]
365  virtual tmp<scalarField> rho(const label patchi) const;
366 
367  //- Heat capacity at constant pressure [J/kg/K]
368  virtual tmp<volScalarField> Cp() const;
369 
370  //- Heat capacity at constant pressure for patch [J/kg/K]
371  virtual tmp<scalarField> Cp
372  (
373  const scalarField& T,
374  const label patchi
375  ) const;
376 
377  //- Heat capacity at constant volume [J/kg/K]
378  virtual tmp<volScalarField> Cv() const;
379 
380  //- Heat capacity at constant volume for patch [J/kg/K]
381  virtual tmp<scalarField> Cv
382  (
383  const scalarField& T,
384  const label patchi
385  ) const;
386 
387  //- Gamma = Cp/Cv []
388  virtual tmp<volScalarField> gamma() const;
389 
390  //- Gamma = Cp/Cv for patch []
391  virtual tmp<scalarField> gamma
392  (
393  const scalarField& T,
394  const label patchi
395  ) const;
396 
397  //- Heat capacity at constant pressure/volume [J/kg/K]
398  virtual tmp<volScalarField> Cpv() const;
399 
400  //- Heat capacity at constant pressure/volume for patch [J/kg/K]
401  virtual tmp<scalarField> Cpv
402  (
403  const scalarField& T,
404  const label patchi
405  ) const;
406 
407  //- Heat capacity ratio []
408  virtual tmp<volScalarField> CpByCpv() const;
409 
410  //- Heat capacity ratio for patch []
411  virtual tmp<scalarField> CpByCpv
412  (
413  const scalarField& T,
414  const label patchi
415  ) const;
416 
417  //- Molecular weight [kg/kmol]
418  virtual tmp<volScalarField> W() const;
419 
420  //- Molecular weight for patch [kg/kmol]
421  virtual tmp<scalarField> W(const label patchi) const;
422 
423 
424  // Fields derived from transport state variables
425 
426  //- Kinematic viscosity of mixture [m^2/s]
427  virtual tmp<volScalarField> nu() const;
428 
429  //- Kinematic viscosity of mixture for patch [m^2/s]
430  virtual tmp<scalarField> nu(const label patchi) const;
431 
432  //- Thermal diffusivity for temperature of mixture [W/m/K]
433  virtual tmp<volScalarField> kappa() const;
434 
435  //- Thermal diffusivity of mixture for patch [W/m/K]
436  virtual tmp<scalarField> kappa
437  (
438  const label patchi
439  ) const;
440 
441  //- Thermal diffusivity for energy of mixture [kg/m/s]
442  virtual tmp<volScalarField> alphahe() const;
443 
444  //- Thermal diffusivity for energy of mixture for patch [kg/m/s]
445  virtual tmp<scalarField> alphahe(const label patchi) const;
446 
447  //- Effective thermal diffusivity of mixture [W/m/K]
449  (
450  const volScalarField& alphat
451  ) const;
452 
453  //- Effective thermal diffusivity of mixture for patch [W/m/K]
454  virtual tmp<scalarField> kappaEff
455  (
456  const scalarField& alphat,
457  const label patchi
458  ) const;
459 
460  //- Effective thermal diffusivity of mixture [W/m/K]
462  (
463  const volScalarField& alphat
464  ) const;
465 
466  //- Effective thermal diffusivity of mixture for patch [W/m/K]
467  virtual tmp<scalarField> alphaEff
468  (
469  const scalarField& alphat,
470  const label patchi
471  ) const;
472 
473 
474  //- Return the phase-averaged reciprocal Cv
475  tmp<volScalarField> rCv() const;
476 
478 
479  //- Indicator of the proximity of the interface
480  // Field values are 1 near and 0 away for the interface.
482 
483  //- Solve for the mixture phase-fractions
484  void solve();
485 };
486 
487 
488 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
489 
490 } // End namespace Foam
491 
492 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
493 
494 #endif
495 
496 // ************************************************************************* //
Template dictionary class which manages the storage associated with it.
Definition: PtrDictionary.H:53
Foam::surfaceFields.
const volVectorField & U() const
Return the velocity.
Hashing function class, shared by all the derived classes.
Definition: string.H:92
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 surfaceScalarField & phi() const
Return the volumetric flux.
virtual tmp< volScalarField > gamma() const
Gamma = Cp/Cv [].
virtual tmp< volScalarField > hc() const
Enthalpy of formation [J/kg].
virtual tmp< volScalarField > ha() const
Absolute enthalpy [J/kg].
virtual volScalarField & p()
Pressure [Pa].
Definition: fluidThermo.C:79
virtual const volScalarField & T() const
Temperature [K].
Definition: basicThermo.C:482
virtual ~multiphaseMixtureThermo()
Destructor.
CGAL::Exact_predicates_exact_constructions_kernel K
virtual tmp< volScalarField > Cp() const
Heat capacity at constant pressure [J/kg/K].
alpha2
Definition: alphaEqn.H:115
virtual tmp< volScalarField > CpByCpv() const
Heat capacity ratio [].
void correctRho(const volScalarField &dp)
Update densities for given pressure change.
void solve()
Solve for the mixture phase-fractions.
Dimension set for the base types.
Definition: dimensionSet.H:120
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
const cellShapeList & cells
virtual tmp< volScalarField > alphahe() const
Thermal diffusivity for energy of mixture [kg/m/s].
virtual tmp< volScalarField > rho() const
Density [kg/m^3].
const dimensionedScalar & b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
tmp< volScalarField > rCv() const
Return the phase-averaged reciprocal Cv.
const PtrDictionary< phaseModel > & phases() const
Return the phases.
A class for handling words, derived from string.
Definition: word.H:59
volScalarField & alpha1(mixture.alpha1())
virtual tmp< scalarField > THE(const scalarField &h, const scalarField &T0, const labelList &cells) const
Temperature from enthalpy/internal energy for cell-set.
const word & name() const
Name function is needed to disambiguate those inherited.
virtual tmp< volScalarField > nu() const
Kinematic viscosity of mixture [m^2/s].
label operator()(const interfacePair &key) const
const Type & second() const
Return second.
Definition: Pair.H:99
TypeName("multiphaseMixtureThermo")
Runtime type information.
virtual tmp< volScalarField > kappaEff(const volScalarField &alphat) const
Effective thermal diffusivity of mixture [W/m/K].
multiphaseMixtureThermo(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
label patchi
virtual tmp< volScalarField > W() const
Molecular weight [kg/kmol].
const surfaceScalarField & rhoPhi() const
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Hash function class for primitives. All non-primitives used to hash entries on hash tables likely nee...
Definition: Hash.H:54
virtual tmp< volScalarField > hs() const
Sensible enthalpy [J/kg].
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
const dimensionedScalar & h
Planck constant.
virtual tmp< volScalarField > Cv() const
Heat capacity at constant volume [J/kg/K].
A class for managing temporary objects.
Definition: PtrList.H:53
tmp< surfaceScalarField > surfaceTensionForce() const
virtual tmp< volScalarField > alphaEff(const volScalarField &alphat) const
Effective thermal diffusivity of mixture [W/m/K].
virtual tmp< volScalarField > Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:366
virtual tmp< volScalarField > kappa() const
Thermal diffusivity for temperature of mixture [W/m/K].
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
const Type & first() const
Return first.
Definition: Pair.H:87
Basic thermodynamic properties based on compressibility.
Definition: psiThermo.H:49
virtual bool isochoric() const
Return true if the equation of state is isochoric.
virtual word thermoName() const
Return the name of the thermo physics.
virtual volScalarField & he()
Enthalpy/Internal energy [J/kg].
Namespace for OpenFOAM.
scalar T0
Definition: createFields.H:22
virtual bool incompressible() const
Return true if the equation of state is incompressible.
virtual void correct()
Update properties.