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 #include "surfaceTensionModel.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 Foam::tmp<Foam::volScalarField> Foam::phasePair::EoH
32 (
33  const volScalarField& d
34 ) const
35 {
36  return
37  mag(dispersed().rho() - continuous().rho())
38  *mag(g())
39  *sqr(d)
40  /sigma();
41 }
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
47 (
48  const phaseModel& phase1,
49  const phaseModel& phase2,
50  const bool ordered
51 )
52 :
53  phasePairKey(phase1.name(), phase2.name(), ordered),
54  phase1_(phase1),
55  phase2_(phase2),
56  g_(phase1.mesh().lookupObject<uniformDimensionedVectorField>("g"))
57 {}
58 
59 
60 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
61 
63 {}
64 
65 
66 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
67 
69 {
71  << "Requested dispersed phase from an unordered pair."
72  << exit(FatalError);
73 
74  return phase1_;
75 }
76 
77 
79 {
81  << "Requested continuous phase from an unordered pair."
82  << exit(FatalError);
83 
84  return phase1_;
85 }
86 
87 
89 {
90  word name2(second());
91  name2[0] = toupper(name2[0]);
92  return first() + "And" + name2;
93 }
94 
95 
97 {
98  return phase1()*phase1().rho() + phase2()*phase2().rho();
99 }
100 
101 
103 {
104  return mag(phase1().U() - phase2().U());
105 }
106 
107 
109 {
110  return dispersed().U() - continuous().U();
111 }
112 
113 
115 {
116  return magUr()*dispersed().d()/continuous().nu();
117 }
118 
119 
121 {
122  return
123  continuous().nu()
124  *continuous().thermo().Cpv()
125  *continuous().rho()
126  /continuous().kappa();
127 }
128 
129 
131 {
132  return EoH(dispersed().d());
133 }
134 
135 
137 {
138  return
139  EoH
140  (
141  dispersed().d()
142  *cbrt(1 + 0.163*pow(Eo(), 0.757))
143  );
144 }
145 
146 
148 {
149  return
150  EoH
151  (
152  dispersed().d()
153  /cbrt(E())
154  );
155 }
156 
157 
159 {
160  return
161  phase1().mesh().lookupObject<surfaceTensionModel>
162  (
164  (
165  surfaceTensionModel::typeName,
167  )
168  ).sigma();
169 }
170 
171 
173 {
174  return
175  mag(g())
176  *continuous().nu()
177  *pow3
178  (
179  continuous().nu()
180  *continuous().rho()
181  /sigma()
182  );
183 }
184 
185 
187 {
188  return Re()*pow(Mo(), 0.23);
189 }
190 
191 
193 {
195  << "Requested aspect ratio of the dispersed phase in an unordered pair"
196  << exit(FatalError);
197 
198  return phase1();
199 }
200 
201 
202 // ************************************************************************* //
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.
UniformDimensionedField< vector > uniformDimensionedVectorField
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.
const word & second() const
Return second.
Definition: Pair.H:99
tmp< volScalarField > sigma() const
Surface tension coefficient.
Definition: phasePairI.H:46
const word & first() const
Return first.
Definition: Pair.H:87
virtual const phaseModel & continuous() const
Continuous phase.
virtual tmp< volScalarField > Cpv() const =0
Heat capacity at constant pressure/volume [J/kg/K].
virtual word name() const
Pair name.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
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.
tmp< volScalarField > d() const
const phaseModel & phase1() const
Definition: phasePairI.H:28
tmp< volScalarField > Ta() const
Takahashi Number.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
const rhoThermo & thermo() const
Return const-access to phase rhoThermo.
Definition: phaseModel.H:120
virtual tmp< volScalarField > E() const
Aspect ratio.
const Mesh & mesh() const
Return mesh.
tmp< volScalarField > rho() const
Average density.
phasePairKey()
Construct null.
virtual ~phasePair()
Destructor.
dimensioned< scalar > mag(const dimensioned< Type > &)
tmp< volScalarField > Eo() const
Eotvos number.
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