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