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-2024 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  Computes the Sauter mean diameter based on a user specified size
29  distribution, defined in terms of size class fractions. Intended for use
30  with a population balance model to account for the evolution of a size
31  distribution by means of coalescence, breakup, drift and nucleation.
32 
33 Usage
34  Excerpt from an exemplary phaseProperties dictionary:
35  \verbatim
36  diameterModel velocityGroup;
37 
38  velocityGroupCoeffs
39  {
40  populationBalance bubbles;
41 
42  shapeModel spherical;
43 
44  sizeGroups
45  (
46  f1 {dSph 1e-3; }
47  f2 {dSph 2e-3; }
48  f3 {dSph 3e-3; }
49  f4 {dSph 4e-3; }
50  f5 {dSph 5e-3; }
51  ...
52  );
53  }
54  \endverbatim
55 
56 See also
57  Foam::diameterModels::sizeGroup
58  Foam::diameterModels::populationBalanceModel
59 
60 SourceFiles
61  velocityGroup.C
62 
63 \*---------------------------------------------------------------------------*/
64 
65 #ifndef velocityGroup_H
66 #define velocityGroup_H
67 
68 #include "diameterModel.H"
69 
70 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
71 
72 namespace Foam
73 {
74 namespace diameterModels
75 {
76 
77 // Forward declaration of classes
78 class sizeGroup;
79 class populationBalanceModel;
80 
81 /*---------------------------------------------------------------------------*\
82  Class velocityGroup Declaration
83 \*---------------------------------------------------------------------------*/
84 
85 class velocityGroup
86 :
87  public diameterModel
88 {
89  // Private Data
90 
91  //- Name of the populationBalance this velocityGroup belongs to
92  word popBalName_;
93 
94  //- Pointer to the populationBalance this velocityGroup belongs to
95  mutable const populationBalanceModel* popBalPtr_;
96 
97  //- sizeGroups belonging to this velocityGroup
98  PtrList<sizeGroup> sizeGroups_;
99 
100  //- Sauter-mean diameter of the phase
101  volScalarField d_;
102 
103 
104  // Private Member Functions
105 
106  tmp<volScalarField> dsm() const;
107 
108  tmp<volScalarField> fSum() const;
109 
110  void scale();
111 
112 
113 public:
114 
115  //- Runtime type information
116  TypeName("velocityGroup");
117 
118 
119  // Constructors
120 
121  //- Construct from dictionary and phase
123  (
125  const phaseModel& phase
126  );
127 
128 
129  //- Destructor
130  virtual ~velocityGroup();
131 
132 
133  // Member Functions
134 
135  //- Return name of populationBalance this velocityGroup belongs to
136  inline const word& popBalName() const;
137 
138  //- Return the populationBalance this velocityGroup belongs to
139  const populationBalanceModel& popBal() const;
140 
141  //- Return sizeGroups belonging to this velocityGroup
142  inline const PtrList<sizeGroup>& sizeGroups() const;
143 
144  //- Return sizeGroups belonging to this velocityGroup
145  inline PtrList<sizeGroup>& sizeGroups();
146 
147  //- Get the diameter field
148  virtual tmp<volScalarField> d() const;
149 
150  //- Get the surface area per unit volume field
151  virtual tmp<volScalarField> Av() const;
152 
153  //- Correct the model
154  virtual void correct();
155 
156  //- Read diameterProperties dictionary
157  virtual bool read(const dictionary& diameterProperties);
158 };
159 
160 
161 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
162 
163 } // End namespace diameterModels
164 } // End namespace Foam
165 
166 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
167 
168 #include "velocityGroupI.H"
169 
170 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171 
172 #endif
173 
174 // ************************************************************************* //
Generic GeometricField class.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
Abstract base-class for dispersed-phase particle diameter models.
Definition: diameterModel.H:52
const phaseModel & phase() const
Return the phase.
const dictionary & diameterProperties() const
Return the phase diameter properties dictionary.
Model for tracking the evolution of a dispersed phase size distribution due to coalescence (synonymou...
Computes the Sauter mean diameter based on a user specified size distribution, defined in terms of si...
Definition: velocityGroup.H:87
virtual void correct()
Correct the model.
const populationBalanceModel & popBal() const
Return the populationBalance this velocityGroup belongs to.
virtual bool read(const dictionary &diameterProperties)
Read diameterProperties dictionary.
virtual tmp< volScalarField > d() const
Get the diameter field.
velocityGroup(const dictionary &diameterProperties, const phaseModel &phase)
Construct from dictionary and phase.
virtual ~velocityGroup()
Destructor.
TypeName("velocityGroup")
Runtime type information.
virtual tmp< volScalarField > Av() const
Get the surface area per unit volume field.
const word & popBalName() const
Return name of populationBalance this velocityGroup belongs to.
const PtrList< sizeGroup > & sizeGroups() const
Return sizeGroups belonging to this velocityGroup.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.