regionSizeDistribution.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) 2013-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::functionObjects::regionSizeDistribution
26 
27 Description
28  Creates a size distribution via interrogating a continuous phase fraction
29  field.
30 
31  Looks up a phase-fraction (alpha) field and splits the mesh into regions
32  based on where the field is below the threshold value. These
33  regions ("droplets") can now be analysed.
34 
35  Regions:
36  - print the regions connected to a user-defined set of patches.
37  (in spray calculation these form the liquid core)
38  - print the regions with too large volume. These are the 'background'
39  regions.
40  - (debug) write regions as a volScalarField
41  - (debug) print for all regions the sum of volume and alpha*volume
42 
43  Output (volume scalar) fields include:
44  - alpha_liquidCore : alpha with outside liquid core set to 0
45  - alpha_background : alpha with outside background set to 0.
46 
47  %Histogram:
48  - determine histogram of diameter (given minDiameter, maxDiameter, nBins)
49  - write graph of number of droplets per bin
50  - write graph of sum, average and deviation of droplet volume per bin
51  - write graph of sum, average and deviation of user-defined fields. For
52  volVectorFields these are those of the 3 components and the magnitude.
53 
54  Example of function object specification:
55  \verbatim
56  regionSizeDistribution1
57  {
58  type regionSizeDistribution;
59  libs ("libfieldFunctionObjects.so");
60  ...
61  alpha alpha;
62  patches (inlet);
63  threshold 0.4;
64  fields (p U);
65  nBins 100;
66  maxDiameter 0.5e-4;
67  minDiameter 0;
68  setFormat gnuplot;
69  }
70  \endverbatim
71 
72 Usage
73  \table
74  Property | Description | Required | Default value
75  type | type name: regionSizeDistribution |yes|
76  alpha | phase field to interrogate | yes |
77  patches | patches from which the liquid core is identified | yes|
78  threshold | phase fraction applied to delimit regions | yes |
79  fields | fields to sample | yes |
80  nBins | number of bins for histogram | yes |
81  maxDiameter | maximum region equivalent diameter | yes |
82  minDiameter | minimum region equivalent diameter | no | 0
83  setFormat | writing format | yes |
84  \endtable
85 
86 See also
87  Foam::functionObject
88  Foam::functionObjects::fvMeshFunctionObject
89  Foam::functionObjects::writeFile
90 
91 SourceFiles
92  regionSizeDistribution.C
93 
94 \*---------------------------------------------------------------------------*/
95 
96 #ifndef functionObjects_regionSizeDistribution_H
97 #define functionObjects_regionSizeDistribution_H
98 
99 #include "fvMeshFunctionObject.H"
100 #include "writeFile.H"
101 #include "setWriter.H"
102 #include "Map.H"
103 #include "volFieldsFwd.H"
104 #include "wordReList.H"
105 #include "coordinateSystem.H"
106 
107 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
108 
109 namespace Foam
110 {
111 
112 // Forward declaration of classes
113 class regionSplit;
114 
115 namespace functionObjects
116 {
117 
118 /*---------------------------------------------------------------------------*\
119  Class regionSizeDistribution Declaration
120 \*---------------------------------------------------------------------------*/
121 
122 class regionSizeDistribution
123 :
124  public fvMeshFunctionObject
125 {
126  // Private Data
127 
128  writeFile file_;
129 
130  //- Name of field
131  word alphaName_;
132 
133  //- Patches to walk from
134  wordReList patchNames_;
135 
136  //- Clip value
137  scalar threshold_;
138 
139  //- Maximum droplet diameter
140  scalar maxDiam_;
141 
142  //- Minimum droplet diameter
143  scalar minDiam_;
144 
145  //- Number of bins
146  label nBins_;
147 
148  //- Names of fields to sample on regions
149  wordList fields_;
150 
151  //- Output formatter to write
152  autoPtr<setWriter> formatterPtr_;
153 
154 
155  // Private Member Functions
156 
157  template<class Type>
158  Map<Type> regionSum(const regionSplit&, const Field<Type>&) const;
159 
160  //- Get data in order
161  template<class Type>
162  List<Type> extractData(const UList<label>& keys, const Map<Type>&)
163  const;
164 
165  //- Write volfields with the parts of alpha which are not
166  // droplets (liquidCore, backGround)
167  void writeAlphaFields
168  (
169  const regionSplit& regions,
170  const Map<label>& keepRegions,
171  const Map<scalar>& regionVolume,
172  const volScalarField& alpha
173  ) const;
174 
175  //- Mark all regions starting at patches
176  Map<label> findPatchRegions(const regionSplit&) const;
177 
178  //- Helper: divide if denom != 0
179  template<class Type>
180  static tmp<Field<Type>> divide(const Field<Type>&, const scalarField&);
181 
182  //- Generate fields for writing
183  template<class Type>
184  void generateFields
185  (
186  const word& fieldName, // name of field
187  const labelList& indices, // index of bin for each region
188  const Field<Type>& sortedField, // per region field data
189  const scalarField& binCount, // per bin number of regions
192  ) const;
193 
194  //- Generate fields for writing. Scalar overload. It is possible to
195  // generate more fields for scalar input.
196  void generateFields
197  (
198  const word& fieldName, // name of field
199  const labelList& indices, // index of bin for each region
200  const scalarField& sortedField, // per region field data
201  const scalarField& binCount, // per bin number of regions
204  ) const;
205 
206  //- Generate fields for writing
207  template<class Type>
208  void generateFields
209  (
210  const word& fieldName, // name of field
211  const Field<Type>& cellField, // per cell field data
212  const regionSplit& regions, // per cell the region(=droplet)
213  const labelList& sortedRegions, // valid regions in sorted order
214  const scalarField& sortedNormalisation,
215  const labelList& indices, // index of bin for each region
216  const scalarField& binCount, // per bin number of regions
219  ) const;
220 
221 
222 public:
223 
224  //- Runtime type information
225  TypeName("regionSizeDistribution");
226 
227 
228  // Constructors
229 
230  //- Construct for given objectRegistry and dictionary.
231  // Allow the possibility to load fields from files
233  (
234  const word& name,
235  const Time& runTime,
236  const dictionary&
237  );
238 
239  //- Disallow default bitwise copy construction
241 
242 
243  // Destructor
244 
245  virtual ~regionSizeDistribution();
246 
247 
248  // Member Functions
249 
250  //- Read the regionSizeDistribution data
251  virtual bool read(const dictionary&);
252 
253  //- Return the list of fields required
254  virtual wordList fields() const;
255 
256  //- Do nothing
257  virtual bool execute();
258 
259  //- Calculate the regionSizeDistribution and write
260  virtual bool write();
261 
262 
263  // Member Operators
264 
265  //- Disallow default bitwise assignment
266  void operator=(const regionSizeDistribution&) = delete;
267 };
268 
269 
270 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
271 
272 } // End namespace functionObjects
273 } // End namespace Foam
274 
275 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
276 
277 #ifdef NoRepository
279 #endif
280 
281 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
282 
283 #endif
284 
285 // ************************************************************************* //
void operator=(const regionSizeDistribution &)=delete
Disallow default bitwise assignment.
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Definition: regionSplit.H:116
virtual bool write()
Calculate the regionSizeDistribution and write.
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const word & name() const
Return the name of this functionObject.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Divide a list of fields.
Definition: divide.H:70
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
regionSizeDistribution(const word &name, const Time &runTime, const dictionary &)
Construct for given objectRegistry and dictionary.
TypeName("regionSizeDistribution")
Runtime type information.
Creates a size distribution via interrogating a continuous phase fraction field.
static List< word > fieldNames
Definition: globalFoam.H:46
virtual wordList fields() const
Return the list of fields required.
Pre-declare SubField and related Field type.
Definition: Field.H:56
A class for handling words, derived from string.
Definition: word.H:59
virtual bool read(const dictionary &)
Read the regionSizeDistribution data.
List< word > wordList
A List of words.
Definition: fileName.H:54
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
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.
A HashTable to objects of type <T> with a label key.
Definition: Map.H:49