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-2025 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"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace diameterModels
35 {
36 namespace breakupModels
37 {
40  (
42  Laakkonen,
44  );
45 }
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
54 (
55  const populationBalanceModel& popBal,
56  const dictionary& dict
57 )
58 :
59  breakupModel(popBal, dict),
60  C1_("C1", dimensionSet(0, -2.0/3.0, 0, 0, 0), dict, 2.25),
61  C2_("C2", dimless, dict, 0.04),
62  C3_("C3", dimless, dict, 0.01)
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
67 
69 (
70  volScalarField::Internal& breakupRate,
71  const label i
72 )
73 {
74  const sizeGroup& fi = popBal_.sizeGroups()[i];
75 
76  const volScalarField::Internal& rhoc = popBal_.continuousPhase().rho();
77 
78  tmp<volScalarField> tsigma(popBal_.sigmaWithContinuousPhase(fi.phase()));
79  const volScalarField::Internal& sigma = tsigma();
80 
81  tmp<volScalarField> tepsilonc(popBal_.continuousTurbulence().epsilon());
82  const volScalarField::Internal& epsilonc = tepsilonc();
83  tmp<volScalarField> tmu(popBal_.continuousPhase().fluidThermo().mu());
84  const volScalarField::Internal muc = tmu();
85 
86  breakupRate =
87  C1_
88  *cbrt(epsilonc)
89  *erfc
90  (
91  sqrt
92  (
93  C2_
94  *sigma
95  /(
96  rhoc*pow(fi.dSph(), 5.0/3.0)
97  *pow(epsilonc, 2.0/3.0)
98  )
99  + C3_
100  *muc
101  /(
102  sqrt(rhoc*fi.phase().rho())
103  *cbrt(epsilonc)
104  *pow(fi.dSph(), 4.0/3.0)
105  )
106  )
107  );
108 }
109 
110 
111 // ************************************************************************* //
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...
Base class for breakup models which provide a total breakup rate and a separate daughter size distrib...
Definition: breakupModel.H:56
Model of Laakkonen et al. (2007). The total breakup rate is calculated by.
Definition: Laakkonen.H:144
virtual void setBreakupRate(volScalarField::Internal &breakupRate, const label i)
Set total breakupRate.
Definition: Laakkonen.C:69
Laakkonen(const populationBalanceModel &popBal, const dictionary &dict)
Definition: Laakkonen.C:54
Model for tracking the evolution of a dispersed phase size distribution due to coalescence (synonymou...
Single size class fraction field representing a fixed particle volume as defined by the user through ...
Definition: sizeGroup.H:96
const dimensionedScalar & dSph() const
Return representative spherical diameter of the sizeGroup.
Definition: sizeGroupI.H:51
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:37
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Dimension set for the base types.
Definition: dimensionSet.H:125
virtual const volScalarField & rho() const =0
Return the density field.
A class for managing temporary objects.
Definition: tmp.H:55
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
addToRunTimeSelectionTable(breakupModel, exponential, dictionary)
Namespace for OpenFOAM.
dimensionedScalar erfc(const dimensionedScalar &ds)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const dimensionSet dimless
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
void cbrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dictionary dict