Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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
27Class
28 Foam::functionObjects::regionSizeDistribution
29
30Group
31 grpFieldFunctionObjects
32
33Description
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 | - | postProcessing/<FO>/<time>/files
68 output field | - | -
69 \endtable
70
71Usage
72 Minimal example by using \c system/controlDict.functions:
73 \verbatim
74 regionSizeDistributionFO
75 {
76 // Mandatory entries
77 type regionSizeDistribution;
78 libs (fieldFunctionObjects);
79 field <word>;
80 patches (<wordList>); // (patch1 patch2 ...);
81 fields (<wordList>); // (field1 field2 ...);
82 threshold <scalar>;
83 maxDiameter <scalar>;
84 nBins <label>;
85 setFormat <word>; // csv, raw, xmgrace, gnuplot
86
87 // Optional entries
88 minDiameter <scalar>;
89 coordinateSystem
90 {
91 origin (0 0 0);
92 rotation
93 {
94 type axes;
95 e3 (0 0 1);
96 e1 (1 0 0);
97 }
98 }
99
100 // Optional downstream iso-plane bins
101 isoPlanes true;
102
103 // Conditional entries
104
105 // when 'isoPlanes' is 'true'
106 // Plane normal and point definition
107 origin (1e-4 0 5e-4);
108 direction (1 0 1);
109
110 // Maximum diameter of the cylinder formed by the origin point
111 // and direction
112 maxD 3e-4;
113
114 // Maximum downstream distance
115 maxDownstream 6e-4;
116
117 // Number of iso-plane bins
118 nDownstreamBins 20;
119
120 // Inherited entries
121 ...
122 }
123 \endverbatim
124
125 where the entries mean:
126 \table
127 Property | Description | Type | Reqd | Deflt
128 type | Type name: regionSizeDistribution | word | yes | -
129 libs | Library name: fieldFunctionObjects | word | yes | -
130 field | Phase field to interrogate | word | yes | -
131 patches | Patches wherefrom the liquid core is identified <!--
132 --> | wordList | yes | -
133 fields | Fields to sample | wordList | yes | -
134 threshold | Phase fraction applied to delimit regions | scalar | yes | -
135 maxDiameter | Maximum droplet diameter | scalar | yes | -
136 minDiameter | Minimum droplet diameter | scalar | no | 0.0
137 nBins | Number of bins for histogram | label | yes | -
138 setFormat | Output format | word | yes | -
139 isoPlanes | Flag for isoPlanes | bool | no | false
140 origin | Origin of the plane when isoPlanes is used <!--
141 --> | vector | yes | -
142 direction <!--
143 --> | Direction of the plane when isoPlanes is used | vector | yes | -
144 maxD | Maximum diameter of the sampling <!--
145 --> cylinder when isoPlanes is used | vector | yes | -
146 nDownstreamBins <!--
147 --> | Number of bins when isoPlanes is used | label | yes | -
148 maxDownstream <!--
149 --> | Maximum distance from origin when isoPlanes is used <!--
150 --> | scalar | yes | -
151 \endtable
152
153 The inherited entries are elaborated in:
154 - \link functionObject.H \endlink
155 - \link writeFile.H \endlink
156
157SourceFiles
158 regionSizeDistribution.C
159 regionSizeDistributionTemplates.C
160
161\*---------------------------------------------------------------------------*/
162
163#ifndef Foam_functionObjects_regionSizeDistribution_H
164#define Foam_functionObjects_regionSizeDistribution_H
165
166#include "fvMeshFunctionObject.H"
167#include "writeFile.H"
168#include "coordSetWriter.H"
169#include "Map.H"
170#include "volFieldsFwd.H"
171#include "wordRes.H"
172#include "coordinateSystem.H"
173
174// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
175
176namespace Foam
177{
178
179// Forward Declarations
180class regionSplit;
181
182namespace functionObjects
183{
184
185/*---------------------------------------------------------------------------*\
186 Class regionSizeDistribution Declaration
187\*---------------------------------------------------------------------------*/
188
190:
191 public functionObjects::fvMeshFunctionObject,
192 public writeFile
193{
194 // Private Data
195
196 //- Number of bins
197 label nBins_;
198
199 //- Name of field
200 word alphaName_;
201
202 //- Clip value
203 scalar threshold_;
204
205 //- Maximum droplet diameter
206 scalar maxDiam_;
207
208 //- Minimum droplet diameter
209 scalar minDiam_;
210
211 //- Patches to walk from
212 wordRes patchNames_;
213
214 //- Names of fields to sample on regions
215 wordRes fields_;
216
217 //- Output formatter to write
218 mutable autoPtr<coordSetWriter> formatterPtr_;
219
220 //- Optional coordinate system
221 autoPtr<coordinateSystem> csysPtr_;
222
223 // Optional extra definition of bins on planes downstream to the origin
224 // point and maximum diameter
225
226 //- Optional origin point
227 vector origin_;
228
229 //- Optional plane normal direction
230 vector direction_;
231
232 //- Optional maximum diameter on plane
233 scalar maxDiameter_;
234
235 //- Optional number of bins for
236 scalar nDownstreamBins_;
237
238 //- Optional maximum downstream coordinate from origin
239 scalar maxDownstream_;
240
241 //- Switch to enable iso-planes sampling
242 bool isoPlanes_;
243
244
245 // Private Member Functions
246
247 //- Write volfields with the parts of alpha which are not
248 //- droplets (liquidCore, backGround)
249 void writeAlphaFields
250 (
251 const regionSplit& regions,
252 const labelHashSet& keepRegions,
253 const Map<scalar>& regionVolume,
254 const volScalarField& alpha
255 ) const;
256
257 //- Mark all regions starting at patches
258 labelHashSet findPatchRegions(const regionSplit&) const;
259
260 //- Helper: divide if denom != 0
261 static tmp<scalarField> divide(const scalarField&, const scalarField&);
262
263 //- Given per-region data calculate per-bin average/deviation and graph
264 void writeGraphs
265 (
266 const word& fieldName, // name of field
267 const scalarField& sortedField, // per region field data
268
269 const labelList& indices, // index of bin for each region
270 const scalarField& binCount, // per bin number of regions
271 const coordSet& coords // graph data for bins
272 ) const;
273
274 //- Given per-cell data calculate per-bin average/deviation and graph
275 void writeGraphs
276 (
277 const word& fieldName, // name of field
278 const scalarField& cellField, // per cell field data
279 const regionSplit& regions, // per cell the region(=droplet)
280 const labelList& sortedRegions, // valid regions in sorted order
281 const scalarField& sortedNormalisation,
282
283 const labelList& indices, // index of bin for each region
284 const scalarField& binCount, // per bin number of regions
285 const coordSet& coords // graph data for bins
286 ) const;
287
288
289public:
290
291 //- Runtime type information
292 TypeName("regionSizeDistribution");
293
294
295 // Constructors
296
297 //- Construct for given objectRegistry and dictionary.
298 // Allow the possibility to load fields from files
300 (
301 const word& name,
302 const Time& runTime,
303 const dictionary&
304 );
305
306 //- No copy construct
308
309 //- No copy assignment
310 void operator=(const regionSizeDistribution&) = delete;
311
312
313 // Destructor
314 virtual ~regionSizeDistribution() = default;
315
317 // Member Functions
318
319 //- Read the function-object dictionary
320 virtual bool read(const dictionary& dict);
321
322 //- Execute the function-object operations
323 virtual bool execute();
324
325 //- Write the function-object results
326 virtual bool write();
327};
328
329
330// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
331
332} // End namespace functionObjects
333} // End namespace Foam
334
335// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
336
337#endif
338
339// ************************************************************************* //
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
Holds list of sampling positions.
Definition coordSet.H:52
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const word & name() const noexcept
Return the name of this functionObject.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Creates a droplet size distribution via interrogating a continuous phase fraction field.
TypeName("regionSizeDistribution")
Runtime type information.
regionSizeDistribution(const word &name, const Time &runTime, const dictionary &)
Construct for given objectRegistry and dictionary.
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
regionSizeDistribution(const regionSizeDistribution &)=delete
No copy construct.
void operator=(const regionSizeDistribution &)=delete
No copy assignment.
virtual bool execute()
Execute the function-object operations.
virtual bool write()
Write the function-object results.
Base class for writing single files from the function objects.
Definition writeFile.H:113
writeFile(const objectRegistry &obr, const fileName &prefix, const word &name="undefined", const bool writeToFile=true, const string &ext=".dat")
Construct from objectRegistry, prefix, fileName.
Definition writeFile.C:200
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
A class for managing temporary objects.
Definition tmp.H:75
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
engineTime & runTime
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition List.H:62
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
Vector< scalar > vector
Definition vector.H:57
volScalarField & alpha
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68
Forwards and collection of common volume field types.