Loading...
Searching...
No Matches
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-2018 OpenFOAM Foundation
9 Copyright (C) 2018-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "phasePair.H"
30#include "phaseSystem.H"
31#include "surfaceTensionModel.H"
32
33// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34
35Foam::tmp<Foam::volScalarField> Foam::phasePair::EoH
36(
37 const volScalarField& d
38) const
39{
40 return
41 mag(dispersed().rho() - continuous().rho())
42 *mag(g())
43 *sqr(d)
44 /sigma();
45}
46
47
48// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49
51(
52 const phaseModel& phase1,
53 const phaseModel& phase2,
54 const bool ordered
55)
56:
57 phasePairKey(phase1.name(), phase2.name(), ordered),
58 phase1_(phase1),
59 phase2_(phase2),
60 g_(phase1.mesh().time().lookupObject<uniformDimensionedVectorField>("g"))
61{}
62
63
64// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
67{}
68
69
70// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
71
73{
75 << "Requested dispersed phase from an unordered pair."
76 << exit(FatalError);
77
78 return phase1_;
79}
80
81
83{
85 << "Requested continuous phase from an unordered pair."
86 << exit(FatalError);
87
88 return phase1_;
89}
90
91
94 word name2(second());
95 name2[0] = toupper(name2[0]);
96 return first() + "And" + name2;
97}
98
99
102 word name1(first());
103 name1[0] = toupper(name1[0]);
104 return second() + "And" + name1;
105}
106
111}
112
117}
118
123}
124
133{
134 return
136 *continuous().thermo().Cpv()
137 *continuous().rho()
138 /continuous().kappa();
139}
140
150 return
151 EoH
153 dispersed().d()
154 *cbrt(1 + 0.163*pow(Eo(), 0.757))
155 );
156}
157
158
160{
161 return
162 EoH
164 dispersed().d()
165 /cbrt(E())
166 );
167}
168
169
171{
172 return
173 phase1().fluid().
174 lookupSubModel<reactingMultiphaseEuler::surfaceTensionModel>
175 (
177 ).sigma();
178}
179
180
182{
183 return
184 mag(g())
185 *continuous().nu()
186 *pow3
187 (
189 *continuous().rho()
190 /sigma()
191 );
192}
193
202{
204 << "Requested aspect ratio of the dispersed phase in an unordered pair"
205 << exit(FatalError);
206
207 return phase1();
208}
209
210
211// ************************************************************************* //
const uniformDimensionedVectorField & g
phaseModel & phase1
phaseModel & phase2
const word & first() const noexcept
Definition Pair.H:137
const word & second() const noexcept
Definition Pair.H:147
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phaseModel.H:56
const phaseSystem & fluid() const
Return the system to which this phase belongs.
Definition phaseModel.C:141
const dimensionedScalar & rho() const
Definition phaseModel.H:193
An ordered or unorder pair of phase names. Typically specified as follows.
bool ordered() const noexcept
Return the ordered flag.
phasePairKey()=default
Default construct.
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
tmp< volScalarField > Eo() const
Eotvos number.
Definition phasePair.C:135
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 ~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 uniformDimensionedVectorField & g() const
Return gravitation acceleration.
Definition phasePairI.H:86
tmp< volScalarField > Re() const
Reynolds number.
Definition phasePair.C:119
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
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 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
tmp< volScalarField > sigma(const phasePairKey &key) const
Return the surface tension coefficient for a pair.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
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)
UniformDimensionedField< vector > uniformDimensionedVectorField
scalarField Re(const UList< complex > &cmplx)
Extract real component.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensionedScalar cbrt(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
volScalarField & nu
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)