Loading...
Searching...
No Matches
IATEsource.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) 2011-2018 OpenFOAM Foundation
9 Copyright (C) 2018-2021 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 "IATEsource.H"
30#include "fvMatrix.H"
31#include "phaseCompressibleTurbulenceModel.H"
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace diameterModels
39{
43}
44
45
46// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
47
50(
51 const word& type,
52 const IATE& iate,
53 const dictionary& dict
54)
55{
56 auto* ctorPtr = dictionaryConstructorTable(type);
57
58 if (!ctorPtr)
59 {
61 (
62 dict,
63 "IATEsource",
64 type,
65 *dictionaryConstructorTablePtr_
66 ) << exit(FatalIOError);
67 }
69 return autoPtr<IATEsource>(ctorPtr(iate, dict));
70}
71
72
73// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74
76{
78 phase().mesh().time().lookupObject<uniformDimensionedVectorField>("g");
79
80 return
81 sqrt(2.0)
82 *pow025
83 (
85 *(otherPhase().rho() - phase().rho())
86 /sqr(otherPhase().rho())
87 )
88 *pow(max(1 - phase(), scalar(0)), 1.75);
93 return sqrt(2*otherPhase().k());
95
97{
98 return max(Ur()*phase().d()/otherPhase().nu(), scalar(1e-3));
99}
100
102{
103 const volScalarField Eo(this->Eo());
104 const volScalarField Re(this->Re());
105
106 return
107 max
108 (
109 min
110 (
111 (16/Re)*(1 + 0.15*pow(Re, 0.687)),
112 48/Re
113 ),
114 8*Eo/(3*(Eo + 4))
115 );
116}
117
119{
123 return
125 *(otherPhase().rho() - phase().rho())
126 /pow3(fluid().sigma());
127}
128
130{
134 return
135 mag(g)*sqr(phase().d())
136 *(otherPhase().rho() - phase().rho())
137 /fluid().sigma();
138}
139
141{
142 return
143 otherPhase().rho()*sqr(Ur())*phase().d()/fluid().sigma();
144}
145
146
147// ************************************************************************* //
label k
const uniformDimensionedVectorField & g
twoPhaseSystem & fluid
const Mesh & mesh() const noexcept
Return const reference to mesh.
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition IOobject.C:450
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
IATE (Interfacial Area Transport Equation) bubble diameter model.
Definition IATE.H:67
IATE (Interfacial Area Transport Equation) bubble diameter model run-time selectable sources.
Definition IATEsource.H:53
tmp< volScalarField > CD() const
Return the bubble drag coefficient.
Definition IATEsource.C:94
tmp< volScalarField > Eo() const
Return the bubble Eotvos number.
Definition IATEsource.C:122
tmp< volScalarField > Mo() const
Return the bubble Morton number.
Definition IATEsource.C:111
const phaseModel & phase() const
Definition IATEsource.H:145
tmp< volScalarField > We() const
Return the bubble Webber number.
Definition IATEsource.C:133
const phaseModel & otherPhase() const
Definition IATEsource.H:155
tmp< volScalarField > Ur() const
Return the bubble relative velocity.
Definition IATEsource.C:68
tmp< volScalarField > Re() const
Return the bubble Reynolds number.
Definition IATEsource.C:89
static autoPtr< IATEsource > New(const word &type, const IATE &iate, const dictionary &dict)
Definition IATEsource.C:43
tmp< volScalarField > Ut() const
Return the bubble turbulent velocity.
Definition IATEsource.C:84
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const Time & time() const noexcept
Return time registry.
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phase.H:53
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition error.H:637
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
UniformDimensionedField< vector > uniformDimensionedVectorField
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
dimensionedScalar pow4(const dimensionedScalar &ds)
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
scalarField Re(const UList< complex > &cmplx)
Extract real component.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensionedScalar pow025(const dimensionedScalar &ds)
volScalarField & nu
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
volScalarField & e
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
Various UniformDimensionedField types.