shortestPathSet.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) 2017-2019 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 
28 #include "shortestPathSet.H"
29 #include "meshSearch.H"
30 #include "DynamicList.H"
31 #include "topoDistanceData.H"
33 #include "FaceCellWave.H"
34 #include "syncTools.H"
35 #include "fvMesh.H"
36 #include "volFields.H"
37 #include "OBJstream.H"
38 #include "PatchTools.H"
39 #include "foamVtkSurfaceWriter.H"
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45  defineTypeNameAndDebug(shortestPathSet, 0);
46  addToRunTimeSelectionTable(sampledSet, shortestPathSet, word);
47 }
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
52 Foam::label Foam::shortestPathSet::findMinFace
53 (
54  const polyMesh& mesh,
55  const label cellI,
56  const List<topoDistanceData>& allFaceInfo,
57  const PackedBoolList& isLeakPoint,
58  const bool distanceMode,
59  const point& origin
60 )
61 {
62  const cell& cFaces2 = mesh.cells()[cellI];
63 
64  // 1. Get topologically nearest face
65 
66  label minDist = labelMax;
67  label minFaceI = -1;
68  label nMin = 0;
69  forAll(cFaces2, i)
70  {
71  label faceI = cFaces2[i];
72  const topoDistanceData& info = allFaceInfo[faceI];
73  if (info.distance() < minDist)
74  {
75  minDist = info.distance();
76  minFaceI = faceI;
77  nMin = 1;
78  }
79  else if (info.distance() == minDist)
80  {
81  nMin++;
82  }
83  }
84 
85  if (nMin > 1)
86  {
87  // 2. Check all faces with minDist for minimum (or maximum)
88  // distance to origin
89  if (distanceMode)
90  {
91  scalar minDist2 = ROOTVGREAT;
92  forAll(cFaces2, i)
93  {
94  label faceI = cFaces2[i];
95  if (allFaceInfo[faceI].distance() == minDist)
96  {
97  scalar d2 = magSqr(mesh.faceCentres()[faceI]-origin);
98  if (d2 < minDist2)
99  {
100  minDist2 = d2;
101  minFaceI = faceI;
102  }
103  }
104  }
105  }
106  else
107  {
108  // Avoid leak points
109  label minLeakPoints = labelMax;
110  forAll(cFaces2, i)
111  {
112  label faceI = cFaces2[i];
113  if (allFaceInfo[faceI].distance() == minDist)
114  {
115  // Count number of leak points
116  label nLeak = 0;
117  {
118  const face& f = mesh.faces()[faceI];
119  forAll(f, fp)
120  {
121  if (isLeakPoint[f[fp]])
122  {
123  nLeak++;
124  }
125  }
126  }
127 
128  if (nLeak < minLeakPoints)
129  {
130  minLeakPoints = nLeak;
131  minFaceI = faceI;
132  }
133  }
134  }
135  }
136  }
137 
138  return minFaceI;
139 }
140 
141 
142 void Foam::shortestPathSet::calculateDistance
143 (
144  const label iter,
145  const polyMesh& mesh,
146  const label cellI,
147 
148  List<topoDistanceData>& allFaceInfo,
149  List<topoDistanceData>& allCellInfo
150 ) const
151 {
152  int dummyTrackData = 0;
153 
154  // Seed faces on cell1
155  DynamicList<topoDistanceData> faceDist;
156  DynamicList<label> cFaces1;
157 
158  if (cellI != -1)
159  {
160  const labelList& cFaces = mesh.cells()[cellI];
161  faceDist.reserve(cFaces.size());
162  cFaces1.reserve(cFaces.size());
163 
164  for (label facei : cFaces)
165  {
166  if (!allFaceInfo[facei].valid(dummyTrackData))
167  {
168  cFaces1.append(facei);
169  faceDist.append(topoDistanceData(123, 0));
170  }
171  }
172  }
173 
174 
175 
176  // Walk through face-cell wave till all cells are reached
177  FaceCellWave
178  <
179  topoDistanceData
180  > wallDistCalc
181  (
182  mesh,
183  cFaces1,
184  faceDist,
185  allFaceInfo,
186  allCellInfo,
187  mesh.globalData().nTotalCells()+1 // max iterations
188  );
189 
190 
191  if (debug)
192  {
193  const fvMesh& fm = refCast<const fvMesh>(mesh);
194 
195  const_cast<fvMesh&>(fm).setInstance(fm.time().timeName());
196  fm.write();
198  (
199  IOobject
200  (
201  "allCellInfo" + Foam::name(iter),
202  fm.time().timeName(),
203  fm,
206  ),
207  fm,
209  );
210  forAll(allCellInfo, celli)
211  {
212  fld[celli] = allCellInfo[celli].distance();
213  }
214  forAll(fld.boundaryField(), patchi)
215  {
216  const polyPatch& pp = mesh.boundaryMesh()[patchi];
217  SubList<topoDistanceData> p(pp.patchSlice(allFaceInfo));
218  scalarField pfld(fld.boundaryField()[patchi].size());
219  forAll(pfld, i)
220  {
221  pfld[i] = 1.0*p[i].distance();
222  }
223  fld.boundaryFieldRef()[patchi] == pfld;
224  }
225  //Note: do not swap cell values so do not do
226  //fld.correctBoundaryConditions();
227  fld.write();
228  }
229 }
230 
231 
232 void Foam::shortestPathSet::sync
233 (
234  const polyMesh& mesh,
235  PackedBoolList& isLeakFace,
236  PackedBoolList& isLeakPoint,
237  const label celli,
238  point& origin,
239  bool& findMinDistance
240 ) const
241 {
243  (
244  mesh,
245  isLeakPoint,
246  orEqOp<unsigned int>(),
247  0u
248  );
250  (
251  mesh,
252  isLeakFace,
253  orEqOp<unsigned int>()
254  );
255  // Sync origin, findMinDistance
256  {
257  typedef Tuple2<label, Tuple2<point, bool>> ProcData;
258 
259  ProcData searchData
260  (
261  celli,
262  Tuple2<point, bool>(origin, findMinDistance)
263  );
265  (
266  searchData,
267  [](ProcData& x, const ProcData& y)
268  {
269  if (y.first() != -1)
270  {
271  x = y;
272  }
273  }
274  );
275  origin = searchData.second().first();
276  findMinDistance = searchData.second().second();
277  }
278 }
279 
280 
281 bool Foam::shortestPathSet::touchesWall
282 (
283  const polyMesh& mesh,
284  const label facei,
285 
286  PackedBoolList& isLeakFace,
287  const PackedBoolList& isLeakPoint
288 ) const
289 {
290  // Check if facei touches leakPoint
291  const face& f = mesh.faces()[facei];
292  forAll(f, fp)
293  {
294  if (isLeakPoint[f[fp]])
295  {
296  if (isLeakFace.set(facei))
297  {
298  return true;
299  }
300  }
301  }
302 
303  return false;
304 }
305 
306 
307 void Foam::shortestPathSet::genSamples
308 (
309  const bool markLeakPath,
310  const label maxIter,
311  const polyMesh& mesh,
312  const boolList& isBlockedFace,
313  const point& insidePoint,
314  const label insideCelli,
315  const point& outsidePoint,
316 
317  DynamicList<point>& samplingPts,
318  DynamicList<label>& samplingCells,
319  DynamicList<label>& samplingFaces,
320  DynamicList<label>& samplingSegments,
321  DynamicList<scalar>& samplingCurveDist,
322  PackedBoolList& isLeakCell,
323  PackedBoolList& isLeakFace,
324  PackedBoolList& isLeakPoint
325 ) const
326 {
327  const topoDistanceData maxData(labelMax, labelMax);
328 
329  // Get the target point
330  label outsideCelli = mesh.findCell(outsidePoint);
331 
332  // Maintain overall track length. Used to make curveDist continuous.
333  label trackLength = 0;
334 
335  for (label iter = 0; iter < maxIter; iter++)
336  {
337  List<topoDistanceData> allFaceInfo(mesh.nFaces());
338  List<topoDistanceData> allCellInfo(mesh.nCells());
339 
340  // Mark blocked faces with high distance
341  forAll(isBlockedFace, facei)
342  {
343  if (isBlockedFace[facei])
344  {
345  allFaceInfo[facei] = maxData;
346  }
347  }
348 
349  if (markLeakPath)
350  {
351  // Mark any previously found leak path. This marks all
352  // faces of all cells on the path. This will make sure that
353  // ultimately the path will go through another gap.
354  forAll(isLeakCell, celli)
355  {
356  if (celli != insideCelli && celli != outsideCelli)
357  {
358  if (isLeakCell[celli])
359  {
360  allCellInfo[celli] = maxData;
361  const cell& cFaces = mesh.cells()[celli];
362  for (auto facei : cFaces)
363  {
364  allFaceInfo[facei] = maxData;
365  }
366  }
367  }
368  }
369  }
370 
371  // Mark any previously found leak faces. These are faces that
372  // are (point-)connected to an existing boundary.
373  forAll(isLeakFace, facei)
374  {
375  if (isLeakFace[facei])
376  {
377  allFaceInfo[facei] = maxData;
378  }
379  }
380 
381 
382  // Pass1: Get distance to insideCelli
383 
384  calculateDistance(iter, mesh, insideCelli, allFaceInfo, allCellInfo);
385 
386 
387 
388  // Pass2: walk from outside points backwards. Note: could be done
389  // using FaceCellWave as well but is overly complex since
390  // does not allow logic comparing all faces of a cell.
391 
392  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
393 
394 
395  bool targetFound = false;
396  if (outsideCelli != -1)
397  {
398  int dummyTrackData;
399  targetFound = allCellInfo[outsideCelli].valid(dummyTrackData);
400  if (!targetFound)
401  {
403  << "Point " << outsidePoint
404  << " not reachable by walk from " << insidePoint
405  << ". Probably mesh has island/regions."
406  << " Skipped route detection." << endl;
407  }
408  }
409  reduce(targetFound, orOp<bool>());
410  if (!targetFound)
411  {
412  break;
413  }
414 
415 
416  // Start with given target cell and walk back
417  // If point happens to be on multiple processors, random pick
418  label frontCellI = outsideCelli;
419  point origin(outsidePoint);
420  bool findMinDistance = true;
421 
422  while (true)
423  {
424  label frontFaceI = -1;
425 
426  // Work within same processor
427  if (frontCellI != -1)
428  {
429  // Find face with lowest distance from seed
430  // x | x 2 1 2 2 | x | x
431  // --- + --- + -1- + -2- + --- + ---
432  // x | 1 1 0 1 1 | x | x
433  // --- + --- + -1- + -2- + --- + ---
434  // x | x 2 1 2 2 3 3 4 4
435  // --- + --- + --- + -3- + -4- + -5-
436  // x | x 3 2 3 3 4 4 5 5
437  // e.g. if we start from cell with value = 4, we have
438  // neighbour faces 4, 4, 5, 5. Choose 4 (least distance
439  // to seed) and continue...
440 
441  frontFaceI = findMinFace
442  (
443  mesh,
444  frontCellI,
445  allFaceInfo,
446  isLeakPoint,
447  findMinDistance, // mode : find min or find max
448  origin
449  );
450 
451 
452  // Loop until we hit a boundary face
453  PackedBoolList isNewLeakPoint(isLeakPoint);
454  while (mesh.isInternalFace(frontFaceI))
455  {
456  if (isBlockedFace.size() && isBlockedFace[frontFaceI])
457  {
458  // Pout<< " Found blocked face" << endl;
459  frontCellI = -1;
460  break;
461  }
462 
463  // Step to neighbouring cell
464  label nbrCellI = mesh.faceOwner()[frontFaceI];
465  if (nbrCellI == frontCellI)
466  {
467  nbrCellI = mesh.faceNeighbour()[frontFaceI];
468  }
469 
470  if (nbrCellI == insideCelli)
471  {
472  // Pout<< " Found connection seed cell!" << endl;
473  frontCellI = -1;
474  break;
475  }
476 
477  frontCellI = nbrCellI;
478 
479  // Pick best face on cell
480  frontFaceI = findMinFace
481  (
482  mesh,
483  frontCellI,
484  allFaceInfo,
485  isLeakPoint,
486  findMinDistance,
487  origin
488  );
489 
490  const topoDistanceData& cInfo = allCellInfo[frontCellI];
491  const topoDistanceData& fInfo = allFaceInfo[frontFaceI];
492 
493  if (fInfo.distance() <= cInfo.distance())
494  {
495  samplingPts.append(mesh.cellCentres()[frontCellI]);
496  samplingCells.append(frontCellI);
497  samplingFaces.append(-1);
498  samplingSegments.append(iter);
499  samplingCurveDist.append
500  (
501  trackLength+cInfo.distance()
502  );
503  isLeakCell.set(frontCellI);
504 
505  if
506  (
507  touchesWall
508  (
509  mesh,
510  frontFaceI,
511  isLeakFace,
512  isLeakPoint
513  )
514  )
515  {
516  isNewLeakPoint.set(mesh.faces()[frontFaceI]);
517  origin = mesh.faceCentres()[frontFaceI];
518  findMinDistance = false;
519  }
520 
521  }
522  }
523  isLeakPoint.transfer(isNewLeakPoint);
524  }
525 
526  // Situation 1: we found the destination cell (do nothing),
527  // frontCellI is -1 on all processors
528  if (returnReduce(frontCellI == -1, andOp<bool>()))
529  {
530  break;
531  }
532 
533  // Situation 2: we're on a coupled patch and might need to
534  // switch processor/cell. We need to transfer:
535  // -frontface -frontdistance -leak points/faces
536  boolList isFront(mesh.nBoundaryFaces(), false);
537 
538  if (frontFaceI != -1)
539  {
540  isFront[frontFaceI-mesh.nInternalFaces()] = true;
541  }
543 
544  label minCellDistance = labelMax;
545  if (frontCellI != -1)
546  {
547  minCellDistance = allCellInfo[frontCellI].distance();
548  }
549  reduce(minCellDistance, minOp<label>());
550 
551  // Sync all leak data
552  sync
553  (
554  mesh,
555  isLeakFace,
556  isLeakPoint,
557  frontCellI,
558  origin,
559  findMinDistance
560  );
561 
562 
563 
564  const label oldFrontFaceI = frontFaceI;
565  frontCellI = -1;
566  frontFaceI = -1;
567  forAll(pbm, patchI)
568  {
569  const polyPatch& pp = pbm[patchI];
570  forAll(pp, i)
571  {
572  label faceI = pp.start()+i;
573  if
574  (
575  oldFrontFaceI == -1
576  && isFront[faceI-mesh.nInternalFaces()]
577  && (isBlockedFace.empty() || !isBlockedFace[faceI])
578  )
579  {
580  frontFaceI = faceI;
581  frontCellI = pp.faceCells()[i];
582  break;
583  }
584  }
585 
586  if
587  (
588  frontCellI != -1
589  && allCellInfo[frontCellI].distance() < minCellDistance
590  )
591  {
592  const topoDistanceData& cInfo = allCellInfo[frontCellI];
593 
594  samplingPts.append(mesh.cellCentres()[frontCellI]);
595  samplingCells.append(frontCellI);
596  samplingFaces.append(-1);
597  samplingSegments.append(iter);
598  samplingCurveDist.append
599  (
600  trackLength+cInfo.distance()
601  );
602  isLeakCell.set(frontCellI);
603 
604  // Check if frontFacei touches leakPoint
605  if
606  (
607  touchesWall
608  (
609  mesh,
610  frontFaceI,
611  isLeakFace,
612  isLeakPoint
613  )
614  )
615  {
616  isLeakPoint.set(mesh.faces()[frontFaceI]);
617  origin = mesh.faceCentres()[frontFaceI];
618  findMinDistance = false;
619  }
620 
621  break;
622  }
623 
624  // When seed cell is isolated by processor boundaries
625  if (insideCelli != -1 && frontCellI == insideCelli)
626  {
627  // Pout<< " Found connection seed cell!" << endl;
628  frontCellI = -1;
629  break;
630  }
631  }
632 
633  // Sync all leak data
634  sync
635  (
636  mesh,
637  isLeakFace,
638  isLeakPoint,
639  frontCellI,
640  origin,
641  findMinDistance
642  );
643  }
644 
645 
646  // Recalculate the overall trackLength
647  trackLength = returnReduce
648  (
649  (
650  samplingCurveDist.size()
651  ? samplingCurveDist.last()
652  : 0
653  ),
654  maxOp<scalar>()
655  );
656 
657 
658  if (debug)
659  {
660  const fvMesh& fm = refCast<const fvMesh>(mesh);
661 
662  const_cast<fvMesh&>(fm).setInstance(fm.time().timeName());
664  (
665  IOobject
666  (
667  "isLeakCell",
668  fm.time().timeName(),
669  fm,
672  ),
673  fm,
675  );
676  forAll(isLeakCell, celli)
677  {
678  fld[celli] = isLeakCell[celli];
679  }
680  fld.correctBoundaryConditions();
681  fld.write();
682  }
683  }
684 }
685 
686 
687 void Foam::shortestPathSet::genSamples
688 (
689  const bool markLeakPath,
690  const label maxIter,
691  const polyMesh& mesh,
692  const labelUList& wallPatches,
693  const boolList& isBlockedFace
694 )
695 {
696  // Storage for sample points
697  DynamicList<point> samplingPts;
698  DynamicList<label> samplingCells;
699  DynamicList<label> samplingFaces;
700  DynamicList<label> samplingSegments;
701  DynamicList<scalar> samplingCurveDist;
702 
703  // Seed faces and points on 'real' boundary
704  PackedBoolList isLeakFace(mesh.nFaces());
705  PackedBoolList isLeakPoint(mesh.nPoints());
706  {
707  // Real boundaries
708  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
709  for (label patchi : wallPatches)
710  {
711  const polyPatch& pp = pbm[patchi];
712  forAll(pp, i)
713  {
714  isLeakPoint.set(pp[i]);
715  }
716  }
717 
718  // Meshed boundaries
719  forAll(isBlockedFace, facei)
720  {
721  if (isBlockedFace[facei])
722  {
723  isLeakPoint.set(mesh.faces()[facei]);
724  }
725  }
726 
728  (
729  mesh,
730  isLeakPoint,
731  orEqOp<unsigned int>(),
732  0u
733  );
734  }
735 
736 
737  // All cells along leak paths
738  PackedBoolList isLeakCell(mesh.nCells());
739 
740  label prevSegmenti = 0;
741  scalar prevDistance = 0.0;
742 
743  for (auto insidePoint : insidePoints_)
744  {
745  const label insideCelli = mesh.findCell(insidePoint);
746 
747  for (auto outsidePoint : outsidePoints_)
748  {
749  const label nOldSamples = samplingSegments.size();
750 
751  // Append any leak path to sampling*
752  genSamples
753  (
754  markLeakPath,
755  maxIter,
756  mesh,
757  isBlockedFace,
758  insidePoint,
759  insideCelli,
760  outsidePoint,
761 
762  samplingPts,
763  samplingCells,
764  samplingFaces,
765  samplingSegments,
766  samplingCurveDist,
767  isLeakCell,
768  isLeakFace,
769  isLeakPoint
770  );
771 
772  // Make segment, distance consecutive
773  label maxSegment = 0;
774  scalar maxDistance = 0.0;
775  for (label i = nOldSamples; i < samplingSegments.size(); ++i)
776  {
777  samplingSegments[i] += prevSegmenti;
778  maxSegment = max(maxSegment, samplingSegments[i]);
779  samplingCurveDist[i] += prevDistance;
780  maxDistance = max(maxDistance, samplingCurveDist[i]);
781  }
782  prevSegmenti = returnReduce(maxSegment, maxOp<label>());
783  prevDistance = returnReduce(maxDistance, maxOp<scalar>());
784  }
785  }
786 
787  if (debug)
788  {
789  const faceList leakFaces(mesh.faces(), isLeakFace.sortedToc());
790 
791  OBJstream str(mesh.time().path()/"isLeakFace.obj");
792  str.write(leakFaces, mesh.points(), false);
793  Pout<< "Writing " << leakFaces.size() << " faces to " << str.name()
794  << endl;
795  }
796 
797 
798  samplingPts.shrink();
799  samplingCells.shrink();
800  samplingFaces.shrink();
801  samplingSegments.shrink();
802  samplingCurveDist.shrink();
803 
804  // Move into *this
805  setSamples
806  (
807  std::move(samplingPts),
808  std::move(samplingCells),
809  std::move(samplingFaces),
810  std::move(samplingSegments),
811  std::move(samplingCurveDist)
812  );
813 
814  if (debug)
815  {
816  write(Info);
817  }
818 }
819 
820 
821 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
822 
824 (
825  const word& name,
826  const polyMesh& mesh,
827  const meshSearch& searchEngine,
828  const word& axis,
829  const bool markLeakPath,
830  const label maxIter,
831  const labelUList& wallPatches,
832  const pointField& insidePoints,
833  const pointField& outsidePoints,
834  const boolList& isBlockedFace
835 )
836 :
837  sampledSet(name, mesh, searchEngine, axis),
838  insidePoints_(insidePoints),
839  outsidePoints_(outsidePoints)
840 {
841  if (debug)
842  {
843  fileName outputDir =
844  (
845  mesh.time().globalPath()
847  / mesh.pointsInstance()
848  );
849  outputDir.clean();
850 
851  Info<< "shortestPathSet : Writing blocked faces to "
852  << outputDir << endl;
853 
854  const indirectPrimitivePatch setPatch
855  (
857  (
858  mesh.faces(),
859  findIndices(isBlockedFace, true)
860  ),
861  mesh.points()
862  );
863 
864  if (Pstream::parRun())
865  {
866  // Topological merge
867  labelList pointToGlobal;
868  labelList uniqueMeshPointLabels;
870  autoPtr<globalIndex> globalFaces;
871  faceList mergedFaces;
872  pointField mergedPoints;
874  (
875  mesh,
876  setPatch.localFaces(),
877  setPatch.meshPoints(),
878  setPatch.meshPointMap(),
879 
880  pointToGlobal,
881  uniqueMeshPointLabels,
882  globalPoints,
883  globalFaces,
884 
885  mergedFaces,
886  mergedPoints
887  );
888 
889  // Write
890  if (Pstream::master())
891  {
893  (
894  mergedPoints,
895  mergedFaces,
896  (outputDir / "blockedFace"),
897  false // serial only - already merged
898  );
899 
900  writer.writeGeometry();
901  }
902  }
903  else
904  {
906  (
907  setPatch.localPoints(),
908  setPatch.localFaces(),
909  (outputDir / "blockedFace"),
910  false // serial only - redundant
911  );
912 
913  writer.writeGeometry();
914  }
915  }
916 
917  genSamples(markLeakPath, maxIter, mesh, wallPatches, isBlockedFace);
918 }
919 
920 
922 (
923  const word& name,
924  const polyMesh& mesh,
925  const meshSearch& searchEngine,
926  const dictionary& dict
927 )
928 :
929  sampledSet(name, mesh, searchEngine, dict),
930  insidePoints_(dict.get<pointField>("insidePoints")),
931  outsidePoints_(dict.get<pointField>("outsidePoints"))
932 {
933  const label maxIter(dict.lookupOrDefault<label>("maxIter", 50));
934  const bool markLeakPath(dict.lookupOrDefault<bool>("markLeakPath", true));
935 
936  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
937 
938  DynamicList<label> wallPatches(pbm.size());
939  forAll(pbm, patchi)
940  {
941  const polyPatch& pp = pbm[patchi];
942  if (!pp.coupled() && !isA<emptyPolyPatch>(pp))
943  {
944  wallPatches.append(patchi);
945  }
946  }
947 
948  genSamples(markLeakPath, maxIter, mesh, wallPatches, boolList());
949 }
950 
951 
952 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::sampledSet
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:83
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
volFields.H
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1038
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
topoDistanceData.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::globalPoints
Calculates points shared by more than two processor patches or cyclic patches.
Definition: globalPoints.H:102
Foam::labelMax
constexpr label labelMax
Definition: label.H:65
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::meshSearch
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition: meshSearch.H:60
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
PatchTools.H
Foam::DynamicList< label >
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const
Number of internal faces.
Definition: primitiveMeshI.H:78
Foam::primitiveMesh::nFaces
label nFaces() const
Number of mesh faces.
Definition: primitiveMeshI.H:90
Foam::polyPatch::coupled
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:327
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:414
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:138
Foam::globalMeshData::nTotalCells
label nTotalCells() const
Return total number of cells in decomposed mesh.
Definition: globalMeshData.H:394
Foam::combineReduce
void combineReduce(const List< UPstream::commsStruct > &comms, T &Value, const CombineOp &cop, const int tag, const label comm)
Definition: PstreamCombineReduceOps.H:54
Foam::syncTools::swapBoundaryFaceList
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues)
Swap coupled boundary face values. Uses eqOp.
Definition: syncTools.H:439
Foam::boolList
List< bool > boolList
A List of bools.
Definition: List.H:72
Foam::List::append
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:182
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:435
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::shortestPathSet::shortestPathSet
shortestPathSet(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const word &axis, const bool markLeakPath, const label maxIter, const labelUList &wallPatches, const pointField &insidePoints, const pointField &outsidePoints, const boolList &blockedFace)
Construct from components. blockedFace is an optional specification.
Definition: shortestPathSet.C:824
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
Foam::dictionary::lookupOrDefault
T lookupOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.H:1241
foamVtkSurfaceWriter.H
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
syncTools.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::syncTools::syncPointList
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
Definition: syncToolsTemplates.C:747
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::polyMesh::pointsInstance
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:815
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::syncTools::syncFaceList
static void syncFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop)
Synchronize values on all mesh faces.
Definition: syncTools.H:390
Foam::primitiveMesh::nCells
label nCells() const
Number of mesh cells.
Definition: primitiveMeshI.H:96
Foam::label
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:62
Foam::Field< vector >
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::OSstream::name
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:91
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1076
Foam::IndirectList
A List with indirect addressing.
Definition: IndirectList.H:56
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::functionObject::outputPrefix
static word outputPrefix
Directory prefix.
Definition: functionObject.H:259
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
fld
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputLagrangian.H:23
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::primitiveMesh::nBoundaryFaces
label nBoundaryFaces() const
Number of boundary faces (== nFaces - nInternalFaces)
Definition: primitiveMeshI.H:84
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::writer
Base class for graphics format writing. Entry points are.
Definition: writer.H:80
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::distance
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:438
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::polyMesh::findCell
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1453
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:176
shortestPathSet.H
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1063
Foam::sampledSet::mesh
const polyMesh & mesh() const
Definition: sampledSet.H:281
meshSearch.H
Foam::Time::path
fileName path() const
Return path.
Definition: Time.H:303
f
labelList f(nPoints)
Foam::vtk::surfaceWriter
Write faces/points (optionally with fields) as a vtp file or a legacy vtk file.
Definition: foamVtkSurfaceWriter.H:68
Foam::faceList
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
Foam::List< bool >
Foam::primitiveMesh::isInternalFace
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
Definition: primitiveMeshI.H:102
Foam::primitiveMesh::faceCentres
const vectorField & faceCentres() const
Definition: primitiveMeshFaceCentresAndAreas.C:145
FaceCellWave.H
Foam::UList< label >
Foam::primitiveMesh::nPoints
label nPoints() const
Number of mesh points.
Definition: primitiveMeshI.H:37
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:35
Foam::TimePaths::globalPath
fileName globalPath() const
Return global path for the case.
Definition: TimePathsI.H:72
insidePoints
insidePoints((1 2 3))
Foam::fileName::clean
static bool clean(std::string &str)
Cleanup filename.
Definition: fileName.C:298
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:246
DynamicList.H
Foam::polyMesh::globalData
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1241
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::labelUList
UList< label > labelUList
A UList of labels.
Definition: UList.H:79
Foam::PatchTools::gatherAndMerge
static void gatherAndMerge(const scalar mergeDist, const PrimitivePatch< Face, FaceList, PointField, PointType > &p, Field< PointType > &mergedPoints, List< Face > &mergedFaces, labelList &pointMergeMap)
Gather points and faces onto master and merge into single patch.
Definition: PatchToolsGatherAndMerge.C:44
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
OBJstream.H
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1082
Foam::PackedBoolList
bitSet PackedBoolList
Compatibility name. Superseded (MAR-2018) by Foam::bitSet.
Definition: PackedBoolList.H:40
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:90