compressibleMultiphaseMixture.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-2022 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::compressibleMultiphaseMixture
26 
27 Description
28 
29 SourceFiles
30  compressibleMultiphaseMixture.C
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #ifndef compressibleMultiphaseMixture_H
35 #define compressibleMultiphaseMixture_H
36 
37 #include "phaseModel.H"
38 #include "PtrDictionary.H"
39 #include "volFields.H"
40 #include "surfaceFields.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 /*---------------------------------------------------------------------------*\
48  Class compressibleMultiphaseMixture Declaration
49 \*---------------------------------------------------------------------------*/
50 
52 :
53  public IOdictionary,
54  public viscosity
55 {
56 public:
57 
58  class interfacePair
59  :
60  public Pair<word>
61  {
62  public:
63 
64  class hash
65  :
66  public Hash<interfacePair>
67  {
68  public:
69 
70  hash()
71  {}
72 
73  label operator()(const interfacePair& key) const
74  {
75  return word::hash()(key.first()) + word::hash()(key.second());
76  }
77  };
78 
79 
80  // Constructors
81 
83  {}
84 
85  interfacePair(const word& alpha1Name, const word& alpha2Name)
86  :
87  Pair<word>(alpha1Name, alpha2Name)
88  {}
89 
91  :
92  Pair<word>(alpha1.name(), alpha2.name())
93  {}
94 
95 
96  // Friend Operators
97 
98  friend bool operator==
99  (
100  const interfacePair& a,
101  const interfacePair& b
102  )
103  {
104  return
105  (
106  ((a.first() == b.first()) && (a.second() == b.second()))
107  || ((a.first() == b.second()) && (a.second() == b.first()))
108  );
109  }
110 
111  friend bool operator!=
112  (
113  const interfacePair& a,
114  const interfacePair& b
115  )
116  {
117  return (!(a == b));
118  }
119  };
120 
121 
122 private:
123 
124  // Private Data
125 
126  const fvMesh& mesh_;
127 
128  //- Pressure
129  volScalarField p_;
130 
131  //- Mixture temperature
132  volScalarField T_;
133 
134  //- Dictionary of phases
136 
137  volScalarField rho_;
138 
139  const volVectorField& U_;
140  const surfaceScalarField& phi_;
141 
142  surfaceScalarField rhoPhi_;
143 
144  volScalarField alphas_;
145 
147  sigmaTable;
148 
149  sigmaTable sigmas_;
150  dimensionSet dimSigma_;
151 
152  //- Stabilisation for normalisation of the interface normal
153  const dimensionedScalar deltaN_;
154 
155 
156  // Private Member Functions
157 
158  void calcAlphas();
159 
160  void solveAlphas(const scalar cAlpha);
161 
163  (
164  const volScalarField& alpha1,
165  const volScalarField& alpha2
166  ) const;
167 
169  (
170  const volScalarField& alpha1,
171  const volScalarField& alpha2
172  ) const;
173 
175  (
176  const phaseModel& alpha1,
177  const phaseModel& alpha2
178  ) const;
179 
180 
181 public:
182 
183  //- Runtime type information
184  TypeName("compressibleMultiphaseMixture");
185 
186 
187  // Constructors
188 
189  //- Construct from components
191  (
192  const volVectorField& U,
193  const surfaceScalarField& phi
194  );
195 
196 
197  //- Destructor
199  {}
200 
201 
202  // Member Functions
203 
204  //- Return the phases
205  const PtrDictionary<phaseModel>& phases() const
206  {
207  return phases_;
208  }
209 
210  //- Return non-const access to the phases
212  {
213  return phases_;
214  }
215 
216  //- Return pressure [Pa]
217  volScalarField& p()
218  {
219  return p_;
220  }
221 
222  //- Return mixture temperature [K]
223  volScalarField& T()
224  {
225  return T_;
226  }
227 
228  //- Return mixture density [kg/m^3]
229  const volScalarField& rho() const
230  {
231  return rho_;
232  }
233 
234  //- Return mixture velocity
235  const volVectorField& U() const
236  {
237  return U_;
238  }
239 
240  //- Return the volumetric flux
241  const surfaceScalarField& phi() const
242  {
243  return phi_;
244  }
246  const surfaceScalarField& rhoPhi() const
247  {
248  return rhoPhi_;
249  }
250 
251  //- Correct the thermodynamics of each phase
252  virtual void correctThermo();
253 
254  //- Update properties
255  virtual void correct();
256 
257  //- Update densities for given pressure change
258  void correctRho(const volScalarField& dp);
259 
260  //- Kinematic viscosity of mixture [m^2/s]
261  virtual tmp<volScalarField> nu() const;
262 
263  //- Kinematic viscosity of mixture for patch [m^2/s]
264  virtual tmp<scalarField> nu(const label patchi) const;
265 
266  //- Effective thermal diffusivity of mixture [W/m/K]
268  (
269  const volScalarField& alphat
270  ) const;
271 
272  //- Return the phase-averaged reciprocal Cv
273  tmp<volScalarField> rCv() const;
274 
276 
277  //- Indicator of the proximity of the interface
278  // Field values are 1 near and 0 away for the interface.
280 
281  //- Solve for the mixture phase-fractions
282  void solve();
283 };
284 
285 
286 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
287 
288 } // End namespace Foam
289 
290 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
291 
292 #endif
293 
294 // ************************************************************************* //
Template dictionary class which manages the storage associated with it.
Definition: PtrDictionary.H:53
Foam::surfaceFields.
const volScalarField & rho() const
Return mixture density [kg/m^3].
virtual tmp< volScalarField > nu() const
Kinematic viscosity of mixture [m^2/s].
Hashing function class, shared by all the derived classes.
Definition: string.H:92
compressibleMultiphaseMixture(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
const surfaceScalarField & phi() const
Return the volumetric flux.
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
volScalarField & alpha1(mixture.alpha1())
volScalarField & p()
Return pressure [Pa].
CGAL::Exact_predicates_exact_constructions_kernel K
alpha2
Definition: alphaEqn.H:115
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
Dimension set for the base types.
Definition: dimensionSet.H:121
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
const fileName & name() const
Return the dictionary name.
Definition: dictionary.H:109
void solve()
Solve for the mixture phase-fractions.
tmp< volScalarField > rCv() const
Return the phase-averaged reciprocal Cv.
const volVectorField & U() const
Return mixture velocity.
A class for handling words, derived from string.
Definition: word.H:59
void correctRho(const volScalarField &dp)
Update densities for given pressure change.
virtual void correctThermo()
Correct the thermodynamics of each phase.
TypeName("compressibleMultiphaseMixture")
Runtime type information.
volScalarField & T()
Return mixture temperature [K].
Abstract base class for all fluid physical properties.
Definition: viscosity.H:49
const Type & second() const
Return second.
Definition: Pair.H:110
virtual tmp< volScalarField > alphaEff(const volScalarField &alphat) const
Effective thermal diffusivity of mixture [W/m/K].
const surfaceScalarField & rhoPhi() const
label patchi
virtual void correct()
Update properties.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
Hash function class for primitives. All non-primitives used to hash entries on hash tables likely nee...
Definition: Hash.H:52
const PtrDictionary< phaseModel > & phases() const
Return the phases.
A class for managing temporary objects.
Definition: PtrList.H:53
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:98
tmp< surfaceScalarField > surfaceTensionForce() const
Namespace for OpenFOAM.