velocityGroup.H
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) 2017-2019 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::diameterModels::velocityGroup
28 
29 Description
30  This diameterModel is intended for use with a populationBalanceModel in
31  order to simulate polydispersed bubbly or particulate flows. It can hold any
32  number of sizeGroups from which the Sauter mean diameter is calculated. It
33  can also be used as a diameterModel without a populationBalance and would
34  then behave like a constantDiameter model. In this case, some arbitrary name
35  must be entered for the populationBalance keyword.
36 
37 Usage
38  \table
39  Property | Description
40  populationBalance | Name of the corresponding populationBalance
41  formFactor | Form factor for converting diameter into volume
42  sizeGroups | List of sizeGroups
43  \endtable
44 
45  Example
46  \verbatim
47  diameterModel velocityGroup;
48  velocityGroupCoeffs
49  {
50  populationBalance bubbles;
51 
52  formFactor 0.5235987756;
53 
54  sizeGroups
55  (
56  f0{d 1.00e-3; value 0;}
57  f1{d 1.08e-3; value 0;}
58  f2{d 1.16e-3; value 0.25;}
59  f3{d 1.25e-3; value 0.5;}
60  f4{d 1.36e-3; value 0.25;}
61  f5{d 1.46e-3; value 0;}
62  ...
63  );
64  }
65  \endverbatim
66 
67 See also
68  Foam::diameterModels::sizeGroup
69  Foam::diameterModels::populationBalanceModel
70 
71 SourceFiles
72  velocityGroup.C
73 
74 \*---------------------------------------------------------------------------*/
75 
76 #ifndef velocityGroup_H
77 #define velocityGroup_H
78 
79 #include "diameterModel.H"
80 #include "multivariateScheme.H"
81 #include "convectionScheme.H"
82 
83 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
84 
85 namespace Foam
86 {
87 namespace diameterModels
88 {
89 
90 class sizeGroup;
91 
92 /*---------------------------------------------------------------------------*\
93  Class velocityGroup Declaration
94 \*---------------------------------------------------------------------------*/
95 
96 class velocityGroup
97 :
98  public diameterModel
99 {
100  // Private data
101 
102  //- Name of the populationBalance this velocityGroup belongs to
103  word popBalName_;
104 
105  //- Reference field from which the sizeGroup fields are derived
106  volScalarField f_;
107 
108  //- Form factor relating diameter and volume
109  dimensionedScalar formFactor_;
110 
111  //- sizeGroups belonging to this velocityGroup
112  PtrList<sizeGroup> sizeGroups_;
113 
114  //- Sum of sizeGroup volume fractions
115  volScalarField fSum_;
116 
117  //- Number-based Sauter-mean diameter of the phase
118  volScalarField d_;
119 
120  //- Multivariate convection scheme
121  tmp<fv::convectionScheme<scalar>> mvConvection_;
122 
123  //- Table of fields for multivariate convection
125 
126  //- Mass transfer rate
127  volScalarField dmdt_;
128 
129 
130  // Private member functions
131 
132  tmp<volScalarField> dsm() const;
133 
134  tmp<volScalarField> fSum() const;
135 
136  void renormalize();
137 
138  tmp<Foam::fv::convectionScheme<Foam::scalar>> mvconvection() const;
139 
140 
141 public:
142 
143  //- Runtime type information
144  TypeName("velocityGroup");
145 
146 
147  // Constructors
148 
149  //- Construct from components
151  (
153  const phaseModel& phase
154  );
155 
156 
157  //- Destructor
158  virtual ~velocityGroup();
159 
160 
161  // Member Functions
162 
163  //- Return name of populationBalance this velocityGroup belongs to
164  inline const word& popBalName() const;
165 
166  //- Return reference field for sizeGroup's
167  inline const volScalarField& f() const;
168 
169  //- Return the form factor
170  inline const dimensionedScalar& formFactor() const;
171 
172  //- Return sizeGroups belonging to this velocityGroup
173  inline const PtrList<sizeGroup>& sizeGroups() const;
174 
175  //- Return const-reference to multivariate convectionScheme
176  inline const tmp<fv::convectionScheme<scalar>>& mvConvection() const;
177 
178  //- Return const-reference to the mass transfer rate
179  inline const volScalarField& dmdt() const;
180 
181  //- Return reference to the mass transfer rate
182  inline volScalarField& dmdtRef();
183 
184  //- Corrections before populationBalanceModel::solve()
185  void preSolve();
186 
187  //- Corrections after populationBalanceModel::solve()
188  void postSolve();
189 
190  //- Read diameterProperties dictionary
191  virtual bool read(const dictionary& diameterProperties);
192 
193  //- Return the Sauter-mean diameter
194  virtual tmp<volScalarField> d() const;
195 };
196 
197 
198 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
199 
200 } // End namespace diameterModels
201 } // End namespace Foam
202 
203 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
204 
205 #include "velocityGroupI.H"
206 
207 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
208 
209 #endif
210 
211 // ************************************************************************* //
Foam::diameterModels::velocityGroup::velocityGroup
velocityGroup(const dictionary &diameterProperties, const phaseModel &phase)
Construct from components.
Definition: velocityGroup.C:147
Foam::phaseModel
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:54
convectionScheme.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
multivariateScheme.H
Foam::phase
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:54
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::diameterModels::velocityGroup::formFactor
const dimensionedScalar & formFactor() const
Return the form factor.
Definition: velocityGroupI.H:45
Foam::diameterModels::velocityGroup::read
virtual bool read(const dictionary &diameterProperties)
Read diameterProperties dictionary.
Definition: velocityGroup.C:323
Foam::diameterModel::diameterProperties
const dictionary & diameterProperties() const
Return the phase diameter properties dictionary.
Definition: diameterModel.H:110
Foam::diameterModel
Abstract base-class for dispersed-phase particle diameter models.
Definition: diameterModel.H:53
Foam::diameterModels::velocityGroup::~velocityGroup
virtual ~velocityGroup()
Destructor.
Definition: velocityGroup.C:277
Foam::diameterModels::velocityGroup::mvConvection
const tmp< fv::convectionScheme< scalar > > & mvConvection() const
Return const-reference to multivariate convectionScheme.
Definition: velocityGroupI.H:59
Foam::diameterModels::velocityGroup::popBalName
const word & popBalName() const
Return name of populationBalance this velocityGroup belongs to.
Definition: velocityGroupI.H:31
velocityGroupI.H
Foam::diameterModels::velocityGroup::postSolve
void postSolve()
Corrections after populationBalanceModel::solve()
Definition: velocityGroup.C:290
Foam::diameterModels::velocityGroup::dmdt
const volScalarField & dmdt() const
Return const-reference to the mass transfer rate.
Definition: velocityGroupI.H:66
Foam::diameterModels::velocityGroup::dmdtRef
volScalarField & dmdtRef()
Return reference to the mass transfer rate.
Definition: velocityGroupI.H:72
Foam::multivariateSurfaceInterpolationScheme
Abstract base class for multi-variate surface interpolation schemes.
Definition: multivariateSurfaceInterpolationScheme.H:52
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
Foam::diameterModels::velocityGroup
This diameterModel is intended for use with a populationBalanceModel in order to simulate polydispers...
Definition: velocityGroup.H:107
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::dimensioned< scalar >
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::diameterModels::velocityGroup::f
const volScalarField & f() const
Return reference field for sizeGroup's.
Definition: velocityGroupI.H:38
Foam::diameterModels::velocityGroup::TypeName
TypeName("velocityGroup")
Runtime type information.
Foam::diameterModels::velocityGroup::preSolve
void preSolve()
Corrections before populationBalanceModel::solve()
Definition: velocityGroup.C:284
Foam::diameterModels::velocityGroup::d
virtual tmp< volScalarField > d() const
Return the Sauter-mean diameter.
Definition: velocityGroup.C:332
Foam::diameterModels::velocityGroup::sizeGroups
const PtrList< sizeGroup > & sizeGroups() const
Return sizeGroups belonging to this velocityGroup.
Definition: velocityGroupI.H:52
Foam::GeometricField< scalar, fvPatchField, volMesh >