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  This function object creates a size distribution via interrogating a
32  continuous phase fraction 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::writeFile
100 
101 SourceFiles
102  regionSizeDistribution.C
103 
104 \*---------------------------------------------------------------------------*/
105 
106 #ifndef functionObjects_regionSizeDistribution_H
107 #define functionObjects_regionSizeDistribution_H
108 
109 #include "writeFile.H"
110 #include "writer.H"
111 #include "Map.H"
112 #include "volFieldsFwd.H"
113 #include "wordReList.H"
114 #include "coordinateSystem.H"
115 
116 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
117 
118 namespace Foam
119 {
120 
121 // Forward declaration of classes
122 class regionSplit;
123 
124 namespace functionObjects
125 {
126 
127 /*---------------------------------------------------------------------------*\
128  Class regionSizeDistribution Declaration
129 \*---------------------------------------------------------------------------*/
130 
131 class regionSizeDistribution
132 :
133  public writeFile
134 {
135  // Private data
136 
137  //- Name of field
138  word alphaName_;
139 
140  //- Patches to walk from
141  wordReList patchNames_;
142 
143  //- Clip value
144  scalar threshold_;
145 
146  //- Maximum droplet diameter
147  scalar maxDiam_;
148 
149  //- Minimum droplet diameter
150  scalar minDiam_;
151 
152  //- Mumber of bins
153  label nBins_;
154 
155  //- Names of fields to sample on regions
156  wordReList fields_;
157 
158  //- Output formatter to write
159  autoPtr<writer<scalar>> formatterPtr_;
160 
161  //- Optional coordinate system
162  autoPtr<coordinateSystem> coordSysPtr_;
163 
164 
165  // Private Member Functions
166 
167  template<class Type>
168  Map<Type> regionSum(const regionSplit&, const Field<Type>&) const;
169 
170  //- Get data in order
171  template<class Type>
172  List<Type> extractData(const UList<label>& keys, const Map<Type>&)
173  const;
174 
175  void writeGraph
176  (
177  const coordSet& coords,
178  const word& valueName,
179  const scalarField& values
180  ) const;
181 
182  //- Write volfields with the parts of alpha which are not
183  // droplets (liquidCore, backGround)
184  void writeAlphaFields
185  (
186  const regionSplit& regions,
187  const Map<label>& keepRegions,
188  const Map<scalar>& regionVolume,
189  const volScalarField& alpha
190  ) const;
191 
192  //- Mark all regions starting at patches
193  Map<label> findPatchRegions(const polyMesh&, const regionSplit&) const;
194 
195  //- Helper: divide if denom != 0
196  static tmp<scalarField> divide(const scalarField&, const scalarField&);
197 
198  //- Given per-region data calculate per-bin average/deviation and graph
199  void writeGraphs
200  (
201  const word& fieldName, // name of field
202  const labelList& indices, // index of bin for each region
203  const scalarField& sortedField, // per region field data
204  const scalarField& binCount, // per bin number of regions
205  const coordSet& coords // graph data for bins
206  ) const;
207 
208  //- Given per-cell data calculate per-bin average/deviation and graph
209  void writeGraphs
210  (
211  const word& fieldName, // name of field
212  const scalarField& cellField, // per cell field data
213 
214  const regionSplit& regions, // per cell the region(=droplet)
215  const labelList& sortedRegions, // valid regions in sorted order
216  const scalarField& sortedNormalisation,
217 
218  const labelList& indices, // index of bin for each region
219  const scalarField& binCount, // per bin number of regions
220  const coordSet& coords // graph data for bins
221  ) const;
222 
223  //- Disallow default bitwise copy construct
225 
226  //- Disallow default bitwise assignment
227  void operator=(const regionSizeDistribution&);
228 
229 
230 public:
231 
232  //- Runtime type information
233  TypeName("regionSizeDistribution");
234 
235 
236  // Constructors
237 
238  //- Construct for given objectRegistry and dictionary.
239  // Allow the possibility to load fields from files
241  (
242  const word& name,
243  const Time& runTime,
244  const dictionary&
245  );
246 
247 
248  // Destructor
249 
250  virtual ~regionSizeDistribution();
251 
252 
253  // Member Functions
254 
255  //- Read the regionSizeDistribution data
256  virtual bool read(const dictionary&);
257 
258  //- Do nothing
259  virtual bool execute();
260 
261  //- Calculate the regionSizeDistribution and write
262  virtual bool write();
263 };
264 
265 
266 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
267 
268 } // End namespace functionObjects
269 } // End namespace Foam
270 
271 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
272 
273 #ifdef NoRepository
275 #endif
276 
277 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
278 
279 #endif
280 
281 // ************************************************************************* //
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.
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
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.
This function object creates a size distribution via interrogating a continuous phase fraction field...
Holds list of sampling positions.
Definition: coordSet.H:49
const word & name() const
Return the name of this functionObject.
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.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
A class for managing temporary objects.
Definition: PtrList.H:54
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Namespace for OpenFOAM.