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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2014-2015 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
26\*---------------------------------------------------------------------------*/
27
28#include "phasePair.H"
29
30// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
31
32Foam::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 dimensionedVector& g,
52 const scalarTable& sigmaTable,
53 const bool ordered
54)
55:
56 phasePairKey(phase1.name(), phase2.name(), ordered),
57 phase1_(phase1),
58 phase2_(phase2),
59 g_(g),
60 sigma_
61 (
62 "sigma",
63 dimensionSet(1, 0, -2, 0, 0),
64 sigmaTable
65 [
67 (
68 phase1.name(),
69 phase2.name(),
70 false
71 )
72 ]
73 )
74{}
75
76
77// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
78
80{}
81
82
83// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
84
86{
88 << "Requested dispersed phase from an unordered pair."
89 << exit(FatalError);
90
91 return phase1_;
92}
93
94
96{
98 << "Requested continuous phase from an unordered pair."
99 << exit(FatalError);
100
101 return phase1_;
102}
103
104
106{
107 word name2(phase2().name());
108 name2[0] = toupper(name2[0]);
109 return phase1().name() + "And" + name2;
110}
111
112
114{
115 return phase1()*phase1().rho() + phase2()*phase2().rho();
116}
117
118
120{
121 return mag(phase1().U() - phase2().U());
122}
123
124
126{
127 return dispersed().U() - continuous().U();
128}
129
130
132{
133 return magUr()*dispersed().d()/continuous().nu();
134}
135
136
138{
139 return
140 continuous().nu()*continuous().Cp()*continuous().rho()
141 /continuous().kappa();
142}
143
144
146{
147 return EoH(dispersed().d());
148}
149
150
152{
153 return
154 EoH
155 (
156 dispersed().d()
157 *cbrt(1 + 0.163*pow(Eo(), 0.757))
158 );
159}
160
161
163{
164 return
165 EoH
166 (
167 dispersed().d()
168 /cbrt(E())
169 );
170}
171
172
174{
175 return
176 mag(g())
177 *continuous().nu()
178 *pow3(continuous().nu()*continuous().rho()/sigma());
179}
180
181
183{
184 return Re()*pow(Mo(), 0.23);
185}
186
187
189{
191 << "Requested aspect ratio of the dispersed phase in an unordered pair"
192 << exit(FatalError);
193
194 return phase1();
195}
196
197
198// ************************************************************************* //
const uniformDimensionedVectorField & g
phaseModel & phase1
phaseModel & phase2
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:109
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:61
const dimensionedScalar & rho() const
Definition: phaseModel.H:171
const word & name() const
Definition: phaseModel.H:144
An ordered or unorder pair of phase names. Typically specified as follows.
Definition: phasePairKey.H:68
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:56
tmp< volScalarField > EoH1() const
Eotvos number based on hydraulic diameter type 1.
Definition: phasePair.C:148
tmp< volScalarField > Eo() const
Eotvos number.
Definition: phasePair.C:142
tmp< volScalarField > Mo() const
Morton Number.
Definition: phasePair.C:181
virtual word name() const
Pair name.
Definition: phasePair.C:69
tmp< volScalarField > magUr() const
Relative velocity magnitude.
Definition: phasePair.C:114
tmp< volScalarField > EoH2() const
Eotvos number based on hydraulic diameter type 2.
Definition: phasePair.C:159
virtual ~phasePair()=default
Destructor.
Definition: phasePair.C:66
tmp< volScalarField > Pr() const
Prandtl number.
Definition: phasePair.C:132
virtual const phaseModel & continuous() const
Continuous phase.
Definition: phasePair.C:82
const uniformDimensionedVectorField & g() const
Return gravitation acceleration.
Definition: phasePairI.H:93
tmp< volScalarField > Re() const
Reynolds number.
Definition: phasePair.C:126
tmp< volScalarField > sigma() const
Surface tension coefficient.
Definition: phasePair.C:170
tmp< volVectorField > Ur() const
Relative velocity.
Definition: phasePair.C:120
tmp< volScalarField > rho() const
Average density.
Definition: phasePair.C:108
virtual const phaseModel & dispersed() const
Dispersed phase.
Definition: phasePair.C:72
virtual tmp< volScalarField > E() const
Aspect ratio.
Definition: phasePair.C:201
tmp< volScalarField > Ta() const
Takahashi Number.
Definition: phasePair.C:195
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
scalarField Re(const UList< complex > &cf)
Extract real component.
Definition: complexField.C:159
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
dimensionedScalar cbrt(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
volScalarField & nu