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