IATEsource.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) 2011-2023 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 "IATEsource.H"
27 #include "fvMatrix.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace diameterModels
36 {
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
44 
47 (
48  const word& type,
49  const IATE& iate,
50  const dictionary& dict
51 )
52 {
53  dictionaryConstructorTable::iterator cstrIter =
54  dictionaryConstructorTablePtr_->find(type);
55 
56  if (cstrIter == dictionaryConstructorTablePtr_->end())
57  {
59  << "Unknown IATE source type "
60  << type << nl << nl
61  << "Valid IATE source types : " << endl
62  << dictionaryConstructorTablePtr_->sortedToc()
63  << exit(FatalError);
64  }
65 
66  return autoPtr<IATEsource>(cstrIter()(iate, dict));
67 }
68 
69 
70 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
71 
73 {
76 
77  return
78  sqrt(2.0)
79  *pow025
80  (
81  sigma()*mag(g)
82  *(otherPhase().rho() - phase().rho())
83  /sqr(otherPhase().rho())
84  )
85  *pow(max(1 - phase(), scalar(0)), 1.75);
86 }
87 
89 {
90  return sqrt(2*otherPhase().k());
91 }
92 
94 {
95  return max(Ur()*phase().d()/otherPhase().fluidThermo().nu(), scalar(1e-3));
96 }
97 
99 {
100  const volScalarField Eo(this->Eo());
101  const volScalarField Re(this->Re());
102 
103  return
104  max
105  (
106  min
107  (
108  (16/Re)*(1 + 0.15*pow(Re, 0.687)),
109  48/Re
110  ),
111  8*Eo/(3*(Eo + 4))
112  );
113 }
114 
116 {
119 
120  return
121  mag(g)*pow4(otherPhase().fluidThermo().nu())*sqr(otherPhase().rho())
122  *(otherPhase().rho() - phase().rho())
123  /pow3(sigma());
124 }
125 
127 {
130 
131  return
132  mag(g)*sqr(phase().d())
133  *(otherPhase().rho() - phase().rho())
134  /sigma();
135 }
136 
138 {
139  return
140  otherPhase().rho()*sqr(Ut())*phase().d()/sigma();
141 }
142 
143 
144 // ************************************************************************* //
label k
Generic GeometricField class.
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:312
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
IATE (Interfacial Area Transport Equation) bubble diameter model.
Definition: IATE.H:69
IATE (Interfacial Area Transport Equation) bubble diameter model run-time selectable sources.
Definition: IATEsource.H:54
tmp< volScalarField > CD() const
Return the bubble drag coefficient.
Definition: IATEsource.C:98
tmp< volScalarField > Eo() const
Return the bubble Eotvos number.
Definition: IATEsource.C:126
tmp< volScalarField > Mo() const
Return the bubble Morton number.
Definition: IATEsource.C:115
tmp< volScalarField > We() const
Return the bubble turbulent Weber number.
Definition: IATEsource.C:137
tmp< volScalarField > Ur() const
Return the bubble relative velocity.
Definition: IATEsource.C:72
tmp< volScalarField > Re() const
Return the bubble Reynolds number.
Definition: IATEsource.C:93
static autoPtr< IATEsource > New(const word &type, const IATE &iate, const dictionary &dict)
Definition: IATEsource.C:47
tmp< volScalarField > Ut() const
Return the bubble turbulent velocity.
Definition: IATEsource.C:88
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:57
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
defineRunTimeSelectionTable(IATEsource, dictionary)
defineTypeNameAndDebug(constant, 0)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const doubleScalar e
Definition: doubleScalar.H:106
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
dimensionedScalar pow3(const dimensionedScalar &ds)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
dimensionedScalar pow4(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
error FatalError
static const char nl
Definition: Ostream.H:266
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dimensionedScalar pow025(const dimensionedScalar &ds)
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
dictionary dict