Loading...
Searching...
No Matches
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{
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());
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// ************************************************************************* //
#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 noexcept
Return const reference to mesh.
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())
@ NO_REGISTER
Do not request registration (bool: false).
const phaseModel & phase() const
Return the phase.
virtual tmp< volScalarField > d() const
Return the Sauter-mean diameter.
Definition IATE.H:155
const IATE & iate_
Reference to the IATE this source applies to.
Definition IATEsource.H:62
tmp< volScalarField > We() const
Return the bubble Webber number.
Definition IATEsource.C:133
tmp< volScalarField > Ut() const
Return the bubble turbulent velocity.
Definition IATEsource.C:84
Turbulence-induced break-up IATE source as defined in paper:
virtual tmp< volScalarField > R() const
turbulentBreakUp(const IATE &iate, const dictionary &dict)
virtual tmp< fvScalarMatrix > R(const volScalarField &alphai, volScalarField &kappai) const
const Type & value() const noexcept
Return const reference to value.
const volVectorField & U() const
Definition phaseModel.H:198
A class for managing temporary objects.
Definition tmp.H:75
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
Namespace for OpenFOAM.
dimensionedScalar exp(const dimensionedScalar &ds)
const dimensionSet dimless
Dimensionless.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar sqrt(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299