regionSizeDistribution.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2013-2016 OpenFOAM Foundation
9  Copyright (C) 2016-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::functionObjects::regionSizeDistribution
29 
30 Group
31  grpFieldFunctionObjects
32 
33 Description
34  Creates a droplet size distribution via interrogating a continuous phase
35  fraction field.
36 
37  Looks up a phase-fraction (alpha) field and splits the mesh into regions
38  based on where the field is below the threshold value. These
39  regions ("droplets") can now be analysed.
40 
41  Regions:
42  - print the regions connected to a user-defined set of patches.
43  (in spray calculation these form the liquid core)
44  - print the regions with too large volume. These are the 'background'
45  regions.
46  - (debug) write regions as a volScalarField
47  - (debug) print for all regions the sum of volume and alpha*volume
48 
49  Output (volume scalar) fields include:
50  - alpha_liquidCore : alpha with outside liquid core set to 0
51  - alpha_background : alpha with outside background set to 0.
52 
53  Histogram:
54  - determine histogram of diameter (given minDiameter, maxDiameter, nBins)
55  - write graph of number of droplets per bin
56  - write graph of sum, average and deviation of droplet volume per bin
57  - write graph of sum, average and deviation of user-defined fields. For
58  volVectorFields these are those of the 3 components and the magnitude.
59  - (optional) write graph of histogram of centroids on iso planes
60  downstream of the injector determined by origin, direction and maxDiameter
61  up to maxDownstream
62 
63  Operands:
64  \table
65  Operand | Type | Location
66  input | - | -
67  output file | - <!--
68  --> | $FOAM_CASE/postProcessing/<FO>/<time>/<files>
69  output field | - | -
70  \endtable
71 
72 Usage
73  Minimal example by using \c system/controlDict.functions:
74  \verbatim
75  regionSizeDistribution1
76  {
77  // Mandatory entries (unmodifiable)
78  type regionSizeDistribution;
79  libs (fieldFunctionObjects);
80 
81  // Mandatory entries (runtime modifiable)
82  field alpha;
83  patches (inlet);
84  fields (p U);
85  threshold 0.4;
86  maxDiameter 5e-5;
87  nBins 100;
88  setFormat gnuplot;
89 
90  // Optional entries (runtime modifiable)
91  minDiameter 0.0;
92  coordinateSystem
93  {
94  origin (0 0 0);
95  rotation
96  {
97  type axes;
98  e3 (0 0 1);
99  e1 (1 0 0);
100  }
101  }
102 
103  // Optional downstream iso-plane bins
104  isoPlanes true;
105 
106  // Mandatory entries if isoPlanes=true (runtime modifiable)
107  // Plane normal and point definition
108  origin (1e-4 0 5e-4);
109  direction (1 0 1);
110 
111  // Maximum diameter of the cylinder formed by the origin point
112  // and direction
113  maxD 3e-4;
114 
115  // Maximum downstream distance
116  maxDownstream 6e-4;
117 
118  // Number of iso-plane bins
119  nDownstreamBins 20;
120 
121  // Optional (inherited) entries
122  ...
123  }
124  \endverbatim
125 
126  where the entries mean:
127  \table
128  Property | Description | Type | Req'd | Dflt
129  type | Type name: regionSizeDistribution | word | yes | -
130  libs | Library name: fieldFunctionObjects | word | yes | -
131  field | Phase field to interrogate | word | yes | -
132  patches | Patches wherefrom the liquid core is identified <!--
133  --> | wordList | yes | -
134  fields | Fields to sample | wordList | yes | -
135  threshold | Phase fraction applied to delimit regions | scalar | yes | -
136  maxDiameter | Maximum droplet diameter | scalar | yes | -
137  minDiameter | Minimum droplet diameter | scalar | no | 0.0
138  nBins | Number of bins for histogram | label | yes | -
139  setFormat | Output format | word | yes | -
140  isoPlanes | Flag for isoPlanes | bool | no | false
141  origin | Origin of the plane when isoPlanes is used <!--
142  --> | vector | yes | -
143  direction <!--
144  --> | Direction of the plane when isoPlanes is used | vector | yes | -
145  maxD | Maximum diameter of the sampling <!--
146  --> cylinder when isoPlanes is used | vector | yes | -
147  nDownstreamBins <!--
148  --> | Number of bins when isoPlanes is used | label | yes | -
149  maxDownstream <!--
150  --> | Maximum distance from origin when isoPlanes is used <!--
151  --> | scalar | yes | -
152  \endtable
153 
154  The inherited entries are elaborated in:
155  - \link functionObject.H \endlink
156  - \link writeFile.H \endlink
157 
158  Usage by the \c postProcess utility is not available.
159 
160 See also
161  - Foam::functionObject
162  - Foam::functionObjects::fvMeshFunctionObject
163  - Foam::functionObjects::writeFile
164  - ExtendedCodeGuide::functionObjects::field::regionSizeDistribution
165 
166 SourceFiles
167  regionSizeDistribution.C
168  regionSizeDistributionTemplates.C
169 
170 \*---------------------------------------------------------------------------*/
171 
172 #ifndef Foam_functionObjects_regionSizeDistribution_H
173 #define Foam_functionObjects_regionSizeDistribution_H
174 
175 #include "fvMeshFunctionObject.H"
176 #include "writeFile.H"
177 #include "coordSetWriter.H"
178 #include "Map.H"
179 #include "volFieldsFwd.H"
180 #include "wordRes.H"
181 #include "coordinateSystem.H"
182 
183 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
184 
185 namespace Foam
186 {
187 
188 class regionSplit;
189 
190 namespace functionObjects
191 {
192 
193 /*---------------------------------------------------------------------------*\
194  Class regionSizeDistribution Declaration
195 \*---------------------------------------------------------------------------*/
196 
197 class regionSizeDistribution
198 :
199  public functionObjects::fvMeshFunctionObject,
200  public writeFile
201 {
202  // Private Data
203 
204  //- Number of bins
205  label nBins_;
206 
207  //- Name of field
208  word alphaName_;
209 
210  //- Clip value
211  scalar threshold_;
212 
213  //- Maximum droplet diameter
214  scalar maxDiam_;
215 
216  //- Minimum droplet diameter
217  scalar minDiam_;
218 
219  //- Patches to walk from
220  wordRes patchNames_;
221 
222  //- Names of fields to sample on regions
223  wordRes fields_;
224 
225  //- Output formatter to write
226  mutable autoPtr<coordSetWriter> formatterPtr_;
227 
228  //- Optional coordinate system
229  autoPtr<coordinateSystem> csysPtr_;
230 
231  // Optional extra definition of bins on planes downstream to the origin
232  // point and maximum diameter
233 
234  //- Optional origin point
235  vector origin_;
236 
237  //- Optional plane normal direction
238  vector direction_;
239 
240  //- Optional maximum diameter on plane
241  scalar maxDiameter_;
242 
243  //- Optional number of bins for
244  scalar nDownstreamBins_;
245 
246  //- Optional maximum downstream coordinate from origin
247  scalar maxDownstream_;
248 
249  //- Switch to enable iso-planes sampling
250  bool isoPlanes_;
251 
252 
253  // Private Member Functions
254 
255  //- Write volfields with the parts of alpha which are not
256  //- droplets (liquidCore, backGround)
257  void writeAlphaFields
258  (
259  const regionSplit& regions,
260  const labelHashSet& keepRegions,
261  const Map<scalar>& regionVolume,
262  const volScalarField& alpha
263  ) const;
264 
265  //- Mark all regions starting at patches
266  labelHashSet findPatchRegions(const regionSplit&) const;
267 
268  //- Helper: divide if denom != 0
269  static tmp<scalarField> divide(const scalarField&, const scalarField&);
270 
271  //- Given per-region data calculate per-bin average/deviation and graph
272  void writeGraphs
273  (
274  const word& fieldName, // name of field
275  const scalarField& sortedField, // per region field data
276 
277  const labelList& indices, // index of bin for each region
278  const scalarField& binCount, // per bin number of regions
279  const coordSet& coords // graph data for bins
280  ) const;
281 
282  //- Given per-cell data calculate per-bin average/deviation and graph
283  void writeGraphs
284  (
285  const word& fieldName, // name of field
286  const scalarField& cellField, // per cell field data
287  const regionSplit& regions, // per cell the region(=droplet)
288  const labelList& sortedRegions, // valid regions in sorted order
289  const scalarField& sortedNormalisation,
290 
291  const labelList& indices, // index of bin for each region
292  const scalarField& binCount, // per bin number of regions
293  const coordSet& coords // graph data for bins
294  ) const;
295 
296 
297 public:
298 
299  //- Runtime type information
300  TypeName("regionSizeDistribution");
301 
302 
303  // Constructors
304 
305  //- Construct for given objectRegistry and dictionary.
306  // Allow the possibility to load fields from files
308  (
309  const word& name,
310  const Time& runTime,
311  const dictionary&
312  );
313 
314  //- No copy construct
316 
317  //- No copy assignment
318  void operator=(const regionSizeDistribution&) = delete;
319 
320 
321  // Destructor
322  virtual ~regionSizeDistribution() = default;
323 
324 
325  // Member Functions
327  //- Read the regionSizeDistribution data
328  virtual bool read(const dictionary&);
329 
330  //- Do nothing
331  virtual bool execute();
332 
333  //- Calculate the regionSizeDistribution and write
334  virtual bool write();
335 };
336 
337 
338 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
339 
340 } // End namespace functionObjects
341 } // End namespace Foam
342 
343 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
344 
345 #endif
346 
347 // ************************************************************************* //
void operator=(const regionSizeDistribution &)=delete
No copy assignment.
virtual bool write()
Calculate the regionSizeDistribution and write.
Forwards and collection of common volume field types.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
engineTime & runTime
regionSizeDistribution(const word &name, const Time &runTime, const dictionary &)
Construct for given objectRegistry and dictionary.
TypeName("regionSizeDistribution")
Runtime type information.
const word & name() const noexcept
Return the name of this functionObject.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Vector< scalar > vector
Definition: vector.H:57
virtual bool read(const dictionary &)
Read the regionSizeDistribution data.
List< label > labelList
A List of labels.
Definition: List.H:62
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Namespace for OpenFOAM.