All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-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 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  sizeGroups | List of sizeGroups
40  \endtable
41 
42  Example
43  \verbatim
44  diameterModel velocityGroup;
45  velocityGroupCoeffs
46  {
47  populationBalance bubbles;
48 
49  shapeModel constant;
50 
51  sizeGroups
52  (
53  f0{dSph 1.00e-3; value 0;}
54  f1{dSph 1.08e-3; value 0;}
55  f2{dSph 1.16e-3; value 0.25;}
56  f3{dSph 1.25e-3; value 0.5;}
57  f4{dSph 1.36e-3; value 0.25;}
58  f5{dSph 1.46e-3; value 0;}
59  ...
60  );
61  }
62  \endverbatim
63 
64 See also
65  Foam::diameterModels::sizeGroup
66  Foam::diameterModels::populationBalanceModel
67 
68 SourceFiles
69  velocityGroup.C
70 
71 \*---------------------------------------------------------------------------*/
72 
73 #ifndef velocityGroup_H
74 #define velocityGroup_H
75 
76 #include "diameterModel.H"
77 #include "multivariateScheme.H"
78 #include "convectionScheme.H"
79 
80 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
81 
82 namespace Foam
83 {
84 namespace diameterModels
85 {
86 
87 // Forward declaration of classes
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_;
103  //- Sum of the sizeGroup volume fractions and reference field from which
104  // the sizeGroup fields are derived
105  volScalarField f_;
106 
107  //- sizeGroups belonging to this velocityGroup
108  PtrList<sizeGroup> sizeGroups_;
109 
110  //- Sauter-mean diameter of the phase
111  volScalarField d_;
112 
113 
114  // Private Member Functions
115 
116  tmp<volScalarField> dsm() const;
117 
118  tmp<volScalarField> fSum() const;
119 
120  tmp<volScalarField> N() const;
121 
122  void scale();
123 
124 
125 public:
126 
127  //- Runtime type information
128  TypeName("velocityGroup");
129 
130 
131  // Constructors
132 
133  //- Construct from dictionary and phase
135  (
137  const phaseModel& phase
138  );
139 
140 
141  //- Destructor
142  virtual ~velocityGroup();
143 
144 
145  // Member Functions
146 
147  //- Return name of populationBalance this velocityGroup belongs to
148  inline const word& popBalName() const;
149 
150  //- Return reference field for sizeGroup's
151  inline const volScalarField& f() const;
152 
153  //- Return sizeGroups belonging to this velocityGroup
154  inline const PtrList<sizeGroup>& sizeGroups() const;
155 
156  //- Get the diameter field
157  virtual tmp<volScalarField> d() const;
158 
159  //- Get the surface area per unit volume field
160  virtual tmp<volScalarField> a() const;
161 
162  //- Correct the model
163  virtual void correct();
164 
165  //- Read diameterProperties dictionary
166  virtual bool read(const dictionary& diameterProperties);
167 };
168 
169 
170 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171 
172 } // End namespace diameterModels
173 } // End namespace Foam
174 
175 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
176 
177 #include "velocityGroupI.H"
178 
179 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
180 
181 #endif
182 
183 // ************************************************************************* //
virtual tmp< volScalarField > d() const
Get the diameter field.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
velocityGroup(const dictionary &diameterProperties, const phaseModel &phase)
Construct from dictionary and phase.
virtual tmp< volScalarField > a() const
Get the surface area per unit volume field.
TypeName("velocityGroup")
Runtime type information.
const word & popBalName() const
Return name of populationBalance this velocityGroup belongs to.
A class for handling words, derived from string.
Definition: word.H:59
virtual ~velocityGroup()
Destructor.
const volScalarField & f() const
Return reference field for sizeGroup&#39;s.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
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.
virtual void correct()
Correct the model.