populationBalanceSizeDistribution.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::functionObjects::populationBalanceSizeDistribution
26 
27 Description
28  Writes out the size distribution determined by a population balance model,
29  either for the entire domain or a cell zone. Requires solver post-
30  processing.
31 
32  The following function object specification for example returns the volume-
33  based number density function:
34 
35  Example of function object specification:
36  \verbatim
37  numberDensity
38  {
39  type populationBalanceSizeDistribution;
40  libs ("libmultiphaseEulerFunctionObjects.so");
41  writeControl writeTime;
42  populationBalance bubbles;
43  functionType numberDensity;
44  coordinateType volume;
45  setFormat raw;
46  }
47  \endverbatim
48 
49 Usage
50  \table
51  Property | Description | Required | Default
52  populationBalance | population balance name | yes |
53  functionType | function type | yes |
54  coordinateType | particle property | yes |
55  allCoordinates | write all coordinate values | no | false
56  normalise | divide by total concentration | no | false
57  logTransform | class width based on log of coordinate\\
58  | no | false
59  weightType | weighting in case of field-dependent particle\\
60  properties | no\\
61  | numberConcentration
62  cellZone | cellZone | yes |
63  name | name of cellZone if required | no |
64  setFormat | output format | yes |
65  \endtable
66 
67 See also
68  Foam::diameterModels::populationBalanceModel
69  Foam::functionObjects::fvMeshFunctionObject
70  Foam::fvCellZone
71  Foam::functionObject
72 
73 SourceFiles
74  populationBalanceSizeDistribution.C
75 
76 \*---------------------------------------------------------------------------*/
77 
78 #ifndef populationBalanceSizeDistribution_H
79 #define populationBalanceSizeDistribution_H
80 
81 #include "fvMeshFunctionObject.H"
82 #include "fvCellZone.H"
83 #include "populationBalanceModel.H"
84 #include "writeFile.H"
85 #include "setWriter.H"
86 
87 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
88 
89 namespace Foam
90 {
91 namespace functionObjects
92 {
93 
94 /*---------------------------------------------------------------------------*\
95  Class populationBalanceSizeDistribution Declaration
96 \*---------------------------------------------------------------------------*/
97 
98 class populationBalanceSizeDistribution
99 :
100  public fvMeshFunctionObject
101 {
102 public:
103 
104  // Public Data Types
105 
106  //- Function type enumeration
107  enum functionType
108  {
115  };
116 
117  //- Function type names
118  static const NamedEnum<functionType, 6> functionTypeNames_;
119 
120  //- Coordinate type enumeration
121  enum coordinateType
122  {
123  volume,
124  area,
125  diameter,
127  };
128 
129  //- Coordinate type names
130  static const NamedEnum<coordinateType, 4> coordinateTypeNames_;
131 
132  //- Enumeration for the weight types
133  enum class weightType
134  {
138  cellVolume
139  };
140 
141  //- Names of the weight types
142  static const NamedEnum<weightType, 4> weightTypeNames_;
143 
144 
145 private:
146 
147  // Private Data
148 
149  //- Write file
150  writeFile file_;
151 
152  //- Reference to mesh
153  const fvMesh& mesh_;
154 
155  //- cellZone
156  fvCellZone zone_;
157 
158  //- Name of the population balance
159  const word popBalName_;
160 
161  //- Function type
162  functionType functionType_;
163 
164  //- Coordinate type
165  coordinateType coordinateType_;
166 
167  //- Add values for all coordinate types to output
168  Switch allCoordinates_;
169 
170  //- Normalise result through division by sum
171  Switch normalise_;
172 
173  //- Log transform
174  Switch logTransform_;
175 
176  //- Weight types, relevant if particle properties are field dependent
177  weightType weightType_;
178 
179  //- Set formatter
180  autoPtr<setWriter> formatterPtr_;
181 
182 
183  // Private Member Functions
184 
185  //- Function type symbolic name for shorter file header
186  word functionTypeSymbolicName();
187 
188  //- Coordinate type symbolic name for shorter file header
189  word coordinateTypeSymbolicName(const coordinateType& cType);
190 
191  //- Filter a field according to cellIds
192  tmp<scalarField> filterField(const scalarField& field) const;
193 
194  //- Field averaged coordinate value
195  scalar averageCoordinateValue
196  (
197  const diameterModels::sizeGroup& fi,
198  const coordinateType& cType
199  );
200 
201  //- Weighted average
202  scalar weightedAverage
203  (
204  const scalarField& fld,
205  const diameterModels::sizeGroup& fi
206  );
207 
208 
209 public:
210 
211  //- Runtime type information
212  TypeName("populationBalanceSizeDistribution");
213 
214 
215  // Constructors
216 
217  //- Construct from Time and dictionary
219  (
220  const word& name,
221  const Time& runTime,
222  const dictionary& dict
223  );
224 
225  //- Disallow default bitwise copy construction
227  (
229  ) = delete;
230 
231 
232  //- Destructor
234 
235 
236  // Member Functions
237 
238  //- Return the list of fields required
239  virtual wordList fields() const
240  {
241  return wordList::null();
242  }
243 
244  //- Read the populationBalanceSizeDistribution data
245  virtual bool read(const dictionary&);
246 
247  //- Execute, currently does nothing
248  virtual bool execute();
249 
250  //- Calculate and write the size distribution
251  virtual bool write();
252 
253 
254  // Mesh changes
255 
256  //- Update for mesh motion
257  virtual void movePoints(const polyMesh&);
258 
259  //- Update topology using the given map
260  virtual void topoChange(const polyTopoChangeMap&);
261 
262  //- Update from another mesh using the given map
263  virtual void mapMesh(const polyMeshMap&);
264 
265  //- Redistribute or update using the given distribution map
266  virtual void distribute(const polyDistributionMap&);
267 
268 
269  // Member Operators
270 
271  //- Disallow default bitwise assignment
272  void operator=(const populationBalanceSizeDistribution&) = delete;
273 };
274 
275 
276 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
277 
278 } // End namespace functionObjects
279 } // End namespace Foam
280 
281 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
282 
283 #endif
284 
285 // ************************************************************************* //
static const List< word > & null()
Return a null List.
Definition: ListI.H:118
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Single size class fraction field representing a fixed particle volume as defined by the user through ...
Definition: sizeGroup.H:96
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 functionObject.
Writes out the size distribution determined by a population balance model, either for the entire doma...
virtual wordList fields() const
Return the list of fields required.
TypeName("populationBalanceSizeDistribution")
Runtime type information.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
void operator=(const populationBalanceSizeDistribution &)=delete
Disallow default bitwise assignment.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
populationBalanceSizeDistribution(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
virtual void movePoints(const polyMesh &)
Update for mesh motion.
static const NamedEnum< weightType, 4 > weightTypeNames_
Names of the weight types.
static const NamedEnum< functionType, 6 > functionTypeNames_
Function type names.
static const NamedEnum< coordinateType, 4 > coordinateTypeNames_
Coordinate type names.
virtual bool write()
Calculate and write the size distribution.
virtual bool read(const dictionary &)
Read the populationBalanceSizeDistribution data.
functionObject base class for writing single files
Definition: writeFile.H:56
cellZone selection or generation class with caching of zone volume
Definition: fvCellZone.H:94
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter) { positions[particlei++]=iter().position(mesh);} for(i=0;i< pTraits< point >::nComponents;i++) { forAll(positions, particlei) { gmvFile<< component(positions[particlei], i)<< ' ';} gmvFile<< nl;}}forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Namespace for OpenFOAM.
dictionary dict