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-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::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"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 
46 /*---------------------------------------------------------------------------*\
47  Class phasePair Declaration
48 \*---------------------------------------------------------------------------*/
49 
50 class phasePair
51 :
52  public phasePairKey
53 {
54 public:
55 
56  // Hash table types
57 
58  //- Dictionary hash table
60  dictTable;
61 
62  //- Scalar hash table
65 
66 
67 private:
68 
69  // Private Data
70 
71  //- Phase 1
72  const phaseModel& phase1_;
73 
74  //- Phase 2
75  const phaseModel& phase2_;
76 
77  //- Gravitational acceleration
79 
80 
81  // Private Member Functions
82 
83  // Etvos number for given diameter
84  tmp<volScalarField> EoH(const volScalarField& d) const;
85 
86 
87 public:
88 
89  // Constructors
90 
91  //- Construct from two phases and gravity
92  phasePair
93  (
94  const phaseModel& phase1,
95  const phaseModel& phase2,
96  const bool ordered = false
97  );
98 
99 
100  //- Destructor
101  virtual ~phasePair();
102 
103 
104  // Member Functions
105 
106  //- Dispersed phase
107  virtual const phaseModel& dispersed() const;
108 
109  //- Continuous phase
110  virtual const phaseModel& continuous() const;
111 
112  //- Pair name
113  virtual word name() const;
114 
115  //- Other pair name
116  virtual word otherName() 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  //- Surface tension coefficient
143  tmp<volScalarField> sigma() const;
144 
145  //- Morton Number
146  tmp<volScalarField> Mo() const;
147 
148  //- Takahashi Number
149  tmp<volScalarField> Ta() const;
150 
151  //- Aspect ratio
152  virtual tmp<volScalarField> E() const;
153 
154  // Access
155 
156  //- Return phase 1
157  inline const phaseModel& phase1() const;
158 
159  //- Return phase 2
160  inline const phaseModel& phase2() const;
161 
162  //- Return true if this phasePair contains the given phase
163  inline bool contains(const phaseModel& phase) const;
164 
165  //- Return the other phase relative to the given phase
166  // Generates a FatalError if this phasePair does not contain
167  // the given phase
168  inline const phaseModel& otherPhase(const phaseModel& phase) const;
169 
170  //- Return the index of the given phase. Generates a FatalError if
171  // this phasePair does not contain the given phase
172  inline label index(const phaseModel& phase) const;
173 
174  //- Return gravitation acceleration
175  inline const uniformDimensionedVectorField& g() const;
176 
177  //- Return the mesh
178  inline const fvMesh& mesh() const;
179 
180 
181  //- STL const_iterator
182  class const_iterator
183  {
184  // Private Data
185 
186  //- Reference to the pair for which this is an iterator
187  const phasePair& pair_;
188 
189  //- Current index
190  label index_;
191 
192  //- Construct an iterator with the given index
193  inline const_iterator(const phasePair&, const label index);
194 
195  public:
197  friend class phasePair;
198 
199  // Constructors
200 
201  //- Construct from pair, moving to its 'begin' position
202  inline explicit const_iterator(const phasePair&);
203 
204 
205  // Access
206 
207  //- Return the current index
208  inline label index() const;
209 
210 
211  // Member Operators
212 
213  inline bool operator==(const const_iterator&) const;
214 
215  inline bool operator!=(const const_iterator&) const;
216 
217  inline const phaseModel& operator*() const;
218  inline const phaseModel& operator()() const;
219 
220  inline const phaseModel& otherPhase() const;
221 
222  inline const_iterator& operator++();
223  inline const_iterator operator++(int);
224  };
225 
226 
227  //- const_iterator set to the beginning of the pair
228  inline const_iterator cbegin() const;
229 
230  //- const_iterator set to beyond the end of the pair
231  inline const_iterator cend() const;
232 
233  //- const_iterator set to the beginning of the pair
234  inline const_iterator begin() const;
235 
236  //- const_iterator set to beyond the end of the pair
237  inline const_iterator end() const;
238 };
239 
240 
241 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
242 
243 } // End namespace Foam
244 
245 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
246 
247 #include "phasePairI.H"
248 
249 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
250 
251 #endif
252 
253 // ************************************************************************* //
STL const_iterator.
Definition: phasePair.H:181
const phaseModel & operator*() const
Definition: phasePairI.H:147
HashTable< dictionary, phasePairKey, phasePairKey::hash > dictTable
Dictionary hash table.
Definition: phasePair.H:59
bool ordered() const
Return the ordered flag.
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 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
bool operator==(const const_iterator &) const
Definition: phasePairI.H:129
tmp< volVectorField > Ur() const
Relative velocity.
virtual word name() const
Pair name.
tmp< volScalarField > magUr() const
Relative velocity magnitude.
const phaseModel & operator()() const
Definition: phasePairI.H:161
bool operator!=(const const_iterator &) const
Definition: phasePairI.H:138
tmp< volScalarField > Pr() const
Prandtl number.
tmp< volScalarField > sigma() const
Surface tension coefficient.
label index(const phaseModel &phase) const
Return the index of the given phase. Generates a FatalError if.
Definition: phasePairI.H:70
tmp< volScalarField > Eo() const
Eotvos number.
const_iterator cbegin() const
const_iterator set to the beginning of the pair
Definition: phasePairI.H:198
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_iterator & operator++()
Definition: phasePairI.H:182
tmp< volScalarField > Re() const
Reynolds number.
const uniformDimensionedVectorField & g() const
Return gravitation acceleration.
Definition: phasePairI.H:91
An STL-conforming hash table.
Definition: HashTable.H:61
virtual const phaseModel & dispersed() const
Dispersed phase.
const_iterator cend() const
const_iterator set to beyond the end of the pair
Definition: phasePairI.H:204
tmp< volScalarField > Mo() const
Morton Number.
label index() const
Return the current index.
Definition: phasePairI.H:122
const_iterator end() const
const_iterator set to beyond the end of the pair
Definition: phasePairI.H:216
const fvMesh & mesh() const
Return the mesh.
Definition: phasePairI.H:97
const phaseModel & otherPhase(const phaseModel &phase) const
Return the other phase relative to the given phase.
Definition: phasePairI.H:47
virtual tmp< volScalarField > E() const
Aspect ratio.
const phaseModel & otherPhase() const
Definition: phasePairI.H:168
const_iterator begin() const
const_iterator set to the beginning of the pair
Definition: phasePairI.H:210
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:52
virtual ~phasePair()
Destructor.
bool contains(const phaseModel &phase) const
Return true if this phasePair contains the given phase.
Definition: phasePairI.H:40
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.
virtual word otherName() const
Other pair name.
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.