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-2019 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  field 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  coordinateSystem
70  {
71  type cartesian;
72  origin (0 0 0);
73  e3 (0 1 1);
74  e1 (1 0 0);
75  }
76  }
77  \endverbatim
78 
79 Usage
80  \table
81  Property | Description | Required | Default value
82  type | type name: regionSizeDistribution |yes|
83  field | phase field to interrogate | yes |
84  patches | patches from which the liquid core is identified | yes|
85  threshold | phase fraction applied to delimit regions | yes |
86  fields | fields to sample | yes |
87  nBins | number of bins for histogram | yes |
88  maxDiameter | maximum region equivalent diameter | yes |
89  minDiameter | minimum region equivalent diameter | no | 0
90  setFormat | writing format | yes |
91  coordinateSystem | transformation for vector fields | no |
92  \endtable
93 
94 See also
95  Foam::functionObject
96  Foam::functionObjects::fvMeshFunctionObject
97  Foam::functionObjects::writeFile
98 
99 SourceFiles
100  regionSizeDistribution.C
101 
102 \*---------------------------------------------------------------------------*/
103 
104 #ifndef functionObjects_regionSizeDistribution_H
105 #define functionObjects_regionSizeDistribution_H
106 
107 #include "fvMeshFunctionObject.H"
108 #include "writeFile.H"
109 #include "writer.H"
110 #include "Map.H"
111 #include "volFieldsFwd.H"
112 #include "wordReList.H"
113 #include "coordinateSystem.H"
114 
115 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
116 
117 namespace Foam
118 {
119 
120 // Forward declaration of classes
121 class regionSplit;
122 
123 namespace functionObjects
124 {
125 
126 /*---------------------------------------------------------------------------*\
127  Class regionSizeDistribution Declaration
128 \*---------------------------------------------------------------------------*/
129 
130 class regionSizeDistribution
131 :
132  public fvMeshFunctionObject
133 {
134  // Private Data
135 
136  writeFile file_;
137 
138  //- Name of field
139  word alphaName_;
140 
141  //- Patches to walk from
142  wordReList patchNames_;
143 
144  //- Clip value
145  scalar threshold_;
146 
147  //- Maximum droplet diameter
148  scalar maxDiam_;
149 
150  //- Minimum droplet diameter
151  scalar minDiam_;
152 
153  //- Number of bins
154  label nBins_;
155 
156  //- Names of fields to sample on regions
157  wordReList fields_;
158 
159  //- Output formatter to write
160  autoPtr<writer<scalar>> formatterPtr_;
161 
162  //- Optional coordinate system
163  autoPtr<coordinateSystem> coordSysPtr_;
164 
165 
166  // Private Member Functions
167 
168  template<class Type>
169  Map<Type> regionSum(const regionSplit&, const Field<Type>&) const;
170 
171  //- Get data in order
172  template<class Type>
173  List<Type> extractData(const UList<label>& keys, const Map<Type>&)
174  const;
175 
176  void writeGraph
177  (
178  const coordSet& coords,
179  const word& valueName,
180  const scalarField& values
181  ) const;
182 
183  //- Write volfields with the parts of alpha which are not
184  // droplets (liquidCore, backGround)
185  void writeAlphaFields
186  (
187  const regionSplit& regions,
188  const Map<label>& keepRegions,
189  const Map<scalar>& regionVolume,
190  const volScalarField& alpha
191  ) const;
192 
193  //- Mark all regions starting at patches
194  Map<label> findPatchRegions(const regionSplit&) const;
195 
196  //- Helper: divide if denom != 0
197  static tmp<scalarField> divide(const scalarField&, const scalarField&);
198 
199  //- Given per-region data calculate per-bin average/deviation and graph
200  void writeGraphs
201  (
202  const word& fieldName, // name of field
203  const labelList& indices, // index of bin for each region
204  const scalarField& sortedField, // per region field data
205  const scalarField& binCount, // per bin number of regions
206  const coordSet& coords // graph data for bins
207  ) const;
208 
209  //- Given per-cell data calculate per-bin average/deviation and graph
210  void writeGraphs
211  (
212  const word& fieldName, // name of field
213  const scalarField& cellField, // per cell field data
214 
215  const regionSplit& regions, // per cell the region(=droplet)
216  const labelList& sortedRegions, // valid regions in sorted order
217  const scalarField& sortedNormalisation,
218 
219  const labelList& indices, // index of bin for each region
220  const scalarField& binCount, // per bin number of regions
221  const coordSet& coords // graph data for bins
222  ) const;
223 
224 
225 public:
226 
227  //- Runtime type information
228  TypeName("regionSizeDistribution");
229 
230 
231  // Constructors
232 
233  //- Construct for given objectRegistry and dictionary.
234  // Allow the possibility to load fields from files
236  (
237  const word& name,
238  const Time& runTime,
239  const dictionary&
240  );
241 
242  //- Disallow default bitwise copy construction
244 
245 
246  // Destructor
247 
248  virtual ~regionSizeDistribution();
249 
250 
251  // Member Functions
252 
253  //- Read the regionSizeDistribution data
254  virtual bool read(const dictionary&);
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: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:158
engineTime & runTime
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
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.
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
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.