phasePair.C
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-2018 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 #include "phaseSystem.H"
28 #include "surfaceTensionModel.H"
29 
30 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
31 
32 Foam::tmp<Foam::volScalarField> Foam::phasePair::EoH
33 (
34  const volScalarField& d
35 ) const
36 {
37  return
38  mag(dispersed().rho() - continuous().rho())
39  *mag(g())
40  *sqr(d)
41  /sigma();
42 }
43 
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
48 (
49  const phaseModel& phase1,
50  const phaseModel& phase2,
51  const bool ordered
52 )
53 :
54  phasePairKey(phase1.name(), phase2.name(), ordered),
55  phase1_(phase1),
56  phase2_(phase2),
57  g_(phase1.mesh().lookupObject<uniformDimensionedVectorField>("g"))
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
62 
64 {}
65 
66 
67 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
68 
70 {
72  << "Requested dispersed phase from an unordered pair."
73  << exit(FatalError);
74 
75  return phase1_;
76 }
77 
78 
80 {
82  << "Requested continuous phase from an unordered pair."
83  << exit(FatalError);
84 
85  return phase1_;
86 }
87 
88 
90 {
91  word name2(second());
92  name2[0] = toupper(name2[0]);
93  return first() + "And" + name2;
94 }
95 
96 
98 {
99  word name1(first());
100  name1[0] = toupper(name1[0]);
101  return second() + "And" + name1;
102 }
103 
104 
106 {
107  return phase1()*phase1().rho() + phase2()*phase2().rho();
108 }
109 
110 
112 {
113  return mag(phase1().U() - phase2().U());
114 }
115 
116 
118 {
119  return dispersed().U() - continuous().U();
120 }
121 
122 
124 {
125  return magUr()*dispersed().d()/continuous().nu();
126 }
127 
128 
130 {
131  return
132  continuous().nu()
133  *continuous().thermo().Cpv()
134  *continuous().rho()
135  /continuous().kappa();
136 }
137 
138 
140 {
141  return EoH(dispersed().d());
142 }
143 
144 
146 {
147  return
148  EoH
149  (
150  dispersed().d()
151  *cbrt(1 + 0.163*pow(Eo(), 0.757))
152  );
153 }
154 
155 
157 {
158  return
159  EoH
160  (
161  dispersed().d()
162  /cbrt(E())
163  );
164 }
165 
166 
168 {
169  return
170  phase1().fluid().lookupSubModel<surfaceTensionModel>
171  (
172  phasePair(phase1(), phase2())
173  ).sigma();
174 }
175 
176 
178 {
179  return
180  mag(g())
181  *continuous().nu()
182  *pow3
183  (
184  continuous().nu()
185  *continuous().rho()
186  /sigma()
187  );
188 }
189 
190 
192 {
193  return Re()*pow(Mo(), 0.23);
194 }
195 
196 
198 {
200  << "Requested aspect ratio of the dispersed phase in an unordered pair"
201  << exit(FatalError);
202 
203  return phase1();
204 }
205 
206 
207 // ************************************************************************* //
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
UniformDimensionedField< vector > uniformDimensionedVectorField
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.
const phaseModel & phase1() const
Return phase 1.
Definition: phasePairI.H:28
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const modelType & lookupSubModel(const phasePair &key) const
Return a sub model between a phase pair.
tmp< volVectorField > Ur() const
Relative velocity.
virtual word name() const
Pair name.
tmp< volScalarField > magUr() const
Relative velocity magnitude.
tmp< volScalarField > d() const
tmp< volScalarField > Pr() const
Prandtl number.
virtual tmp< volScalarField > Cpv() const =0
Heat capacity at constant pressure/volume [J/kg/K].
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
tmp< volScalarField > Eo() const
Eotvos number.
const volVectorField & U() const
Definition: phaseModel.H:170
tmp< volScalarField > EoH1() const
Eotvos number based on hydraulic diameter type 1.
dynamicFvMesh & mesh
virtual const phaseModel & continuous() const
Continuous phase.
A class for handling words, derived from string.
Definition: word.H:59
const rhoThermo & thermo() const
Return const-access to phase rhoThermo.
Definition: phaseModel.H:120
dimensionedScalar cbrt(const dimensionedScalar &ds)
const uniformDimensionedVectorField & g() const
Return gravitation acceleration.
Definition: phasePairI.H:91
tmp< volScalarField > Re() const
Reynolds number.
virtual const phaseModel & dispersed() const
Dispersed phase.
const dimensionedScalar & kappa() const
Definition: phaseModel.H:155
const dimensionedScalar & nu() const
Return the laminar viscosity.
Definition: phaseModel.H:150
tmp< volScalarField > Mo() const
Morton Number.
virtual const phaseModel & continuous() const
Continuous phase.
virtual tmp< volScalarField > E() const
Aspect ratio.
const word & second() const
Return second.
Definition: Pair.H:99
virtual tmp< volScalarField > E() const
Aspect ratio.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
U
Definition: pEqn.H:72
const dimensionedScalar & rho() const
Definition: phaseModel.H:165
const phaseSystem & fluid() const
Return the system to which this phase belongs.
phasePairKey()
Construct null.
virtual ~phasePair()
Destructor.
tmp< volScalarField > Ta() const
Takahashi Number.
dimensioned< scalar > mag(const dimensioned< Type > &)
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
const word & first() const
Return first.
Definition: Pair.H:87
volScalarField & nu
virtual const phaseModel & dispersed() const
Dispersed phase.
tmp< volScalarField > rho() const
Average density.