sizeGroup.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::sizeGroup
26 
27 Description
28  Single size class fraction field representing a fixed particle volume as
29  defined by the user through the corresponding sphere equivalent diameter.
30 
31  If present, the field is read from the start time directory, e.g. during
32  restart. Otherwise, it is constructed as a uniform field using the size
33  class fraction value provided by the user. In the latter case, the boundary
34  condition types are taken from the 'f.<phaseName>' field, which represents
35  the sum of all size group fractions of a phase. The user specified value is
36  also applied at fixed value boundary conditions and used as inlet value for
37  mixed boundary conditions.
38 
39  An alternative diameter field is provided by the selected shape model, e.g.
40  a collisional diameter, which is then utilised in selected population
41  balance submodels, e.g. for modelling fractal aggregation.
42 
43 Usage
44  \table
45  Property | Description | Required | Default value
46  dSph | Sphere equivalent diameter | yes | none
47  value | Initial and boundary condition value of\\
48  size class fraction | yes | none
49  shapeModel | Shape model providing an alternative diameter field\\
50  | yes | none
51  \endtable
52 
53 See also
54  Foam::diameterModels::populationBalanceModel
55  Foam::diameterModels::velocityGroup
56  Foam::diameterModels::shapeModel
57 
58 SourceFiles
59  sizeGroup.C
60 
61 \*---------------------------------------------------------------------------*/
62 
63 #ifndef sizeGroup_H
64 #define sizeGroup_H
65 
66 #include "dictionaryEntry.H"
67 #include "velocityGroup.H"
68 
69 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
70 
71 namespace Foam
72 {
73 namespace diameterModels
74 {
75 
76 class shapeModel;
77 
78 /*---------------------------------------------------------------------------*\
79  Class sizeGroup Declaration
80 \*---------------------------------------------------------------------------*/
81 
82 class sizeGroup
83 :
84  public volScalarField
85 {
86 private:
87 
88  // Private Data
89 
90  dictionary dict_;
91 
92  //- Phase this sizeGroup belongs to
93  const phaseModel& phase_;
94 
95  //- VelocityGroup this sizeGroup belongs to
96  const velocityGroup& velocityGroup_;
97 
98  //- Sphere equivalent diameter of the sizeGroup
99  const dimensionedScalar dSph_;
100 
101  //- Representative volume of the sizeGroup
102  const dimensionedScalar x_;
103 
104  //- Initial value and value at boundaries
105  const scalar value_;
106 
107  //- Label of this sizeGroup within the corresponding populationBalance
108  mutable autoPtr<label> i_;
109 
110  //- Model for describing the representative shape of the elements in the
111  // sizeGroup
112  autoPtr<shapeModel> shapeModel_;
113 
114 
115 public:
116 
117  // Constructors
118 
119  sizeGroup
120  (
121  const word& name,
122  const dictionary& dict,
123  const phaseModel& phase,
125  const fvMesh& mesh
126  );
127 
128  //- Return clone
129  autoPtr<sizeGroup> clone() const;
130 
131  //- Return a pointer to a new sizeGroup created on freestore
132  // from Istream
133  class iNew
134  {
135  const phaseModel& phase_;
136  const velocityGroup& velocityGroup_;
137 
138  public:
139 
140  iNew
141  (
142  const phaseModel& phase,
144  )
145  :
146  phase_(phase),
147  velocityGroup_(velocityGroup)
148  {}
149 
151  {
153  return autoPtr<sizeGroup>
154  (
155  new sizeGroup
156  (
157  ent.keyword(),
158  ent,
159  phase_,
160  velocityGroup_,
161  phase_.mesh()
162  )
163  );
164  }
165  };
166 
167 
168  //- Destructor
169  virtual ~sizeGroup();
170 
171 
172  // Member Functions
173 
174  //- Return const-reference to the dictionary
175  inline const dictionary& dict() const;
176 
177  //- Return const-reference to the phase
178  inline const phaseModel& phase() const;
179 
180  //- Return const-reference to the velocityGroup
181  inline const velocityGroup& VelocityGroup() const;
182 
183  //- Return representative spherical diameter of the sizeGroup
184  inline const dimensionedScalar& dSph() const;
185 
186  //- Return representative volume of the sizeGroup
187  inline const dimensionedScalar& x() const;
188 
189  //- Return reference to diameterModel of the phase
191 
192  //- Return const-reference to diameterModel of the phase
193  inline const autoPtr<shapeModel>& shapeModelPtr() const;
194 
195  //- Return index of the size group within the population balance
196  const label& i() const;
197 
198  //- Return representative surface area of the sizeGroup
199  const tmp<volScalarField> a() const;
200 
201  //- Return representative diameter of the sizeGroup
202  const tmp<volScalarField> d() const;
203 
204  //- Correct secondary properties
205  void correct();
206 };
207 
208 
209 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
210 
211 } // End namespace diameterModels
212 } // End namespace Foam
213 
214 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
215 
216 #include "sizeGroupI.H"
217 
218 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
219 
220 #endif
221 
222 // ************************************************************************* //
const Mesh & mesh() const
Return mesh.
const word & name() const
Return name.
Definition: IOobject.H:310
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
Return a pointer to a new sizeGroup created on freestore.
Definition: sizeGroup.H:151
autoPtr< sizeGroup > operator()(Istream &is) const
Definition: sizeGroup.H:167
iNew(const phaseModel &phase, const velocityGroup &velocityGroup)
Definition: sizeGroup.H:158
Single size class fraction field representing a fixed particle volume as defined by the user through ...
Definition: sizeGroup.H:102
const dimensionedScalar & dSph() const
Return representative spherical diameter of the sizeGroup.
Definition: sizeGroupI.H:50
void correct()
Correct secondary properties.
Definition: sizeGroup.C:148
const label & i() const
Return index of the size group within the population balance.
Definition: sizeGroup.C:111
const dimensionedScalar & x() const
Return representative volume of the sizeGroup.
Definition: sizeGroupI.H:57
const velocityGroup & VelocityGroup() const
Return const-reference to the velocityGroup.
Definition: sizeGroupI.H:43
sizeGroup(const word &name, const dictionary &dict, const phaseModel &phase, const velocityGroup &velocityGroup, const fvMesh &mesh)
Definition: sizeGroup.C:38
virtual ~sizeGroup()
Destructor.
Definition: sizeGroup.C:97
const tmp< volScalarField > a() const
Return representative surface area of the sizeGroup.
Definition: sizeGroup.C:135
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:36
autoPtr< shapeModel > & shapeModelPtr()
Return reference to diameterModel of the phase.
Definition: sizeGroupI.H:64
autoPtr< sizeGroup > clone() const
Return clone.
Definition: sizeGroup.C:104
const dictionary & dict() const
Return const-reference to the dictionary.
Definition: sizeGroupI.H:29
const tmp< volScalarField > d() const
Return representative diameter of the sizeGroup.
Definition: sizeGroup.C:142
Computes the Sauter mean diameter based on a user specified size distribution, defined in terms of si...
Definition: velocityGroup.H:86
A keyword and a list of tokens is a 'dictionaryEntry'.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
static const dictionary null
Null dictionary.
Definition: dictionary.H:242
const keyType & keyword() const
Return keyword.
Definition: entry.H:123
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
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.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61