extractEulerianParticles.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) 2015-2020 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::functionObjects::extractEulerianParticles
28 
29 Group
30  grpFieldFunctionObjects
31 
32 Description
33  Generates particle size information from Eulerian calculations, e.g. \c VoF.
34 
35  Operands:
36  \table
37  Operand | Type | Location
38  input | - | -
39  output file | - | -
40  output field 1 | scalarField | $OUTPUT/d
41  output field 2 | scalarField | $OUTPUT/soi
42  output field 3 | labelField | $OUTPUT/tag
43  output field 4 | vectorField | $OUTPUT/U
44  \endtable
45 
46  where \c $OUTPUT=$FOAM_CASE/<time>/lagrangian/eulerianParticleCloud.
47 
48 Usage
49  Minimal example by using \c system/controlDict.functions:
50  \verbatim
51  extractEulerianParticles1
52  {
53  // Mandatory entries (unmodifiable)
54  type extractEulerianParticles;
55  libs (fieldFunctionObjects);
56 
57  // Mandatory entries (runtime modifiable)
58  faceZone f0;
59  alpha alpha.water;
60 
61  // Optional entries (runtime modifiable)
62  alphaThreshold 0.1;
63  nLocations 0;
64  U U;
65  rho rho;
66  phi phi;
67  minDiameter 1e-30;
68  maxDiameter 1e30;
69 
70  // Optional (inherited) entries
71  ...
72  }
73  \endverbatim
74 
75  where the entries mean:
76  \table
77  Property | Description | Type | Req'd | Dflt
78  type | Type name: extractEulerianParticles | word | yes | -
79  libs | Library name: fieldFunctionObjects | word | yes | -
80  faceZone | Name of faceZone used as collection surface | word | yes | -
81  alpha | Name of phase indicator field | word | yes | -
82  alphaThreshold | Threshold for alpha field | scalar | no | 0.1
83  nLocations | Number of injection bins to generate | label | no | 0
84  U | Name of velocity field | word | no | U
85  rho | Name of density field | word | no | rho
86  phi | Name of flux field | word | no | phi
87  minDiameter | Minimum diameter | scalar | no | SMALL
88  maxDiameter | Maximum diameter | scalar | no | GREAT
89  \endtable
90 
91  The inherited entries are elaborated in:
92  - \link functionObject.H \endlink
93  - \link writeFile.H \endlink
94 
95  Usage by the \c postProcess utility is not available.
96 
97 See also
98  - Foam::functionObject
99  - Foam::functionObjects::fvMeshFunctionObject
100  - Foam::functionObjects::writeFile
101  - Foam::eulerianParticle
102  - ExtendedCodeGuide::functionObjects::field::extractEulerianParticles
103 
104 SourceFiles
105  extractEulerianParticles.C
106  extractEulerianParticlesTemplates.C
107 
108 \*---------------------------------------------------------------------------*/
109 
110 #ifndef functionObjects_extractEulerianParticles_H
111 #define functionObjects_extractEulerianParticles_H
112 
113 #include "fvMeshFunctionObject.H"
114 #include "writeFile.H"
115 #include "runTimeSelectionTables.H"
116 #include "polyMesh.H"
117 #include "surfaceFieldsFwd.H"
118 #include "globalIndex.H"
119 #include "eulerianParticle.H"
120 #include "injectedParticleCloud.H"
121 
122 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
123 
124 namespace Foam
125 {
126 namespace functionObjects
127 {
128 /*---------------------------------------------------------------------------*\
129  Class extractEulerianParticlesFunctionObject Declaration
130 \*---------------------------------------------------------------------------*/
131 
132 class extractEulerianParticles
133 :
134  public fvMeshFunctionObject,
135  public writeFile
136 {
137 protected:
138 
139  // Protected Data
140 
141  //- Storage for collected particles
142  injectedParticleCloud cloud_;
143 
144 
145  // faceZone info
146 
147  //- Name of faceZone to sample
148  word faceZoneName_;
149 
150  //- Index of the faceZone
151  label zoneID_;
152 
153  //- Patch indices where faceZone face intersect patch
155 
156  //- Patch face indices where faceZone face intersect patch
158 
159 
160  // Field names
161 
162  //- Name of phase fraction field
163  word alphaName_;
164 
165  //- Value of phase fraction used to identify particle boundaries
166  scalar alphaThreshold_;
167 
168  //- Name of the velocity field, default = U
169  word UName_;
170 
171  //- Name of the density field, default = rho
172  word rhoName_;
173 
174  //- Name of the flux field, default ="rho"
175  word phiName_;
176 
177 
178  // Agglomeration
179 
180  //- Number of sample locations to generate
181  label nInjectorLocations_;
182 
183  //- Agglomeration addressing from fine to coarse (local proc only)
185 
186  //- Global coarse face addressing
187  globalIndex globalCoarseFaces_;
188 
189 
190  // Particle collection info
191 
192  //- Region indices in faceZone faces from last iteration
194 
195  //- Particle properties (partial, being accumulated)
196  List<eulerianParticle> particles_;
197 
198  //- Map from region to index in particles_ list
199  Map<label> regionToParticleMap_;
200 
201  //- Minimum diameter (optional)
202  // Can be used to filter out 'small' particles
203  scalar minDiameter_;
204 
205  //- Maximum diameter (optional)
206  // Can be used to filter out 'large' particles
207  scalar maxDiameter_;
208 
209 
210  // Statistics
211 
212  //- Total number of collected particles
213  label nCollectedParticles_;
214 
215  //- Total collected volume [m3]
216  scalar collectedVolume_;
217 
218  //- Total number of discarded particles
219  label nDiscardedParticles_;
220 
221  //- Total discarded volume [m3]
222  scalar discardedVolume_;
223 
224 
225  // Protected Member Functions
226 
227  //- Check that the faceZone is valid
228  virtual void checkFaceZone();
229 
230  //- Initialise the particle collection bins
231  virtual void initialiseBins();
232 
233  //- Return the volumetric flux
234  virtual tmp<surfaceScalarField> phiU() const;
235 
236  //- Set the blocked faces, i.e. where alpha > alpha threshold value
237  virtual void setBlockedFaces
238  (
239  const surfaceScalarField& alphaf,
240  const faceZone& fz,
241  boolList& blockedFaces
242  );
243 
244  //- Calculate the addressing between regions between iterations
245  //- Returns the number of active regions (particles)
246  virtual void calculateAddressing
247  (
248  const label nRegionsNew,
249  const scalar time,
250  labelList& regionFaceIDs
251  );
253  //- Collect particles that have passed through the faceZone
254  virtual void collectParticle
255  (
256  const scalar time,
257  const label regioni
258  );
259 
260  //- Process latest region information
261  virtual void accumulateParticleInfo
262  (
263  const surfaceScalarField& alphaf,
264  const surfaceScalarField& phi,
265  const labelList& regionFaceIDs,
266  const faceZone& fz
267  );
268 
269  template<class Type>
270  inline Type faceValue
271  (
273  const label localFaceI,
274  const label globalFaceI
275  ) const;
276 
277 
278 public:
279 
280  // Static data members
281 
282  //- Runtime type information
283  TypeName("extractEulerianParticles");
284 
286  // Constructors
287 
288  //- Construct from components
290  (
291  const word& name,
292  const Time& runTime,
293  const dictionary& dict
294  );
295 
296  //- No copy construct
299  //- No copy assignment
300  void operator=(const extractEulerianParticles&) = delete;
301 
302 
303  //- Destructor
304  virtual ~extractEulerianParticles() = default;
305 
306 
307  // Member Functions
309  //- Read the field min/max data
310  virtual bool read(const dictionary&);
311 
312  //- Execute
313  virtual bool execute();
314 
315  //- Write
316  virtual bool write();
317 };
318 
319 
320 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
322 } // End namespace functionObjects
323 } // End namespace Foam
324 
325 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
327 #ifdef NoRepository
329 #endif
330 
331 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
332 
333 #endif
334 
335 // ************************************************************************* //
label nDiscardedParticles_
Total number of discarded particles.
dictionary dict
virtual void checkFaceZone()
Check that the faceZone is valid.
virtual ~extractEulerianParticles()=default
Destructor.
virtual void calculateAddressing(const label nRegionsNew, const scalar time, labelList &regionFaceIDs)
Calculate the addressing between regions between iterations Returns the number of active regions (par...
rDeltaTY field()
virtual void collectParticle(const scalar time, const label regioni)
Collect particles that have passed through the faceZone.
labelList fineToCoarseAddr_
Agglomeration addressing from fine to coarse (local proc only)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Map< label > regionToParticleMap_
Map from region to index in particles_ list.
scalar discardedVolume_
Total discarded volume [m3].
word rhoName_
Name of the density field, default = rho.
labelList regions0_
Region indices in faceZone faces from last iteration.
engineTime & runTime
extractEulerianParticles(const word &name, const Time &runTime, const dictionary &dict)
Construct from components.
globalIndex globalCoarseFaces_
Global coarse face addressing.
Generates particle size information from Eulerian calculations, e.g. VoF.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
const word & name() const noexcept
Return the name of this functionObject.
List< eulerianParticle > particles_
Particle properties (partial, being accumulated)
labelList patchIDs_
Patch indices where faceZone face intersect patch.
virtual bool read(const dictionary &)
Read the field min/max data.
virtual void initialiseBins()
Initialise the particle collection bins.
A class for handling words, derived from Foam::string.
Definition: word.H:63
word phiName_
Name of the flux field, default ="rho".
void operator=(const extractEulerianParticles &)=delete
No copy assignment.
virtual void accumulateParticleInfo(const surfaceScalarField &alphaf, const surfaceScalarField &phi, const labelList &regionFaceIDs, const faceZone &fz)
Process latest region information.
TypeName("extractEulerianParticles")
Runtime type information.
scalar collectedVolume_
Total collected volume [m3].
virtual void setBlockedFaces(const surfaceScalarField &alphaf, const faceZone &fz, boolList &blockedFaces)
Set the blocked faces, i.e. where alpha > alpha threshold value.
label nInjectorLocations_
Number of sample locations to generate.
word UName_
Name of the velocity field, default = U.
scalar alphaThreshold_
Value of phase fraction used to identify particle boundaries.
Type faceValue(const GeometricField< Type, fvsPatchField, surfaceMesh > &field, const label localFaceI, const label globalFaceI) const
virtual tmp< surfaceScalarField > phiU() const
Return the volumetric flux.
injectedParticleCloud cloud_
Storage for collected particles.
label nCollectedParticles_
Total number of collected particles.
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:60
Macros to ease declaration of run-time selection tables.
const Time & time() const
Return time database.
labelList patchFaceIDs_
Patch face indices where faceZone face intersect patch.
List< label > labelList
A List of labels.
Definition: List.H:62
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Namespace for OpenFOAM.