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 | Copyright (C) 2011-2015 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 {
37  defineTypeNameAndDebug(IATEsource, 0);
38  defineRunTimeSelectionTable(IATEsource, dictionary);
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  (
60  "IATEsource::New"
61  "(const word& type, const IATE&, const dictionary&)"
62  ) << "Unknown IATE source type "
63  << type << nl << nl
64  << "Valid IATE source types : " << endl
65  << dictionaryConstructorTablePtr_->sortedToc()
66  << exit(FatalError);
67  }
68 
69  return autoPtr<IATEsource>(cstrIter()(iate, dict));
70 }
71 
72 
73 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74 
76 {
78  phase().db().lookupObject<uniformDimensionedVectorField>("g");
79 
80  return
81  sqrt(2.0)
82  *pow025
83  (
84  fluid().sigma()*mag(g)
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().turbulence().k());
94 }
95 
97 {
98  return max(Ur()*phase().d()/otherPhase().nu(), scalar(1.0e-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 {
121  phase().db().lookupObject<uniformDimensionedVectorField>("g");
122 
123  return
124  mag(g)*pow4(otherPhase().nu())*sqr(otherPhase().rho())
125  *(otherPhase().rho() - phase().rho())
126  /pow3(fluid().sigma());
127 }
128 
130 {
132  phase().db().lookupObject<uniformDimensionedVectorField>("g");
133 
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 // ************************************************************************* //
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar pow3(const dimensionedScalar &ds)
multiphaseSystem & fluid
Definition: createFields.H:10
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
tmp< volScalarField > We() const
Return the bubble Webber number.
dimensioned< scalar > mag(const dimensioned< Type > &)
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
tmp< volScalarField > CD() const
Return the bubble drag coefficient.
Namespace for OpenFOAM.
tmp< volScalarField > Eo() const
Return the bubble Eotvos number.
UniformDimensionedField< vector > uniformDimensionedVectorField
dictionary dict
static const char nl
Definition: Ostream.H:260
const double e
Elementary charge.
Definition: doubleFloat.H:78
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
label k
Boltzmann constant.
dimensionedScalar pow025(const dimensionedScalar &ds)
tmp< volScalarField > Ut() const
Return the bubble turbulent velocity.
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
error FatalError
dimensionedScalar pow4(const dimensionedScalar &ds)
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
volScalarField & nu
tmp< volScalarField > Re() const
Return the bubble Reynolds number.
static autoPtr< IATEsource > New(const word &type, const IATE &iate, const dictionary &dict)
tmp< volScalarField > Mo() const
Return the bubble Morton number.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
tmp< volScalarField > Ur() const
Return the bubble relative velocity.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
A class for managing temporary objects.
Definition: PtrList.H:118
defineTypeNameAndDebug(combustionModel, 0)
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97