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-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
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"
41
42// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43
44namespace Foam
45{
48}
49
50
51// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52
53Foam::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
143void 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
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
284bool 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
313Foam::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
376bool 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.
780Foam::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
911void 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 forAll(leakCcs, i)
1068 {
1069 str.write(leakCcs[i]);
1070 }
1071 }
1072 if (debug & 2)
1073 {
1074 OBJstream str
1075 (
1076 mesh.time().timePath()
1077 /"leakPath" + Foam::name(iter) + ".obj"
1078 );
1079 Pout<< "Writing new leak-path to " << str.name() << endl;
1080 for
1081 (
1082 label samplei = oldNSamplingPts+1;
1083 samplei < samplingPts.size();
1084 samplei++
1085 )
1086 {
1087 Pout<< " passing through cell "
1088 << samplingCells[samplei]
1089 << " at:" << mesh.cellCentres()[samplingCells[samplei]]
1090 << " distance:" << samplingCurveDist[samplei]
1091 << endl;
1092
1093 str.write
1094 (
1096 (
1097 samplingPts[samplei-1],
1098 samplingPts[samplei]
1099 )
1100 );
1101 }
1102 }
1103 }
1104 }
1105
1106 if (debug)
1107 {
1108 const fvMesh& fm = refCast<const fvMesh>(mesh);
1109
1110 const_cast<fvMesh&>(fm).setInstance(fm.time().timeName());
1112 (
1113 IOobject
1114 (
1115 "isLeakCell",
1116 fm.time().timeName(),
1117 fm,
1120 ),
1121 fm,
1123 );
1124 forAll(isLeakCell, celli)
1125 {
1126 fld[celli] = isLeakCell[celli];
1127 }
1128 fld.correctBoundaryConditions();
1129 fld.write();
1130 }
1131
1132 if (maxIter > 1 && iter == maxIter)
1133 {
1134 WarningInFunction << "Did not manage to close gap using " << iter
1135 << " leak paths" << nl << "This can cause problems when using the"
1136 << " paths to close leaks" << endl;
1137 }
1138}
1139
1140
1141void Foam::shortestPathSet::genSamples
1142(
1143 const bool markLeakPath,
1144 const label maxIter,
1145 const polyMesh& mesh,
1146 const labelUList& wallPatches,
1147 const bitSet& isBlockedFace
1148)
1149{
1150 // Storage for sample points
1151 DynamicList<point> samplingPts;
1152 DynamicList<label> samplingCells;
1153 DynamicList<label> samplingFaces;
1154 DynamicList<label> samplingSegments;
1155 DynamicList<scalar> samplingCurveDist;
1156
1157 // Seed faces and points on 'real' boundary
1158 bitSet isBlockedPoint(mesh.nPoints());
1159 {
1160 // Real boundaries
1161 const polyBoundaryMesh& pbm = mesh.boundaryMesh();
1162 for (label patchi : wallPatches)
1163 {
1164 const polyPatch& pp = pbm[patchi];
1165 forAll(pp, i)
1166 {
1167 isBlockedPoint.set(pp[i]);
1168 }
1169 }
1170
1171 // Meshed boundaries
1172 forAll(isBlockedFace, facei)
1173 {
1174 if (isBlockedFace[facei])
1175 {
1176 isBlockedPoint.set(mesh.faces()[facei]);
1177 }
1178 }
1179
1181 (
1182 mesh,
1183 isBlockedPoint,
1184 orEqOp<unsigned int>(),
1185 0u
1186 );
1187
1188 if (debug)
1189 {
1190 mkDir(mesh.time().timePath());
1191 OBJstream str(mesh.time().timePath()/"isBlockedPoint.obj");
1192 for (const auto& pointi : isBlockedPoint)
1193 {
1194 str.write(mesh.points()[pointi]);
1195 }
1196 Pout<< "Writing " << str.nVertices() << " points to " << str.name()
1197 << endl;
1198 }
1199 }
1200
1201
1202 bitSet isLeakPoint(isBlockedPoint);
1203 // Newly closed faces
1204 bitSet isLeakFace(mesh.nFaces());
1205 // All cells along leak paths
1206 bitSet isLeakCell(mesh.nCells());
1207
1208 label prevSegmenti = 0;
1209 scalar prevDistance = 0.0;
1210
1211 for (auto insidePoint : insidePoints_)
1212 {
1213 const label insideCelli = mesh.findCell(insidePoint);
1214
1215 for (auto outsidePoint : outsidePoints_)
1216 {
1217 const label nOldSamples = samplingSegments.size();
1218
1219 // Append any leak path to sampling*
1220 genSamples
1221 (
1222 markLeakPath,
1223 maxIter,
1224 mesh,
1225 isBlockedFace,
1226 insidePoint,
1227 insideCelli,
1228 outsidePoint,
1229
1230 samplingPts,
1231 samplingCells,
1232 samplingFaces,
1233 samplingSegments,
1234 samplingCurveDist,
1235 isLeakCell,
1236 isLeakFace,
1237 isLeakPoint
1238 );
1239
1240 // Make segment, distance consecutive
1241 label maxSegment = 0;
1242 scalar maxDistance = 0.0;
1243 for (label i = nOldSamples; i < samplingSegments.size(); ++i)
1244 {
1245 samplingSegments[i] += prevSegmenti;
1246 maxSegment = max(maxSegment, samplingSegments[i]);
1247 samplingCurveDist[i] += prevDistance;
1248 maxDistance = max(maxDistance, samplingCurveDist[i]);
1249 }
1250 prevSegmenti = returnReduce(maxSegment, maxOp<label>());
1251 prevDistance = returnReduce(maxDistance, maxOp<scalar>());
1252 }
1253 }
1254
1255 // Clean up leak faces (erode open edges). These are leak faces which are
1256 // not connected to another leak face or leak point
1257 erodeFaceSet(mesh, isBlockedPoint, isLeakFace);
1258
1259 leakFaces_ = isLeakFace.sortedToc();
1260
1261
1262 if (debug)
1263 {
1264 const faceList leakFaces(mesh.faces(), leakFaces_);
1265
1266 mkDir(mesh.time().timePath());
1267 OBJstream str(mesh.time().timePath()/"isLeakFace.obj");
1268 str.write(leakFaces, mesh.points(), false);
1269 Pout<< "Writing " << leakFaces.size() << " faces to " << str.name()
1270 << endl;
1271 }
1272
1273
1274 samplingPts.shrink();
1275 samplingCells.shrink();
1276 samplingFaces.shrink();
1277 samplingSegments.shrink();
1278 samplingCurveDist.shrink();
1279
1280 // Move into *this
1281 setSamples
1282 (
1283 std::move(samplingPts),
1284 std::move(samplingCells),
1285 std::move(samplingFaces),
1286 std::move(samplingSegments),
1287 std::move(samplingCurveDist)
1288 );
1289
1290 if (debug)
1291 {
1292 write(Info);
1293 }
1294}
1295
1296
1297// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1298
1300(
1301 const word& name,
1302 const polyMesh& mesh,
1303 const meshSearch& searchEngine,
1304 const word& axis,
1305 const bool markLeakPath,
1306 const label maxIter,
1307 const labelUList& wallPatches,
1308 const pointField& insidePoints,
1309 const pointField& outsidePoints,
1310 const boolList& isBlockedFace
1311)
1312:
1313 sampledSet(name, mesh, searchEngine, axis),
1314 insidePoints_(insidePoints),
1315 outsidePoints_(outsidePoints)
1316{
1317 if (debug)
1318 {
1319 fileName outputDir
1320 (
1321 mesh.time().globalPath()
1324 );
1325 outputDir.clean(); // Remove unneeded ".."
1326
1327 Info<< "shortestPathSet : Writing blocked faces to "
1328 << outputDir << endl;
1329
1330 const indirectPrimitivePatch setPatch
1331 (
1333 (
1334 mesh.faces(),
1335 findIndices(isBlockedFace, true)
1336 ),
1337 mesh.points()
1338 );
1339
1340 if (Pstream::parRun())
1341 {
1342 // Topological merge
1343 labelList pointToGlobal;
1344 labelList uniqueMeshPointLabels;
1346 autoPtr<globalIndex> globalFaces;
1347 faceList mergedFaces;
1348 pointField mergedPoints;
1350 (
1351 mesh,
1352 setPatch.localFaces(),
1353 setPatch.meshPoints(),
1354 setPatch.meshPointMap(),
1355
1356 pointToGlobal,
1357 uniqueMeshPointLabels,
1359 globalFaces,
1360
1361 mergedFaces,
1362 mergedPoints
1363 );
1364
1365 // Write
1366 if (Pstream::master())
1367 {
1369 (
1370 mergedPoints,
1371 mergedFaces,
1372 (outputDir / "blockedFace"),
1373 false // serial only - already merged
1374 );
1375
1376 writer.writeGeometry();
1377 }
1378 }
1379 else
1380 {
1382 (
1383 setPatch.localPoints(),
1384 setPatch.localFaces(),
1385 (outputDir / "blockedFace"),
1386 false // serial only - redundant
1387 );
1388
1389 writer.writeGeometry();
1390 }
1391 }
1392
1393 genSamples
1394 (
1395 markLeakPath,
1396 maxIter,
1397 mesh,
1398 wallPatches,
1399 bitSet(isBlockedFace)
1400 );
1401}
1402
1403
1405(
1406 const word& name,
1407 const polyMesh& mesh,
1408 const meshSearch& searchEngine,
1409 const dictionary& dict
1410)
1411:
1412 sampledSet(name, mesh, searchEngine, dict),
1413 insidePoints_(dict.get<pointField>("insidePoints")),
1414 outsidePoints_(dict.get<pointField>("outsidePoints"))
1415{
1416 const label maxIter(dict.getOrDefault<label>("maxIter", 50));
1417 const bool markLeakPath(dict.getOrDefault("markLeakPath", true));
1418
1419 const polyBoundaryMesh& pbm = mesh.boundaryMesh();
1420
1421 DynamicList<label> wallPatches(pbm.size());
1422 forAll(pbm, patchi)
1423 {
1424 const polyPatch& pp = pbm[patchi];
1425 if (!pp.coupled() && !isA<emptyPolyPatch>(pp))
1426 {
1427 wallPatches.append(patchi);
1428 }
1429 }
1430
1431 genSamples(markLeakPath, maxIter, mesh, wallPatches, bitSet());
1432}
1433
1434
1435// ************************************************************************* //
scalar y
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
vtk::internalMeshWriter writer(topoMesh, topoCells, vtk::formatType::INLINE_ASCII, runTime.path()/"blockTopology")
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:503
A List with indirect addressing.
Definition: IndirectList.H:119
void transfer(List< T > &list)
Definition: List.C:447
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:175
virtual const fileName & name() const
Get the name of the stream.
Definition: OSstream.H:107
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
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, labelList &pointMergeMap, const bool useLocal=false)
Gather points and faces onto master and merge into single patch.
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.
A list of faces which address into the list of points.
const Map< label > & meshPointMap() const
Mesh point map.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
const Field< point_type > & localPoints() const
Return pointField of points in patch.
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
static void combineAllGather(const List< commsStruct > &comms, T &value, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
fileName globalPath() const
Return global path for the case.
Definition: TimePathsI.H:80
fileName timePath() const
Return current time path.
Definition: Time.H:375
label fcIndex(const label i) const noexcept
Definition: UListI.H:60
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:66
const scalarList & distance() const noexcept
Return the cumulative distance.
Definition: coordSet.H:152
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
A class for handling file names.
Definition: fileName.H:76
static bool clean(std::string &str)
Definition: fileName.C:199
static word outputPrefix
Directory prefix.
void sync()
Do all: synchronise all IOFields and objectRegistry.
Definition: syncObjects.C:70
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
const mapDistribute & globalEdgeSlavesMap() const
const indirectPrimitivePatch & coupledPatch() const
Return patch of all coupled faces.
label nTotalCells() const noexcept
Return total number of cells in decomposed mesh.
Calculates points shared by more than two processor patches or cyclic patches.
Definition: globalPoints.H:103
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition: meshSearch.H:61
const Time & time() const noexcept
Return time registry.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1108
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1522
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1121
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1310
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:860
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1127
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1083
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:380
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
label nBoundaryFaces() const noexcept
Number of boundary faces (== nFaces - nInternalFaces)
const vectorField & faceCentres() const
label nInternalFaces() const noexcept
Number of internal faces.
const vectorField & cellCentres() const
label nPoints() const noexcept
Number of mesh points.
label nCells() const noexcept
Number of mesh cells.
label nFaces() const noexcept
Number of mesh faces.
const cellList & cells() const
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:86
const polyMesh & mesh() const noexcept
Definition: sampledSet.H:319
Finds shortest path (in terms of cell centres) to walk on mesh from any point in insidePoints to any ...
splitCell * master() const
Definition: splitCell.H:113
static void syncFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop)
Synchronize values on all mesh faces.
Definition: syncTools.H:396
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues)
Swap coupled boundary face values. Uses eqOp.
Definition: syncTools.H:445
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
Write faces/points (optionally with fields) as a vtp file or a legacy vtk file.
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
const dimensionSet dimless
Dimensionless.
List< label > labelList
A List of labels.
Definition: List.H:66
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:515
line< point, const point & > linePointRef
A line using referred points.
Definition: linePointRef.H:47
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
messageStream Info
Information stream (stdout output on master, null elsewhere)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
vector point
Point is a vector.
Definition: point.H:43
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
constexpr label labelMax
Definition: label.H:61
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
A PrimitivePatch with UIndirectList for the faces, const reference for the point field.
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
List< bool > boolList
A List of bools.
Definition: List.H:64
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
labelList findIndices(const ListType &input, typename ListType::const_reference val, label start=0)
Linear search to find all occurrences of given element.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
UList< label > labelUList
A UList of labels.
Definition: UList.H:85
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
runTime write()
labelList f(nPoints)
insidePoints((1 2 3))
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333