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 -------------------------------------------------------------------------------
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"
29 #include "fvmSup.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace diameterModels
37 {
38 namespace IATEsources
39 {
40  defineTypeNameAndDebug(turbulentBreakUp, 0);
41  addToRunTimeSelectionTable(IATEsource, turbulentBreakUp, dictionary);
42 }
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
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  (
73  IOobject
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 // ************************************************************************* //
const Type & value() const noexcept
Return const reference to value.
dictionary dict
IATE (Interfacial Area Transport Equation) bubble diameter model.
Definition: IATE.H:64
virtual tmp< volScalarField > R() const
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
turbulentBreakUp(const IATE &iate, const dictionary &dict)
dimensionedScalar sqrt(const dimensionedScalar &ds)
IATE (Interfacial Area Transport Equation) bubble diameter model run-time selectable sources...
Definition: IATEsource.H:52
const dimensionSet dimless
Dimensionless.
zeroField Su(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
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 phaseModel & phase() const
Return the phase.
addToRunTimeSelectionTable(IATEsource, dummy, word)
dimensionedScalar exp(const dimensionedScalar &ds)
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition: IOobject.C:456
tmp< volScalarField > Ut() const
Return the bubble turbulent velocity.
Definition: IATEsource.C:84
const Mesh & mesh() const noexcept
Return mesh.
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
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.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
Calculate the finiteVolume matrix for implicit and explicit sources.
Namespace for OpenFOAM.
const dimensionSet & dimensions() const noexcept
Return dimensions.