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-------------------------------------------------------------------------------
10License
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
33namespace Foam
34{
35namespace diameterModels
36{
37namespace IATEsources
38{
39 defineTypeNameAndDebug(turbulentBreakUp, 0);
40 addToRunTimeSelectionTable(IATEsource, turbulentBreakUp, dictionary);
41}
42}
43}
44
45
46// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47
48Foam::diameterModels::IATEsources::turbulentBreakUp::
49turbulentBreakUp
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{
67 (
69 (
71 (
72 "R",
73 iate_.phase().U().time().timeName(),
74 iate_.phase().mesh()
75 ),
76 iate_.phase().U().mesh(),
78 )
79 );
80
81 volScalarField R = tR();
82
83 scalar Cti = Cti_.value();
84 scalar WeCr = WeCr_.value();
85 volScalarField Ut(this->Ut());
86 volScalarField We(this->We());
87 const volScalarField& d(iate_.d()());
88
89 forAll(R, celli)
90 {
91 if (We[celli] > WeCr)
92 {
93 R[celli] =
94 (1.0/3.0)
95 *Cti/d[celli]
96 *Ut[celli]
97 *sqrt(1 - WeCr/We[celli])
98 *exp(-WeCr/We[celli]);
99 }
100 }
101
102 return tR;
103}
104
105
106// ************************************************************************* //
#define R(A, B, C, D, E, F, K, M)
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const Mesh & mesh() const
Return mesh.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const Time & time() const
Return Time associated with the objectRegistry.
Definition: IOobject.C:506
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
const phaseModel & phase() const
Return the phase.
virtual tmp< volScalarField > d() const
Return the Sauter-mean diameter.
Definition: IATE.H:135
const IATE & iate_
Reference to the IATE this source applies to.
Definition: IATEsource.H:64
tmp< volScalarField > We() const
Return the bubble Webber number.
Definition: IATEsource.C:140
tmp< volScalarField > Ut() const
Return the bubble turbulent velocity.
Definition: IATEsource.C:91
virtual tmp< volScalarField > R() const
const Type & value() const
Return const reference to value.
const volVectorField & U() const
Definition: phaseModel.H:176
A class for managing temporary objects.
Definition: tmp.H:65
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
Namespace for OpenFOAM.
dimensionedScalar exp(const dimensionedScalar &ds)
const dimensionSet dimless
Dimensionless.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
dimensionedScalar sqrt(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333