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-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 
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"
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  reduce(targetFound, orOp<bool>());
501  if (!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 (returnReduce(frontCellI == -1, andOp<bool>()))
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());
813  const uindirectPrimitivePatch pp
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 curveDist 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.valid() ? 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  if (isBlockedFace.valid())
1027  {
1028  isBlockedFace.clear();
1029  }
1030  markLeakPath = false;
1031  }
1032  else
1033  {
1034  // Did not mark any additional faces/points. Save current state
1035  // and add faces/points on leakpath to force next walk
1036  // to pass outside of leakpath. This is done until the leakpath
1037  // 'touchesWall' (i.e. edge connected to an original boundary face)
1038 
1039  if (markLeakPath && !foundPath)
1040  {
1041  // Is marking all faces on all paths and no additional path
1042  // found. Also nothing new marked (in isLeakFace) since
1043  // nLeakFaces == nOldLeakFaces) so we're
1044  // giving up. Convert all path faces into leak faces
1045  //Pout<< "** giving up" << endl;
1046  break;
1047  }
1048 
1049 
1050  // Revert to boundaryFaces only
1051  if (!isBlockedFace.valid())
1052  {
1053  //Pout<< "** Starting from original boundary faces." << endl;
1054  isBlockedFace.reset(new bitSet(isBoundaryFace));
1055  }
1056 
1057  markLeakPath = true;
1058 
1059 
1060  if (debug & 2)
1061  {
1062  const pointField leakCcs(mesh.cellCentres(), isLeakCell.toc());
1063  mkDir(mesh.time().timePath());
1064  OBJstream str
1065  (
1066  mesh.time().timePath()
1067  /"isLeakCell" + Foam::name(iter) + ".obj"
1068  );
1069  Pout<< "Writing new isLeakCell to " << str.name() << endl;
1070  forAll(leakCcs, i)
1071  {
1072  str.write(leakCcs[i]);
1073  }
1074  }
1075  if (debug & 2)
1076  {
1077  OBJstream str
1078  (
1079  mesh.time().timePath()
1080  /"leakPath" + Foam::name(iter) + ".obj"
1081  );
1082  Pout<< "Writing new leak-path to " << str.name() << endl;
1083  for
1084  (
1085  label samplei = oldNSamplingPts+1;
1086  samplei < samplingPts.size();
1087  samplei++
1088  )
1089  {
1090  Pout<< " passing through cell "
1091  << samplingCells[samplei]
1092  << " at:" << mesh.cellCentres()[samplingCells[samplei]]
1093  << " distance:" << samplingCurveDist[samplei]
1094  << endl;
1095 
1096  str.write
1097  (
1098  linePointRef
1099  (
1100  samplingPts[samplei-1],
1101  samplingPts[samplei]
1102  )
1103  );
1104  }
1105  }
1106  }
1107  }
1108 
1109  if (debug)
1110  {
1111  const fvMesh& fm = refCast<const fvMesh>(mesh);
1112 
1113  const_cast<fvMesh&>(fm).setInstance(fm.time().timeName());
1115  (
1116  IOobject
1117  (
1118  "isLeakCell",
1119  fm.time().timeName(),
1120  fm,
1123  ),
1124  fm,
1126  );
1127  forAll(isLeakCell, celli)
1128  {
1129  fld[celli] = isLeakCell[celli];
1130  }
1131  fld.correctBoundaryConditions();
1132  fld.write();
1133  }
1134 
1135  if (maxIter > 1 && iter == maxIter)
1136  {
1137  WarningInFunction << "Did not manage to close gap using " << iter
1138  << " leak paths" << nl << "This can cause problems when using the"
1139  << " paths to close leaks" << endl;
1140  }
1141 }
1142 
1143 
1144 void Foam::shortestPathSet::genSamples
1145 (
1146  const bool markLeakPath,
1147  const label maxIter,
1148  const polyMesh& mesh,
1149  const labelUList& wallPatches,
1150  const bitSet& isBlockedFace
1151 )
1152 {
1153  // Storage for sample points
1154  DynamicList<point> samplingPts;
1155  DynamicList<label> samplingCells;
1156  DynamicList<label> samplingFaces;
1157  DynamicList<label> samplingSegments;
1158  DynamicList<scalar> samplingCurveDist;
1159 
1160  // Seed faces and points on 'real' boundary
1161  bitSet isBlockedPoint(mesh.nPoints());
1162  {
1163  // Real boundaries
1164  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
1165  for (label patchi : wallPatches)
1166  {
1167  const polyPatch& pp = pbm[patchi];
1168  forAll(pp, i)
1169  {
1170  isBlockedPoint.set(pp[i]);
1171  }
1172  }
1173 
1174  // Meshed boundaries
1175  forAll(isBlockedFace, facei)
1176  {
1177  if (isBlockedFace[facei])
1178  {
1179  isBlockedPoint.set(mesh.faces()[facei]);
1180  }
1181  }
1182 
1184  (
1185  mesh,
1186  isBlockedPoint,
1187  orEqOp<unsigned int>(),
1188  0u
1189  );
1190 
1191  if (debug)
1192  {
1193  mkDir(mesh.time().timePath());
1194  OBJstream str(mesh.time().timePath()/"isBlockedPoint.obj");
1195  for (const auto& pointi : isBlockedPoint)
1196  {
1197  str.write(mesh.points()[pointi]);
1198  }
1199  Pout<< "Writing " << str.nVertices() << " points to " << str.name()
1200  << endl;
1201  }
1202  }
1203 
1204 
1205  bitSet isLeakPoint(isBlockedPoint);
1206  // Newly closed faces
1207  bitSet isLeakFace(mesh.nFaces());
1208  // All cells along leak paths
1209  bitSet isLeakCell(mesh.nCells());
1210 
1211  label prevSegmenti = 0;
1212  scalar prevDistance = 0.0;
1213 
1214  for (auto insidePoint : insidePoints_)
1215  {
1216  const label insideCelli = mesh.findCell(insidePoint);
1217 
1218  for (auto outsidePoint : outsidePoints_)
1219  {
1220  const label nOldSamples = samplingSegments.size();
1221 
1222  // Append any leak path to sampling*
1223  genSamples
1224  (
1225  markLeakPath,
1226  maxIter,
1227  mesh,
1228  isBlockedFace,
1229  insidePoint,
1230  insideCelli,
1231  outsidePoint,
1232 
1233  samplingPts,
1234  samplingCells,
1235  samplingFaces,
1236  samplingSegments,
1237  samplingCurveDist,
1238  isLeakCell,
1239  isLeakFace,
1240  isLeakPoint
1241  );
1242 
1243  // Make segment, distance consecutive
1244  label maxSegment = 0;
1245  scalar maxDistance = 0.0;
1246  for (label i = nOldSamples; i < samplingSegments.size(); ++i)
1247  {
1248  samplingSegments[i] += prevSegmenti;
1249  maxSegment = max(maxSegment, samplingSegments[i]);
1250  samplingCurveDist[i] += prevDistance;
1251  maxDistance = max(maxDistance, samplingCurveDist[i]);
1252  }
1253  prevSegmenti = returnReduce(maxSegment, maxOp<label>());
1254  prevDistance = returnReduce(maxDistance, maxOp<scalar>());
1255  }
1256  }
1257 
1258  // Clean up leak faces (erode open edges). These are leak faces which are
1259  // not connected to another leak face or leak point
1260  erodeFaceSet(mesh, isBlockedPoint, isLeakFace);
1261 
1262  leakFaces_ = isLeakFace.sortedToc();
1263 
1264 
1265  if (debug)
1266  {
1267  const faceList leakFaces(mesh.faces(), leakFaces_);
1268 
1269  mkDir(mesh.time().timePath());
1270  OBJstream str(mesh.time().timePath()/"isLeakFace.obj");
1271  str.write(leakFaces, mesh.points(), false);
1272  Pout<< "Writing " << leakFaces.size() << " faces to " << str.name()
1273  << endl;
1274  }
1275 
1276 
1277  samplingPts.shrink();
1278  samplingCells.shrink();
1279  samplingFaces.shrink();
1280  samplingSegments.shrink();
1281  samplingCurveDist.shrink();
1282 
1283  // Move into *this
1284  setSamples
1285  (
1286  std::move(samplingPts),
1287  std::move(samplingCells),
1288  std::move(samplingFaces),
1289  std::move(samplingSegments),
1290  std::move(samplingCurveDist)
1291  );
1292 
1293  if (debug)
1294  {
1295  write(Info);
1296  }
1297 }
1298 
1299 
1300 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1301 
1304  const word& name,
1305  const polyMesh& mesh,
1306  const meshSearch& searchEngine,
1307  const word& axis,
1308  const bool markLeakPath,
1309  const label maxIter,
1310  const labelUList& wallPatches,
1311  const pointField& insidePoints,
1312  const pointField& outsidePoints,
1313  const boolList& isBlockedFace
1314 )
1315 :
1316  sampledSet(name, mesh, searchEngine, axis),
1317  insidePoints_(insidePoints),
1318  outsidePoints_(outsidePoints)
1319 {
1320  if (debug)
1321  {
1322  fileName outputDir =
1323  (
1324  mesh.time().globalPath()
1326  / mesh.pointsInstance()
1327  );
1328  outputDir.clean();
1329 
1330  Info<< "shortestPathSet : Writing blocked faces to "
1331  << outputDir << endl;
1332 
1333  const indirectPrimitivePatch setPatch
1334  (
1336  (
1337  mesh.faces(),
1338  findIndices(isBlockedFace, true)
1339  ),
1340  mesh.points()
1341  );
1342 
1343  if (Pstream::parRun())
1344  {
1345  // Topological merge
1346  labelList pointToGlobal;
1347  labelList uniqueMeshPointLabels;
1349  autoPtr<globalIndex> globalFaces;
1350  faceList mergedFaces;
1351  pointField mergedPoints;
1353  (
1354  mesh,
1355  setPatch.localFaces(),
1356  setPatch.meshPoints(),
1357  setPatch.meshPointMap(),
1358 
1359  pointToGlobal,
1360  uniqueMeshPointLabels,
1361  globalPoints,
1362  globalFaces,
1363 
1364  mergedFaces,
1365  mergedPoints
1366  );
1367 
1368  // Write
1369  if (Pstream::master())
1370  {
1372  (
1373  mergedPoints,
1374  mergedFaces,
1375  (outputDir / "blockedFace"),
1376  false // serial only - already merged
1377  );
1378 
1379  writer.writeGeometry();
1380  }
1381  }
1382  else
1383  {
1385  (
1386  setPatch.localPoints(),
1387  setPatch.localFaces(),
1388  (outputDir / "blockedFace"),
1389  false // serial only - redundant
1390  );
1391 
1392  writer.writeGeometry();
1393  }
1394  }
1395 
1396  genSamples
1397  (
1398  markLeakPath,
1399  maxIter,
1400  mesh,
1401  wallPatches,
1402  bitSet(isBlockedFace)
1403  );
1404 }
1405 
1406 
1409  const word& name,
1410  const polyMesh& mesh,
1411  const meshSearch& searchEngine,
1412  const dictionary& dict
1413 )
1414 :
1415  sampledSet(name, mesh, searchEngine, dict),
1416  insidePoints_(dict.get<pointField>("insidePoints")),
1417  outsidePoints_(dict.get<pointField>("outsidePoints"))
1418 {
1419  const label maxIter(dict.getOrDefault<label>("maxIter", 50));
1420  const bool markLeakPath(dict.getOrDefault("markLeakPath", true));
1421 
1422  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
1423 
1424  DynamicList<label> wallPatches(pbm.size());
1425  forAll(pbm, patchi)
1426  {
1427  const polyPatch& pp = pbm[patchi];
1428  if (!pp.coupled() && !isA<emptyPolyPatch>(pp))
1429  {
1430  wallPatches.append(patchi);
1431  }
1432  }
1433 
1434  genSamples(markLeakPath, maxIter, mesh, wallPatches, bitSet());
1435 }
1436 
1437 
1438 // ************************************************************************* //
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:71
volFields.H
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
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:61
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::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:64
uindirectPrimitivePatch.H
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
Foam::globalMeshData::globalEdgeSlavesMap
const mapDistribute & globalEdgeSlavesMap() const
Definition: globalMeshData.C:2265
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
PatchTools.H
Foam::DynamicList< label >
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const
Number of internal faces.
Definition: primitiveMeshI.H:78
Foam::PatchTools::gatherAndMerge
static void gatherAndMerge(const scalar mergeDist, const PrimitivePatch< FaceList, PointField > &p, Field< typename PrimitivePatch< FaceList, PointField >::point_type > &mergedPoints, List< typename PrimitivePatch< FaceList, PointField >::face_type > &mergedFaces, labelList &pointMergeMap)
Gather points and faces onto master and merge into single patch.
Definition: PatchToolsGatherAndMerge.C:38
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:328
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:415
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:69
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:350
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:1303
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
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:296
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::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:67
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::Time::timePath
fileName timePath() const
Return current time path.
Definition: Time.H:375
Foam::OSstream::name
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:107
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:340
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
Foam::PatchTools::matchEdges
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.
Definition: PatchToolsMatch.C:76
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
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
Foam::Ostream::write
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
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::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:439
Foam::indirectPrimitivePatch
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
Definition: indirectPrimitivePatch.H:49
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
Foam::polyMesh::findCell
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1446
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:175
Foam::linePointRef
line< point, const point & > linePointRef
Line using referred points.
Definition: linePointRef.H:47
shortestPathSet.H
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1063
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:372
Foam::sampledSet::mesh
const polyMesh & mesh() const
Definition: sampledSet.H:281
meshSearch.H
Foam::nl
constexpr char nl
Definition: Ostream.H:385
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::globalMeshData::coupledPatch
const indirectPrimitivePatch & coupledPatch() const
Return patch of all coupled faces.
Definition: globalMeshData.C:2066
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:144
FaceCellWave.H
Foam::UList< label >
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::primitiveMesh::nPoints
label nPoints() const
Number of mesh points.
Definition: primitiveMeshI.H:37
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
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:248
Foam::findIndices
labelList findIndices(const ListType &input, typename ListType::const_reference val, label start=0)
Linear search to find all occurrences of given element.
DynamicList.H
Foam::uindirectPrimitivePatch
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
A PrimitivePatch with UIndirectList for the faces, const reference for the point field.
Definition: uindirectPrimitivePatch.H:49
Foam::polyMesh::globalData
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1234
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::labelUList
UList< label > labelUList
A UList of labels.
Definition: UList.H:80
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:122
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::mkDir
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: MSwindows.C:507
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:298
OBJstream.H
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1082
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:85