turbulentShear.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) 2019-2021 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 "turbulentShear.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace diameterModels
35 {
36 namespace coalescenceModels
37 {
38  defineTypeNameAndDebug(turbulentShear, 0);
40  (
41  coalescenceModel,
42  turbulentShear,
43  dictionary
44  );
45 }
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
54 (
55  const populationBalanceModel& popBal,
56  const dictionary& dict
57 )
58 :
59  coalescenceModel(popBal, dict),
60  C_("C", dimless, dict)
61 {}
62 
63 
64 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
65 
66 void
69 (
70  volScalarField& coalescenceRate,
71  const label i,
72  const label j
73 )
74 {
75  const sizeGroup& fi = popBal_.sizeGroups()[i];
76  const sizeGroup& fj = popBal_.sizeGroups()[j];
77 
79  const volScalarField& rho = popBal_.continuousPhase().rho();
81 
82  coalescenceRate += C_*sqrt(epsilon*rho/mu)*pow3(fi.d() + fj.d());
83 }
84 
85 
86 // ************************************************************************* //
turbulentShear(const populationBalanceModel &popBal, const dictionary &dict)
dictionary dict
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 phaseCompressible::momentumTransportModel & continuousTurbulence() const
Return reference to turbulence model of the continuous phase.
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual tmp< volScalarField > rho() const =0
Return the density field.
const dimensionSet dimless
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
const phaseModel & continuousPhase() const
Return continuous phase.
const rhoThermo & thermo() const
Return const-access to phase rhoThermo.
Definition: phaseModel.H:121
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
defineTypeNameAndDebug(combustionModel, 0)
const UPtrList< sizeGroup > & sizeGroups() const
Return the sizeGroups belonging to this populationBalance.
dimensionedScalar pow3(const dimensionedScalar &ds)
virtual tmp< volScalarField > mu() const =0
Dynamic viscosity of mixture [kg/m/s].
virtual void addToCoalescenceRate(volScalarField &coalescenceRate, const label i, const label j)
Add to coalescenceRate.
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
Namespace for OpenFOAM.