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-2015 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 "twoPhaseSystem.H"
31 #include "fvMatrix.H"
33 #include "gravityMeshObject.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace diameterModels
40 {
41  defineTypeNameAndDebug(IATEsource, 0);
42  defineRunTimeSelectionTable(IATEsource, dictionary);
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
48 
51 (
52  const word& type,
53  const IATE& iate,
54  const dictionary& dict
55 )
56 {
57  auto* ctorPtr = dictionaryConstructorTable(type);
58 
59  if (!ctorPtr)
60  {
62  (
63  dict,
64  "IATEsource",
65  type,
66  *dictionaryConstructorTablePtr_
67  ) << exit(FatalIOError);
68  }
69 
70  return autoPtr<IATEsource>(ctorPtr(iate, dict));
71 }
72 
73 
74 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
75 
77 {
79  meshObjects::gravity::New(phase().U().time());
80 
81  return
82  sqrt(2.0)
83  *pow025
84  (
85  fluid().sigma()*mag(g)
86  *(otherPhase().rho() - phase().rho())
87  /sqr(otherPhase().rho())
88  )
89  *pow(max(1 - phase(), scalar(0)), 1.75);
90 }
91 
92 
94 {
95  return sqrt(2*otherPhase().turbulence().k());
96 }
97 
98 
100 {
101  return max(Ur()*phase().d()/otherPhase().nu(), scalar(1.0e-3));
102 }
103 
104 
106 {
107  const volScalarField Eo(this->Eo());
108  const volScalarField Re(this->Re());
109 
110  return
111  max
112  (
113  min
114  (
115  (16/Re)*(1 + 0.15*pow(Re, 0.687)),
116  48/Re
117  ),
118  8*Eo/(3*(Eo + 4))
119  );
120 }
121 
122 
124 {
126  meshObjects::gravity::New(phase().U().time());
127 
128  return
129  mag(g)*pow4(otherPhase().nu())*sqr(otherPhase().rho())
130  *(otherPhase().rho() - phase().rho())
131  /pow3(fluid().sigma());
132 }
133 
134 
136 {
138  meshObjects::gravity::New(phase().U().time());
139 
140  return
141  mag(g)*sqr(phase().d())
142  *(otherPhase().rho() - phase().rho())
143  /fluid().sigma();
144 }
145 
146 
148 {
149  return otherPhase().rho()*sqr(Ur())*phase().d()/fluid().sigma();
150 }
151 
152 
153 // ************************************************************************* //
twoPhaseSystem & fluid
dictionary dict
type
Types of root.
Definition: Roots.H:52
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
UniformDimensionedField< vector > uniformDimensionedVectorField
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
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)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:84
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:752
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
static const gravity & New(const Time &runTime)
Return cached object or construct on Time.
scalarField Re(const UList< complex > &cf)
Extract real component.
Definition: complexField.C:152
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 uniformDimensionedVectorField & g
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
U
Definition: pEqn.H:72
dimensionedScalar pow4(const dimensionedScalar &ds)
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:615
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 ...