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-2018 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"
29#include "fvmSup.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace diameterModels
37{
38namespace IATEsources
39{
42}
43}
44}
45
46
47// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48
49Foam::diameterModels::IATEsources::turbulentBreakUp::
50turbulentBreakUp
51(
52 const IATE& iate,
53 const dictionary& dict
54)
55:
56 IATEsource(iate),
57 Cti_("Cti", dimless, dict),
58 WeCr_("WeCr", dimless, dict)
59{}
60
61
62// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
63
66(
67 const volScalarField& alphai,
68 volScalarField& kappai
69) const
70{
72 (
74 (
75 "turbulentBreakUp:R",
76 iate_.phase().time().timeName(),
77 iate_.phase().mesh()
78 ),
79 iate_.phase().mesh(),
81 );
82
83 const scalar Cti = Cti_.value();
84 const scalar WeCr = WeCr_.value();
85 const volScalarField Ut(this->Ut());
86 const volScalarField We(this->We());
87
88 forAll(R, celli)
89 {
90 if (We[celli] > WeCr)
91 {
92 R[celli] =
93 2*Cti*Ut[celli]*sqrt(1 - WeCr/We[celli])*exp(-WeCr/We[celli]);
94 }
95 }
96
97 return fvm::Su(R, kappai);
98}
99
100
101// ************************************************************************* //
#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 dimensionSet & dimensions() const
Return dimensions.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
IATE (Interfacial Area Transport Equation) bubble diameter model.
Definition: IATE.H:71
IATE (Interfacial Area Transport Equation) bubble diameter model run-time selectable sources.
Definition: IATEsource.H:57
Turbulence-induced break-up IATE source as defined in paper:
virtual tmp< volScalarField > R() const
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A class for managing temporary objects.
Definition: tmp.H:65
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
Calculate the finiteVolume matrix for implicit and explicit sources.
zeroField Su(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
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)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333