phasePair.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) 2014-2019 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::phasePair
26 
27 Description
28 
29 SourceFiles
30  phasePair.C
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #ifndef phasePair_H
35 #define phasePair_H
36 
37 #include "phaseModel.H"
38 #include "phasePairKey.H"
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 
45 /*---------------------------------------------------------------------------*\
46  Class phasePair Declaration
47 \*---------------------------------------------------------------------------*/
48 
49 class phasePair
50 :
51  public phasePairKey
52 {
53 public:
54 
55  // Hash table types
56 
57  //- Dictionary hash table
59  dictTable;
60 
61  //- Scalar hash table
64 
65 
66 private:
67 
68  // Private Data
69 
70  //- Phase 1
71  const phaseModel& phase1_;
72 
73  //- Phase 2
74  const phaseModel& phase2_;
75 
76  //- Gravitational acceleration
77  const dimensionedVector& g_;
78 
79  //- Surface tension coefficient
80  const dimensionedScalar sigma_;
81 
82 
83  // Private Member Functions
84 
85  // Etvos number for given diameter
86  tmp<volScalarField> EoH(const volScalarField& d) const;
87 
88 public:
89 
90  // Constructors
91 
92  //- Construct from two phases, gravity and surface tension table
93  phasePair
94  (
95  const phaseModel& phase1,
96  const phaseModel& phase2,
97  const dimensionedVector& g,
98  const scalarTable& sigmaTable,
99  const bool ordered = false
100  );
101 
102 
103  //- Destructor
104  virtual ~phasePair();
105 
106 
107  // Member Functions
108 
109  //- Dispersed phase
110  virtual const phaseModel& dispersed() const;
111 
112  //- Continuous phase
113  virtual const phaseModel& continuous() const;
114 
115  //- Pair name
116  virtual word name() const;
117 
118  //- Average density
119  tmp<volScalarField> rho() const;
120 
121  //- Relative velocity magnitude
122  tmp<volScalarField> magUr() const;
123 
124  //- Relative velocity
125  tmp<volVectorField> Ur() const;
126 
127  //- Reynolds number
128  tmp<volScalarField> Re() const;
129 
130  //- Prandtl number
131  tmp<volScalarField> Pr() const;
132 
133  //- Eotvos number
134  tmp<volScalarField> Eo() const;
135 
136  //- Eotvos number based on hydraulic diameter type 1
137  tmp<volScalarField> EoH1() const;
138 
139  //- Eotvos number based on hydraulic diameter type 2
140  tmp<volScalarField> EoH2() const;
141 
142  //- Morton Number
143  tmp<volScalarField> Mo() const;
144 
145  //- Takahashi Number
146  tmp<volScalarField> Ta() const;
147 
148  //- Aspect ratio
149  virtual tmp<volScalarField> E() const;
150 
151  // Access
152 
153  // Phase 1
154  inline const phaseModel& phase1() const;
155 
156  // Phase 2
157  inline const phaseModel& phase2() const;
158 
159  // Gravitational acceleration
160  inline const dimensionedVector& g() const;
161 
162  // Surface tension coefficient
163  inline const dimensionedScalar& sigma() const;
164 };
165 
166 
167 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
168 
169 } // End namespace Foam
170 
171 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172 
173 #include "phasePairI.H"
174 
175 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
176 
177 #endif
178 
179 // ************************************************************************* //
HashTable< dictionary, phasePairKey, phasePairKey::hash > dictTable
Dictionary hash table.
Definition: phasePair.H:59
bool ordered() const
Return the ordered flag.
tmp< volScalarField > sigma() const
Surface tension coefficient.
Definition: phasePairI.H:46
const phaseModel & phase2() const
Return phase 2.
Definition: phasePairI.H:34
phasePair(const phaseModel &phase1, const phaseModel &phase2, const bool ordered=false)
Construct from two phases and gravity.
const phaseModel & phase1() const
Return phase 1.
Definition: phasePairI.H:28
tmp< volVectorField > Ur() const
Relative velocity.
virtual word name() const
Pair name.
tmp< volScalarField > magUr() const
Relative velocity magnitude.
Generic dimensioned Type class.
tmp< volScalarField > Pr() const
Prandtl number.
tmp< volScalarField > Eo() const
Eotvos number.
tmp< volScalarField > EoH1() const
Eotvos number based on hydraulic diameter type 1.
virtual const phaseModel & continuous() const
Continuous phase.
A class for handling words, derived from string.
Definition: word.H:59
const uniformDimensionedVectorField & g() const
Return gravitation acceleration.
Definition: phasePairI.H:91
tmp< volScalarField > Re() const
Reynolds number.
An STL-conforming hash table.
Definition: HashTable.H:61
virtual const phaseModel & dispersed() const
Dispersed phase.
tmp< volScalarField > Mo() const
Morton Number.
virtual tmp< volScalarField > E() const
Aspect ratio.
virtual ~phasePair()
Destructor.
tmp< volScalarField > Ta() const
Takahashi Number.
HashTable< scalar, phasePairKey, phasePairKey::hash > scalarTable
Scalar hash table.
Definition: phasePair.H:63
A class for managing temporary objects.
Definition: PtrList.H:53
tmp< volScalarField > EoH2() const
Eotvos number based on hydraulic diameter type 2.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
Namespace for OpenFOAM.
tmp< volScalarField > rho() const
Average density.