IATEturbulentBreakUp.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2013-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "IATEturbulentBreakUp.H"
27 #include "fvmSup.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace diameterModels
35 {
36 namespace IATEsources
37 {
40 }
41 }
42 }
43 
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
49 (
50  const IATE& iate,
51  const dictionary& dict
52 )
53 :
54  IATEsource(iate),
55  Cti_("Cti", dimless, dict),
56  WeCr_("WeCr", dimless, dict)
57 {}
58 
59 
60 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
61 
64 (
65  const volScalarField& alphai,
66  volScalarField& kappai
67 ) const
68 {
70  (
71  IOobject
72  (
73  typedName("R"),
74  iate_.phase().time().name(),
75  iate_.phase().mesh()
76  ),
77  iate_.phase().mesh(),
79  );
80 
81  const scalar Cti = Cti_.value();
82  const scalar WeCr = WeCr_.value();
83  const volScalarField Ut(this->Ut());
84  const volScalarField We(this->We());
85 
86  forAll(R, celli)
87  {
88  if (We[celli] > WeCr)
89  {
90  R[celli] =
91  (Cti/18)*Ut[celli]*sqr(kappai[celli])
92  *sqrt(1 - WeCr/We[celli])*exp(-WeCr/We[celli]);
93  }
94  }
95 
96  return fvm::Su(R, kappai);
97 }
98 
99 
100 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
IATE (Interfacial Area Transport Equation) bubble diameter model.
Definition: IATE.H:69
IATE (Interfacial Area Transport Equation) bubble diameter model run-time selectable sources.
Definition: IATEsource.H:54
Turbulence-induced break-up IATE source as defined in paper:
turbulentBreakUp(const IATE &iate, const dictionary &dict)
virtual tmp< fvScalarMatrix > R(const volScalarField &alphai, volScalarField &kappai) const
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
A class for managing temporary objects.
Definition: tmp.H:55
Calculate the matrix for implicit and explicit sources.
defineTypeNameAndDebug(wallBoiling, 0)
addToRunTimeSelectionTable(IATEsource, wallBoiling, dictionary)
tmp< fvMatrix< Type > > Su(const DimensionedField< Type, volMesh > &, const VolField< Type > &)
Namespace for OpenFOAM.
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimless
const dimensionSet dimTime
dimensionedScalar sqrt(const dimensionedScalar &ds)
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:149
static scalar R(const scalar a, const scalar x)
Definition: invIncGamma.C:102
dictionary dict