phasePair.C
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) 2014-2015 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 \*---------------------------------------------------------------------------*/
25 
26 #include "phasePair.H"
27 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 Foam::tmp<Foam::volScalarField> Foam::phasePair::EoH
31 (
32  const volScalarField& d
33 ) const
34 {
35  return
36  mag(dispersed().rho() - continuous().rho())
37  *mag(g())
38  *sqr(d)
39  /sigma();
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
46 (
47  const phaseModel& phase1,
48  const phaseModel& phase2,
49  const dimensionedVector& g,
50  const scalarTable& sigmaTable,
51  const bool ordered
52 )
53 :
54  phasePairKey(phase1.name(), phase2.name(), ordered),
55  phase1_(phase1),
56  phase2_(phase2),
57  g_(g),
58  sigma_
59  (
60  "sigma",
61  dimensionSet(1, 0, -2, 0, 0),
62  sigmaTable
63  [
65  (
66  phase1.name(),
67  phase2.name(),
68  false
69  )
70  ]
71  )
72 {}
73 
74 
75 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
76 
78 {}
79 
80 
81 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
82 
84 {
86  << "Requested dispersed phase from an unordered pair."
87  << exit(FatalError);
88 
89  return phase1_;
90 }
91 
92 
94 {
96  << "Requested continuous phase from an unordered pair."
97  << exit(FatalError);
98 
99  return phase1_;
100 }
101 
102 
104 {
105  word name2(phase2().name());
106  name2[0] = toupper(name2[0]);
107  return phase1().name() + "And" + name2;
108 }
109 
110 
112 {
113  return phase1()*phase1().rho() + phase2()*phase2().rho();
114 }
115 
116 
118 {
119  return mag(phase1().U() - phase2().U());
120 }
121 
122 
124 {
125  return dispersed().U() - continuous().U();
126 }
127 
128 
130 {
131  return magUr()*dispersed().d()/continuous().nu();
132 }
133 
134 
136 {
137  return
139  /continuous().kappa();
140 }
141 
142 
144 {
145  return EoH(dispersed().d());
146 }
147 
148 
150 {
151  return
152  EoH
153  (
154  dispersed().d()
155  *cbrt(1 + 0.163*pow(Eo(), 0.757))
156  );
157 }
158 
159 
161 {
162  return
163  EoH
164  (
165  dispersed().d()
166  /cbrt(E())
167  );
168 }
169 
170 
172 {
173  return
174  mag(g())
175  *continuous().nu()
176  *pow3(continuous().nu()*continuous().rho()/sigma());
177 }
178 
179 
181 {
182  return Re()*pow(Mo(), 0.23);
183 }
184 
185 
187 {
189  << "Requested aspect ratio of the dispersed phase in an unordered pair"
190  << exit(FatalError);
191 
192  return phase1();
193 }
194 
195 
196 // ************************************************************************* //
const uniformDimensionedVectorField & g() const
Definition: phasePairI.H:40
U
Definition: pEqn.H:83
tmp< volScalarField > Mo() const
Morton Number.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
tmp< volScalarField > EoH2() const
Eotvos number based on hydraulic diameter type 2.
tmp< volVectorField > Ur() const
Relative velocity.
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
phasePair(const phaseModel &phase1, const phaseModel &phase2, const bool ordered=false)
Construct from two phases and gravity.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual const phaseModel & continuous() const
Continuous phase.
tmp< volScalarField > Pr() const
Prandtl number.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
tmp< volScalarField > sigma() const
Surface tension coefficient.
Definition: phasePairI.H:46
virtual const phaseModel & continuous() const
Continuous phase.
virtual word name() const
Pair name.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
A class for handling words, derived from string.
Definition: word.H:59
dimensionedScalar cbrt(const dimensionedScalar &ds)
const dimensionedScalar & rho() const
Definition: phaseModel.H:165
tmp< volScalarField > EoH1() const
Eotvos number based on hydraulic diameter type 1.
const dimensionedScalar & nu() const
Return the laminar viscosity.
Definition: phaseModel.H:150
bool ordered() const
Return the ordered flag.
tmp< volScalarField > Re() const
Reynolds number.
const volVectorField & U() const
Definition: phaseModel.H:170
virtual const phaseModel & dispersed() const
Dispersed phase.
virtual const phaseModel & dispersed() const
Dispersed phase.
virtual word name() const
Pair name.
tmp< volScalarField > d() const
const phaseModel & phase1() const
Definition: phasePairI.H:28
const word & name() const
Definition: phaseModel.H:109
tmp< volScalarField > Ta() const
Takahashi Number.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
virtual tmp< volScalarField > E() const
Aspect ratio.
tmp< volScalarField > rho() const
Average density.
phasePairKey()
Construct null.
const dimensionedScalar & Cp() const
Definition: phaseModel.H:160
virtual ~phasePair()
Destructor.
dimensioned< scalar > mag(const dimensioned< Type > &)
tmp< volScalarField > Eo() const
Eotvos number.
HashTable< scalar, phasePairKey, phasePairKey::hash > scalarTable
Scalar hash table.
Definition: phasePair.H:63
virtual tmp< volScalarField > E() const
Aspect ratio.
A class for managing temporary objects.
Definition: PtrList.H:54
tmp< volScalarField > magUr() const
Relative velocity magnitude.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
volScalarField & nu
const phaseModel & phase2() const
Definition: phasePairI.H:34
const dimensionedScalar & kappa() const
Definition: phaseModel.H:155