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-2025 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  shapeModel | Shape model providing an alternative diameter field\\
48  | yes | none
49  \endtable
50 
51 See also
52  Foam::diameterModels::populationBalanceModel
53  Foam::diameterModels::velocityGroup
54  Foam::diameterModels::shapeModel
55 
56 SourceFiles
57  sizeGroup.C
58 
59 \*---------------------------------------------------------------------------*/
60 
61 #ifndef sizeGroup_H
62 #define sizeGroup_H
63 
64 #include "velocityGroup.H"
65 
66 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67 
68 namespace Foam
69 {
70 namespace diameterModels
71 {
72 
73 // Forward declaration of classes
74 class shapeModel;
75 
76 /*---------------------------------------------------------------------------*\
77  Class sizeGroup Declaration
78 \*---------------------------------------------------------------------------*/
79 
80 class sizeGroup
81 :
82  public volScalarField
83 {
84 private:
85 
86  // Private Data
87 
88  //- Label of this sizeGroup within the corresponding populationBalance
89  const label i_;
90 
91  //- VelocityGroup this sizeGroup belongs to
92  const velocityGroup& group_;
93 
94  //- Sphere equivalent diameter of the sizeGroup
95  const dimensionedScalar dSph_;
96 
97  //- Representative volume of the sizeGroup
98  const dimensionedScalar x_;
99 
100  //- Model for describing the representative shape of the elements in the
101  // sizeGroup
102  autoPtr<shapeModel> shapeModel_;
103 
104 
105 public:
106 
107  // Static Member Functions
108 
109  //- Return IO for a size-group field. Construction helper.
110  static IOobject fieldIo
111  (
112  const word& name,
113  const label i,
114  const velocityGroup& group,
116  const bool registerObject = true
117  );
118 
119  //- Return IO a size-group field. Construction helper.
121  (
122  const word& name,
123  const label i,
124  const velocityGroup& group
125  );
126 
127 
128  // Constructors
129 
130  //- Construct from index, dictionary and velocity group
131  sizeGroup
132  (
133  const label i,
134  const dictionary& dict,
135  const velocityGroup& group
136  );
137 
138  //- Return clone
139  autoPtr<sizeGroup> clone() const
140  {
142  return autoPtr<sizeGroup>(nullptr);
143  }
144 
145  //- Return a pointer to a new sizeGroup created from Istream
146  class iNew
147  {
148  const velocityGroup& group_;
149 
150  mutable label i_;
151 
152 
153  public:
154 
155  iNew(const velocityGroup& group, const label i0)
156  :
157  group_(group),
158  i_(i0)
159  {}
160 
162  {
163  token t(is);
164 
165  // Check and filter for old syntax (remove in due course)
166  if (t.isWord())
167  {
168  const word fName = "f" + Foam::name(i_);
169 
170  if (t.wordToken() != fName)
171  {
173  << "The name '" << t.wordToken() << "' should not "
174  << "have been provided for the fraction field "
175  << IOobject::groupName(fName, group_.phase().name())
176  << " of size-group #" << i_
177  << " of population balance " << group_.popBalName()
178  << ". Only the coefficients are required."
179  << exit(FatalError);
180  }
181  }
182  else
183  {
184  is.putBack(t);
185  }
186 
187  autoPtr<sizeGroup> result(new sizeGroup(i_, is, group_));
188 
189  i_ ++;
190 
191  return result;
192  }
193  };
194 
195 
196  //- Destructor
197  virtual ~sizeGroup();
198 
199 
200  // Member Functions
201 
202  //- Return index of the size group within the population balance
203  inline label i() const;
204 
205  //- Return const-reference to the phase
206  inline const phaseModel& phase() const;
207 
208  //- Return const-reference to the velocityGroup
209  inline const velocityGroup& group() const;
210 
211  //- Return representative spherical diameter of the sizeGroup
212  inline const dimensionedScalar& dSph() const;
213 
214  //- Return representative volume of the sizeGroup
215  inline const dimensionedScalar& x() const;
216 
217  //- Return reference to the shape model
218  shapeModel& shape();
219 
220  //- Return const-reference to the shape model
221  const shapeModel& shape() const;
222 
223  //- Return representative surface area of the sizeGroup
224  const tmp<volScalarField> a() const;
225 
226  //- Return representative diameter of the sizeGroup
227  const tmp<volScalarField> d() const;
228 
229  //- Correct secondary properties
230  void correct();
231 };
232 
233 
234 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
235 
236 } // End namespace diameterModels
237 } // End namespace Foam
238 
239 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
240 
241 #include "sizeGroupI.H"
242 
243 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
244 
245 #endif
246 
247 // ************************************************************************* //
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
readOption
Enumeration defining the read options.
Definition: IOobject.H:117
bool & registerObject()
Register object created from this IOobject with registry if true.
Definition: IOobject.H:343
const word & name() const
Return name.
Definition: IOobject.H:307
static word groupName(Name name, const word &group)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
void putBack(const token &)
Put back token.
Definition: Istream.C:30
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
const phaseModel & phase() const
Return the phase.
iNew(const velocityGroup &group, const label i0)
Definition: sizeGroup.H:168
autoPtr< sizeGroup > operator()(Istream &is) const
Definition: sizeGroup.H:174
autoPtr< sizeGroup > clone() const
Return clone.
Definition: sizeGroup.H:152
const dimensionedScalar & dSph() const
Return representative spherical diameter of the sizeGroup.
Definition: sizeGroupI.H:51
void correct()
Correct secondary properties.
Definition: sizeGroup.C:153
shapeModel & shape()
Return reference to the shape model.
Definition: sizeGroup.C:127
const dimensionedScalar & x() const
Return representative volume of the sizeGroup.
Definition: sizeGroupI.H:58
static tmp< volScalarField > field(const word &name, const label i, const velocityGroup &group)
Return IO a size-group field. Construction helper.
Definition: sizeGroup.C:61
const velocityGroup & group() const
Return const-reference to the velocityGroup.
Definition: sizeGroupI.H:44
static IOobject fieldIo(const word &name, const label i, const velocityGroup &group, const IOobject::readOption r=IOobject::NO_READ, const bool registerObject=true)
Return IO for a size-group field. Construction helper.
Definition: sizeGroup.C:35
virtual ~sizeGroup()
Destructor.
Definition: sizeGroup.C:120
const tmp< volScalarField > a() const
Return representative surface area of the sizeGroup.
Definition: sizeGroup.C:140
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:37
sizeGroup(const label i, const dictionary &dict, const velocityGroup &group)
Construct from index, dictionary and velocity group.
Definition: sizeGroup.C:89
label i() const
Return index of the size group within the population balance.
Definition: sizeGroupI.H:30
const tmp< volScalarField > d() const
Return representative diameter of the sizeGroup.
Definition: sizeGroup.C:147
Computes the Sauter mean diameter based on a user specified size distribution, defined in terms of si...
Definition: velocityGroup.H:87
const word & popBalName() const
Return name of populationBalance this velocityGroup belongs to.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:145
A class for managing temporary objects.
Definition: tmp.H:55
A token holds items read from Istream.
Definition: token.H:73
A class for handling words, derived from string.
Definition: word.H:62
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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:62
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict