LaakkonenAlopaeusAittamaa.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) 2018-2019 OpenFOAM Foundation
9 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
31#include "phaseCompressibleTurbulenceModel.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
37namespace diameterModels
38{
39namespace breakupModels
40{
43 (
47 );
48}
49}
50}
51
52
53// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54
55Foam::diameterModels::breakupModels::LaakkonenAlopaeusAittamaa::
56LaakkonenAlopaeusAittamaa
57(
58 const populationBalanceModel& popBal,
59 const dictionary& dict
60)
61:
62 breakupModel(popBal, dict),
63 C1_
64 (
65 dimensionedScalar::getOrDefault
66 (
67 "C1",
68 dict,
69 dimensionSet(0, -2.0/3.0, 0, 0, 0),
70 6.0
71 )
72 ),
73 C2_(dimensionedScalar::getOrDefault("C2", dict, dimless, 0.04)),
74 C3_(dimensionedScalar::getOrDefault("C3", dict, dimless, 0.01))
75{}
76
77
78// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79
80void
82(
83 volScalarField& breakupRate,
84 const label i
85)
86{
87 const phaseModel& continuousPhase = popBal_.continuousPhase();
88 const sizeGroup& fi = popBal_.sizeGroups()[i];
89
90 breakupRate =
91 C1_*cbrt(popBal_.continuousTurbulence().epsilon())
92 *erfc
93 (
94 sqrt
95 (
96 C2_*popBal_.sigmaWithContinuousPhase(fi.phase())
97 /(
98 continuousPhase.rho()*pow(fi.d(), 5.0/3.0)
99 *pow(popBal_.continuousTurbulence().epsilon(), 2.0/3.0)
100 )
101 + C3_*continuousPhase.mu()
102 /(
103 sqrt(continuousPhase.rho()*fi.phase().rho())
104 *cbrt(popBal_.continuousTurbulence().epsilon())
105 *pow(fi.d(), 4.0/3.0)
106 )
107 )
108 );
109}
110
111
112// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Base class for breakup models which give a total breakup rate and a separate daughter size distributi...
Definition: breakupModel.H:56
Model of Laakkonen et al. (2006). The total breakup rate is calculated by.
virtual void setBreakupRate(volScalarField &breakupRate, const label i)
Set total breakupRate.
Class that solves the univariate population balance equation by means of a class method (also called ...
This class represents a single sizeGroup belonging to a velocityGroup. The main property of a sizeGro...
Definition: sizeGroup.H:99
const dimensionedScalar & d() const
Return representative diameter of the sizeGroup.
Definition: sizeGroupI.H:52
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:38
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:109
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:61
const dimensionedScalar & rho() const
Definition: phaseModel.H:171
virtual tmp< volScalarField > mu() const =0
Return the laminar dynamic viscosity.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
Namespace for OpenFOAM.
dimensionedScalar erfc(const dimensionedScalar &ds)
const dimensionSet dimless
Dimensionless.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar cbrt(const dimensionedScalar &ds)
dictionary dict