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 -------------------------------------------------------------------------------
11 License
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"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace diameterModels
39 {
40  defineTypeNameAndDebug(IATEsource, 0);
41  defineRunTimeSelectionTable(IATEsource, dictionary);
42 }
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  }
68 
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);
89 }
90 
92 {
93  return sqrt(2*otherPhase().k());
94 }
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
124  mag(g)*pow4(otherPhase().nu())*sqr(otherPhase().rho())
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 // ************************************************************************* //
twoPhaseSystem & fluid
dictionary dict
IATE (Interfacial Area Transport Equation) bubble diameter model.
Definition: IATE.H:64
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...
type
Types of root.
Definition: Roots.H:52
const dimensionedScalar & rho() const
Return const-access to phase1 density.
Definition: phase.H:151
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
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)
static autoPtr< IATEsource > New(const word &type, const IATE &iate, const dictionary &dict)
Definition: IATEsource.C:43
Various UniformDimensionedField types.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar pow025(const dimensionedScalar &ds)
tmp< volScalarField > Ur() const
Return the bubble relative velocity.
Definition: IATEsource.C:68
defineTypeNameAndDebug(constant, 0)
label k
Boltzmann constant.
tmp< volScalarField > We() const
Return the bubble Webber number.
Definition: IATEsource.C:133
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
A class for handling words, derived from Foam::string.
Definition: word.H:63
const Time & time() const noexcept
Return time registry.
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:450
tmp< volScalarField > CD() const
Return the bubble drag coefficient.
Definition: IATEsource.C:94
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
tmp< volScalarField > Ut() const
Return the bubble turbulent velocity.
Definition: IATEsource.C:84
const Mesh & mesh() const noexcept
Return mesh.
const uniformDimensionedVectorField & g
scalarField Re(const UList< complex > &cmplx)
Extract real component.
Definition: complexField.C:207
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar pow4(const dimensionedScalar &ds)
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:50
tmp< volScalarField > Eo() const
Return the bubble Eotvos number.
Definition: IATEsource.C:122
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
tmp< volScalarField > Mo() const
Return the bubble Morton number.
Definition: IATEsource.C:111
A class for managing temporary objects.
Definition: HashPtrTable.H:50
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:635
volScalarField & nu
Namespace for OpenFOAM.
defineRunTimeSelectionTable(binaryBreakupModel, dictionary)
tmp< volScalarField > Re() const
Return the bubble Reynolds number.
Definition: IATEsource.C:89
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...