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 -------------------------------------------------------------------------------
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 \*---------------------------------------------------------------------------*/
27 
29 #include "regionSplit2D.H"
30 #include "mathematicalConstants.H"
31 #include "volFields.H"
32 #include "surfaceFields.H"
33 #include "surfaceInterpolate.H"
34 #include "pairPatchAgglomeration.H"
35 #include "emptyPolyPatch.H"
36 #include "coupledPolyPatch.H"
37 #include "binned.H"
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 namespace functionObjects
45 {
46  defineTypeNameAndDebug(extractEulerianParticles, 0);
48  (
49  functionObject,
50  extractEulerianParticles,
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
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 
132  pairPatchAgglomeration ppa
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);
154  Info<< "Created " << returnReduce(nCoarseFaces, sumOp<label>())
155  << " coarse faces" << endl;
156 }
157 
158 
161 {
163 
164  const auto& phi = mesh_.lookupObject<surfaceScalarField>(phiName_);
165 
166  if (phi.dimensions() == dimMass/dimTime)
167  {
168  const auto& rho = mesh_.lookupObject<volScalarField>(rhoName_);
169 
170  return phi/fvc::interpolate(rho);
171  }
172 
173  return phi;
174 }
175 
176 
178 (
179  const surfaceScalarField& alphaf,
180  const faceZone& fz,
181  boolList& blockedFaces
182 )
183 {
185 
186  // Initialise storage for patch and patch-face indices where faceZone
187  // intersects mesh patch(es)
188  patchIDs_.setSize(fz.size(), -1);
189  patchFaceIDs_.setSize(fz.size(), -1);
190 
191  label nBlockedFaces = 0;
192  forAll(fz, localFacei)
193  {
194  const label facei = fz[localFacei];
195 
196  if (mesh_.isInternalFace(facei))
197  {
198  if (alphaf[facei] > alphaThreshold_)
199  {
200  blockedFaces[localFacei] = true;
201  }
202  }
203  else
204  {
205  label patchi = mesh_.boundaryMesh().whichPatch(facei);
206  label patchFacei = -1;
207 
208  const polyPatch& pp = mesh_.boundaryMesh()[patchi];
209  const scalarField& alphafp = alphaf.boundaryField()[patchi];
210  const auto* cpp = isA<coupledPolyPatch>(pp);
211 
212  if (cpp)
213  {
214  patchFacei = (cpp->owner() ? pp.whichFace(facei) : -1);
215  }
216  else if (!isA<emptyPolyPatch>(pp))
217  {
218  patchFacei = pp.whichFace(facei);
219  }
220 
221  if (patchFacei == -1)
222  {
223  patchi = -1;
224  }
225  else if (alphafp[patchFacei] > alphaThreshold_)
226  {
227  blockedFaces[localFacei] = true;
228  }
229 
230  patchIDs_[localFacei] = patchi;
231  patchFaceIDs_[localFacei] = patchFacei;
232  }
233  }
234 
235  DebugInFunction << "Number of blocked faces: " << nBlockedFaces << endl;
236 }
237 
238 
240 (
241  const scalar time,
242  const label regioni
243 )
244 {
245  DebugInFunction << "collectParticle: " << regioni << endl;
246 
247  const label particlei = regionToParticleMap_[regioni];
248  eulerianParticle p = particles_[particlei];
249 
250  if (p.faceIHit != -1 && nInjectorLocations_)
251  {
252  // Use coarse face index for tag output
253  label coarseFacei = fineToCoarseAddr_[p.faceIHit];
254  p.faceIHit = globalCoarseFaces_.toGlobal(coarseFacei);
255  }
256 
257  reduce(p, sumParticleOp<eulerianParticle>());
258 
259  const scalar pDiameter = cbrt(6.0*p.V/constant::mathematical::pi);
260 
261  if ((pDiameter > minDiameter_) && (pDiameter < maxDiameter_))
262  {
263  if (Pstream::master())
264  {
265  const scalar d = cbrt(6.0*p.V/constant::mathematical::pi);
266  const point position = p.VC/(p.V + ROOTVSMALL);
267  const vector U = p.VU/(p.V + ROOTVSMALL);
268  label tag = -1;
269  if (nInjectorLocations_)
270  {
271  tag = p.faceIHit;
272  }
273 
274  injectedParticle* ip = new injectedParticle
275  (
276  mesh_,
277  position,
278  tag,
279  time,
280  d,
281  U,
282  false // not looking to set cell owner etc.
283  );
284 
285  cloud_.addParticle(ip);
286 
287  collectedVolume_ += p.V;
288  }
289 
290  ++nCollectedParticles_;
291  }
292  else
293  {
294  // Discard particles over/under diameter threshold
295  ++nDiscardedParticles_;
296  discardedVolume_ += p.V;
297  }
298 }
299 
300 
302 (
303  const label nNewRegions,
304  const scalar time,
305  labelList& regionFaceIDs
306 )
307 {
309 
310  // Determine mapping between old and new regions so that we can
311  // accumulate particle info
312  labelList oldToNewRegion(particles_.size(), -1);
313  labelList newToNewRegion(identity(nNewRegions));
314 
315  forAll(regionFaceIDs, facei)
316  {
317  label newRegioni = regionFaceIDs[facei];
318  label oldRegioni = regions0_[facei];
319 
320  if (newRegioni != -1 && oldRegioni != -1)
321  {
322  // If old region has split into multiple regions we need to
323  // renumber new regions to maintain connectivity with old regions
324  newToNewRegion[newRegioni] =
325  max(newRegioni, oldToNewRegion[oldRegioni]);
326  oldToNewRegion[oldRegioni] = newRegioni;
327  }
328  }
329 
330  // Create map from new regions to slots in particles list
331  // - filter through new-to-new addressing to identify new particles
332  Pstream::listCombineReduce(newToNewRegion, maxEqOp<label>());
333 
334  label nParticle = -1;
335  labelHashSet newRegions;
336  Map<label> newRegionToParticleMap;
337  forAll(newToNewRegion, newRegioni0)
338  {
339  label newRegioni = newToNewRegion[newRegioni0];
340  if (newRegions.insert(newRegioni))
341  {
342  ++nParticle;
343  }
344 
345  // New particle slot
346  newRegionToParticleMap.insert(newRegioni0, nParticle);
347  }
348 
349  // Accumulate old region data or create a new particle if there is no
350  // mapping from the old-to-new region
351  Pstream::listCombineReduce(oldToNewRegion, maxEqOp<label>());
352 
353  List<eulerianParticle> newParticles(newRegionToParticleMap.size());
354  forAll(oldToNewRegion, oldRegioni)
355  {
356  label newRegioni = oldToNewRegion[oldRegioni];
357  if (newRegioni == -1)
358  {
359  // No mapping from old-to-new - collect new particle
360  DebugInfo
361  << "Collecting particle from oldRegion:" << oldRegioni
362  << endl;
363 
364  collectParticle(time, oldRegioni);
365  }
366  else
367  {
368  // Combine existing particle into new particle
369  label newParticlei = newRegionToParticleMap[newRegioni];
370  label oldParticlei = regionToParticleMap_[oldRegioni];
371 
372  DebugInfo
373  << "Combining newRegioni: " << newRegioni
374  << "(p:" << newParticlei << ") and "
375  << "oldRegioni: " << oldRegioni
376  << "(p:" << oldParticlei << ")"
377  << endl;
378 
379  newParticles[newParticlei] =
380  sumParticleOp<eulerianParticle>()
381  (
382  newParticles[newParticlei],
383  particles_[oldParticlei]
384  );
385  }
386  }
387 
388  // Reset the particles list and addressing for latest available info
389  particles_.transfer(newParticles);
390  regionToParticleMap_ = newRegionToParticleMap;
391 
392  // Reset the region IDs for the next integration step
393  // - these become the oldRegioni's
394  regions0_ = regionFaceIDs;
395 }
396 
397 
399 (
400  const surfaceScalarField& alphaf,
401  const surfaceScalarField& phi,
402  const labelList& regionFaceIDs,
403  const faceZone& fz
404 )
405 {
407 
408  const volVectorField& U = mesh_.lookupObject<volVectorField>(UName_);
410 
411  const scalar deltaT = mesh_.time().deltaTValue();
412  const pointField& faceCentres = mesh_.faceCentres();
413 
414  forAll(regionFaceIDs, localFacei)
415  {
416  const label newRegioni = regionFaceIDs[localFacei];
417  if (newRegioni != -1)
418  {
419  const label particlei = regionToParticleMap_[newRegioni];
420  const label meshFacei = fz[localFacei];
421  eulerianParticle& p = particles_[particlei];
422 
423  if (p.faceIHit < 0)
424  {
425  // New particle - does not exist in particles_ list
426  p.faceIHit = localFacei;
427  p.V = 0;
428  p.VC = vector::zero;
429  p.VU = vector::zero;
430  }
431 
432  // Accumulate particle properties
433  scalar magPhii = mag(faceValue(phi, localFacei, meshFacei));
434  vector Ufi = faceValue(Uf, localFacei, meshFacei);
435  scalar dV = magPhii*deltaT;
436  p.V += dV;
437  p.VC += dV*faceCentres[meshFacei];
438  p.VU += dV*Ufi;
439  }
440  }
441 }
442 
443 
444 // * * * * * * * * * * * * * * * * Constructor * * * * * * * * * * * * * * * //
445 
447 (
448  const word& name,
449  const Time& runTime,
450  const dictionary& dict
451 )
452 :
453  fvMeshFunctionObject(name, runTime, dict),
454  writeFile(runTime, name),
455  cloud_(mesh_, "eulerianParticleCloud"),
456  faceZoneName_(),
457  zoneID_(-1),
458  patchIDs_(),
459  patchFaceIDs_(),
460  alphaName_("alpha"),
461  alphaThreshold_(0.1),
462  UName_("U"),
463  rhoName_("rho"),
464  phiName_("phi"),
465  nInjectorLocations_(0),
466  fineToCoarseAddr_(),
467  globalCoarseFaces_(),
468  regions0_(),
469  particles_(),
470  regionToParticleMap_(),
471  minDiameter_(ROOTVSMALL),
472  maxDiameter_(GREAT),
473  nCollectedParticles_(getProperty<label>("nCollectedParticles", 0)),
474  collectedVolume_(getProperty<scalar>("collectedVolume", 0)),
475  nDiscardedParticles_(getProperty<label>("nDiscardedParticles", 0)),
476  discardedVolume_(getProperty<scalar>("discardedVolume", 0))
477 {
478  if (mesh_.nSolutionD() != 3)
479  {
481  << name << " function object only applicable to 3-D cases"
482  << exit(FatalError);
483  }
484 
486 }
487 
488 
489 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
490 
492 (
493  const dictionary& dict
494 )
495 {
497 
499  {
500  dict.readEntry("faceZone", faceZoneName_);
501  dict.readEntry("alpha", alphaName_);
502 
503  dict.readIfPresent("alphaThreshold", alphaThreshold_);
504  dict.readIfPresent("U", UName_);
505  dict.readIfPresent("rho", rhoName_);
506  dict.readIfPresent("phi", phiName_);
507  dict.readIfPresent("nLocations", nInjectorLocations_);
508  dict.readIfPresent("minDiameter", minDiameter_);
509  dict.readIfPresent("maxDiameter", maxDiameter_);
510 
511  checkFaceZone();
512 
513  if (nInjectorLocations_)
514  {
515  initialiseBins();
516  }
517 
518  return true;
519  }
520 
521  return false;
522 }
523 
524 
526 {
528 
529  Log << type() << " " << name() << " output:" << nl;
530 
531  const volScalarField& alpha =
532  mesh_.lookupObject<volScalarField>(alphaName_);
533 
534  const surfaceScalarField alphaf
535  (
536  typeName + ":alphaf",
538  );
539 
540  const faceZone& fz = mesh_.faceZones()[zoneID_];
542  (
543  IndirectList<face>(mesh_.faces(), fz),
544  mesh_.points()
545  );
546 
547  // Set the blocked faces, i.e. where alpha > alpha threshold value
548  boolList blockedFaces(fz.size(), false);
549  setBlockedFaces(alphaf, fz, blockedFaces);
550 
551  // Split the faceZone according to the blockedFaces
552  // - Returns a list of (disconnected) region index per face zone face
553  regionSplit2D regionFaceIDs(mesh_, patch, blockedFaces);
554 
555  // Global number of regions
556  const label nRegionsNew = regionFaceIDs.nRegions();
557 
558  // Calculate the addressing between the old and new region information
559  // Also collects particles that have traversed the faceZone
560  // - Note: may also update regionFaceIDs
561  calculateAddressing
562  (
563  nRegionsNew,
564  mesh_.time().value(),
565  regionFaceIDs
566  );
567 
568  // Process latest region information
569  tmp<surfaceScalarField> tphi = phiU();
570  accumulateParticleInfo(alphaf, tphi(), regionFaceIDs, fz);
571 
572  Log << " Collected particles : " << nCollectedParticles_ << nl
573  << " Collected volume : " << collectedVolume_ << nl
574  << " Discarded particles : " << nDiscardedParticles_ << nl
575  << " Discarded volume : " << discardedVolume_ << nl
576  << " Particles in progress : " << particles_.size() << nl
577  << endl;
578 
579  return true;
580 }
581 
582 
584 {
586 
587  cloud_.write();
588 
589  setProperty("nCollectedParticles", nCollectedParticles_);
590  setProperty("collectedVolume", collectedVolume_);
591  setProperty("nDiscardedParticles", nDiscardedParticles_);
592  setProperty("discardedVolume", discardedVolume_);
593 
594  return true;
595 }
596 
597 
598 // ************************************************************************* //
Foam::surfaceFields.
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
virtual void checkFaceZone()
Check that the faceZone is valid.
defineTypeNameAndDebug(ObukhovLength, 0)
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...
virtual void collectParticle(const scalar time, const label regioni)
Collect particles that have passed through the faceZone.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
Lightweight class to store particle data derived from VOF calculations, with special handling for inp...
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
Splits a patch into regions based on a mask field. Result is a globally consistent label list of regi...
Definition: regionSplit2D.H:55
labelList regions0_
Region indices in faceZone faces from last iteration.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1049
extractEulerianParticles(const word &name, const Time &runTime, const dictionary &dict)
Construct from components.
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Macros for easy insertion into run-time selection tables.
const word & name() const noexcept
Return the name of this functionObject.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
A list of faces which address into the list of points.
void setSize(const label n)
Alias for resize()
Definition: List.H:316
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
virtual bool read(const dictionary &)
Read the field min/max data.
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Definition: ZoneMesh.C:629
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
Definition: labelLists.C:44
virtual void initialiseBins()
Initialise the particle collection bins.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
#define DebugInFunction
Report an information message using Foam::Info.
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
dimensionedScalar cbrt(const dimensionedScalar &ds)
virtual const word & type() const =0
Runtime type information.
virtual void accumulateParticleInfo(const surfaceScalarField &alphaf, const surfaceScalarField &phi, const labelList &regionFaceIDs, const faceZone &fz)
Process latest region information.
constexpr scalar pi(M_PI)
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
Vector< scalar > vector
Definition: vector.H:57
autoPtr< surfaceVectorField > Uf
#define DebugInfo
Report an information message using Foam::Info.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
int debug
Static debugging option.
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:670
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
virtual bool read(const dictionary &dict)
Read.
Definition: writeFile.C:241
label nRegions() const noexcept
Return the global number of regions.
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.
virtual void setBlockedFaces(const surfaceScalarField &alphaf, const faceZone &fz, boolList &blockedFaces)
Set the blocked faces, i.e. where alpha > alpha threshold value.
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
Definition: polyMesh.C:893
label nInjectorLocations_
Number of sample locations to generate.
U
Definition: pEqn.H:72
vector point
Point is a vector.
Definition: point.H:37
virtual tmp< surfaceScalarField > phiU() const
Return the volumetric flux.
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
#define Log
Definition: PDRblock.C:28
const std::string patch
OpenFOAM patch number as a std::string.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
wordList names() const
A list of the zone names.
Definition: ZoneMesh.C:362
messageStream Info
Information stream (stdout output on master, null elsewhere)
constexpr label labelMax
Definition: label.H:55
virtual bool read(const dictionary &dict)
Read optional controls.
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:60
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
A List with indirect addressing.
Definition: IndirectList.H:60
const fvMesh & mesh_
Reference to the fvMesh.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
static void listCombineReduce(List< T > &values, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
After completion all processors have the same data.