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-2020 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 functionObjects_regionSizeDistribution_H
173 #define functionObjects_regionSizeDistribution_H
174 
175 #include "fvMeshFunctionObject.H"
176 #include "writeFile.H"
177 #include "writer.H"
178 #include "Map.H"
179 #include "volFieldsFwd.H"
180 #include "wordRes.H"
181 #include "coordinateSystem.H"
182 #include "Switch.H"
183 
184 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185 
186 namespace Foam
187 {
188 
189 class regionSplit;
190 
191 namespace functionObjects
192 {
193 
194 /*---------------------------------------------------------------------------*\
195  Class regionSizeDistribution Declaration
196 \*---------------------------------------------------------------------------*/
197 
198 class regionSizeDistribution
199 :
200  public fvMeshFunctionObject,
201  public writeFile
202 {
203  // Private Data
204 
205  //- Number of bins
206  label nBins_;
207 
208  //- Name of field
209  word alphaName_;
210 
211  //- Clip value
212  scalar threshold_;
213 
214  //- Maximum droplet diameter
215  scalar maxDiam_;
216 
217  //- Minimum droplet diameter
218  scalar minDiam_;
219 
220  //- Patches to walk from
221  wordRes patchNames_;
222 
223  //- Names of fields to sample on regions
224  wordRes fields_;
225 
226  //- Output formatter to write
227  autoPtr<writer<scalar>> formatterPtr_;
228 
229  //- Optional coordinate system
230  autoPtr<coordinateSystem> csysPtr_;
231 
232  // Optional extra definition of bins on planes downstream to the origin
233  // point and maximum diameter
234 
235  //- Optional origin point
236  vector origin_;
237 
238  //- Optional plane normal direction
239  vector direction_;
240 
241  //- Optional maximum diameter on plane
242  scalar maxDiameter_;
243 
244  //- Optional number of bins for
245  scalar nDownstreamBins_;
246 
247  //- Optional maximum downstream coordinate from origin
248  scalar maxDownstream_;
249 
250  //- Switch to enable iso-planes sampling
251  bool isoPlanes_;
252 
253 
254  // Private Member Functions
255 
256  template<class Type>
257  Map<Type> regionSum(const regionSplit&, const Field<Type>&) const;
258 
259  //- Get data in order
260  template<class Type>
261  List<Type> extractData(const labelUList& keys, const Map<Type>&) const;
262 
263  void writeGraph
264  (
265  const coordSet& coords,
266  const word& valueName,
267  const scalarField& values
268  ) const;
269 
270  //- Write volfields with the parts of alpha which are not
271  //- droplets (liquidCore, backGround)
272  void writeAlphaFields
273  (
274  const regionSplit& regions,
275  const Map<label>& keepRegions,
276  const Map<scalar>& regionVolume,
277  const volScalarField& alpha
278  ) const;
279 
280  //- Mark all regions starting at patches
281  Map<label> findPatchRegions(const regionSplit&) const;
282 
283  //- Helper: divide if denom != 0
284  static tmp<scalarField> divide(const scalarField&, const scalarField&);
285 
286  //- Given per-region data calculate per-bin average/deviation and graph
287  void writeGraphs
288  (
289  const word& fieldName, // name of field
290  const labelList& indices, // index of bin for each region
291  const scalarField& sortedField, // per region field data
292  const scalarField& binCount, // per bin number of regions
293  const coordSet& coords // graph data for bins
294  ) const;
295 
296  //- Given per-cell data calculate per-bin average/deviation and graph
297  void writeGraphs
298  (
299  const word& fieldName, // name of field
300  const scalarField& cellField, // per cell field data
301  const regionSplit& regions, // per cell the region(=droplet)
302  const labelList& sortedRegions, // valid regions in sorted order
303  const scalarField& sortedNormalisation,
304  const labelList& indices, // index of bin for each region
305  const scalarField& binCount, // per bin number of regions
306  const coordSet& coords // graph data for bins
307  ) const;
308 
309 
310 public:
311 
312  //- Runtime type information
313  TypeName("regionSizeDistribution");
314 
315 
316  // Constructors
317 
318  //- Construct for given objectRegistry and dictionary.
319  // Allow the possibility to load fields from files
321  (
322  const word& name,
323  const Time& runTime,
324  const dictionary&
325  );
326 
327  //- No copy construct
329 
330  //- No copy assignment
331  void operator=(const regionSizeDistribution&) = delete;
332 
333 
334  // Destructor
335  virtual ~regionSizeDistribution() = default;
336 
337 
338  // Member Functions
339 
340  //- Read the regionSizeDistribution data
341  virtual bool read(const dictionary&);
342 
343  //- Do nothing
344  virtual bool execute();
345 
346  //- Calculate the regionSizeDistribution and write
347  virtual bool write();
348 };
349 
350 
351 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
352 
353 } // End namespace functionObjects
354 } // End namespace Foam
355 
356 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
357 
358 #ifdef NoRepository
360 #endif
361 
362 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
363 
364 #endif
365 
366 // ************************************************************************* //
Foam::functionObjects::regionSizeDistribution::write
virtual bool write()
Calculate the regionSizeDistribution and write.
Definition: regionSizeDistribution.C:384
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
runTime
engineTime & runTime
Definition: createEngineTime.H:13
volFieldsFwd.H
wordRes.H
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
writeFile.H
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
fvMeshFunctionObject.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::HashTableOps::values
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:149
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::Map
A HashTable to objects of type <T> with a label key.
Definition: lumpedPointController.H:69
Foam::functionObjects::regionSizeDistribution::TypeName
TypeName("regionSizeDistribution")
Runtime type information.
Foam::functionObjects::regionSizeDistribution::regionSizeDistribution
regionSizeDistribution(const word &name, const Time &runTime, const dictionary &)
Construct for given objectRegistry and dictionary.
Definition: regionSizeDistribution.C:314
Foam::functionObjects::fvMeshFunctionObject
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Definition: fvMeshFunctionObject.H:64
coordinateSystem.H
Foam::functionObjects::regionSizeDistribution::read
virtual bool read(const dictionary &)
Read the regionSizeDistribution data.
Definition: regionSizeDistribution.C:332
Map.H
Foam::functionObjects::regionSizeDistribution
Creates a droplet size distribution via interrogating a continuous phase fraction field.
Definition: regionSizeDistribution.H:331
Foam::Field
Generic templated field type.
Definition: Field.H:63
regionSizeDistributionTemplates.C
Switch.H
Foam::regionSplit
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Definition: regionSplit.H:140
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::coordSet
Holds list of sampling positions.
Definition: coordSet.H:53
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::functionObject::name
const word & name() const noexcept
Return the name of this functionObject.
Definition: functionObject.C:143
Foam::Vector< scalar >
Foam::List< Type >
Foam::UList< label >
Foam::wordRes
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:51
Foam::functionObjects::regionSizeDistribution::execute
virtual bool execute()
Do nothing.
Definition: regionSizeDistribution.C:378
Foam::functionObjects::writeFile
Base class for writing single files from the function objects.
Definition: writeFile.H:119
Foam::labelUList
UList< label > labelUList
A UList of labels.
Definition: UList.H:85
writer.H
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::functionObjects::regionSizeDistribution::operator=
void operator=(const regionSizeDistribution &)=delete
No copy assignment.
Foam::functionObjects::regionSizeDistribution::~regionSizeDistribution
virtual ~regionSizeDistribution()=default