turbulentBreakUp.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) 2013-2015 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "turbulentBreakUp.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace diameterModels
36 {
37 namespace IATEsources
38 {
39  defineTypeNameAndDebug(turbulentBreakUp, 0);
40  addToRunTimeSelectionTable(IATEsource, turbulentBreakUp, dictionary);
41 }
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
50 (
51  const IATE& iate,
52  const dictionary& dict
53 )
54 :
55  IATEsource(iate),
56  Cti_("Cti", dimless, dict),
57  WeCr_("WeCr", dimless, dict)
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
62 
65 {
66  auto tR = volScalarField::New
67  (
68  "R",
70  iate_.phase().U().mesh(),
72  );
73  auto R = tR();
74 
75  const scalar Cti = Cti_.value();
76  const scalar WeCr = WeCr_.value();
77  volScalarField Ut(this->Ut());
78  volScalarField We(this->We());
79  const tmp<volScalarField> td(iate_.d());
80  const volScalarField& d = td();
81 
82  forAll(R, celli)
83  {
84  if (We[celli] > WeCr)
85  {
86  R[celli] =
87  (1.0/3.0)
88  *Cti/d[celli]
89  *Ut[celli]
90  *sqrt(1 - WeCr/We[celli])
91  *exp(-WeCr/We[celli]);
92  }
93  }
94 
95  return tR;
96 }
97 
98 
99 // ************************************************************************* //
const Type & value() const noexcept
Return const reference to value.
dictionary dict
virtual tmp< volScalarField > R() const
turbulentBreakUp(const IATE &iate, const dictionary &dict)
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual tmp< volScalarField > d() const
Return the Sauter-mean diameter.
Definition: IATE.H:155
const dimensionSet dimless
Dimensionless.
tmp< volScalarField > We() const
Return the bubble Webber number.
Definition: IATEsource.C:133
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
const volVectorField & U() const
Definition: phaseModel.H:198
const phaseModel & phase() const
Return the phase.
addToRunTimeSelectionTable(IATEsource, dummy, word)
dimensionedScalar exp(const dimensionedScalar &ds)
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
tmp< volScalarField > Ut() const
Return the bubble turbulent velocity.
Definition: IATEsource.C:84
const Mesh & mesh() const noexcept
Return mesh.
const IATE & iate_
Reference to the IATE this source applies to.
Definition: IATEsource.H:62
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Do not request registration (bool: false)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127