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-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::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  //- Label of this sizeGroup within the corresponding populationBalance
91  const label i_;
92 
93  //- VelocityGroup this sizeGroup belongs to
94  const velocityGroup& group_;
95 
96  //- Sphere equivalent diameter of the sizeGroup
97  const dimensionedScalar dSph_;
98 
99  //- Representative volume of the sizeGroup
100  const dimensionedScalar x_;
101 
102  //- Model for describing the representative shape of the elements in the
103  // sizeGroup
104  autoPtr<shapeModel> shapeModel_;
105 
106 
107 public:
108 
109  // Static Member Functions
110 
111  //- Return IO for a size-group field. Construction helper.
112  static IOobject fieldIo
113  (
114  const word& name,
115  const label i,
116  const velocityGroup& group,
118  const bool registerObject = true
119  );
120 
121  //- Return IO a size-group field. Construction helper.
123  (
124  const word& name,
125  const label i,
126  const velocityGroup& group
127  );
128 
129 
130  // Constructors
131 
132  //- Construct from index, dictionary and velocity group
133  sizeGroup
134  (
135  const label i,
136  const dictionary& dict,
137  const velocityGroup& group
138  );
139 
140  //- Return clone
141  autoPtr<sizeGroup> clone() const;
142 
143  //- Return a pointer to a new sizeGroup created from Istream
144  class iNew
145  {
146  const velocityGroup& group_;
147 
148  mutable label i_;
149 
150 
151  public:
152 
153  iNew(const velocityGroup& group, const label i0)
154  :
155  group_(group),
156  i_(i0)
157  {}
158 
160  {
161  token t(is);
162 
163  // Check and filter for old syntax (remove in due course)
164  if (t.isWord())
165  {
166  const word fName = "f" + Foam::name(i_);
167 
168  if (t.wordToken() != fName)
169  {
171  << "The name '" << t.wordToken() << "' should not "
172  << "have been provided for the fraction field "
173  << IOobject::groupName(fName, group_.phase().name())
174  << " of size-group #" << i_
175  << " of population balance " << group_.popBalName()
176  << ". Only the coefficients are required."
177  << exit(FatalError);
178  }
179  }
180  else
181  {
182  is.putBack(t);
183  }
184 
185  autoPtr<sizeGroup> result(new sizeGroup(i_, is, group_));
186 
187  i_ ++;
188 
189  return result;
190  }
191  };
192 
193 
194  //- Destructor
195  virtual ~sizeGroup();
196 
197 
198  // Member Functions
199 
200  //- Return index of the size group within the population balance
201  label i() const;
202 
203  //- Return const-reference to the phase
204  inline const phaseModel& phase() const;
205 
206  //- Return const-reference to the velocityGroup
207  inline const velocityGroup& group() const;
208 
209  //- Return representative spherical diameter of the sizeGroup
210  inline const dimensionedScalar& dSph() const;
211 
212  //- Return representative volume of the sizeGroup
213  inline const dimensionedScalar& x() const;
214 
215  //- Return reference to diameterModel of the phase
216  inline autoPtr<shapeModel>& shapeModelPtr();
217 
218  //- Return const-reference to diameterModel of the phase
219  inline const autoPtr<shapeModel>& shapeModelPtr() const;
220 
221  //- Return representative surface area of the sizeGroup
222  const tmp<volScalarField> a() const;
223 
224  //- Return representative diameter of the sizeGroup
225  const tmp<volScalarField> d() const;
226 
227  //- Correct secondary properties
228  void correct();
229 };
230 
231 
232 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
233 
234 } // End namespace diameterModels
235 } // End namespace Foam
236 
237 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
238 
239 #include "sizeGroupI.H"
240 
241 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
242 
243 #endif
244 
245 // ************************************************************************* //
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:346
const word & name() const
Return name.
Definition: IOobject.H:310
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.
Return a pointer to a new sizeGroup created from Istream.
Definition: sizeGroup.H:162
iNew(const velocityGroup &group, const label i0)
Definition: sizeGroup.H:170
autoPtr< sizeGroup > operator()(Istream &is) const
Definition: sizeGroup.H:176
const dimensionedScalar & dSph() const
Return representative spherical diameter of the sizeGroup.
Definition: sizeGroupI.H:51
void correct()
Correct secondary properties.
Definition: sizeGroup.C:148
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:135
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:37
autoPtr< shapeModel > & shapeModelPtr()
Return reference to diameterModel of the phase.
Definition: sizeGroupI.H:65
autoPtr< sizeGroup > clone() const
Return clone.
Definition: sizeGroup.C:127
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:142
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 keyword definitions, which are a keyword followed by any number of values (e....
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 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
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
error FatalError
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict