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-2023 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; value 1.0;}
47  f2 {dSph 2e-3; value 0.0;}
48  f3 {dSph 3e-3; value 0.0;}
49  f4 {dSph 4e-3; value 0.0;}
50  f5 {dSph 5e-3; value 0.0;}
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  tmp<volScalarField> N() const;
111 
112  void scale();
113 
114 
115 public:
116 
117  //- Runtime type information
118  TypeName("velocityGroup");
119 
120 
121  // Constructors
122 
123  //- Construct from dictionary and phase
125  (
127  const phaseModel& phase
128  );
129 
130 
131  //- Destructor
132  virtual ~velocityGroup();
133 
134 
135  // Member Functions
136 
137  //- Return name of populationBalance this velocityGroup belongs to
138  inline const word& popBalName() const;
139 
140  //- Return the populationBalance this velocityGroup belongs to
141  const populationBalanceModel& popBal() const;
142 
143  //- Return sizeGroups belonging to this velocityGroup
144  inline const PtrList<sizeGroup>& sizeGroups() const;
145 
146  //- Return sizeGroups belonging to this velocityGroup
147  inline PtrList<sizeGroup>& sizeGroups();
148 
149  //- Get the diameter field
150  virtual tmp<volScalarField> d() const;
151 
152  //- Get the surface area per unit volume field
153  virtual tmp<volScalarField> Av() const;
154 
155  //- Correct the model
156  virtual void correct();
157 
158  //- Read diameterProperties dictionary
159  virtual bool read(const dictionary& diameterProperties);
160 };
161 
162 
163 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
164 
165 } // End namespace diameterModels
166 } // End namespace Foam
167 
168 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169 
170 #include "velocityGroupI.H"
171 
172 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
173 
174 #endif
175 
176 // ************************************************************************* //
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 keyword definitions, which are a keyword followed by any number of values (e....
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.