extractEulerianParticles.C
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-------------------------------------------------------------------------------
10License
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\*---------------------------------------------------------------------------*/
27
29#include "regionSplit2D.H"
31#include "volFields.H"
32#include "surfaceFields.H"
33#include "surfaceInterpolate.H"
35#include "emptyPolyPatch.H"
36#include "coupledPolyPatch.H"
37#include "binned.H"
39
40// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41
42namespace Foam
43{
44namespace functionObjects
45{
48 (
52 );
53}
54}
55
56// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
57
59{
61
63 if (zoneID_ == -1)
64 {
66 << "Unable to find faceZone " << faceZoneName_
67 << ". Available faceZones are: " << mesh_.faceZones().names()
68 << exit(FatalError);
69 }
70
71 const faceZone& fz = mesh_.faceZones()[zoneID_];
72 const label nFaces = fz.size();
73 const label allFaces = returnReduce(nFaces, sumOp<label>());
74
75 if (allFaces < nInjectorLocations_)
76 {
78 << "faceZone " << faceZoneName_
79 << ": Number of faceZone faces (" << allFaces
80 << ") is less than the number of requested locations ("
81 << nInjectorLocations_ << ")."
82 << exit(FatalError);
83 }
84
85 Info<< type() << " " << name() << " output:" << nl
86 << " faceZone : " << faceZoneName_ << nl
87 << " faces : " << allFaces << nl
88 << endl;
89
90 // Initialise old iteration blocked faces
91 // Note: for restart, this info needs to be written/read
92 regions0_.setSize(fz.size(), -1);
93}
94
95
97{
99
100 if (!nInjectorLocations_)
101 {
102 return;
103 }
104
105 const faceZone& fz = mesh_.faceZones()[zoneID_];
106
107 // Agglomerate faceZone faces into nInjectorLocations_ global locations
108 const indirectPrimitivePatch patch
109 (
110 IndirectList<face>(mesh_.faces(), fz),
111 mesh_.points()
112 );
113
114 const label nFaces = fz.size();
115 label nLocations = nInjectorLocations_;
116
117 if (Pstream::parRun())
118 {
119 label nGlobalFaces = returnReduce(nFaces, sumOp<label>());
120 scalar fraction = scalar(nFaces)/scalar(nGlobalFaces);
121 nLocations = ceil(fraction*nInjectorLocations_);
122 if (debug)
123 {
124 Pout<< "nFaces:" << nFaces
125 << ", nGlobalFaces:" << nGlobalFaces
126 << ", fraction:" << fraction
127 << ", nLocations:" << nLocations
128 << endl;
129 }
130 }
131
133 (
134 patch.localFaces(),
135 patch.localPoints(),
136 10,
137 50,
138 nLocations,
139 labelMax,
140 180
141 );
142
143 ppa.agglomerate();
144
145 label nCoarseFaces = 0;
146 if (nFaces != 0)
147 {
148 fineToCoarseAddr_ = ppa.restrictTopBottomAddressing();
149 nCoarseFaces = max(fineToCoarseAddr_) + 1;
150 }
151
152 globalCoarseFaces_ = globalIndex(nCoarseFaces);
153
154 Info<< "Created " << returnReduce(nCoarseFaces, sumOp<label>())
155 << " coarse faces" << endl;
156}
157
158
161{
163
165 (
166 mesh_.lookupObject<surfaceScalarField>(phiName_)
167 );
168
169 if (phi.dimensions() == dimMass/dimTime)
170 {
171 const volScalarField& rho =
172 mesh_.lookupObject<volScalarField>(rhoName_);
173
174 return phi/fvc::interpolate(rho);
175 }
176
177 return phi;
178}
179
180
182(
183 const surfaceScalarField& alphaf,
184 const faceZone& fz,
185 boolList& blockedFaces
186)
187{
189
190 // Initialise storage for patch and patch-face indices where faceZone
191 // intersects mesh patch(es)
192 patchIDs_.setSize(fz.size(), -1);
193 patchFaceIDs_.setSize(fz.size(), -1);
194
195 label nBlockedFaces = 0;
196 forAll(fz, localFacei)
197 {
198 const label facei = fz[localFacei];
199
200 if (mesh_.isInternalFace(facei))
201 {
202 if (alphaf[facei] > alphaThreshold_)
203 {
204 blockedFaces[localFacei] = true;
205 }
206 }
207 else
208 {
209 label patchi = mesh_.boundaryMesh().whichPatch(facei);
210 label patchFacei = -1;
211
212 const polyPatch& pp = mesh_.boundaryMesh()[patchi];
213 const scalarField& alphafp = alphaf.boundaryField()[patchi];
214 const auto* cpp = isA<coupledPolyPatch>(pp);
215
216 if (cpp)
217 {
218 patchFacei = (cpp->owner() ? pp.whichFace(facei) : -1);
219 }
220 else if (!isA<emptyPolyPatch>(pp))
221 {
222 patchFacei = pp.whichFace(facei);
223 }
224
225 if (patchFacei == -1)
226 {
227 patchi = -1;
228 }
229 else if (alphafp[patchFacei] > alphaThreshold_)
230 {
231 blockedFaces[localFacei] = true;
232 }
233
234 patchIDs_[localFacei] = patchi;
235 patchFaceIDs_[localFacei] = patchFacei;
236 }
237 }
238
239 DebugInFunction << "Number of blocked faces: " << nBlockedFaces << endl;
240}
241
242
244(
245 const scalar time,
246 const label regioni
247)
248{
249 DebugInFunction << "collectParticle: " << regioni << endl;
250
251 const label particlei = regionToParticleMap_[regioni];
252 eulerianParticle p = particles_[particlei];
253
254 if (p.faceIHit != -1 && nInjectorLocations_)
255 {
256 // Use coarse face index for tag output
257 label coarseFacei = fineToCoarseAddr_[p.faceIHit];
258 p.faceIHit = globalCoarseFaces_.toGlobal(coarseFacei);
259 }
260
262
263 const scalar pDiameter = cbrt(6.0*p.V/constant::mathematical::pi);
264
265 if ((pDiameter > minDiameter_) && (pDiameter < maxDiameter_))
266 {
267 if (Pstream::master())
268 {
269 const scalar d = cbrt(6.0*p.V/constant::mathematical::pi);
270 const point position = p.VC/(p.V + ROOTVSMALL);
271 const vector U = p.VU/(p.V + ROOTVSMALL);
272 label tag = -1;
273 if (nInjectorLocations_)
274 {
275 tag = p.faceIHit;
276 }
277
279 (
280 mesh_,
281 position,
282 tag,
283 time,
284 d,
285 U,
286 false // not looking to set cell owner etc.
287 );
288
289 cloud_.addParticle(ip);
290
291 collectedVolume_ += p.V;
292 }
293
294 ++nCollectedParticles_;
295 }
296 else
297 {
298 // Discard particles over/under diameter threshold
299 ++nDiscardedParticles_;
300 discardedVolume_ += p.V;
301 }
302}
303
304
306(
307 const label nNewRegions,
308 const scalar time,
309 labelList& regionFaceIDs
310)
311{
313
314 // Determine mapping between old and new regions so that we can
315 // accumulate particle info
316 labelList oldToNewRegion(particles_.size(), -1);
317 labelList newToNewRegion(identity(nNewRegions));
318
319 forAll(regionFaceIDs, facei)
320 {
321 label newRegioni = regionFaceIDs[facei];
322 label oldRegioni = regions0_[facei];
323
324 if (newRegioni != -1 && oldRegioni != -1)
325 {
326 // If old region has split into multiple regions we need to
327 // renumber new regions to maintain connectivity with old regions
328 newToNewRegion[newRegioni] =
329 max(newRegioni, oldToNewRegion[oldRegioni]);
330 oldToNewRegion[oldRegioni] = newRegioni;
331 }
332 }
333
334 // Create map from new regions to slots in particles list
335 // - filter through new-to-new addressing to identify new particles
337
338 label nParticle = -1;
339 labelHashSet newRegions;
340 Map<label> newRegionToParticleMap;
341 forAll(newToNewRegion, newRegioni0)
342 {
343 label newRegioni = newToNewRegion[newRegioni0];
344 if (newRegions.insert(newRegioni))
345 {
346 ++nParticle;
347 }
348
349 // New particle slot
350 newRegionToParticleMap.insert(newRegioni0, nParticle);
351 }
352
353 // Accumulate old region data or create a new particle if there is no
354 // mapping from the old-to-new region
356
357 List<eulerianParticle> newParticles(newRegionToParticleMap.size());
358 forAll(oldToNewRegion, oldRegioni)
359 {
360 label newRegioni = oldToNewRegion[oldRegioni];
361 if (newRegioni == -1)
362 {
363 // No mapping from old-to-new - collect new particle
365 << "Collecting particle from oldRegion:" << oldRegioni
366 << endl;
367
368 collectParticle(time, oldRegioni);
369 }
370 else
371 {
372 // Combine existing particle into new particle
373 label newParticlei = newRegionToParticleMap[newRegioni];
374 label oldParticlei = regionToParticleMap_[oldRegioni];
375
377 << "Combining newRegioni: " << newRegioni
378 << "(p:" << newParticlei << ") and "
379 << "oldRegioni: " << oldRegioni
380 << "(p:" << oldParticlei << ")"
381 << endl;
382
383 newParticles[newParticlei] =
385 (
386 newParticles[newParticlei],
387 particles_[oldParticlei]
388 );
389 }
390 }
391
392 // Reset the particles list and addressing for latest available info
393 particles_.transfer(newParticles);
394 regionToParticleMap_ = newRegionToParticleMap;
395
396 // Reset the region IDs for the next integration step
397 // - these become the oldRegioni's
398 regions0_ = regionFaceIDs;
399}
400
401
403(
404 const surfaceScalarField& alphaf,
405 const surfaceScalarField& phi,
406 const labelList& regionFaceIDs,
407 const faceZone& fz
408)
409{
411
412 const volVectorField& U = mesh_.lookupObject<volVectorField>(UName_);
414
415 const scalar deltaT = mesh_.time().deltaTValue();
416 const pointField& faceCentres = mesh_.faceCentres();
417
418 forAll(regionFaceIDs, localFacei)
419 {
420 const label newRegioni = regionFaceIDs[localFacei];
421 if (newRegioni != -1)
422 {
423 const label particlei = regionToParticleMap_[newRegioni];
424 const label meshFacei = fz[localFacei];
425 eulerianParticle& p = particles_[particlei];
426
427 if (p.faceIHit < 0)
428 {
429 // New particle - does not exist in particles_ list
430 p.faceIHit = localFacei;
431 p.V = 0;
432 p.VC = vector::zero;
433 p.VU = vector::zero;
434 }
435
436 // Accumulate particle properties
437 scalar magPhii = mag(faceValue(phi, localFacei, meshFacei));
438 vector Ufi = faceValue(Uf, localFacei, meshFacei);
439 scalar dV = magPhii*deltaT;
440 p.V += dV;
441 p.VC += dV*faceCentres[meshFacei];
442 p.VU += dV*Ufi;
443 }
444 }
445}
446
447
448// * * * * * * * * * * * * * * * * Constructor * * * * * * * * * * * * * * * //
449
451(
452 const word& name,
453 const Time& runTime,
454 const dictionary& dict
455)
456:
459 cloud_(mesh_, "eulerianParticleCloud"),
460 faceZoneName_(word::null),
461 zoneID_(-1),
462 patchIDs_(),
463 patchFaceIDs_(),
464 alphaName_("alpha"),
465 alphaThreshold_(0.1),
466 UName_("U"),
467 rhoName_("rho"),
468 phiName_("phi"),
469 nInjectorLocations_(0),
470 fineToCoarseAddr_(),
471 globalCoarseFaces_(),
472 regions0_(),
473 particles_(),
474 regionToParticleMap_(),
475 minDiameter_(ROOTVSMALL),
476 maxDiameter_(GREAT),
477 nCollectedParticles_(getProperty<label>("nCollectedParticles", 0)),
478 collectedVolume_(getProperty<scalar>("collectedVolume", 0)),
479 nDiscardedParticles_(getProperty<label>("nDiscardedParticles", 0)),
480 discardedVolume_(getProperty<scalar>("discardedVolume", 0))
481{
482 if (mesh_.nSolutionD() != 3)
483 {
485 << name << " function object only applicable to 3-D cases"
486 << exit(FatalError);
487 }
488
489 read(dict);
490}
491
492
493// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
494
496(
497 const dictionary& dict
498)
499{
501
503 {
504 dict.readEntry("faceZone", faceZoneName_);
505 dict.readEntry("alpha", alphaName_);
506
507 dict.readIfPresent("alphaThreshold", alphaThreshold_);
508 dict.readIfPresent("U", UName_);
509 dict.readIfPresent("rho", rhoName_);
510 dict.readIfPresent("phi", phiName_);
511 dict.readIfPresent("nLocations", nInjectorLocations_);
512 dict.readIfPresent("minDiameter", minDiameter_);
513 dict.readIfPresent("maxDiameter", maxDiameter_);
514
515 checkFaceZone();
516
517 if (nInjectorLocations_)
518 {
519 initialiseBins();
520 }
521
522 return true;
523 }
524
525 return false;
526}
527
528
530{
532
533 Log << type() << " " << name() << " output:" << nl;
534
535 const volScalarField& alpha =
536 mesh_.lookupObject<volScalarField>(alphaName_);
537
538 const surfaceScalarField alphaf
539 (
540 typeName + ":alphaf",
542 );
543
544 const faceZone& fz = mesh_.faceZones()[zoneID_];
545 const indirectPrimitivePatch patch
546 (
547 IndirectList<face>(mesh_.faces(), fz),
548 mesh_.points()
549 );
550
551 // Set the blocked faces, i.e. where alpha > alpha threshold value
552 boolList blockedFaces(fz.size(), false);
553 setBlockedFaces(alphaf, fz, blockedFaces);
554
555 // Split the faceZone according to the blockedFaces
556 // - Returns a list of (disconnected) region index per face zone face
557 regionSplit2D regionFaceIDs(mesh_, patch, blockedFaces);
558
559 // Global number of regions
560 const label nRegionsNew = regionFaceIDs.nRegions();
561
562 // Calculate the addressing between the old and new region information
563 // Also collects particles that have traversed the faceZone
564 // - Note: may also update regionFaceIDs
565 calculateAddressing
566 (
567 nRegionsNew,
568 mesh_.time().value(),
569 regionFaceIDs
570 );
571
572 // Process latest region information
573 tmp<surfaceScalarField> tphi = phiU();
574 accumulateParticleInfo(alphaf, tphi(), regionFaceIDs, fz);
575
576 Log << " Collected particles : " << nCollectedParticles_ << nl
577 << " Collected volume : " << collectedVolume_ << nl
578 << " Discarded particles : " << nDiscardedParticles_ << nl
579 << " Discarded volume : " << discardedVolume_ << nl
580 << " Particles in progress : " << particles_.size() << nl
581 << endl;
582
583 return true;
584}
585
586
588{
590
591 cloud_.write();
592
593 setProperty("nCollectedParticles", nCollectedParticles_);
594 setProperty("collectedVolume", collectedVolume_);
595 setProperty("nDiscardedParticles", nDiscardedParticles_);
596 setProperty("discardedVolume", discardedVolume_);
597
598 return true;
599}
600
601
602// ************************************************************************* //
#define Log
Definition: PDRblock.C:35
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
surfaceScalarField & phi
const Boundary & boundaryField() const
Return const-reference to the boundary field.
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:180
label size() const noexcept
The number of elements in table.
Definition: HashTableI.H:52
A List with indirect addressing.
Definition: IndirectList.H:119
void setSize(const label n)
Alias for resize()
Definition: List.H:218
A HashTable to objects of type <T> with a label key.
Definition: Map.H:60
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
A list of faces which address into the list of points.
static void listCombineAllGather(const List< commsStruct > &comms, List< T > &values, const CombineOp &cop, const int tag, const label comm)
After completion all processors have the same data.
virtual bool read()
Re-read model coefficients if they have changed.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Definition: ZoneMesh.C:525
wordList names() const
A list of the zone names.
Definition: ZoneMesh.C:304
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:67
Abstract base-class for Time/database function objects.
const word & name() const noexcept
Return the name of this functionObject.
virtual const word & type() const =0
Runtime type information.
Generates particle size information from Eulerian calculations, e.g. VoF.
virtual void accumulateParticleInfo(const surfaceScalarField &alphaf, const surfaceScalarField &phi, const labelList &regionFaceIDs, const faceZone &fz)
Process latest region information.
virtual tmp< surfaceScalarField > phiU() const
Return the volumetric flux.
label nInjectorLocations_
Number of sample locations to generate.
virtual void setBlockedFaces(const surfaceScalarField &alphaf, const faceZone &fz, boolList &blockedFaces)
Set the blocked faces, i.e. where alpha > alpha threshold value.
virtual void collectParticle(const scalar time, const label regioni)
Collect particles that have passed through the faceZone.
labelList regions0_
Region indices in faceZone faces from last iteration.
virtual void initialiseBins()
Initialise the particle collection bins.
virtual void checkFaceZone()
Check that the faceZone is valid.
virtual bool read(const dictionary &)
Read the field min/max data.
virtual void calculateAddressing(const label nRegionsNew, const scalar time, labelList &regionFaceIDs)
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
Computes the magnitude of an input field.
Definition: mag.H:153
Base class for writing single files from the function objects.
Definition: writeFile.H:120
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
Primarily stores particle properties so that it can be injected at a later time. Note that this store...
Primitive patch pair agglomerate method.
const labelList & restrictTopBottomAddressing() const
Return restriction from top level to bottom level.
void agglomerate()
Agglomerate patch.
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:498
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
Definition: polyMesh.C:900
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:451
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:315
Splits a patch into regions based on a mask field. Result is a globally consistent label list of regi...
Definition: regionSplit2D.H:62
label nRegions() const noexcept
Return the global number of regions.
Definition: regionSplit2D.H:99
splitCell * master() const
Definition: splitCell.H:113
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
volScalarField & p
engineTime & runTime
autoPtr< surfaceVectorField > Uf
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const labelList nFaces(UPstream::listGatherValues< label >(aMesh.nFaces()))
#define DebugInfo
Report an information message using Foam::Info.
#define DebugInFunction
Report an information message using Foam::Info.
constexpr scalar pi(M_PI)
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:38
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
constexpr label labelMax
Definition: label.H:61
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
error FatalError
dimensionedScalar cbrt(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
volScalarField & alpha
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
A non-counting (dummy) refCount.
Definition: refCount.H:59
Foam::surfaceFields.