All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Laakkonen.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) 2018-2022 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 "Laakkonen.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace diameterModels
34 {
35 namespace breakupModels
36 {
37  defineTypeNameAndDebug(Laakkonen, 0);
39  (
40  breakupModel,
41  Laakkonen,
42  dictionary
43  );
44 }
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
53 (
54  const populationBalanceModel& popBal,
55  const dictionary& dict
56 )
57 :
58  breakupModel(popBal, dict),
59  C1_
60  (
61  dimensionedScalar::lookupOrDefault
62  (
63  "C1",
64  dict,
65  dimensionSet(0, -2.0/3.0, 0, 0, 0),
66  2.25
67  )
68  ),
69  C2_(dimensionedScalar::lookupOrDefault("C2", dict, dimless, 0.04)),
70  C3_(dimensionedScalar::lookupOrDefault("C3", dict, dimless, 0.01))
71 {}
72 
73 
74 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
75 
76 void
78 (
79  volScalarField& breakupRate,
80  const label i
81 )
82 {
83  const phaseModel& continuousPhase = popBal_.continuousPhase();
84  const sizeGroup& fi = popBal_.sizeGroups()[i];
85 
86  breakupRate =
88  *erfc
89  (
90  sqrt
91  (
92  C2_*popBal_.sigmaWithContinuousPhase(fi.phase())
93  /(
94  continuousPhase.rho()*pow(fi.dSph(), 5.0/3.0)
96  )
97  + C3_*continuousPhase.thermo().mu()
98  /(
99  sqrt(continuousPhase.rho()*fi.phase().rho())
101  *pow(fi.dSph(), 4.0/3.0)
102  )
103  )
104  );
105 }
106 
107 
108 // ************************************************************************* //
dictionary dict
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
virtual void setBreakupRate(volScalarField &breakupRate, const label i)
Set total breakupRate.
Laakkonen(const populationBalanceModel &popBal, const dictionary &dict)
const phaseCompressible::momentumTransportModel & continuousTurbulence() const
Return reference to turbulence model of the continuous phase.
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionSet dimless
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
const phaseModel & continuousPhase() const
Return continuous phase.
dimensionedScalar cbrt(const dimensionedScalar &ds)
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
Definition: breakupModel.H:62
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
const UPtrList< sizeGroup > & sizeGroups() const
Return the size groups belonging to this populationBalance.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar erfc(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const tmp< volScalarField > sigmaWithContinuousPhase(const phaseModel &dispersedPhase) const
Return the surface tension coefficient between a given dispersed.
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
Namespace for OpenFOAM.