multiphaseSystem.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 | Copyright (C) 2011-2016 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::multiphaseSystem
26 
27 Description
28  Incompressible multi-phase mixture with built in solution for the
29  phase fractions with interface compression for interface-capturing.
30 
31  Derived from transportModel so that it can be unsed in conjunction with
32  the incompressible turbulence models.
33 
34  Surface tension and contact-angle is handled for the interface
35  between each phase-pair.
36 
37 SourceFiles
38  multiphaseSystem.C
39 
40 \*---------------------------------------------------------------------------*/
41 
42 #ifndef multiphaseSystem_H
43 #define multiphaseSystem_H
44 
46 #include "IOdictionary.H"
47 #include "phaseModel.H"
48 #include "PtrDictionary.H"
49 #include "volFields.H"
50 #include "surfaceFields.H"
51 #include "dragModel.H"
52 #include "HashPtrTable.H"
53 
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
56 namespace Foam
57 {
58 
59 /*---------------------------------------------------------------------------*\
60  Class multiphaseSystem Declaration
61 \*---------------------------------------------------------------------------*/
62 
63 class multiphaseSystem
64 :
65  public IOdictionary,
66  public transportModel
67 {
68 
69 public:
70 
71  class interfacePair
72  :
73  public Pair<word>
74  {
75  public:
76 
77  class symmHash
78  :
79  public Hash<interfacePair>
80  {
81  public:
82 
83  symmHash()
84  {}
85 
86  label operator()(const interfacePair& key) const
87  {
88  return word::hash()(key.first()) + word::hash()(key.second());
89  }
90  };
91 
92  class hash
93  :
94  public Hash<interfacePair>
95  {
96  public:
97 
98  hash()
99  {}
101  label operator()(const interfacePair& key) const
102  {
103  return word::hash()(key.first(), word::hash()(key.second()));
104  }
105  };
106 
107 
108  // Constructors
110  interfacePair()
111  {}
113  interfacePair(const word& alpha1Name, const word& alpha2Name)
114  :
115  Pair<word>(alpha1Name, alpha2Name)
116  {}
119  :
120  Pair<word>(alpha1.name(), alpha2.name())
121  {}
122 
123 
124  // Friend Operators
125 
126  friend bool operator==
127  (
128  const interfacePair& a,
129  const interfacePair& b
130  )
131  {
132  return
133  (
134  ((a.first() == b.first()) && (a.second() == b.second()))
135  || ((a.first() == b.second()) && (a.second() == b.first()))
136  );
137  }
138 
139  friend bool operator!=
140  (
141  const interfacePair& a,
142  const interfacePair& b
143  )
144  {
145  return (!(a == b));
146  }
147  };
148 
149 
152 
155 
156 
157 private:
158 
159  // Private data
160 
161  //- Dictionary of phases
163 
164  const fvMesh& mesh_;
165  const surfaceScalarField& phi_;
166 
167  volScalarField alphas_;
168 
171 
174 
175  scalarCoeffSymmTable sigmas_;
176  dimensionSet dimSigma_;
177 
178  scalarCoeffSymmTable cAlphas_;
179 
180  scalarCoeffTable Cvms_;
181 
184 
185  dragModelTable dragModels_;
186 
187  //- Stabilisation for normalisation of the interface normal
188  const dimensionedScalar deltaN_;
189 
190  //- Conversion factor for degrees into radians
191  static const scalar convertToRad;
192 
193 
194  // Private member functions
195 
196  void calcAlphas();
197 
198  void solveAlphas();
199 
201  (
202  const volScalarField& alpha1,
203  const volScalarField& alpha2
204  ) const;
205 
207  (
208  const volScalarField& alpha1,
209  const volScalarField& alpha2
210  ) const;
211 
212  void correctContactAngle
213  (
214  const phaseModel& alpha1,
215  const phaseModel& alpha2,
216  surfaceVectorField::Boundary& nHatb
217  ) const;
218 
220  (
221  const phaseModel& alpha1,
222  const phaseModel& alpha2
223  ) const;
224 
225 
226 public:
227 
228  // Constructors
229 
230  //- Construct from components
232  (
233  const volVectorField& U,
234  const surfaceScalarField& phi
235  );
236 
237 
238  //- Destructor
239  virtual ~multiphaseSystem()
240  {}
241 
242 
243  // Member Functions
244 
245  //- Return the phases
246  const PtrDictionary<phaseModel>& phases() const
247  {
248  return phases_;
249  }
250 
251  //- Return the phases
253  {
254  return phases_;
255  }
256 
257  //- Return the mixture density
258  tmp<volScalarField> rho() const;
259 
260  //- Return the mixture density for patch
261  tmp<scalarField> rho(const label patchi) const;
262 
263  //- Return the mixture laminar viscosity
264  tmp<volScalarField> nu() const;
265 
266  //- Return the laminar viscosity for patch
267  tmp<scalarField> nu(const label patchi) const;
268 
269  //- Return the virtual-mass coefficient for the given phase
270  tmp<volScalarField> Cvm(const phaseModel& phase) const;
271 
272  //- Return the virtual-mass source for the given phase
273  tmp<volVectorField> Svm(const phaseModel& phase) const;
274 
275  //- Return the table of drag models
276  const dragModelTable& dragModels() const
277  {
278  return dragModels_;
279  }
280 
281  //- Return the drag coefficients for all of the interfaces
283 
284  //- Return the sum of the drag coefficients for the given phase
286  (
287  const phaseModel& phase,
288  const dragCoeffFields& dragCoeffs
289  ) const;
290 
292 
293  //- Indicator of the proximity of the interface
294  // Field values are 1 near and 0 away for the interface.
296 
297  //- Solve for the mixture phase-fractions
298  void solve();
299 
300  //- Dummy correct
301  void correct()
302  {}
303 
304  //- Read base transportProperties dictionary
305  bool read();
306 };
307 
308 
309 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
310 
311 } // End namespace Foam
312 
313 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
314 
315 #endif
316 
317 // ************************************************************************* //
tmp< volScalarField > nu() const
Return the mixture laminar viscosity.
Template dictionary class which manages the storage associated with it.
Definition: PtrDictionary.H:53
Foam::surfaceFields.
Hashing function class, shared by all the derived classes.
Definition: string.H:90
const PtrDictionary< phaseModel > & phases() const
Return the phases.
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
tmp< surfaceScalarField > surfaceTension(const phaseModel &phase) const
tmp< volVectorField > U() const
Return the mixture velocity.
const Type & second() const
Return second.
Definition: Pair.H:99
autoPtr< dragCoeffFields > dragCoeffs() const
Return the drag coefficients for all of the interfaces.
const Type & first() const
Return first.
Definition: Pair.H:87
tmp< volScalarField > Cvm(const phaseModel &phase) const
Return the virtual-mass coefficient for the given phase.
CGAL::Exact_predicates_exact_constructions_kernel K
HashPtrTable< dragModel, interfacePair, interfacePair::symmHash > dragModelTable
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
const surfaceScalarField & phi() const
Constant access the mixture flux.
Definition: phaseSystemI.H:55
tmp< volScalarField > dragCoeff(const phaseModel &phase, const dragCoeffFields &dragCoeffs) const
Return the sum of the drag coefficients for the given phase.
Dimension set for the base types.
Definition: dimensionSet.H:118
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
A class for handling words, derived from string.
Definition: word.H:59
const word & name() const
Name function is needed to disambiguate those inherited.
Definition: IOdictionary.C:181
HashPtrTable< volScalarField, interfacePair, interfacePair::symmHash > dragCoeffFields
alpha2
Definition: alphaEqn.H:112
void solve()
Solve for the mixture phase-fractions.
tmp< volScalarField > rho() const
Return the mixture density.
tmp< volVectorField > Svm(const phaseModel &phase) const
Return the virtual-mass source for the given phase.
bool read()
Read base transportProperties dictionary.
const dragModelTable & dragModels() const
Return the table of drag models.
volScalarField & alpha1
void correct()
Dummy correct.
label patchi
Base-class for all transport models used by the incompressible turbulence models. ...
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
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:52
label operator()(const interfacePair &key) const
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
A class for managing temporary objects.
Definition: PtrList.H:54
multiphaseSystem(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
Incompressible multi-phase mixture with built in solution for the phase fractions with interface comp...
virtual ~multiphaseSystem()
Destructor.
Namespace for OpenFOAM.