Loading...
Searching...
No Matches
phasePair.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2014-2018 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Class
27 Foam::phasePair
28
29Description
30
31SourceFiles
32 phasePair.C
33
34\*---------------------------------------------------------------------------*/
35
36#ifndef phasePair_H
37#define phasePair_H
38
39#include "phaseModel.H"
40#include "phasePairKey.H"
42
43// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44
45namespace Foam
46{
47
48/*---------------------------------------------------------------------------*\
49 Class phasePair Declaration
50\*---------------------------------------------------------------------------*/
51
52class phasePair
53:
54 public phasePairKey
55{
56public:
57
58 // Hash table types
60 //- Dictionary hash table
63
64 //- Scalar hash table
67
68
69private:
70
71 // Private data
72
73 //- Phase 1
74 const phaseModel& phase1_;
75
76 //- Phase 2
77 const phaseModel& phase2_;
78
79 //- Gravitational acceleration
81
82
83 // Private member functions
84
85 // Etvos number for given diameter
86 tmp<volScalarField> EoH(const volScalarField& d) const;
87
88
89public:
90
91 // Constructors
92
93 //- Construct from two phases and gravity
95 (
96 const phaseModel& phase1,
97 const phaseModel& phase2,
98 const bool ordered = false
99 );
100
101
102 //- Destructor
103 virtual ~phasePair();
104
105
106 // Member Functions
107
108 //- Dispersed phase
109 virtual const phaseModel& dispersed() const;
110
111 //- Continuous phase
112 virtual const phaseModel& continuous() const;
113
114 //- Pair name
115 virtual word name() const;
116
117 //- Other pair name
118 virtual word otherName() const;
119
120 //- Average density
121 tmp<volScalarField> rho() const;
122
123 //- Relative velocity magnitude
125
126 //- Relative velocity
127 tmp<volVectorField> Ur() const;
128
129 //- Reynolds number
131
132 //- Prandtl number
133 tmp<volScalarField> Pr() const;
134
135 //- Eotvos number
136 tmp<volScalarField> Eo() const;
137
138 //- Eotvos number based on hydraulic diameter type 1
140
141 //- Eotvos number based on hydraulic diameter type 2
143
144 //- Surface tension coefficient
146
147 //- Morton Number
148 tmp<volScalarField> Mo() const;
149
150 //- Takahashi Number
151 tmp<volScalarField> Ta() const;
152
153 //- Aspect ratio
154 virtual tmp<volScalarField> E() const;
155
156 // Access
157
158 //- Return phase 1
159 inline const phaseModel& phase1() const;
160
161 //- Return phase 2
162 inline const phaseModel& phase2() const;
163
164 //- Return true if this phasePair contains the given phase
165 inline bool contains(const phaseModel& phase) const;
166
167 //- Return the other phase relative to the given phase
168 // Generates a FatalError if this phasePair does not contain
169 // the given phase
170 inline const phaseModel& otherPhase(const phaseModel& phase) const;
171
172 //- Return the index of the given phase. Generates a FatalError if
173 // this phasePair does not contain the given phase
174 inline label index(const phaseModel& phase) const;
175
176 //- Return gravitation acceleration
177 inline const uniformDimensionedVectorField& g() const;
178
179
180 //- STL const_iterator
181 class const_iterator
182 {
183 // Private data
184
185 //- Reference to the pair for which this is an iterator
186 const phasePair& pair_;
187
188 //- Current index
189 label index_;
190
191 //- Construct an iterator with the given index
192 inline const_iterator(const phasePair&, const label index);
193
194 public:
195
196 friend class phasePair;
197
198 // Constructors
199
200 //- Construct from pair, moving to its 'begin' position
201 inline explicit const_iterator(const phasePair&);
203
204 // Access
205
206 //- Return the current index
207 inline label index() const;
208
209
210 // Member operators
211
212 inline bool operator==(const const_iterator&) const;
213
214 inline bool operator!=(const const_iterator&) const;
215
216 inline const phaseModel& operator*() const;
217 inline const phaseModel& operator()() const;
218
219 inline const phaseModel& otherPhase() const;
220
221 inline const_iterator& operator++();
222 inline const_iterator operator++(int);
223 };
224
225
226 //- const_iterator set to the beginning of the pair
227 inline const_iterator cbegin() const;
228
229 //- const_iterator set to beyond the end of the pair
230 inline const_iterator cend() const;
231
232 //- const_iterator set to the beginning of the pair
233 inline const_iterator begin() const;
234
235 //- const_iterator set to beyond the end of the pair
236 inline const_iterator end() const;
237};
239
240// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
241
242} // End namespace Foam
243
244// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
245
246#include "phasePairI.H"
247
248// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
249
250#endif
251
252// ************************************************************************* //
const T * const_iterator
Random access iterator for traversing FixedList.
Definition FixedList.H:135
A HashTable similar to std::unordered_map.
Definition HashTable.H:124
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phaseModel.H:56
An ordered or unorder pair of phase names. Typically specified as follows.
bool ordered() const noexcept
Return the ordered flag.
const phaseModel & operator*() const
Definition phasePairI.H:136
bool operator!=(const const_iterator &) const
Definition phasePairI.H:127
const_iterator & operator++()
Definition phasePairI.H:171
const phaseModel & otherPhase() const
Definition phasePairI.H:157
const phaseModel & operator()() const
Definition phasePairI.H:150
bool operator==(const const_iterator &) const
Definition phasePairI.H:118
label index() const
Return the current index.
Definition phasePairI.H:111
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition phasePair.H:52
const_iterator cbegin() const
const_iterator set to the beginning of the pair
Definition phasePairI.H:187
tmp< volScalarField > EoH1() const
Eotvos number based on hydraulic diameter type 1.
Definition phasePair.C:141
virtual word otherName() const
Other pair name.
Definition phasePair.C:93
const phaseModel & phase2() const
Return phase 2.
tmp< volScalarField > Eo() const
Eotvos number.
Definition phasePair.C:135
const phaseModel & otherPhase(const phaseModel &phase) const
Return the other phase relative to the given phase.
Definition phasePairI.H:42
tmp< volScalarField > Mo() const
Morton Number.
Definition phasePair.C:174
virtual word name() const
Pair name.
Definition phasePair.C:62
tmp< volScalarField > magUr() const
Relative velocity magnitude.
Definition phasePair.C:107
tmp< volScalarField > EoH2() const
Eotvos number based on hydraulic diameter type 2.
Definition phasePair.C:152
virtual word name() const
Pair name.
virtual ~phasePair()=default
Destructor.
Definition phasePair.C:59
tmp< volScalarField > Pr() const
Prandtl number.
Definition phasePair.C:125
virtual const phaseModel & continuous() const
Continuous phase.
Definition phasePair.C:75
const_iterator begin() const
const_iterator set to the beginning of the pair
Definition phasePairI.H:199
label index(const phaseModel &phase) const
Return the index of the given phase. Generates a FatalError if.
Definition phasePairI.H:65
const uniformDimensionedVectorField & g() const
Return gravitation acceleration.
Definition phasePairI.H:86
tmp< volScalarField > Re() const
Reynolds number.
Definition phasePair.C:119
const phaseModel & phase1() const
Return phase 1.
HashTable< scalar, phasePairKey, phasePairKey::hash > scalarTable
Scalar hash table.
Definition phasePair.H:65
bool contains(const phaseModel &phase) const
Return true if this phasePair contains the given phase.
Definition phasePairI.H:35
virtual ~phasePair()
Destructor.
tmp< volScalarField > sigma() const
Surface tension coefficient.
Definition phasePair.C:163
tmp< volVectorField > Ur() const
Relative velocity.
Definition phasePair.C:113
tmp< volScalarField > rho() const
Average density.
Definition phasePair.C:101
virtual const phaseModel & dispersed() const
Dispersed phase.
Definition phasePair.C:65
const_iterator end() const
const_iterator set to beyond the end of the pair
Definition phasePairI.H:205
phasePair(const multiphaseInter::phaseModel &phase1, const multiphaseInter::phaseModel &phase2, const bool ordered=false)
Construct from two phases.
Definition phasePair.C:28
virtual tmp< volScalarField > E() const
Aspect ratio.
Definition phasePair.C:194
const_iterator cend() const
const_iterator set to beyond the end of the pair
Definition phasePairI.H:193
const multiphaseInter::phaseModel & phase1() const
Definition phasePairI.H:23
const multiphaseInter::phaseModel & phase2() const
Definition phasePairI.H:29
tmp< volScalarField > Ta() const
Takahashi Number.
Definition phasePair.C:188
HashTable< dictionary, phasePairKey, phasePairKey::hash > dictTable
Dictionary hash table.
Definition phasePair.H:59
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phase.H:53
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
UniformDimensionedField< vector > uniformDimensionedVectorField
Various UniformDimensionedField types.