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-2022 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 #include "indirectPrimitivePatch.H"
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46  defineTypeNameAndDebug(shortestPathSet, 0);
47  addToRunTimeSelectionTable(sampledSet, shortestPathSet, word);
48 }
49 
50 
51 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52 
53 Foam::label Foam::shortestPathSet::findMinFace
54 (
55  const polyMesh& mesh,
56  const label cellI,
57  const List<topoDistanceData<label>>& allFaceInfo,
58  const bitSet& isLeakPoint,
59  const bool distanceMode,
60  const point& origin
61 )
62 {
63  const cell& cFaces2 = mesh.cells()[cellI];
64 
65  // 1. Get topologically nearest face
66 
67  label minDist = labelMax;
68  label minFaceI = -1;
69  label nMin = 0;
70  forAll(cFaces2, i)
71  {
72  label faceI = cFaces2[i];
73  const topoDistanceData<label>& info = allFaceInfo[faceI];
74  if (info.distance() < minDist)
75  {
76  minDist = info.distance();
77  minFaceI = faceI;
78  nMin = 1;
79  }
80  else if (info.distance() == minDist)
81  {
82  nMin++;
83  }
84  }
85 
86  if (nMin > 1)
87  {
88  // 2. Check all faces with minDist for minimum (or maximum)
89  // distance to origin
90  if (distanceMode)
91  {
92  scalar minDist2 = ROOTVGREAT;
93  forAll(cFaces2, i)
94  {
95  label faceI = cFaces2[i];
96  if (allFaceInfo[faceI].distance() == minDist)
97  {
98  scalar d2 = magSqr(mesh.faceCentres()[faceI]-origin);
99  if (d2 < minDist2)
100  {
101  minDist2 = d2;
102  minFaceI = faceI;
103  }
104  }
105  }
106  }
107  else
108  {
109  // Avoid leak points i.e. preferentially stay away from the wall
110  label minLeakPoints = labelMax;
111  forAll(cFaces2, i)
112  {
113  label faceI = cFaces2[i];
114  if (allFaceInfo[faceI].distance() == minDist)
115  {
116  // Count number of leak points
117  label nLeak = 0;
118  {
119  const face& f = mesh.faces()[faceI];
120  forAll(f, fp)
121  {
122  if (isLeakPoint[f[fp]])
123  {
124  nLeak++;
125  }
126  }
127  }
128 
129  if (nLeak < minLeakPoints)
130  {
131  minLeakPoints = nLeak;
132  minFaceI = faceI;
133  }
134  }
135  }
136  }
137  }
138 
139  return minFaceI;
140 }
141 
142 
143 void Foam::shortestPathSet::calculateDistance
144 (
145  const label iter,
146  const polyMesh& mesh,
147  const label cellI,
148 
149  List<topoDistanceData<label>>& allFaceInfo,
150  List<topoDistanceData<label>>& allCellInfo
151 ) const
152 {
153  int dummyTrackData = 0;
154 
155  // Seed faces on cell1
156  DynamicList<topoDistanceData<label>> faceDist;
157  DynamicList<label> cFaces1;
158 
159  if (cellI != -1)
160  {
161  const labelList& cFaces = mesh.cells()[cellI];
162  faceDist.reserve(cFaces.size());
163  cFaces1.reserve(cFaces.size());
164 
165  for (label facei : cFaces)
166  {
167  if (!allFaceInfo[facei].valid(dummyTrackData))
168  {
169  cFaces1.append(facei);
170  faceDist.append(topoDistanceData<label>(0, 123));
171  }
172  }
173  }
174 
175 
176 
177  // Walk through face-cell wave till all cells are reached
178  FaceCellWave
179  <
180  topoDistanceData<label>
181  > wallDistCalc
182  (
183  mesh,
184  cFaces1,
185  faceDist,
186  allFaceInfo,
187  allCellInfo,
188  mesh.globalData().nTotalCells()+1 // max iterations
189  );
190 
191 
192  if (debug)
193  {
194  const fvMesh& fm = refCast<const fvMesh>(mesh);
195 
196  //const_cast<fvMesh&>(fm).setInstance(fm.time().timeName());
197  //fm.write();
199  (
200  IOobject
201  (
202  "allCellInfo" + Foam::name(iter),
203  fm.time().timeName(),
204  fm,
207  ),
208  fm,
210  );
211  forAll(allCellInfo, celli)
212  {
213  fld[celli] = allCellInfo[celli].distance();
214  }
215  forAll(fld.boundaryField(), patchi)
216  {
217  const polyPatch& pp = mesh.boundaryMesh()[patchi];
218  SubList<topoDistanceData<label>> p(pp.patchSlice(allFaceInfo));
219  scalarField pfld(fld.boundaryField()[patchi].size());
220  forAll(pfld, i)
221  {
222  pfld[i] = 1.0*p[i].distance();
223  }
224  fld.boundaryFieldRef()[patchi] == pfld;
225  }
226  //Note: do not swap cell values so do not do
227  //fld.correctBoundaryConditions();
228  Pout<< "Writing distance field for initial cell " << cellI
229  << " to " << fld.objectPath() << endl;
230  fld.write();
231  }
232 }
233 
234 
235 void Foam::shortestPathSet::sync
236 (
237  const polyMesh& mesh,
238  bitSet& isLeakFace,
239  bitSet& isLeakPoint,
240  const label celli,
241  point& origin,
242  bool& findMinDistance
243 ) const
244 {
246  (
247  mesh,
248  isLeakPoint,
249  orEqOp<unsigned int>(),
250  0u
251  );
253  (
254  mesh,
255  isLeakFace,
256  orEqOp<unsigned int>()
257  );
258  // Sync origin, findMinDistance
259  {
260  typedef Tuple2<label, Tuple2<point, bool>> ProcData;
261 
262  ProcData searchData
263  (
264  celli,
265  Tuple2<point, bool>(origin, findMinDistance)
266  );
268  (
269  searchData,
270  [](ProcData& x, const ProcData& y)
271  {
272  if (y.first() != -1)
273  {
274  x = y;
275  }
276  }
277  );
278  origin = searchData.second().first();
279  findMinDistance = searchData.second().second();
280  }
281 }
282 
283 
284 bool Foam::shortestPathSet::touchesWall
285 (
286  const polyMesh& mesh,
287  const label facei,
288 
289  bitSet& isLeakFace,
290  const bitSet& isLeakPoint
291 ) const
292 {
293  // Check if facei touches leakPoint
294  const face& f = mesh.faces()[facei];
295  forAll(f, fp)
296  {
297  const label nextFp = f.fcIndex(fp);
298 
299  if (isLeakPoint[f[fp]] && isLeakPoint[f[nextFp]]) // edge on boundary
300  //if (isLeakPoint[f[fp]]) // point on boundary
301  {
302  if (isLeakFace.set(facei))
303  {
304  return true;
305  }
306  }
307  }
308 
309  return false;
310 }
311 
312 
313 Foam::bitSet Foam::shortestPathSet::pathFaces
314 (
315  const polyMesh& mesh,
316  const bitSet& isLeakCell
317 ) const
318 {
319  // Calculates the list of faces inbetween leak cells.
320  // Does not account for the fact that the cells might be from different
321  // paths...
322 
323  const polyBoundaryMesh& patches = mesh.boundaryMesh();
324  const labelList& own = mesh.faceOwner();
325  const labelList& nbr = mesh.faceNeighbour();
326 
327  // Get remote value of bitCell. Note: keep uncoupled boundary faces false.
328  boolList nbrLeakCell(mesh.nBoundaryFaces(), false);
329  {
330  for (const polyPatch& pp : patches)
331  {
332  if (pp.coupled())
333  {
334  label bFacei = pp.start()-mesh.nInternalFaces();
335 
336  const labelUList& faceCells = pp.faceCells();
337 
338  for (const label celli : faceCells)
339  {
340  nbrLeakCell[bFacei] = isLeakCell[celli];
341  ++bFacei;
342  }
343  }
344  }
345 
347  (
348  mesh,
349  nbrLeakCell
350  );
351  }
352 
353 
354  bitSet isLeakFace(mesh.nFaces());
355 
356  // Internal faces
357  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
358  {
359  if (isLeakCell[own[facei]] && isLeakCell[nbr[facei]])
360  {
361  isLeakFace.set(facei);
362  }
363  }
364  // Boundary faces
365  for (label facei = mesh.nInternalFaces(); facei < mesh.nFaces(); facei++)
366  {
367  if (isLeakCell[own[facei]] && nbrLeakCell[facei-mesh.nInternalFaces()])
368  {
369  isLeakFace.set(facei);
370  }
371  }
372  return isLeakFace;
373 }
374 
375 
376 bool Foam::shortestPathSet::genSingleLeakPath
377 (
378  const bool markLeakPath,
379  const label iter,
380  const polyMesh& mesh,
381  const bitSet& isBlockedFace,
382  const point& insidePoint,
383  const label insideCelli,
384  const point& outsidePoint,
385  const label outsideCelli,
386 
387  // Generated sampling points
388  const scalar trackLength,
389  DynamicList<point>& samplingPts,
390  DynamicList<label>& samplingCells,
391  DynamicList<label>& samplingFaces,
392  DynamicList<label>& samplingSegments,
393  DynamicList<scalar>& samplingCurveDist,
394 
395  // State of current leak paths
396  bitSet& isLeakCell,
397  bitSet& isLeakFace,
398  bitSet& isLeakPoint,
399 
400  // Work storage
401  List<topoDistanceData<label>>& allFaceInfo,
402  List<topoDistanceData<label>>& allCellInfo
403 ) const
404 {
405  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
406  const topoDistanceData<label> maxData(labelMax, labelMax);
407 
408 
409  allFaceInfo.setSize(mesh.nFaces());
410  allFaceInfo = topoDistanceData<label>();
411  allCellInfo.setSize(mesh.nCells());
412  allCellInfo = topoDistanceData<label>();
413 
414  // Mark blocked faces with high distance
415  forAll(isBlockedFace, facei)
416  {
417  if (isBlockedFace[facei])
418  {
419  allFaceInfo[facei] = maxData;
420  }
421  }
422 
423  if (markLeakPath)
424  {
425  // Mark any previously found leak path. This marks all
426  // faces of all cells on the path. This will make sure that
427  // ultimately the path will go through another gap.
428  forAll(isLeakCell, celli)
429  {
430  if (celli != insideCelli && celli != outsideCelli)
431  {
432  if (isLeakCell[celli])
433  {
434  allCellInfo[celli] = maxData;
435  //- Mark all faces of the cell
436  //const cell& cFaces = mesh.cells()[celli];
437  //for (auto facei : cFaces)
438  //{
439  // allFaceInfo[facei] = maxData;
440  //}
441  }
442  }
443  }
444 
445  //- Mark only faces inbetween two leak cells
446  bitSet isLeakCellWithout(isLeakCell);
447  if (insideCelli != -1)
448  {
449  isLeakCellWithout.unset(insideCelli);
450  }
451  if (outsideCelli != -1)
452  {
453  isLeakCellWithout.unset(outsideCelli);
454  }
455  const bitSet isPathFace(pathFaces(mesh, isLeakCellWithout));
456  forAll(isPathFace, facei)
457  {
458  if (isPathFace[facei])
459  {
460  allFaceInfo[facei] = maxData;
461  }
462  }
463  }
464 
465  // Mark any previously found leak faces. These are faces that
466  // are (point- or edge-)connected to an existing boundary.
467  forAll(isLeakFace, facei)
468  {
469  if (isLeakFace[facei])
470  {
471  allFaceInfo[facei] = maxData;
472  }
473  }
474 
475 
476  // Pass1: Get distance to insideCelli
477 
478  calculateDistance(iter, mesh, insideCelli, allFaceInfo, allCellInfo);
479 
480 
481 
482  // Pass2: walk from outside points backwards. Note: could be done
483  // using FaceCellWave as well but is overly complex since
484  // does not allow logic comparing all faces of a cell.
485 
486  bool targetFound = false;
487  if (outsideCelli != -1)
488  {
489  int dummyTrackData;
490  targetFound = allCellInfo[outsideCelli].valid(dummyTrackData);
491  if (iter == 0 && !targetFound)
492  {
494  << "Point " << outsidePoint
495  << " not reachable by walk from " << insidePoint
496  << ". Probably mesh has island/regions."
497  << " Skipped route detection." << endl;
498  }
499  }
500 
501  if (!returnReduceOr(targetFound))
502  {
503  //Pout<< "now :"
504  // << " nLeakCell:"
505  // << returnReduce(isLeakCell.count(), sumOp<label>())
506  // << " nLeakFace:"
507  // << returnReduce(isLeakFace.count(), sumOp<label>())
508  // << " nLeakPoint:"
509  // << returnReduce(isLeakPoint.count(), sumOp<label>())
510  // << endl;
511 
512  return false;
513  }
514 
515 
516  // Start with given target cell and walk back
517  // If point happens to be on multiple processors, random pick
518  label frontCellI = outsideCelli;
519  point origin(outsidePoint);
520  bool findMinDistance = true;
521 
522  while (true)
523  {
524  // We have a single face front (frontFaceI). Walk from face to cell
525  // to face etc until we reach the destination cell.
526 
527  label frontFaceI = -1;
528 
529  // Work within same processor
530  if (frontCellI != -1)
531  {
532  // Find face with lowest distance from seed. In below figure the
533  // seed cell is marked with distance 0. It is surrounded by faces
534  // and cells with distance 1. The layer outside is marked with
535  // distance 2 etc etc.
536  //
537  // x | x 2 1 2 2 | x | x
538  // --- + --- + -1- + -2- + --- + ---
539  // x | 1 1 0 1 1 | x | x
540  // --- + --- + -1- + -2- + --- + ---
541  // x | x 2 1 2 2 3 3 4 4
542  // --- + --- + --- + -3- + -4- + -5-
543  // x | x 3 2 3 3 4 4 5 5
544  //
545  // e.g. if we start backwards search from cell with value = 4,
546  // we have neighbour faces 4, 4, 5, 5. Choose 4 (least distance
547  // to seed) and continue...
548 
549  frontFaceI = findMinFace
550  (
551  mesh,
552  frontCellI,
553  allFaceInfo,
554  isLeakPoint,
555  findMinDistance, // mode : find min or find max
556  origin
557  );
558 
559 
560  // Loop until we hit a boundary face
561  bitSet isNewLeakPoint(isLeakPoint);
562  while (mesh.isInternalFace(frontFaceI))
563  {
564  if (isBlockedFace.size() && isBlockedFace[frontFaceI])
565  {
566  // This should not happen since we never walk through
567  // a blocked face. However ...
568  // Pout<< " Found blocked face" << endl;
569  frontCellI = -1;
570  break;
571  }
572 
573  // Step to neighbouring cell
574  label nbrCellI = mesh.faceOwner()[frontFaceI];
575  if (nbrCellI == frontCellI)
576  {
577  nbrCellI = mesh.faceNeighbour()[frontFaceI];
578  }
579 
580  if (nbrCellI == insideCelli)
581  {
582  // Reached destination. This is the normal exit.
583  frontCellI = -1;
584  break;
585  }
586 
587  frontCellI = nbrCellI;
588 
589  // Pick best face on cell
590  frontFaceI = findMinFace
591  (
592  mesh,
593  frontCellI,
594  allFaceInfo,
595  isLeakPoint,
596  findMinDistance,
597  origin
598  );
599 
600  const topoDistanceData<label>& cInfo = allCellInfo[frontCellI];
601  const topoDistanceData<label>& fInfo = allFaceInfo[frontFaceI];
602 
603  if (fInfo.distance() <= cInfo.distance())
604  {
605  // Found valid next cell,face. Mark on path
606  samplingPts.append(mesh.cellCentres()[frontCellI]);
607  samplingCells.append(frontCellI);
608  samplingFaces.append(-1);
609  samplingSegments.append(iter);
610  samplingCurveDist.append
611  (
612  trackLength+cInfo.distance()
613  );
614  isLeakCell.set(frontCellI);
615 
616  // Check if front face uses a boundary point.
617  if
618  (
619  touchesWall
620  (
621  mesh,
622  frontFaceI,
623  isLeakFace,
624  isLeakPoint
625  )
626  )
627  {
628  //Pout<< "** added leak face:" << frontFaceI << " at:"
629  //<< mesh.faceCentres()[frontFaceI] << endl;
630  isNewLeakPoint.set(mesh.faces()[frontFaceI]);
631  origin = mesh.faceCentres()[frontFaceI];
632  findMinDistance = false;
633  }
634  }
635  }
636  isLeakPoint.transfer(isNewLeakPoint);
637  }
638 
639  // Situation 1: we found the destination cell (do nothing),
640  // frontCellI is -1 on all processors
641  if (returnReduceAnd(frontCellI == -1))
642  {
643  //Pout<< "now :"
644  // << " nLeakCell:"
645  // << returnReduce(isLeakCell.count(), sumOp<label>())
646  // << " nLeakFace:"
647  // << returnReduce(isLeakFace.count(), sumOp<label>())
648  // << " nLeakPoint:"
649  // << returnReduce(isLeakPoint.count(), sumOp<label>())
650  // << endl;
651 
652  break;
653  }
654 
655  // Situation 2: we're on a coupled patch and might need to
656  // switch processor/cell. We need to transfer:
657  // -frontface -frontdistance -leak points/faces
658  boolList isFront(mesh.nBoundaryFaces(), false);
659 
660  if (frontFaceI != -1)
661  {
662  isFront[frontFaceI-mesh.nInternalFaces()] = true;
663  }
665 
666  label minCellDistance = labelMax;
667  if (frontCellI != -1)
668  {
669  minCellDistance = allCellInfo[frontCellI].distance();
670  }
671  reduce(minCellDistance, minOp<label>());
672 
673  // Sync all leak data
674  sync
675  (
676  mesh,
677  isLeakFace,
678  isLeakPoint,
679  frontCellI,
680  origin,
681  findMinDistance
682  );
683 
684 
685  // Bit tricky:
686  // - processor might get frontFaceI/frontCellI in through multiple faces
687  // (even through different patches?)
688  // - so loop over all (coupled) patches and pick the best frontCellI
689 
690  const label oldFrontFaceI = frontFaceI;
691  frontCellI = -1;
692  frontFaceI = -1;
693  forAll(pbm, patchI)
694  {
695  const polyPatch& pp = pbm[patchI];
696  forAll(pp, i)
697  {
698  label faceI = pp.start()+i;
699  if
700  (
701  oldFrontFaceI == -1
702  && isFront[faceI-mesh.nInternalFaces()]
703  && (isBlockedFace.empty() || !isBlockedFace[faceI])
704  )
705  {
706  frontFaceI = faceI;
707  frontCellI = pp.faceCells()[i];
708  break;
709  }
710  }
711 
712  if
713  (
714  frontCellI != -1
715  && allCellInfo[frontCellI].distance() < minCellDistance
716  )
717  {
718  const topoDistanceData<label>& cInfo = allCellInfo[frontCellI];
719 
720  samplingPts.append(mesh.cellCentres()[frontCellI]);
721  samplingCells.append(frontCellI);
722  samplingFaces.append(-1);
723  samplingSegments.append(iter);
724  samplingCurveDist.append
725  (
726  trackLength+cInfo.distance()
727  );
728  isLeakCell.set(frontCellI);
729 
730  // Check if frontFacei touches leakPoint
731  if
732  (
733  touchesWall
734  (
735  mesh,
736  frontFaceI,
737  isLeakFace,
738  isLeakPoint
739  )
740  )
741  {
742  //Pout<< "** added leak BOUNDARY face:" << frontFaceI
743  // << " at:" << mesh.faceCentres()[frontFaceI] << endl;
744  isLeakPoint.set(mesh.faces()[frontFaceI]);
745  origin = mesh.faceCentres()[frontFaceI];
746  findMinDistance = false;
747  }
748 
749  break;
750  }
751 
752  // When seed cell is isolated by processor boundaries
753  if (insideCelli != -1 && frontCellI == insideCelli)
754  {
755  // Pout<< " Found connection seed cell!" << endl;
756  frontCellI = -1;
757  break;
758  }
759  }
760 
761  // Sync all leak data
762  sync
763  (
764  mesh,
765  isLeakFace,
766  isLeakPoint,
767  frontCellI,
768  origin,
769  findMinDistance
770  );
771  }
772 
773  return true;
774 }
775 
776 
777 // Clean up leak faces (erode open edges). These are leak faces which are
778 // not connected to another leak face or leak point. Parallel consistent.
779 // Returns overall number of faces deselected.
780 Foam::label Foam::shortestPathSet::erodeFaceSet
781 (
782  const polyMesh& mesh,
783  const bitSet& isBlockedPoint,
784  bitSet& isLeakFace
785 ) const
786 {
787  if
788  (
789  (isBlockedPoint.size() != mesh.nPoints())
790  || (isLeakFace.size() != mesh.nFaces())
791  )
792  {
793  FatalErrorInFunction << "Problem :"
794  << " isBlockedPoint:" << isBlockedPoint.size()
795  << " isLeakFace:" << isLeakFace.size()
796  << exit(FatalError);
797  }
798 
799  const globalMeshData& globalData = mesh.globalData();
800  const mapDistribute& map = globalData.globalEdgeSlavesMap();
802 
803 
804  label nTotalEroded = 0;
805 
806  while (true)
807  {
808  bitSet newIsLeakFace(isLeakFace);
809 
810  // Get number of edges
811 
812  const labelList meshFaceIDs(isLeakFace.toc());
814  (
815  UIndirectList<face>(mesh.faces(), meshFaceIDs),
816  mesh.points()
817  );
818 
819  // Count number of faces per edge
820  const labelListList& edgeFaces = pp.edgeFaces();
821  labelList nEdgeFaces(edgeFaces.size());
822  forAll(edgeFaces, edgei)
823  {
824  nEdgeFaces[edgei] = edgeFaces[edgei].size();
825  }
826 
827  // Match pp edges to coupled edges
828  labelList patchEdges;
829  labelList coupledEdges;
830  bitSet sameEdgeOrientation;
832  (
833  pp,
834  cpp,
835  patchEdges,
836  coupledEdges,
837  sameEdgeOrientation
838  );
839 
840 
841  // Convert patch-edge data into cpp-edge data
842  labelList coupledNEdgeFaces(map.constructSize(), Zero);
843  UIndirectList<label>(coupledNEdgeFaces, coupledEdges) =
844  UIndirectList<label>(nEdgeFaces, patchEdges);
845 
846  // Synchronise
847  globalData.syncData
848  (
849  coupledNEdgeFaces,
850  globalData.globalEdgeSlaves(),
851  globalData.globalEdgeTransformedSlaves(),
852  map,
853  plusEqOp<label>()
854  );
855 
856  // Convert back from cpp-edge to patch-edge
857  UIndirectList<label>(nEdgeFaces, patchEdges) =
858  UIndirectList<label>(coupledNEdgeFaces, coupledEdges);
859 
860  // Remove any open edges
861  label nEroded = 0;
862  forAll(nEdgeFaces, edgei)
863  {
864  if (nEdgeFaces[edgei] == 1)
865  {
866  const edge& e = pp.edges()[edgei];
867  const label mp0 = pp.meshPoints()[e[0]];
868  const label mp1 = pp.meshPoints()[e[1]];
869 
870  if (!isBlockedPoint[mp0] || !isBlockedPoint[mp1])
871  {
872  // Edge is not on wall so is open
873  const label patchFacei = edgeFaces[edgei][0];
874  const label meshFacei = pp.addressing()[patchFacei];
875  //Pout<< "Eroding face:" << meshFacei
876  // << " at:" << mesh.faceCentres()[meshFacei]
877  // << " since has open edge:" << mesh.points()[mp0]
878  // << mesh.points()[mp1] << endl;
879 
880  if (newIsLeakFace.unset(meshFacei))
881  {
882  nEroded++;
883  }
884  }
885  }
886  }
887 
888  reduce(nEroded, sumOp<label>());
889  nTotalEroded += nEroded;
890 
891  if (nEroded == 0)
892  {
893  break;
894  }
895 
896  // Make sure we make the same decision on both sides
898  (
899  mesh,
900  newIsLeakFace,
901  orEqOp<unsigned int>()
902  );
903 
904  isLeakFace = std::move(newIsLeakFace);
905  }
906 
907  return nTotalEroded;
908 }
909 
910 
911 void Foam::shortestPathSet::genSamples
912 (
913  const bool addLeakPath,
914  const label maxIter,
915  const polyMesh& mesh,
916  const bitSet& isBoundaryFace,
917  const point& insidePoint,
918  const label insideCelli,
919  const point& outsidePoint,
920 
921  DynamicList<point>& samplingPts,
922  DynamicList<label>& samplingCells,
923  DynamicList<label>& samplingFaces,
924  DynamicList<label>& samplingSegments,
925  DynamicList<scalar>& samplingCurveDist,
926  bitSet& isLeakCell,
927  bitSet& isLeakFace,
928  bitSet& isLeakPoint
929 ) const
930 {
931  // Mark all paths needed to close a single combination of insidePoint,
932  // outsidePoint. The output is
933  // - isLeakCell : is cell along a leak path
934  // - isLeakFace : is face along a leak path (so inbetween two leakCells)
935  // - isLeakPoint : is point on a leakFace
936 
937 
938  const topoDistanceData<label> maxData(labelMax, labelMax);
939 
940  // Get the target point
941  const label outsideCelli = mesh.findCell(outsidePoint);
942 
943  // Maintain overall track length. Used to make curve distance continuous.
944  scalar trackLength = 0;
945 
946  List<topoDistanceData<label>> allFaceInfo(mesh.nFaces());
947  List<topoDistanceData<label>> allCellInfo(mesh.nCells());
948 
949 
950  // Boundary face + additional temporary blocks (to force leakpath to
951  // outside)
952  autoPtr<bitSet> isBlockedFace;
953 
954  label iter;
955  bool markLeakPath = false;
956 
957 
958  for (iter = 0; iter < maxIter; iter++)
959  {
960  const label nOldLeakFaces = returnReduce
961  (
962  isLeakFace.count(),
963  sumOp<label>()
964  );
965  const label oldNSamplingPts(samplingPts.size());
966 
967  bool foundPath = genSingleLeakPath
968  (
969  markLeakPath,
970  iter,
971  mesh,
972  (isBlockedFace ? *isBlockedFace : isBoundaryFace),
973  insidePoint,
974  insideCelli,
975  outsidePoint,
976  outsideCelli,
977 
978  // Generated sampling points
979  trackLength,
980  samplingPts,
981  samplingCells,
982  samplingFaces,
983  samplingSegments,
984  samplingCurveDist,
985 
986  // State of current leak paths
987  isLeakCell,
988  isLeakFace,
989  isLeakPoint,
990 
991  // Work storage
992  allFaceInfo,
993  allCellInfo
994  );
995 
996  // Recalculate the overall trackLength
997  trackLength = returnReduce
998  (
999  (
1000  samplingCurveDist.size()
1001  ? samplingCurveDist.last()
1002  : 0
1003  ),
1004  maxOp<scalar>()
1005  );
1006 
1007  const label nLeakFaces = returnReduce
1008  (
1009  isLeakFace.count(),
1010  sumOp<label>()
1011  );
1012 
1013  if (!foundPath && !markLeakPath)
1014  {
1015  // In mark-boundary-face-only mode and already fully blocked the
1016  // path to outsideCell so we're good
1017  break;
1018  }
1019 
1020 
1021  if (nLeakFaces > nOldLeakFaces)
1022  {
1023  // Normal operation: walking has closed some wall-connected faces
1024  // If previous iteration was markLeakPath-mode make sure to revert
1025  // to normal operation (i.e. face marked in isLeakFace)
1026  isBlockedFace.reset(nullptr);
1027  markLeakPath = false;
1028  }
1029  else
1030  {
1031  // Did not mark any additional faces/points. Save current state
1032  // and add faces/points on leakpath to force next walk
1033  // to pass outside of leakpath. This is done until the leakpath
1034  // 'touchesWall' (i.e. edge connected to an original boundary face)
1035 
1036  if (markLeakPath && !foundPath)
1037  {
1038  // Is marking all faces on all paths and no additional path
1039  // found. Also nothing new marked (in isLeakFace) since
1040  // nLeakFaces == nOldLeakFaces) so we're
1041  // giving up. Convert all path faces into leak faces
1042  //Pout<< "** giving up" << endl;
1043  break;
1044  }
1045 
1046 
1047  // Revert to boundaryFaces only
1048  if (!isBlockedFace)
1049  {
1050  //Pout<< "** Starting from original boundary faces." << endl;
1051  isBlockedFace.reset(new bitSet(isBoundaryFace));
1052  }
1053 
1054  markLeakPath = true;
1055 
1056 
1057  if (debug & 2)
1058  {
1059  const pointField leakCcs(mesh.cellCentres(), isLeakCell.toc());
1060  mkDir(mesh.time().timePath());
1061  OBJstream str
1062  (
1063  mesh.time().timePath()
1064  /"isLeakCell" + Foam::name(iter) + ".obj"
1065  );
1066  Pout<< "Writing new isLeakCell to " << str.name() << endl;
1067  str.write(leakCcs);
1068  }
1069  if (debug & 2)
1070  {
1071  OBJstream str
1072  (
1073  mesh.time().timePath()
1074  /"leakPath" + Foam::name(iter) + ".obj"
1075  );
1076  Pout<< "Writing new leak-path to " << str.name() << endl;
1077  for
1078  (
1079  label samplei = oldNSamplingPts+1;
1080  samplei < samplingPts.size();
1081  samplei++
1082  )
1083  {
1084  Pout<< " passing through cell "
1085  << samplingCells[samplei]
1086  << " at:" << mesh.cellCentres()[samplingCells[samplei]]
1087  << " distance:" << samplingCurveDist[samplei]
1088  << endl;
1089 
1090  str.writeLine
1091  (
1092  samplingPts[samplei-1],
1093  samplingPts[samplei]
1094  );
1095  }
1096  }
1097  }
1098  }
1099 
1100  if (debug)
1101  {
1102  const fvMesh& fm = refCast<const fvMesh>(mesh);
1103 
1104  const_cast<fvMesh&>(fm).setInstance(fm.time().timeName());
1106  (
1107  IOobject
1108  (
1109  "isLeakCell",
1110  fm.time().timeName(),
1111  fm,
1114  ),
1115  fm,
1117  );
1118  for (const label celli : isLeakCell)
1119  {
1120  fld[celli] = scalar(1);
1121  }
1122  fld.correctBoundaryConditions();
1123  fld.write();
1124  }
1125 
1126  if (maxIter > 1 && iter == maxIter)
1127  {
1128  WarningInFunction << "Did not manage to close gap using " << iter
1129  << " leak paths" << nl << "This can cause problems when using the"
1130  << " paths to close leaks" << endl;
1131  }
1132 }
1133 
1134 
1135 void Foam::shortestPathSet::genSamples
1136 (
1137  const bool markLeakPath,
1138  const label maxIter,
1139  const polyMesh& mesh,
1140  const labelUList& wallPatches,
1141  const bitSet& isBlockedFace
1142 )
1143 {
1144  // Storage for sample points
1145  DynamicList<point> samplingPts;
1146  DynamicList<label> samplingCells;
1147  DynamicList<label> samplingFaces;
1148  DynamicList<label> samplingSegments;
1149  DynamicList<scalar> samplingCurveDist;
1150 
1151  // Seed faces and points on 'real' boundary
1152  bitSet isBlockedPoint(mesh.nPoints());
1153  {
1154  // Real boundaries
1155  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
1156  for (label patchi : wallPatches)
1157  {
1158  const polyPatch& pp = pbm[patchi];
1159  forAll(pp, i)
1160  {
1161  isBlockedPoint.set(pp[i]);
1162  }
1163  }
1164 
1165  // Meshed boundaries
1166  forAll(isBlockedFace, facei)
1167  {
1168  if (isBlockedFace[facei])
1169  {
1170  isBlockedPoint.set(mesh.faces()[facei]);
1171  }
1172  }
1173 
1175  (
1176  mesh,
1177  isBlockedPoint,
1178  orEqOp<unsigned int>(),
1179  0u
1180  );
1181 
1182  if (debug)
1183  {
1184  mkDir(mesh.time().timePath());
1185  OBJstream str(mesh.time().timePath()/"isBlockedPoint.obj");
1186  for (const auto& pointi : isBlockedPoint)
1187  {
1188  str.write(mesh.points()[pointi]);
1189  }
1190  Pout<< "Writing " << str.nVertices()
1191  << " points to " << str.name() << endl;
1192  }
1193  }
1194 
1195 
1196  bitSet isLeakPoint(isBlockedPoint);
1197  // Newly closed faces
1198  bitSet isLeakFace(mesh.nFaces());
1199  // All cells along leak paths
1200  bitSet isLeakCell(mesh.nCells());
1201 
1202  label prevSegmenti = 0;
1203  scalar prevDistance = 0.0;
1204 
1205  for (auto insidePoint : insidePoints_)
1206  {
1207  const label insideCelli = mesh.findCell(insidePoint);
1208 
1209  for (auto outsidePoint : outsidePoints_)
1210  {
1211  const label nOldSamples = samplingSegments.size();
1212 
1213  // Append any leak path to sampling*
1214  genSamples
1215  (
1216  markLeakPath,
1217  maxIter,
1218  mesh,
1219  isBlockedFace,
1220  insidePoint,
1221  insideCelli,
1222  outsidePoint,
1223 
1224  samplingPts,
1225  samplingCells,
1226  samplingFaces,
1227  samplingSegments,
1228  samplingCurveDist,
1229  isLeakCell,
1230  isLeakFace,
1231  isLeakPoint
1232  );
1233 
1234  // Make segment, distance consecutive
1235  label maxSegment = 0;
1236  scalar maxDistance = 0.0;
1237  for (label i = nOldSamples; i < samplingSegments.size(); ++i)
1238  {
1239  samplingSegments[i] += prevSegmenti;
1240  maxSegment = max(maxSegment, samplingSegments[i]);
1241  samplingCurveDist[i] += prevDistance;
1242  maxDistance = max(maxDistance, samplingCurveDist[i]);
1243  }
1244  prevSegmenti = returnReduce(maxSegment, maxOp<label>());
1245  prevDistance = returnReduce(maxDistance, maxOp<scalar>());
1246  }
1247  }
1248 
1249  // Clean up leak faces (erode open edges). These are leak faces which are
1250  // not connected to another leak face or leak point
1251  erodeFaceSet(mesh, isBlockedPoint, isLeakFace);
1252 
1253  leakFaces_ = isLeakFace.sortedToc();
1254 
1255 
1256  if (debug)
1257  {
1258  const faceList leakFaces(mesh.faces(), leakFaces_);
1259 
1260  mkDir(mesh.time().timePath());
1261  OBJstream str(mesh.time().timePath()/"isLeakFace.obj");
1262  str.write(leakFaces, mesh.points(), false);
1263  Pout<< "Writing " << leakFaces.size()
1264  << " faces to " << str.name() << endl;
1265  }
1266 
1267 
1268  samplingPts.shrink();
1269  samplingCells.shrink();
1270  samplingFaces.shrink();
1271  samplingSegments.shrink();
1272  samplingCurveDist.shrink();
1273 
1274  // Move into *this
1275  setSamples
1276  (
1277  std::move(samplingPts),
1278  std::move(samplingCells),
1279  std::move(samplingFaces),
1280  std::move(samplingSegments),
1281  std::move(samplingCurveDist)
1282  );
1283 
1284  if (debug)
1285  {
1286  write(Info);
1287  }
1288 }
1289 
1290 
1291 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1294 (
1295  const word& name,
1296  const polyMesh& mesh,
1297  const meshSearch& searchEngine,
1298  const word& axis,
1299  const bool markLeakPath,
1300  const label maxIter,
1301  const labelUList& wallPatches,
1302  const pointField& insidePoints,
1303  const pointField& outsidePoints,
1304  const boolList& isBlockedFace
1305 )
1306 :
1307  sampledSet(name, mesh, searchEngine, axis),
1308  insidePoints_(insidePoints),
1309  outsidePoints_(outsidePoints)
1310 {
1311  if (debug)
1312  {
1313  fileName outputDir
1314  (
1315  mesh.time().globalPath()
1317  / mesh.pointsInstance()
1318  );
1319  outputDir.clean(); // Remove unneeded ".."
1320 
1321  Info<< "shortestPathSet : Writing blocked faces to "
1322  << outputDir << endl;
1323 
1324  const indirectPrimitivePatch setPatch
1325  (
1327  (
1328  mesh.faces(),
1329  findIndices(isBlockedFace, true)
1330  ),
1331  mesh.points()
1332  );
1333 
1334  if (Pstream::parRun())
1335  {
1336  // Topological merge
1337  labelList pointToGlobal;
1338  labelList uniqueMeshPointLabels;
1340  autoPtr<globalIndex> globalFaces;
1341  faceList mergedFaces;
1342  pointField mergedPoints;
1344  (
1345  mesh,
1346  setPatch.localFaces(),
1347  setPatch.meshPoints(),
1348  setPatch.meshPointMap(),
1349 
1350  pointToGlobal,
1351  uniqueMeshPointLabels,
1352  globalPoints,
1353  globalFaces,
1354 
1355  mergedFaces,
1356  mergedPoints
1357  );
1358 
1359  // Write
1360  if (Pstream::master())
1361  {
1363  (
1364  mergedPoints,
1365  mergedFaces,
1366  (outputDir / "blockedFace"),
1367  false // serial only - already merged
1368  );
1369 
1370  writer.writeGeometry();
1371  }
1372  }
1373  else
1374  {
1376  (
1377  setPatch.localPoints(),
1378  setPatch.localFaces(),
1379  (outputDir / "blockedFace"),
1380  false // serial only - redundant
1381  );
1382 
1383  writer.writeGeometry();
1384  }
1385  }
1386 
1387  genSamples
1388  (
1389  markLeakPath,
1390  maxIter,
1391  mesh,
1392  wallPatches,
1393  bitSet(isBlockedFace)
1394  );
1395 }
1396 
1399 (
1400  const word& name,
1401  const polyMesh& mesh,
1402  const meshSearch& searchEngine,
1403  const dictionary& dict
1404 )
1405 :
1406  sampledSet(name, mesh, searchEngine, dict),
1407  insidePoints_(dict.get<pointField>("insidePoints")),
1408  outsidePoints_(dict.get<pointField>("outsidePoints"))
1409 {
1410  const label maxIter(dict.getOrDefault<label>("maxIter", 50));
1411  const bool markLeakPath(dict.getOrDefault("markLeakPath", true));
1412 
1414 
1415  DynamicList<label> wallPatches(pbm.size());
1416  forAll(pbm, patchi)
1417  {
1418  const polyPatch& pp = pbm[patchi];
1419  if (!pp.coupled() && !isA<emptyPolyPatch>(pp))
1420  {
1421  wallPatches.append(patchi);
1422  }
1423  }
1424 
1425  genSamples(markLeakPath, maxIter, mesh, wallPatches, bitSet());
1426 }
1427 
1428 
1429 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
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.
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search...
Definition: meshSearch.H:56
const polyBoundaryMesh & pbm
dictionary dict
A class for handling file names.
Definition: fileName.H:72
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
vtk::lineWriter writer(edgeCentres, edgeList::null(), fileName(aMesh.time().globalPath()/"finiteArea-edgesCentres"))
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition: List.C:326
label nPoints() const noexcept
Number of mesh points.
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
fileName timePath() const
Return current time path = path/timeName.
Definition: Time.H:520
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:517
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1122
const mapDistribute & globalEdgeSlavesMap() const
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
static void matchEdges(const PrimitivePatch< FaceList1, PointField1 > &p1, const PrimitivePatch< FaceList2, PointField2 > &p2, labelList &p1EdgeLabels, labelList &p2EdgeLabels, bitSet &sameOrientation)
Find corresponding edges on patches sharing the same points.
Write faces/points (optionally with fields) as a vtp file or a legacy vtk file.
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
const scalarList & distance() const noexcept
Return the cumulative distance.
Definition: coordSet.H:176
static void syncFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop)
Synchronize values on all mesh faces.
Definition: syncTools.H:432
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
labelList findIndices(const ListType &input, typename ListType::const_reference val, label start=0)
Linear search to find all occurrences of given element.
const cellList & cells() const
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
label nFaces() const noexcept
Number of mesh faces.
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field...
virtual const fileName & name() const override
Get the name of the output serial stream. (eg, the name of the Fstream file name) ...
Definition: OSstream.H:128
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.
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
Macros for easy insertion into run-time selection tables.
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
const labelList & meshPoints() const
Return labelList of mesh points in patch.
static void gatherAndMerge(const scalar mergeDist, const PrimitivePatch< FaceList, PointField > &pp, Field< typename PrimitivePatch< FaceList, PointField >::point_type > &mergedPoints, List< typename PrimitivePatch< FaceList, PointField >::face_type > &mergedFaces, globalIndex &pointAddr, globalIndex &faceAddr, labelList &pointMergeMap=const_cast< labelList &>(labelList::null()), const bool useLocal=false)
Gather points and faces onto master and merge into single patch.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
label fcIndex(const label i) const noexcept
The forward circular index. The next index in the list which returns to the first at the end of the l...
Definition: UListI.H:97
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
scalar y
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
A list of faces which address into the list of points.
bool returnReduceAnd(const bool value, const label comm=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:853
dynamicFvMesh & mesh
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:614
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:79
fileName globalPath() const
Return global path for the case = rootPath/globalCaseName. Same as TimePaths::globalPath() ...
Definition: Time.H:512
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const Time & time() const noexcept
Return time registry.
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1116
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
const labelListList & edgeFaces() const
Return edge-face addressing.
const globalMeshData & globalData() const
Return parallel info (demand-driven)
Definition: polyMesh.C:1311
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
label nInternalFaces() const noexcept
Number of internal faces.
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
const vectorField & cellCentres() const
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
static void combineReduce(const List< commsStruct > &comms, T &value, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) applying cop to inplace combine value from different processors...
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
const polyMesh & mesh() const noexcept
Definition: sampledSet.H:373
int debug
Static debugging option.
defineTypeNameAndDebug(combustionModel, 0)
static bool clean(std::string &str)
Cleanup filename string, possibly applies other transformations such as changing the path separator e...
Definition: fileName.C:192
labelList f(nPoints)
label nTotalCells() const noexcept
Return total number of cells in decomposed mesh.
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))
const indirectPrimitivePatch & coupledPatch() const
Return patch of all coupled faces.
const vectorField & faceCentres() const
Calculates points shared by more than two processor patches or cyclic patches.
Definition: globalPoints.H:98
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
A PrimitivePatch with UIndirectList for the faces, const reference for the point field.
vector point
Point is a vector.
Definition: point.H:37
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
#define WarningInFunction
Report a warning using Foam::Warning.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label nCells() const noexcept
Number of mesh cells.
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
const polyBoundaryMesh & patches
Nothing to be read.
static word outputPrefix
Directory prefix.
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.
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1511
messageStream Info
Information stream (stdout output on master, null elsewhere)
constexpr label labelMax
Definition: label.H:55
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
insidePoints((1 2 3))
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
label nBoundaryFaces() const noexcept
Number of boundary faces (== nFaces - nInternalFaces)
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
List< bool > boolList
A List of bools.
Definition: List.H:60
A List with indirect addressing.
Definition: IndirectList.H:60
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues)
Swap coupled boundary face values. Uses eqOp.
Definition: syncTools.H:485
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127