polyTopoChange.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2015-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "polyTopoChange.H"
30#include "polyMesh.H"
31#include "polyAddPoint.H"
32#include "polyModifyPoint.H"
33#include "polyRemovePoint.H"
34#include "polyAddFace.H"
35#include "polyModifyFace.H"
36#include "polyRemoveFace.H"
37#include "polyAddCell.H"
38#include "polyModifyCell.H"
39#include "polyRemoveCell.H"
40#include "CircularBuffer.H"
41#include "CompactListList.H"
42#include "objectMap.H"
43#include "processorPolyPatch.H"
44#include "mapPolyMesh.H"
45
46// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
47
48namespace Foam
49{
51}
52
53
54// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
55
56// Renumber with special handling for merged items (marked with <-1)
57void Foam::polyTopoChange::renumberReverseMap
58(
59 const labelUList& oldToNew,
60 DynamicList<label>& elems
61)
62{
63 forAll(elems, elemI)
64 {
65 const label val = elems[elemI];
66
67 if (val >= 0)
68 {
69 elems[elemI] = oldToNew[val];
70 }
71 else if (val < -1)
72 {
73 const label mergedVal = -val-2;
74 elems[elemI] = -oldToNew[mergedVal]-2;
75 }
76 }
77}
78
79
80void Foam::polyTopoChange::renumber
81(
82 const labelUList& oldToNew,
83 labelHashSet& labels
84)
85{
86 labelHashSet newSet(labels.capacity());
87
88 for (const label val : labels)
89 {
90 const label newVal = oldToNew[val];
91
92 if (newVal >= 0)
93 {
94 newSet.insert(newVal);
95 }
96 }
97
98 labels.transfer(newSet);
99}
100
101
102// Renumber and remove -1 elements.
103void Foam::polyTopoChange::renumberCompact
104(
105 const labelUList& oldToNew,
106 labelList& elems
107)
108{
109 label nElem = 0;
110
111 for (const label val : elems)
112 {
113 const label newVal = oldToNew[val];
114
115 if (newVal != -1)
116 {
117 elems[nElem++] = newVal;
118 }
119 }
120 elems.setSize(nElem);
121}
122
123
124void Foam::polyTopoChange::countMap
125(
126 const labelUList& map,
127 const labelUList& reverseMap,
128 label& nAdd,
129 label& nInflate,
130 label& nMerge,
131 label& nRemove
132)
133{
134 nAdd = 0;
135 nInflate = 0;
136 nMerge = 0;
137 nRemove = 0;
138
139 forAll(map, newCelli)
140 {
141 const label oldCelli = map[newCelli];
142
143 if (oldCelli >= 0)
144 {
145 if
146 (
147 oldCelli < reverseMap.size()
148 && reverseMap[oldCelli] == newCelli
149 )
150 {
151 // unchanged
152 }
153 else
154 {
155 // Added (from another cell v.s. inflated from face/point)
156 nAdd++;
157 }
158 }
159 else if (oldCelli == -1)
160 {
161 // Created from nothing
162 nInflate++;
163 }
164 else
165 {
167 << " new:" << newCelli << abort(FatalError);
168 }
169 }
170
171 forAll(reverseMap, oldCelli)
172 {
173 const label newCelli = reverseMap[oldCelli];
174
175 if (newCelli >= 0)
176 {
177 // unchanged
178 }
179 else if (newCelli == -1)
180 {
181 // removed
182 nRemove++;
183 }
184 else
185 {
186 // merged into -newCelli-2
187 nMerge++;
188 }
189 }
190}
191
192
193void Foam::polyTopoChange::writeMeshStats(const polyMesh& mesh, Ostream& os)
194{
195 const polyBoundaryMesh& patches = mesh.boundaryMesh();
196
197 labelList patchSizes(patches.size());
198 labelList patchStarts(patches.size());
199 forAll(patches, patchi)
200 {
201 patchSizes[patchi] = patches[patchi].size();
202 patchStarts[patchi] = patches[patchi].start();
203 }
204
205 os << " Points : " << mesh.nPoints() << nl
206 << " Faces : " << mesh.nFaces() << nl
207 << " Cells : " << mesh.nCells() << nl
208 << " PatchSizes : " << patchSizes << nl
209 << " PatchStarts : " << patchStarts << nl
210 << endl;
211}
212
213
214void Foam::polyTopoChange::getMergeSets
215(
216 const labelUList& reverseCellMap,
217 const labelUList& cellMap,
218 List<objectMap>& cellsFromCells
219)
220{
221 // Per new cell the number of old cells that have been merged into it
222 labelList nMerged(cellMap.size(), 1);
223
224 forAll(reverseCellMap, oldCelli)
225 {
226 const label newCelli = reverseCellMap[oldCelli];
227
228 if (newCelli < -1)
229 {
230 label mergeCelli = -newCelli-2;
231
232 nMerged[mergeCelli]++;
233 }
234 }
235
236 // From merged cell to set index
237 labelList cellToMergeSet(cellMap.size(), -1);
238
239 label nSets = 0;
240
241 forAll(nMerged, celli)
242 {
243 if (nMerged[celli] > 1)
244 {
245 cellToMergeSet[celli] = nSets++;
246 }
247 }
248
249 // Collect cell labels.
250 // Each objectMap will have
251 // - index : new mesh cell label
252 // - masterObjects : list of old cells that have been merged. Element 0
253 // will be the original destination cell label.
254
255 cellsFromCells.setSize(nSets);
256
257 forAll(reverseCellMap, oldCelli)
258 {
259 const label newCelli = reverseCellMap[oldCelli];
260
261 if (newCelli < -1)
262 {
263 const label mergeCelli = -newCelli-2;
264
265 // oldCelli was merged into mergeCelli
266
267 const label setI = cellToMergeSet[mergeCelli];
268
269 objectMap& mergeSet = cellsFromCells[setI];
270
271 if (mergeSet.masterObjects().empty())
272 {
273 // First occurrence of master cell mergeCelli
274
275 mergeSet.index() = mergeCelli;
276 mergeSet.masterObjects().setSize(nMerged[mergeCelli]);
277
278 // old master label
279 mergeSet.masterObjects()[0] = cellMap[mergeCelli];
280
281 // old slave label
282 mergeSet.masterObjects()[1] = oldCelli;
283
284 nMerged[mergeCelli] = 2;
285 }
286 else
287 {
288 mergeSet.masterObjects()[nMerged[mergeCelli]++] = oldCelli;
289 }
290 }
291 }
292}
293
294
295bool Foam::polyTopoChange::hasValidPoints(const face& f) const
296{
297 for (const label fp : f)
298 {
299 if (fp < 0 || fp >= points_.size())
300 {
301 return false;
302 }
303 }
304 return true;
305}
306
307
309{
311 forAll(f, fp)
312 {
313 if (f[fp] < 0 && f[fp] >= points_.size())
314 {
316 << "Problem." << abort(FatalError);
317 }
318 points[fp] = points_[f[fp]];
319 }
320 return points;
321}
322
323
324void Foam::polyTopoChange::checkFace
325(
326 const face& f,
327 const label facei,
328 const label own,
329 const label nei,
330 const label patchi,
331 const label zoneI
332) const
333{
334 if (nei == -1)
335 {
336 if (own == -1 && zoneI != -1)
337 {
338 // retired face
339 }
340 else if (patchi == -1 || patchi >= nPatches_)
341 {
343 << "Face has no neighbour (so external) but does not have"
344 << " a valid patch" << nl
345 << "f:" << f
346 << " facei(-1 if added face):" << facei
347 << " own:" << own << " nei:" << nei
348 << " patchi:" << patchi << nl;
349 if (hasValidPoints(f))
350 {
352 << "points (removed points marked with "
353 << vector::max << ") " << facePoints(f);
354 }
356 }
357 }
358 else
359 {
360 if (patchi != -1)
361 {
363 << "Cannot both have valid patchi and neighbour" << nl
364 << "f:" << f
365 << " facei(-1 if added face):" << facei
366 << " own:" << own << " nei:" << nei
367 << " patchi:" << patchi << nl;
368 if (hasValidPoints(f))
369 {
371 << "points (removed points marked with "
372 << vector::max << ") : " << facePoints(f);
373 }
375 }
376
377 if (nei <= own)
378 {
380 << "Owner cell label should be less than neighbour cell label"
381 << nl
382 << "f:" << f
383 << " facei(-1 if added face):" << facei
384 << " own:" << own << " nei:" << nei
385 << " patchi:" << patchi << nl;
386 if (hasValidPoints(f))
387 {
389 << "points (removed points marked with "
390 << vector::max << ") : " << facePoints(f);
391 }
393 }
394 }
395
396 if (f.size() < 3 || f.found(-1))
397 {
399 << "Illegal vertices in face"
400 << nl
401 << "f:" << f
402 << " facei(-1 if added face):" << facei
403 << " own:" << own << " nei:" << nei
404 << " patchi:" << patchi << nl;
405 if (hasValidPoints(f))
406 {
408 << "points (removed points marked with "
409 << vector::max << ") : " << facePoints(f);
410 }
412 }
413 if (facei >= 0 && facei < faces_.size() && faceRemoved(facei))
414 {
416 << "Face already marked for removal"
417 << nl
418 << "f:" << f
419 << " facei(-1 if added face):" << facei
420 << " own:" << own << " nei:" << nei
421 << " patchi:" << patchi << nl;
422 if (hasValidPoints(f))
423 {
425 << "points (removed points marked with "
426 << vector::max << ") : " << facePoints(f);
427 }
429 }
430 forAll(f, fp)
431 {
432 if (f[fp] < points_.size() && pointRemoved(f[fp]))
433 {
435 << "Face uses removed vertices"
436 << nl
437 << "f:" << f
438 << " facei(-1 if added face):" << facei
439 << " own:" << own << " nei:" << nei
440 << " patchi:" << patchi << nl;
441 if (hasValidPoints(f))
442 {
444 << "points (removed points marked with "
445 << vector::max << ") : " << facePoints(f);
446 }
448 }
449 }
450}
451
452
453void Foam::polyTopoChange::makeCells
454(
455 const label nActiveFaces,
456 labelList& cellFaces,
457 labelList& cellFaceOffsets
458) const
459{
460 cellFaces.setSize(2*nActiveFaces);
461 cellFaceOffsets.setSize(cellMap_.size() + 1);
462
463 // Faces per cell
464 labelList nNbrs(cellMap_.size(), Zero);
465
466 // 1. Count faces per cell
467
468 for (label facei = 0; facei < nActiveFaces; facei++)
469 {
470 if (faceOwner_[facei] < 0)
471 {
472 pointField newPoints;
473 if (facei < faces_.size())
474 {
475 const face& f = faces_[facei];
476 newPoints.setSize(f.size(), vector::max);
477 forAll(f, fp)
478 {
479 if (f[fp] < points_.size())
480 {
481 newPoints[fp] = points_[f[fp]];
482 }
483 }
484 }
485
486
488 << "Face " << facei << " is active but its owner has"
489 << " been deleted. This is usually due to deleting cells"
490 << " without modifying exposed faces to be boundary faces."
491 << exit(FatalError);
492 }
493 nNbrs[faceOwner_[facei]]++;
494 }
495 for (label facei = 0; facei < nActiveFaces; facei++)
496 {
497 if (faceNeighbour_[facei] >= 0)
498 {
499 nNbrs[faceNeighbour_[facei]]++;
500 }
501 }
502
503 // 2. Calculate offsets
504
505 cellFaceOffsets[0] = 0;
506 forAll(nNbrs, celli)
507 {
508 cellFaceOffsets[celli+1] = cellFaceOffsets[celli] + nNbrs[celli];
509 }
510
511 // 3. Fill faces per cell
512
513 // reset the whole list to use as counter
514 nNbrs = 0;
515
516 for (label facei = 0; facei < nActiveFaces; facei++)
517 {
518 label celli = faceOwner_[facei];
519
520 cellFaces[cellFaceOffsets[celli] + nNbrs[celli]++] = facei;
521 }
522
523 for (label facei = 0; facei < nActiveFaces; facei++)
524 {
525 label celli = faceNeighbour_[facei];
526
527 if (celli >= 0)
528 {
529 cellFaces[cellFaceOffsets[celli] + nNbrs[celli]++] = facei;
530 }
531 }
532
533 // Last offset points to beyond end of cellFaces.
534 cellFaces.setSize(cellFaceOffsets[cellMap_.size()]);
535}
536
537
538// Create cell-cell addressing. Called after compaction (but before ordering)
539// of faces
540void Foam::polyTopoChange::makeCellCells
541(
542 const label nActiveFaces,
543 CompactListList<label>& cellCells
544) const
545{
546 // Neighbours per cell
547 labelList nNbrs(cellMap_.size(), Zero);
548
549 // 1. Count neighbours (through internal faces) per cell
550
551 for (label facei = 0; facei < nActiveFaces; facei++)
552 {
553 if (faceNeighbour_[facei] >= 0)
554 {
555 nNbrs[faceOwner_[facei]]++;
556 nNbrs[faceNeighbour_[facei]]++;
557 }
558 }
559
560 // 2. Construct csr
561 cellCells.setSize(nNbrs);
562
563
564 // 3. Fill faces per cell
565
566 // reset the whole list to use as counter
567 nNbrs = 0;
568
569 for (label facei = 0; facei < nActiveFaces; facei++)
570 {
571 label nei = faceNeighbour_[facei];
572
573 if (nei >= 0)
574 {
575 label own = faceOwner_[facei];
576 cellCells(own, nNbrs[own]++) = nei;
577 cellCells(nei, nNbrs[nei]++) = own;
578 }
579 }
580}
581
582
583// Cell ordering (based on bandCompression).
584// Handles removed cells. Returns number of remaining cells.
585Foam::label Foam::polyTopoChange::getCellOrder
586(
587 const CompactListList<label>& cellCellAddressing,
588 labelList& oldToNew
589) const
590{
591 const label nOldCells(cellCellAddressing.size());
592
593 // Which cells are visited/unvisited
594 bitSet unvisited(nOldCells, true);
595
596 // Exclude removed cells
597 for (label celli = 0; celli < nOldCells; ++celli)
598 {
599 if (cellRemoved(celli))
600 {
601 unvisited.unset(celli);
602 }
603 }
604
605 // The new output order
606 labelList newOrder(unvisited.count());
607
608
609 // Various work arrays
610 // ~~~~~~~~~~~~~~~~~~~
611
612 //- Neighbour cells
613 DynamicList<label> nbrCells;
614
615 //- Neighbour ordering
616 DynamicList<label> nbrOrder;
617
618 //- Corresponding weights for neighbour cells
619 DynamicList<label> weights;
620
621 // FIFO buffer for renumbering.
622 CircularBuffer<label> queuedCells(1024);
623
624
625 label cellInOrder = 0;
626
627 while (true)
628 {
629 // For a disconnected region find the lowest connected cell.
630 label currCelli = -1;
631 label minCount = labelMax;
632
633 for (const label celli : unvisited)
634 {
635 const label nbrCount = cellCellAddressing[celli].size();
636
637 if (minCount > nbrCount)
638 {
639 minCount = nbrCount;
640 currCelli = celli;
641 }
642 }
643
644 if (currCelli == -1)
645 {
646 break;
647 }
648
649
650 // Starting from currCelli - walk breadth-first
651
652 queuedCells.append(currCelli);
653
654 // Loop through queuedCells list. Add the first cell into the
655 // cell order if it has not already been visited and ask for its
656 // neighbours. If the neighbour in question has not been visited,
657 // add it to the end of the queuedCells list
658
659 while (!queuedCells.empty())
660 {
661 // Process as FIFO
662 currCelli = queuedCells.first();
663 queuedCells.pop_front();
664
665 if (unvisited.test(currCelli))
666 {
667 // First visit...
668 unvisited.unset(currCelli);
669
670 // Add into cellOrder
671 newOrder[cellInOrder] = currCelli;
672 ++cellInOrder;
673
674 // Find if the neighbours have been visited
675 const auto& neighbours = cellCellAddressing[currCelli];
676
677 // Add in increasing order of connectivity
678
679 // 1. Count neighbours of unvisited neighbours
680 nbrCells.clear();
681 weights.clear();
682
683 for (const label nbr : neighbours)
684 {
685 const label nbrCount = cellCellAddressing[nbr].size();
686
687 if (unvisited.test(nbr))
688 {
689 // Not visited (or removed), add to the list
690 nbrCells.append(nbr);
691 weights.append(nbrCount);
692 }
693 }
694
695 // Resize DynamicList prior to sortedOrder
696 nbrOrder.resize_nocopy(weights.size());
697
698 // 2. Ascending order
699 Foam::sortedOrder(weights, nbrOrder);
700
701 // 3. Add to FIFO in sorted order
702 for (const label nbrIdx : nbrOrder)
703 {
704 queuedCells.append(nbrCells[nbrIdx]);
705 }
706 }
707 }
708 }
709
710 // Now we have new-to-old in newOrder.
711 newOrder.resize(cellInOrder); // Extra safety, but should be a no-op
712
713 // Invert to get old-to-new. Make sure removed (i.e. unmapped) cells are -1.
714 oldToNew = invert(nOldCells, newOrder);
715
716 return cellInOrder;
717}
718
719
720// Determine order for faces:
721// - upper-triangular order for internal faces
722// - external faces after internal faces and in patch order.
723void Foam::polyTopoChange::getFaceOrder
724(
725 const label nActiveFaces,
726 const labelUList& cellFaces,
727 const labelUList& cellFaceOffsets,
728
729 labelList& oldToNew,
730 labelList& patchSizes,
731 labelList& patchStarts
732) const
733{
734 oldToNew.setSize(faceOwner_.size());
735 oldToNew = -1;
736
737 // First unassigned face
738 label newFacei = 0;
739
740 labelList nbr;
741 labelList order;
742
743 forAll(cellMap_, celli)
744 {
745 label startOfCell = cellFaceOffsets[celli];
746 label nFaces = cellFaceOffsets[celli+1] - startOfCell;
747
748 // Neighbouring cells
749 nbr.setSize(nFaces);
750
751 for (label i = 0; i < nFaces; i++)
752 {
753 label facei = cellFaces[startOfCell + i];
754
755 label nbrCelli = faceNeighbour_[facei];
756
757 if (facei >= nActiveFaces)
758 {
759 // Retired face.
760 nbr[i] = -1;
761 }
762 else if (nbrCelli != -1)
763 {
764 // Internal face. Get cell on other side.
765 if (nbrCelli == celli)
766 {
767 nbrCelli = faceOwner_[facei];
768 }
769
770 if (celli < nbrCelli)
771 {
772 // Celli is master
773 nbr[i] = nbrCelli;
774 }
775 else
776 {
777 // nbrCell is master. Let it handle this face.
778 nbr[i] = -1;
779 }
780 }
781 else
782 {
783 // External face. Do later.
784 nbr[i] = -1;
785 }
786 }
787
788 sortedOrder(nbr, order);
789
790 for (const label index : order)
791 {
792 if (nbr[index] != -1)
793 {
794 oldToNew[cellFaces[startOfCell + index]] = newFacei++;
795 }
796 }
797 }
798
799
800 // Pick up all patch faces in patch face order.
801 patchStarts.setSize(nPatches_);
802 patchStarts = 0;
803 patchSizes.setSize(nPatches_);
804 patchSizes = 0;
805
806 if (nPatches_ > 0)
807 {
808 patchStarts[0] = newFacei;
809
810 for (label facei = 0; facei < nActiveFaces; facei++)
811 {
812 if (region_[facei] >= 0)
813 {
814 patchSizes[region_[facei]]++;
815 }
816 }
817
818 label facei = patchStarts[0];
819
820 forAll(patchStarts, patchi)
821 {
822 patchStarts[patchi] = facei;
823 facei += patchSizes[patchi];
824 }
825 }
826
827 //if (debug)
828 //{
829 // Pout<< "patchSizes:" << patchSizes << nl
830 // << "patchStarts:" << patchStarts << endl;
831 //}
832
833 labelList workPatchStarts(patchStarts);
834
835 for (label facei = 0; facei < nActiveFaces; facei++)
836 {
837 if (region_[facei] >= 0)
838 {
839 oldToNew[facei] = workPatchStarts[region_[facei]]++;
840 }
841 }
842
843 // Retired faces.
844 for (label facei = nActiveFaces; facei < oldToNew.size(); facei++)
845 {
846 oldToNew[facei] = facei;
847 }
848
849 // Check done all faces.
850 forAll(oldToNew, facei)
851 {
852 if (oldToNew[facei] == -1)
853 {
855 << "Did not determine new position"
856 << " for face " << facei
857 << " owner " << faceOwner_[facei]
858 << " neighbour " << faceNeighbour_[facei]
859 << " region " << region_[facei] << endl
860 << "This is usually caused by not specifying a patch for"
861 << " a boundary face." << nl
862 << "Switch on the polyTopoChange::debug flag to catch"
863 << " this error earlier." << nl;
864 if (hasValidPoints(faces_[facei]))
865 {
867 << "points (removed points marked with "
868 << vector::max << ") " << facePoints(faces_[facei]);
869 }
871 }
872 }
873}
874
875
876// Reorder and compact faces according to map.
877void Foam::polyTopoChange::reorderCompactFaces
878(
879 const label newSize,
880 const labelUList& oldToNew
881)
882{
883 reorder(oldToNew, faces_);
884 faces_.setCapacity(newSize);
885
886 reorder(oldToNew, region_);
887 region_.setCapacity(newSize);
888
889 reorder(oldToNew, faceOwner_);
890 faceOwner_.setCapacity(newSize);
891
892 reorder(oldToNew, faceNeighbour_);
893 faceNeighbour_.setCapacity(newSize);
894
895 // Update faceMaps.
896 reorder(oldToNew, faceMap_);
897 faceMap_.setCapacity(newSize);
898
899 renumberReverseMap(oldToNew, reverseFaceMap_);
900
901 renumberKey(oldToNew, faceFromPoint_);
902 renumberKey(oldToNew, faceFromEdge_);
903 inplaceReorder(oldToNew, flipFaceFlux_);
904 flipFaceFlux_.setCapacity(newSize);
905 renumberKey(oldToNew, faceZone_);
906 inplaceReorder(oldToNew, faceZoneFlip_);
907 faceZoneFlip_.setCapacity(newSize);
908}
909
910
912{
913 points_.shrink();
914 pointMap_.shrink();
915 reversePointMap_.shrink();
916
917 faces_.shrink();
918 region_.shrink();
919 faceOwner_.shrink();
920 faceNeighbour_.shrink();
921 faceMap_.shrink();
922 reverseFaceMap_.shrink();
923
924 cellMap_.shrink();
925 reverseCellMap_.shrink();
926 cellZone_.shrink();
927}
928
929
930// Compact all and orders points and faces:
931// - points into internal followed by external points
932// - internalfaces upper-triangular
933// - externalfaces after internal ones.
935(
936 const bool orderCells,
937 const bool orderPoints,
938 label& nInternalPoints,
939 labelList& patchSizes,
940 labelList& patchStarts
941)
942{
943 shrink();
944
945 // Compact points
946 label nActivePoints = 0;
947 {
948 labelList localPointMap(points_.size(), -1);
949 label newPointi = 0;
950
951 if (!orderPoints)
952 {
953 nInternalPoints = -1;
954
955 forAll(points_, pointi)
956 {
957 if (!pointRemoved(pointi) && !retiredPoints_.found(pointi))
958 {
959 localPointMap[pointi] = newPointi++;
960 }
961 }
962 nActivePoints = newPointi;
963 }
964 else
965 {
966 forAll(points_, pointi)
967 {
968 if (!pointRemoved(pointi) && !retiredPoints_.found(pointi))
969 {
970 nActivePoints++;
971 }
972 }
973
974 // Mark boundary points
975 forAll(faceOwner_, facei)
976 {
977 if
978 (
979 !faceRemoved(facei)
980 && faceOwner_[facei] >= 0
981 && faceNeighbour_[facei] < 0
982 )
983 {
984 // Valid boundary face
985 const face& f = faces_[facei];
986
987 forAll(f, fp)
988 {
989 label pointi = f[fp];
990
991 if (localPointMap[pointi] == -1)
992 {
993 if
994 (
995 pointRemoved(pointi)
996 || retiredPoints_.found(pointi)
997 )
998 {
1000 << "Removed or retired point " << pointi
1001 << " in face " << f
1002 << " at position " << facei << endl
1003 << "Probably face has not been adapted for"
1004 << " removed points." << abort(FatalError);
1005 }
1006 localPointMap[pointi] = newPointi++;
1007 }
1008 }
1009 }
1010 }
1011
1012 label nBoundaryPoints = newPointi;
1013 nInternalPoints = nActivePoints - nBoundaryPoints;
1014
1015 // Move the boundary addressing up
1016 forAll(localPointMap, pointi)
1017 {
1018 if (localPointMap[pointi] != -1)
1019 {
1020 localPointMap[pointi] += nInternalPoints;
1021 }
1022 }
1023
1024 newPointi = 0;
1025
1026 // Mark internal points
1027 forAll(faceOwner_, facei)
1028 {
1029 if
1030 (
1031 !faceRemoved(facei)
1032 && faceOwner_[facei] >= 0
1033 && faceNeighbour_[facei] >= 0
1034 )
1035 {
1036 // Valid internal face
1037 const face& f = faces_[facei];
1038
1039 forAll(f, fp)
1040 {
1041 label pointi = f[fp];
1042
1043 if (localPointMap[pointi] == -1)
1044 {
1045 if
1046 (
1047 pointRemoved(pointi)
1048 || retiredPoints_.found(pointi)
1049 )
1050 {
1052 << "Removed or retired point " << pointi
1053 << " in face " << f
1054 << " at position " << facei << endl
1055 << "Probably face has not been adapted for"
1056 << " removed points." << abort(FatalError);
1057 }
1058 localPointMap[pointi] = newPointi++;
1059 }
1060 }
1061 }
1062 }
1063
1064 if (newPointi != nInternalPoints)
1065 {
1067 << "Problem." << abort(FatalError);
1068 }
1069 newPointi = nActivePoints;
1070 }
1071
1072 for (const label pointi : retiredPoints_)
1073 {
1074 localPointMap[pointi] = newPointi++;
1075 }
1076
1077
1078 if (debug)
1079 {
1080 Pout<< "Points : active:" << nActivePoints
1081 << " removed:" << points_.size()-newPointi << endl;
1082 }
1083
1084 reorder(localPointMap, points_);
1085 points_.setCapacity(newPointi);
1086
1087 // Update pointMaps
1088 reorder(localPointMap, pointMap_);
1089 pointMap_.setCapacity(newPointi);
1090 renumberReverseMap(localPointMap, reversePointMap_);
1091
1092 renumberKey(localPointMap, pointZone_);
1093 renumber(localPointMap, retiredPoints_);
1094
1095 // Use map to relabel face vertices
1096 forAll(faces_, facei)
1097 {
1098 face& f = faces_[facei];
1099
1100 //labelList oldF(f);
1101 renumberCompact(localPointMap, f);
1102
1103 if (!faceRemoved(facei) && f.size() < 3)
1104 {
1106 << "Created illegal face " << f
1107 //<< " from face " << oldF
1108 << " at position:" << facei
1109 << " when filtering removed points"
1110 << abort(FatalError);
1111 }
1112 }
1113 }
1114
1115
1116 // Compact faces.
1117 {
1118 labelList localFaceMap(faces_.size(), -1);
1119 label newFacei = 0;
1120
1121 forAll(faces_, facei)
1122 {
1123 if (!faceRemoved(facei) && faceOwner_[facei] >= 0)
1124 {
1125 localFaceMap[facei] = newFacei++;
1126 }
1127 }
1128 nActiveFaces_ = newFacei;
1129
1130 forAll(faces_, facei)
1131 {
1132 if (!faceRemoved(facei) && faceOwner_[facei] < 0)
1133 {
1134 // Retired face
1135 localFaceMap[facei] = newFacei++;
1136 }
1137 }
1138
1139 if (debug)
1140 {
1141 Pout<< "Faces : active:" << nActiveFaces_
1142 << " removed:" << faces_.size()-newFacei << endl;
1143 }
1144
1145 // Reorder faces.
1146 reorderCompactFaces(newFacei, localFaceMap);
1147 }
1148
1149 // Compact cells.
1150 {
1151 labelList localCellMap;
1152 label newCelli;
1153
1154 if (orderCells)
1155 {
1156 // Construct cellCell addressing
1157 CompactListList<label> cellCells;
1158 makeCellCells(nActiveFaces_, cellCells);
1159
1160 // Cell ordering (based on bandCompression). Handles removed cells.
1161 newCelli = getCellOrder(cellCells, localCellMap);
1162 }
1163 else
1164 {
1165 // Compact out removed cells
1166 localCellMap.setSize(cellMap_.size());
1167 localCellMap = -1;
1168
1169 newCelli = 0;
1170 forAll(cellMap_, celli)
1171 {
1172 if (!cellRemoved(celli))
1173 {
1174 localCellMap[celli] = newCelli++;
1175 }
1176 }
1177 }
1178
1179 if (debug)
1180 {
1181 Pout<< "Cells : active:" << newCelli
1182 << " removed:" << cellMap_.size()-newCelli << endl;
1183 }
1184
1185 // Renumber -if cells reordered or -if cells removed
1186 if (orderCells || (newCelli != cellMap_.size()))
1187 {
1188 reorder(localCellMap, cellMap_);
1189 cellMap_.setCapacity(newCelli);
1190 renumberReverseMap(localCellMap, reverseCellMap_);
1191
1192 reorder(localCellMap, cellZone_);
1193 cellZone_.setCapacity(newCelli);
1194
1195 renumberKey(localCellMap, cellFromPoint_);
1196 renumberKey(localCellMap, cellFromEdge_);
1197 renumberKey(localCellMap, cellFromFace_);
1198
1199 // Renumber owner/neighbour. Take into account if neighbour suddenly
1200 // gets lower cell than owner.
1201 forAll(faceOwner_, facei)
1202 {
1203 label own = faceOwner_[facei];
1204 label nei = faceNeighbour_[facei];
1205
1206 if (own >= 0)
1207 {
1208 // Update owner
1209 faceOwner_[facei] = localCellMap[own];
1210
1211 if (nei >= 0)
1212 {
1213 // Update neighbour.
1214 faceNeighbour_[facei] = localCellMap[nei];
1215
1216 // Check if face needs reversing.
1217 if
1218 (
1219 faceNeighbour_[facei] >= 0
1220 && faceNeighbour_[facei] < faceOwner_[facei]
1221 )
1222 {
1223 faces_[facei].flip();
1224 std::swap(faceOwner_[facei], faceNeighbour_[facei]);
1225 flipFaceFlux_.flip(facei);
1226 faceZoneFlip_.flip(facei);
1227 }
1228 }
1229 }
1230 else if (nei >= 0)
1231 {
1232 // Update neighbour.
1233 faceNeighbour_[facei] = localCellMap[nei];
1234 }
1235 }
1236 }
1237 }
1238
1239 // Reorder faces into upper-triangular and patch ordering
1240 {
1241 // Create cells (packed storage)
1242 labelList cellFaces;
1243 labelList cellFaceOffsets;
1244 makeCells(nActiveFaces_, cellFaces, cellFaceOffsets);
1245
1246 // Do upper triangular order and patch sorting
1247 labelList localFaceMap;
1248 getFaceOrder
1249 (
1250 nActiveFaces_,
1251 cellFaces,
1252 cellFaceOffsets,
1253
1254 localFaceMap,
1255 patchSizes,
1256 patchStarts
1257 );
1258
1259 // Reorder faces.
1260 reorderCompactFaces(localFaceMap.size(), localFaceMap);
1261 }
1262}
1263
1264
1265// Find faces to interpolate to create value for new face. Only used if
1266// face was inflated from edge or point. Internal faces should only be
1267// created from internal faces, external faces only from external faces
1268// (and ideally the same patch)
1269// Is bit problematic if there are no faces to select, i.e. in polyDualMesh
1270// an internal face can be created from a boundary edge with no internal
1271// faces connected to it.
1272Foam::labelList Foam::polyTopoChange::selectFaces
1273(
1274 const primitiveMesh& mesh,
1275 const labelUList& faceLabels,
1276 const bool internalFacesOnly
1277)
1278{
1279 label nFaces = 0;
1280
1281 forAll(faceLabels, i)
1282 {
1283 label facei = faceLabels[i];
1284
1285 if (internalFacesOnly == mesh.isInternalFace(facei))
1286 {
1287 nFaces++;
1288 }
1289 }
1290
1291 labelList collectedFaces;
1292
1293 if (nFaces == 0)
1294 {
1295 // Did not find any faces of the correct type so just use any old
1296 // face.
1297 collectedFaces = faceLabels;
1298 }
1299 else
1300 {
1301 collectedFaces.setSize(nFaces);
1302
1303 nFaces = 0;
1304
1305 forAll(faceLabels, i)
1306 {
1307 label facei = faceLabels[i];
1308
1309 if (internalFacesOnly == mesh.isInternalFace(facei))
1310 {
1311 collectedFaces[nFaces++] = facei;
1312 }
1313 }
1314 }
1315
1316 return collectedFaces;
1317}
1318
1319
1320// Calculate pointMap per patch (so from patch point label to old patch point
1321// label)
1322void Foam::polyTopoChange::calcPatchPointMap
1323(
1324 const UList<Map<label>>& oldPatchMeshPointMaps,
1325 const labelUList& patchMap,
1326 const polyBoundaryMesh& boundary,
1327 labelListList& patchPointMap
1328) const
1329{
1330 patchPointMap.setSize(patchMap.size());
1331
1332 forAll(patchMap, patchi)
1333 {
1334 const label oldPatchi = patchMap[patchi];
1335
1336 if (oldPatchi != -1)
1337 {
1338 const labelList& meshPoints = boundary[patchi].meshPoints();
1339
1340 const Map<label>& oldMeshPointMap =
1341 oldPatchMeshPointMaps[oldPatchi];
1342
1343 labelList& curPatchPointRnb = patchPointMap[patchi];
1344
1345 curPatchPointRnb.setSize(meshPoints.size());
1346
1347 forAll(meshPoints, i)
1348 {
1349 if (meshPoints[i] < pointMap_.size())
1350 {
1351 // Check if old point was part of same patch
1352 curPatchPointRnb[i] = oldMeshPointMap.lookup
1353 (
1354 pointMap_[meshPoints[i]],
1355 -1
1356 );
1357 }
1358 else
1359 {
1360 curPatchPointRnb[i] = -1;
1361 }
1362 }
1363 }
1364 }
1365}
1366
1367
1368void Foam::polyTopoChange::calcFaceInflationMaps
1369(
1370 const polyMesh& mesh,
1371 List<objectMap>& facesFromPoints,
1372 List<objectMap>& facesFromEdges,
1373 List<objectMap>& facesFromFaces
1374) const
1375{
1376 // Faces inflated from points
1377 // ~~~~~~~~~~~~~~~~~~~~~~~~~~
1378
1379 facesFromPoints.setSize(faceFromPoint_.size());
1380
1381 if (faceFromPoint_.size())
1382 {
1383 label nFacesFromPoints = 0;
1384
1385 // Collect all still existing faces connected to this point.
1386 forAllConstIters(faceFromPoint_, iter)
1387 {
1388 const label facei = iter.key();
1389 const label pointi = iter.val();
1390
1391 // Get internal or patch faces using point on old mesh
1392 const bool internal = (region_[facei] == -1);
1393
1394 facesFromPoints[nFacesFromPoints++] = objectMap
1395 (
1396 facei,
1397 selectFaces
1398 (
1399 mesh,
1400 mesh.pointFaces()[pointi],
1401 internal
1402 )
1403 );
1404 }
1405 }
1406
1407
1408 // Faces inflated from edges
1409 // ~~~~~~~~~~~~~~~~~~~~~~~~~
1410
1411 facesFromEdges.setSize(faceFromEdge_.size());
1412
1413 if (faceFromEdge_.size())
1414 {
1415 label nFacesFromEdges = 0;
1416
1417 // Collect all still existing faces connected to this edge.
1418 forAllConstIters(faceFromEdge_, iter)
1419 {
1420 const label facei = iter.key();
1421 const label edgei = iter.val();
1422
1423 // Get internal or patch faces using edge on old mesh
1424 const bool internal = (region_[facei] == -1);
1425
1426 facesFromEdges[nFacesFromEdges++] = objectMap
1427 (
1428 facei,
1429 selectFaces
1430 (
1431 mesh,
1432 mesh.edgeFaces(edgei),
1433 internal
1434 )
1435 );
1436 }
1437 }
1438
1439
1440 // Faces from face merging
1441 // ~~~~~~~~~~~~~~~~~~~~~~~
1442
1443 getMergeSets
1444 (
1445 reverseFaceMap_,
1446 faceMap_,
1447 facesFromFaces
1448 );
1449}
1450
1451
1452void Foam::polyTopoChange::calcCellInflationMaps
1453(
1454 const polyMesh& mesh,
1455 List<objectMap>& cellsFromPoints,
1456 List<objectMap>& cellsFromEdges,
1457 List<objectMap>& cellsFromFaces,
1458 List<objectMap>& cellsFromCells
1459) const
1460{
1461 cellsFromPoints.setSize(cellFromPoint_.size());
1462
1463 if (cellFromPoint_.size())
1464 {
1465 label nCellsFromPoints = 0;
1466
1467 // Collect all still existing cells connected to this point.
1468 forAllConstIters(cellFromPoint_, iter)
1469 {
1470 const label celli = iter.key();
1471 const label pointi = iter.val();
1472
1473 cellsFromPoints[nCellsFromPoints++] = objectMap
1474 (
1475 celli,
1476 mesh.pointCells()[pointi]
1477 );
1478 }
1479 }
1480
1481
1482 cellsFromEdges.setSize(cellFromEdge_.size());
1483
1484 if (cellFromEdge_.size())
1485 {
1486 label nCellsFromEdges = 0;
1487
1488 // Collect all still existing cells connected to this edge.
1489 forAllConstIters(cellFromEdge_, iter)
1490 {
1491 const label celli = iter.key();
1492 const label edgei = iter.val();
1493
1494 cellsFromEdges[nCellsFromEdges++] = objectMap
1495 (
1496 celli,
1497 mesh.edgeCells()[edgei]
1498 );
1499 }
1500 }
1501
1502
1503 cellsFromFaces.setSize(cellFromFace_.size());
1504
1505 if (cellFromFace_.size())
1506 {
1507 label nCellsFromFaces = 0;
1508
1509 labelList twoCells(2);
1510
1511 // Collect all still existing faces connected to this point.
1512 forAllConstIters(cellFromFace_, iter)
1513 {
1514 const label celli = iter.key();
1515 const label oldFacei = iter.val();
1516
1517 if (mesh.isInternalFace(oldFacei))
1518 {
1519 twoCells[0] = mesh.faceOwner()[oldFacei];
1520 twoCells[1] = mesh.faceNeighbour()[oldFacei];
1521 cellsFromFaces[nCellsFromFaces++] = objectMap
1522 (
1523 celli,
1524 twoCells
1525 );
1526 }
1527 else
1528 {
1529 cellsFromFaces[nCellsFromFaces++] = objectMap
1530 (
1531 celli,
1532 labelList(1, mesh.faceOwner()[oldFacei])
1533 );
1534 }
1535 }
1536 }
1537
1538
1539 // Cells from cell merging
1540 // ~~~~~~~~~~~~~~~~~~~~~~~
1541
1542 getMergeSets
1543 (
1544 reverseCellMap_,
1545 cellMap_,
1546 cellsFromCells
1547 );
1548}
1549
1550
1551void Foam::polyTopoChange::resetZones
1552(
1553 const polyMesh& mesh,
1554 polyMesh& newMesh,
1555 labelListList& pointZoneMap,
1556 labelListList& faceZoneFaceMap,
1557 labelListList& cellZoneMap
1558) const
1559{
1560 // pointZones
1561 // ~~~~~~~~~~
1562
1563 pointZoneMap.setSize(mesh.pointZones().size());
1564 {
1565 const pointZoneMesh& pointZones = mesh.pointZones();
1566
1567 // Count points per zone
1568
1569 labelList nPoints(pointZones.size(), Zero);
1570
1571 forAllConstIters(pointZone_, iter)
1572 {
1573 const label pointi = iter.key();
1574 const label zonei = iter.val();
1575
1576 if (zonei < 0 || zonei >= pointZones.size())
1577 {
1579 << "Illegal zoneID " << zonei << " for point "
1580 << pointi << " coord " << mesh.points()[pointi]
1581 << abort(FatalError);
1582 }
1583 nPoints[zonei]++;
1584 }
1585
1586 // Distribute points per zone
1587
1588 labelListList addressing(pointZones.size());
1589 forAll(addressing, zonei)
1590 {
1591 addressing[zonei].setSize(nPoints[zonei]);
1592 }
1593 nPoints = 0;
1594
1595 forAllConstIters(pointZone_, iter)
1596 {
1597 const label pointi = iter.key();
1598 const label zonei = iter.val();
1599
1600 addressing[zonei][nPoints[zonei]++] = pointi;
1601 }
1602 // Sort the addressing
1603 forAll(addressing, zonei)
1604 {
1605 stableSort(addressing[zonei]);
1606 }
1607
1608 // So now we both have old zones and the new addressing.
1609 // Invert the addressing to get pointZoneMap.
1610 forAll(addressing, zonei)
1611 {
1612 const pointZone& oldZone = pointZones[zonei];
1613 const labelList& newZoneAddr = addressing[zonei];
1614
1615 labelList& curPzRnb = pointZoneMap[zonei];
1616 curPzRnb.setSize(newZoneAddr.size());
1617
1618 forAll(newZoneAddr, i)
1619 {
1620 if (newZoneAddr[i] < pointMap_.size())
1621 {
1622 curPzRnb[i] = oldZone.whichPoint(pointMap_[newZoneAddr[i]]);
1623 }
1624 else
1625 {
1626 curPzRnb[i] = -1;
1627 }
1628 }
1629 }
1630
1631 // Reset the addressing on the zone
1632 newMesh.pointZones().clearAddressing();
1633 forAll(newMesh.pointZones(), zonei)
1634 {
1635 if (debug)
1636 {
1637 Pout<< "pointZone:" << zonei
1638 << " name:" << newMesh.pointZones()[zonei].name()
1639 << " size:" << addressing[zonei].size()
1640 << endl;
1641 }
1642
1643 newMesh.pointZones()[zonei] = addressing[zonei];
1644 }
1645 }
1646
1647
1648 // faceZones
1649 // ~~~~~~~~~
1650
1651 faceZoneFaceMap.setSize(mesh.faceZones().size());
1652 {
1653 const faceZoneMesh& faceZones = mesh.faceZones();
1654
1655 labelList nFaces(faceZones.size(), Zero);
1656
1657 forAllConstIters(faceZone_, iter)
1658 {
1659 const label facei = iter.key();
1660 const label zonei = iter.val();
1661
1662 if (zonei < 0 || zonei >= faceZones.size())
1663 {
1665 << "Illegal zoneID " << zonei << " for face "
1666 << facei
1667 << abort(FatalError);
1668 }
1669 nFaces[zonei]++;
1670 }
1671
1672 labelListList addressing(faceZones.size());
1673 boolListList flipMode(faceZones.size());
1674
1675 forAll(addressing, zonei)
1676 {
1677 addressing[zonei].setSize(nFaces[zonei]);
1678 flipMode[zonei].setSize(nFaces[zonei]);
1679 }
1680 nFaces = 0;
1681
1682 forAllConstIters(faceZone_, iter)
1683 {
1684 const label facei = iter.key();
1685 const label zonei = iter.val();
1686
1687 const label index = nFaces[zonei]++;
1688
1689 addressing[zonei][index] = facei;
1690 flipMode[zonei][index] = faceZoneFlip_.test(facei);
1691 }
1692
1693 // Sort the addressing
1694 forAll(addressing, zonei)
1695 {
1696 labelList newToOld(sortedOrder(addressing[zonei]));
1697 {
1698 labelList newAddressing(addressing[zonei].size());
1699 forAll(newAddressing, i)
1700 {
1701 newAddressing[i] = addressing[zonei][newToOld[i]];
1702 }
1703 addressing[zonei].transfer(newAddressing);
1704 }
1705 {
1706 boolList newFlipMode(flipMode[zonei].size());
1707 forAll(newFlipMode, i)
1708 {
1709 newFlipMode[i] = flipMode[zonei][newToOld[i]];
1710 }
1711 flipMode[zonei].transfer(newFlipMode);
1712 }
1713 }
1714
1715 // So now we both have old zones and the new addressing.
1716 // Invert the addressing to get faceZoneFaceMap.
1717 forAll(addressing, zonei)
1718 {
1719 const faceZone& oldZone = faceZones[zonei];
1720 const labelList& newZoneAddr = addressing[zonei];
1721
1722 labelList& curFzFaceRnb = faceZoneFaceMap[zonei];
1723
1724 curFzFaceRnb.setSize(newZoneAddr.size());
1725
1726 forAll(newZoneAddr, i)
1727 {
1728 if (newZoneAddr[i] < faceMap_.size())
1729 {
1730 curFzFaceRnb[i] =
1731 oldZone.whichFace(faceMap_[newZoneAddr[i]]);
1732 }
1733 else
1734 {
1735 curFzFaceRnb[i] = -1;
1736 }
1737 }
1738 }
1739
1740
1741 // Reset the addressing on the zone
1742 newMesh.faceZones().clearAddressing();
1743 forAll(newMesh.faceZones(), zonei)
1744 {
1745 if (debug)
1746 {
1747 Pout<< "faceZone:" << zonei
1748 << " name:" << newMesh.faceZones()[zonei].name()
1749 << " size:" << addressing[zonei].size()
1750 << endl;
1751 }
1752
1753 newMesh.faceZones()[zonei].resetAddressing
1754 (
1755 addressing[zonei],
1756 flipMode[zonei]
1757 );
1758 }
1759 }
1760
1761
1762 // cellZones
1763 // ~~~~~~~~~
1764
1765 cellZoneMap.setSize(mesh.cellZones().size());
1766 {
1767 const cellZoneMesh& cellZones = mesh.cellZones();
1768
1769 labelList nCells(cellZones.size(), Zero);
1770
1771 forAll(cellZone_, celli)
1772 {
1773 const label zonei = cellZone_[celli];
1774
1775 if (zonei >= cellZones.size())
1776 {
1778 << "Illegal zoneID " << zonei << " for cell "
1779 << celli << abort(FatalError);
1780 }
1781
1782 if (zonei >= 0)
1783 {
1784 nCells[zonei]++;
1785 }
1786 }
1787
1788 labelListList addressing(cellZones.size());
1789 forAll(addressing, zonei)
1790 {
1791 addressing[zonei].setSize(nCells[zonei]);
1792 }
1793 nCells = 0;
1794
1795 forAll(cellZone_, celli)
1796 {
1797 const label zonei = cellZone_[celli];
1798
1799 if (zonei >= 0)
1800 {
1801 addressing[zonei][nCells[zonei]++] = celli;
1802 }
1803 }
1804 // Sort the addressing
1805 forAll(addressing, zonei)
1806 {
1807 stableSort(addressing[zonei]);
1808 }
1809
1810 // So now we both have old zones and the new addressing.
1811 // Invert the addressing to get cellZoneMap.
1812 forAll(addressing, zonei)
1813 {
1814 const cellZone& oldZone = cellZones[zonei];
1815 const labelList& newZoneAddr = addressing[zonei];
1816
1817 labelList& curCellRnb = cellZoneMap[zonei];
1818
1819 curCellRnb.setSize(newZoneAddr.size());
1820
1821 forAll(newZoneAddr, i)
1822 {
1823 if (newZoneAddr[i] < cellMap_.size())
1824 {
1825 curCellRnb[i] =
1826 oldZone.whichCell(cellMap_[newZoneAddr[i]]);
1827 }
1828 else
1829 {
1830 curCellRnb[i] = -1;
1831 }
1832 }
1833 }
1834
1835 // Reset the addressing on the zone
1836 newMesh.cellZones().clearAddressing();
1837 forAll(newMesh.cellZones(), zonei)
1838 {
1839 if (debug)
1840 {
1841 Pout<< "cellZone:" << zonei
1842 << " name:" << newMesh.cellZones()[zonei].name()
1843 << " size:" << addressing[zonei].size()
1844 << endl;
1845 }
1846
1847 newMesh.cellZones()[zonei] = addressing[zonei];
1848 }
1849 }
1850}
1851
1852
1853void Foam::polyTopoChange::calcFaceZonePointMap
1854(
1855 const polyMesh& mesh,
1856 const UList<Map<label>>& oldFaceZoneMeshPointMaps,
1857 labelListList& faceZonePointMap
1858) const
1859{
1860 const faceZoneMesh& faceZones = mesh.faceZones();
1861
1862 faceZonePointMap.setSize(faceZones.size());
1863
1864 forAll(faceZones, zonei)
1865 {
1866 const faceZone& newZone = faceZones[zonei];
1867
1868 const labelList& newZoneMeshPoints = newZone().meshPoints();
1869
1870 const Map<label>& oldZoneMeshPointMap = oldFaceZoneMeshPointMaps[zonei];
1871
1872 labelList& curFzPointRnb = faceZonePointMap[zonei];
1873
1874 curFzPointRnb.setSize(newZoneMeshPoints.size());
1875
1876 forAll(newZoneMeshPoints, pointi)
1877 {
1878 if (newZoneMeshPoints[pointi] < pointMap_.size())
1879 {
1880 curFzPointRnb[pointi] =
1881 oldZoneMeshPointMap.lookup
1882 (
1883 pointMap_[newZoneMeshPoints[pointi]],
1884 -1
1885 );
1886 }
1887 else
1888 {
1889 curFzPointRnb[pointi] = -1;
1890 }
1891 }
1892 }
1893}
1894
1895
1896void Foam::polyTopoChange::reorderCoupledFaces
1897(
1898 const bool syncParallel,
1899 const polyBoundaryMesh& oldBoundary,
1900 const labelUList& patchMap, // new to old patches
1901 const labelUList& patchStarts,
1902 const labelUList& patchSizes,
1903 const pointField& points
1904)
1905{
1906 // Mapping for faces (old to new). Extends over all mesh faces for
1907 // convenience (could be just the external faces)
1908 labelList oldToNew(identity(faces_.size()));
1909
1910 // Rotation on new faces.
1911 labelList rotation(faces_.size(), Zero);
1912
1913 PstreamBuffers pBufs(Pstream::commsTypes::nonBlocking);
1914
1915 // Send ordering
1916 forAll(patchMap, patchi)
1917 {
1918 const label oldPatchi = patchMap[patchi];
1919
1920 if
1921 (
1922 oldPatchi != -1
1923 && (syncParallel || !isA<processorPolyPatch>(oldBoundary[oldPatchi]))
1924 )
1925 {
1926 oldBoundary[oldPatchi].initOrder
1927 (
1928 pBufs,
1930 (
1931 SubList<face>
1932 (
1933 faces_,
1934 patchSizes[patchi],
1935 patchStarts[patchi]
1936 ),
1937 points
1938 )
1939 );
1940 }
1941 }
1942
1943 if (syncParallel)
1944 {
1945 pBufs.finishedSends();
1946 }
1947
1948 // Receive and calculate ordering
1949
1950 bool anyChanged = false;
1951
1952 forAll(patchMap, patchi)
1953 {
1954 const label oldPatchi = patchMap[patchi];
1955
1956 if
1957 (
1958 oldPatchi != -1
1959 && (syncParallel || !isA<processorPolyPatch>(oldBoundary[oldPatchi]))
1960 )
1961 {
1962 labelList patchFaceMap(patchSizes[patchi], -1);
1963 labelList patchFaceRotation(patchSizes[patchi], Zero);
1964
1965 const bool changed = oldBoundary[oldPatchi].order
1966 (
1967 pBufs,
1969 (
1970 SubList<face>
1971 (
1972 faces_,
1973 patchSizes[patchi],
1974 patchStarts[patchi]
1975 ),
1976 points
1977 ),
1978 patchFaceMap,
1979 patchFaceRotation
1980 );
1981
1982 if (changed)
1983 {
1984 // Merge patch face reordering into mesh face reordering table
1985 const label start = patchStarts[patchi];
1986
1987 forAll(patchFaceMap, patchFacei)
1988 {
1989 oldToNew[patchFacei + start] =
1990 start + patchFaceMap[patchFacei];
1991 }
1992
1993 forAll(patchFaceRotation, patchFacei)
1994 {
1995 rotation[patchFacei + start] =
1996 patchFaceRotation[patchFacei];
1997 }
1998
1999 anyChanged = true;
2000 }
2001 }
2002 }
2003
2004 if (syncParallel)
2005 {
2006 reduce(anyChanged, orOp<bool>());
2007 }
2008
2009 if (anyChanged)
2010 {
2011 // Reorder faces according to oldToNew.
2012 reorderCompactFaces(oldToNew.size(), oldToNew);
2013
2014 // Rotate faces (rotation is already in new face indices).
2015 forAll(rotation, facei)
2016 {
2017 if (rotation[facei] != 0)
2018 {
2019 inplaceRotateList<List, label>(faces_[facei], rotation[facei]);
2020 }
2021 }
2022 }
2023}
2024
2025
2026void Foam::polyTopoChange::compactAndReorder
2027(
2028 const polyMesh& mesh,
2029 const labelUList& patchMap, // from new to old patch
2030 const bool syncParallel,
2031 const bool orderCells,
2032 const bool orderPoints,
2033
2034 label& nInternalPoints,
2035 pointField& newPoints,
2036 labelList& patchSizes,
2037 labelList& patchStarts,
2038 List<objectMap>& pointsFromPoints,
2039 List<objectMap>& facesFromPoints,
2040 List<objectMap>& facesFromEdges,
2041 List<objectMap>& facesFromFaces,
2042 List<objectMap>& cellsFromPoints,
2043 List<objectMap>& cellsFromEdges,
2044 List<objectMap>& cellsFromFaces,
2045 List<objectMap>& cellsFromCells,
2046 List<Map<label>>& oldPatchMeshPointMaps,
2047 labelList& oldPatchNMeshPoints,
2048 labelList& oldPatchStarts,
2049 List<Map<label>>& oldFaceZoneMeshPointMaps
2050)
2051{
2052 if (patchMap.size() != nPatches_)
2053 {
2055 << "polyTopoChange was constructed with a mesh with "
2056 << nPatches_ << " patches." << endl
2057 << "The mesh now provided has a different number of patches "
2058 << patchMap.size()
2059 << " which is illegal" << endl
2060 << abort(FatalError);
2061 }
2062
2063 // Remove any holes from points/faces/cells and sort faces.
2064 // Sets nActiveFaces_.
2065 compact(orderCells, orderPoints, nInternalPoints, patchSizes, patchStarts);
2066
2067 // Transfer points to pointField. points_ are now cleared!
2068 // Only done since e.g. reorderCoupledFaces requires pointField.
2069 newPoints.transfer(points_);
2070
2071 // Reorder any coupled faces
2072 reorderCoupledFaces
2073 (
2074 syncParallel,
2076 patchMap,
2077 patchStarts,
2078 patchSizes,
2079 newPoints
2080 );
2081
2082
2083 // Calculate inflation/merging maps
2084 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2085 // These are for the new face(/point/cell) the old faces whose value
2086 // needs to be
2087 // averaged/summed to get the new value. These old faces come either from
2088 // merged old faces (face remove into other face),
2089 // the old edgeFaces (inflate from edge) or the old pointFaces (inflate
2090 // from point). As an additional complexity will use only internal faces
2091 // to create new value for internal face and vice versa only patch
2092 // faces to to create patch face value.
2093
2094 // For point only point merging
2095 getMergeSets
2096 (
2097 reversePointMap_,
2098 pointMap_,
2099 pointsFromPoints
2100 );
2101
2102 calcFaceInflationMaps
2103 (
2104 mesh,
2105 facesFromPoints,
2106 facesFromEdges,
2107 facesFromFaces
2108 );
2109
2110 calcCellInflationMaps
2111 (
2112 mesh,
2113 cellsFromPoints,
2114 cellsFromEdges,
2115 cellsFromFaces,
2116 cellsFromCells
2117 );
2118
2119 // Clear inflation info
2120 {
2121 faceFromPoint_.clearStorage();
2122 faceFromEdge_.clearStorage();
2123
2124 cellFromPoint_.clearStorage();
2125 cellFromEdge_.clearStorage();
2126 cellFromFace_.clearStorage();
2127 }
2128
2129
2130 const polyBoundaryMesh& boundary = mesh.boundaryMesh();
2131
2132 // Grab patch mesh point maps
2133 oldPatchMeshPointMaps.setSize(boundary.size());
2134 oldPatchNMeshPoints.setSize(boundary.size());
2135 oldPatchStarts.setSize(boundary.size());
2136
2137 forAll(boundary, patchi)
2138 {
2139 // Copy old face zone mesh point maps
2140 oldPatchMeshPointMaps[patchi] = boundary[patchi].meshPointMap();
2141 oldPatchNMeshPoints[patchi] = boundary[patchi].meshPoints().size();
2142 oldPatchStarts[patchi] = boundary[patchi].start();
2143 }
2144
2145 // Grab old face zone mesh point maps.
2146 // These need to be saved before resetting the mesh and are used
2147 // later on to calculate the faceZone pointMaps.
2148 oldFaceZoneMeshPointMaps.setSize(mesh.faceZones().size());
2149
2150 forAll(mesh.faceZones(), zonei)
2151 {
2152 const faceZone& oldZone = mesh.faceZones()[zonei];
2153
2154 oldFaceZoneMeshPointMaps[zonei] = oldZone().meshPointMap();
2155 }
2156}
2157
2158
2159// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
2160
2161// Construct from components
2162Foam::polyTopoChange::polyTopoChange(const label nPatches, const bool strict)
2163:
2164 strict_(strict),
2165 nPatches_(nPatches),
2166 points_(0),
2167 pointMap_(0),
2168 reversePointMap_(0),
2169 pointZone_(0),
2170 retiredPoints_(0),
2171 faces_(0),
2172 region_(0),
2173 faceOwner_(0),
2174 faceNeighbour_(0),
2175 faceMap_(0),
2176 reverseFaceMap_(0),
2177 faceFromPoint_(0),
2178 faceFromEdge_(0),
2179 flipFaceFlux_(0),
2180 faceZone_(0),
2181 faceZoneFlip_(0),
2182 nActiveFaces_(0),
2183 cellMap_(0),
2184 reverseCellMap_(0),
2185 cellFromPoint_(0),
2186 cellFromEdge_(0),
2187 cellFromFace_(0),
2188 cellZone_(0)
2189{}
2190
2191
2192// Construct from components
2194(
2195 const polyMesh& mesh,
2196 const bool strict
2197)
2198:
2199 strict_(strict),
2200 nPatches_(0),
2201 points_(0),
2202 pointMap_(0),
2203 reversePointMap_(0),
2204 pointZone_(0),
2205 retiredPoints_(0),
2206 faces_(0),
2207 region_(0),
2208 faceOwner_(0),
2209 faceNeighbour_(0),
2210 faceMap_(0),
2211 reverseFaceMap_(0),
2212 faceFromPoint_(0),
2213 faceFromEdge_(0),
2214 flipFaceFlux_(0),
2215 faceZone_(0),
2216 faceZoneFlip_(0),
2217 nActiveFaces_(0),
2218 cellMap_(0),
2219 reverseCellMap_(0),
2220 cellFromPoint_(0),
2221 cellFromEdge_(0),
2222 cellFromFace_(0),
2223 cellZone_(0)
2224{
2225 addMesh
2226 (
2227 mesh,
2232 );
2233}
2234
2235
2236// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
2237
2239{
2240 points_.clearStorage();
2241 pointMap_.clearStorage();
2242 reversePointMap_.clearStorage();
2243 pointZone_.clearStorage();
2244 retiredPoints_.clearStorage();
2245
2246 faces_.clearStorage();
2247 region_.clearStorage();
2248 faceOwner_.clearStorage();
2249 faceNeighbour_.clearStorage();
2250 faceMap_.clearStorage();
2251 reverseFaceMap_.clearStorage();
2252 faceFromPoint_.clearStorage();
2253 faceFromEdge_.clearStorage();
2254 flipFaceFlux_.clearStorage();
2255 faceZone_.clearStorage();
2256 faceZoneFlip_.clearStorage();
2257 nActiveFaces_ = 0;
2258
2259 cellMap_.clearStorage();
2260 reverseCellMap_.clearStorage();
2261 cellZone_.clearStorage();
2262 cellFromPoint_.clearStorage();
2263 cellFromEdge_.clearStorage();
2264 cellFromFace_.clearStorage();
2265}
2266
2267
2269(
2270 const polyMesh& mesh,
2271 const labelUList& patchMap,
2272 const labelUList& pointZoneMap,
2273 const labelUList& faceZoneMap,
2274 const labelUList& cellZoneMap
2275)
2276{
2277 label maxRegion = nPatches_ - 1;
2278 for (const label regioni : patchMap)
2279 {
2280 maxRegion = max(maxRegion, regioni);
2281 }
2282 nPatches_ = maxRegion + 1;
2283
2284
2285 // Add points
2286 {
2287 const pointField& points = mesh.points();
2288 const pointZoneMesh& pointZones = mesh.pointZones();
2289
2290 // Extend
2291 points_.setCapacity(points_.size() + points.size());
2292 pointMap_.setCapacity(pointMap_.size() + points.size());
2293 reversePointMap_.setCapacity(reversePointMap_.size() + points.size());
2294 pointZone_.resize(pointZone_.size() + points.size()/128);
2295
2296 // Precalc offset zones
2297 labelList newZoneID(points.size(), -1);
2298
2299 forAll(pointZones, zonei)
2300 {
2301 const labelList& pointLabels = pointZones[zonei];
2302
2303 for (const label pointi : pointLabels)
2304 {
2305 newZoneID[pointi] = pointZoneMap[zonei];
2306 }
2307 }
2308
2309 // Add points in mesh order
2310 for (label pointi = 0; pointi < mesh.nPoints(); pointi++)
2311 {
2312 addPoint
2313 (
2314 points[pointi],
2315 pointi,
2316 newZoneID[pointi],
2317 true
2318 );
2319 }
2320 }
2321
2322 // Add cells
2323 {
2324 const cellZoneMesh& cellZones = mesh.cellZones();
2325
2326 // Resize
2327
2328 // Note: polyMesh does not allow retired cells anymore. So allCells
2329 // always equals nCells
2330 label nAllCells = mesh.nCells();
2331
2332 cellMap_.setCapacity(cellMap_.size() + nAllCells);
2333 reverseCellMap_.setCapacity(reverseCellMap_.size() + nAllCells);
2334 cellFromPoint_.resize(cellFromPoint_.size() + nAllCells/128);
2335 cellFromEdge_.resize(cellFromEdge_.size() + nAllCells/128);
2336 cellFromFace_.resize(cellFromFace_.size() + nAllCells/128);
2337 cellZone_.setCapacity(cellZone_.size() + nAllCells);
2338
2339
2340 // Precalc offset zones
2341 labelList newZoneID(nAllCells, -1);
2342
2343 forAll(cellZones, zonei)
2344 {
2345 const labelList& cellLabels = cellZones[zonei];
2346
2347 for (const label celli : cellLabels)
2348 {
2349 if (newZoneID[celli] != -1)
2350 {
2352 << "Cell:" << celli
2353 << " centre:" << mesh.cellCentres()[celli]
2354 << " is in two zones:"
2355 << cellZones[newZoneID[celli]].name()
2356 << " and " << cellZones[zonei].name() << endl
2357 << " This is not supported."
2358 << " Continuing with first zone only." << endl;
2359 }
2360 else
2361 {
2362 newZoneID[celli] = cellZoneMap[zonei];
2363 }
2364 }
2365 }
2366
2367 // Add cells in mesh order
2368 for (label celli = 0; celli < nAllCells; celli++)
2369 {
2370 // Add cell from cell
2371 addCell(-1, -1, -1, celli, newZoneID[celli]);
2372 }
2373 }
2374
2375 // Add faces
2376 {
2378 const faceList& faces = mesh.faces();
2379 const labelList& faceOwner = mesh.faceOwner();
2380 const labelList& faceNeighbour = mesh.faceNeighbour();
2381 const faceZoneMesh& faceZones = mesh.faceZones();
2382
2383 // Resize
2384 label nAllFaces = mesh.faces().size();
2385
2386 faces_.setCapacity(faces_.size() + nAllFaces);
2387 region_.setCapacity(region_.size() + nAllFaces);
2388 faceOwner_.setCapacity(faceOwner_.size() + nAllFaces);
2389 faceNeighbour_.setCapacity(faceNeighbour_.size() + nAllFaces);
2390 faceMap_.setCapacity(faceMap_.size() + nAllFaces);
2391 reverseFaceMap_.setCapacity(reverseFaceMap_.size() + nAllFaces);
2392 faceFromPoint_.resize(faceFromPoint_.size() + nAllFaces/128);
2393 faceFromEdge_.resize(faceFromEdge_.size() + nAllFaces/128);
2394 flipFaceFlux_.setCapacity(faces_.size() + nAllFaces);
2395 faceZone_.resize(faceZone_.size() + nAllFaces/128);
2396 faceZoneFlip_.setCapacity(faces_.size() + nAllFaces);
2397
2398
2399 // Precalc offset zones
2400 labelList newZoneID(nAllFaces, -1);
2401 boolList zoneFlip(nAllFaces, false);
2402
2403 forAll(faceZones, zonei)
2404 {
2405 const labelList& faceLabels = faceZones[zonei];
2406 const boolList& flipMap = faceZones[zonei].flipMap();
2407
2408 forAll(faceLabels, facei)
2409 {
2410 newZoneID[faceLabels[facei]] = faceZoneMap[zonei];
2411 zoneFlip[faceLabels[facei]] = flipMap[facei];
2412 }
2413 }
2414
2415 // Add faces in mesh order
2416
2417 // 1. Internal faces
2418 for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
2419 {
2420 addFace
2421 (
2422 faces[facei],
2423 faceOwner[facei],
2424 faceNeighbour[facei],
2425 -1, // masterPointID
2426 -1, // masterEdgeID
2427 facei, // masterFaceID
2428 false, // flipFaceFlux
2429 -1, // patchID
2430 newZoneID[facei], // zoneID
2431 zoneFlip[facei] // zoneFlip
2432 );
2433 }
2434
2435 // 2. Patch faces
2436 forAll(patches, patchi)
2437 {
2438 const polyPatch& pp = patches[patchi];
2439
2440 if (pp.start() != faces_.size())
2441 {
2443 << "Problem : "
2444 << "Patch " << pp.name() << " starts at " << pp.start()
2445 << endl
2446 << "Current face counter at " << faces_.size() << endl
2447 << "Are patches in incremental order?"
2448 << abort(FatalError);
2449 }
2450 forAll(pp, patchFacei)
2451 {
2452 const label facei = pp.start() + patchFacei;
2453
2454 addFace
2455 (
2456 faces[facei],
2457 faceOwner[facei],
2458 -1, // neighbour
2459 -1, // masterPointID
2460 -1, // masterEdgeID
2461 facei, // masterFaceID
2462 false, // flipFaceFlux
2463 patchMap[patchi], // patchID
2464 newZoneID[facei], // zoneID
2465 zoneFlip[facei] // zoneFlip
2466 );
2467 }
2468 }
2469 }
2470}
2471
2472
2474(
2475 const label nPoints,
2476 const label nFaces,
2477 const label nCells
2478)
2479{
2480 points_.setCapacity(nPoints);
2481 pointMap_.setCapacity(nPoints);
2482 reversePointMap_.setCapacity(nPoints);
2483 pointZone_.resize(pointZone_.size() + nPoints/128);
2484
2485 faces_.setCapacity(nFaces);
2486 region_.setCapacity(nFaces);
2487 faceOwner_.setCapacity(nFaces);
2488 faceNeighbour_.setCapacity(nFaces);
2489 faceMap_.setCapacity(nFaces);
2490 reverseFaceMap_.setCapacity(nFaces);
2491 faceFromPoint_.resize(faceFromPoint_.size() + nFaces/128);
2492 faceFromEdge_.resize(faceFromEdge_.size() + nFaces/128);
2493 flipFaceFlux_.setCapacity(nFaces);
2494 faceZone_.resize(faceZone_.size() + nFaces/128);
2495 faceZoneFlip_.setCapacity(nFaces);
2496
2497 cellMap_.setCapacity(nCells);
2498 reverseCellMap_.setCapacity(nCells);
2499 cellFromPoint_.resize(cellFromPoint_.size() + nCells/128);
2500 cellFromEdge_.resize(cellFromEdge_.size() + nCells/128);
2501 cellFromFace_.resize(cellFromFace_.size() + nCells/128);
2502 cellZone_.setCapacity(nCells);
2503}
2504
2505
2507{
2508 if (isType<polyAddPoint>(action))
2509 {
2510 const polyAddPoint& pap = refCast<const polyAddPoint>(action);
2511
2512 return addPoint
2513 (
2514 pap.newPoint(),
2515 pap.masterPointID(),
2516 pap.zoneID(),
2517 pap.inCell()
2518 );
2519 }
2520 else if (isType<polyModifyPoint>(action))
2521 {
2522 const polyModifyPoint& pmp = refCast<const polyModifyPoint>(action);
2523
2524 modifyPoint
2525 (
2526 pmp.pointID(),
2527 pmp.newPoint(),
2528 pmp.zoneID(),
2529 pmp.inCell()
2530 );
2531
2532 return -1;
2533 }
2534 else if (isType<polyRemovePoint>(action))
2535 {
2536 const polyRemovePoint& prp = refCast<const polyRemovePoint>(action);
2537
2538 removePoint(prp.pointID(), prp.mergePointID());
2539
2540 return -1;
2541 }
2542 else if (isType<polyAddFace>(action))
2543 {
2544 const polyAddFace& paf = refCast<const polyAddFace>(action);
2545
2546 return addFace
2547 (
2548 paf.newFace(),
2549 paf.owner(),
2550 paf.neighbour(),
2551 paf.masterPointID(),
2552 paf.masterEdgeID(),
2553 paf.masterFaceID(),
2554 paf.flipFaceFlux(),
2555 paf.patchID(),
2556 paf.zoneID(),
2557 paf.zoneFlip()
2558 );
2559 }
2560 else if (isType<polyModifyFace>(action))
2561 {
2562 const polyModifyFace& pmf = refCast<const polyModifyFace>(action);
2563
2564 modifyFace
2565 (
2566 pmf.newFace(),
2567 pmf.faceID(),
2568 pmf.owner(),
2569 pmf.neighbour(),
2570 pmf.flipFaceFlux(),
2571 pmf.patchID(),
2572 pmf.zoneID(),
2573 pmf.zoneFlip()
2574 );
2575
2576 return -1;
2577 }
2578 else if (isType<polyRemoveFace>(action))
2579 {
2580 const polyRemoveFace& prf = refCast<const polyRemoveFace>(action);
2581
2582 removeFace(prf.faceID(), prf.mergeFaceID());
2583
2584 return -1;
2585 }
2586 else if (isType<polyAddCell>(action))
2587 {
2588 const polyAddCell& pac = refCast<const polyAddCell>(action);
2589
2590 return addCell
2591 (
2592 pac.masterPointID(),
2593 pac.masterEdgeID(),
2594 pac.masterFaceID(),
2595 pac.masterCellID(),
2596 pac.zoneID()
2597 );
2598 }
2599 else if (isType<polyModifyCell>(action))
2600 {
2601 const polyModifyCell& pmc = refCast<const polyModifyCell>(action);
2602
2603 if (pmc.removeFromZone())
2604 {
2605 modifyCell(pmc.cellID(), -1);
2606 }
2607 else
2608 {
2609 modifyCell(pmc.cellID(), pmc.zoneID());
2610 }
2611
2612 return -1;
2613 }
2614 else if (isType<polyRemoveCell>(action))
2615 {
2616 const polyRemoveCell& prc = refCast<const polyRemoveCell>(action);
2617
2618 removeCell(prc.cellID(), prc.mergeCellID());
2619
2620 return -1;
2621 }
2622 else
2623 {
2625 << "Unknown type of topoChange: " << action.type()
2626 << abort(FatalError);
2627
2628 // Dummy return to keep compiler happy
2629 return -1;
2630 }
2631}
2632
2633
2635(
2636 const point& pt,
2637 const label masterPointID,
2638 const label zoneID,
2639 const bool inCell
2640)
2641{
2642 const label pointi = points_.size();
2643
2644 points_.append(pt);
2645 pointMap_.append(masterPointID);
2646 reversePointMap_.append(pointi);
2647
2648 if (zoneID >= 0)
2649 {
2650 pointZone_.insert(pointi, zoneID);
2651 }
2652
2653 if (!inCell)
2654 {
2655 retiredPoints_.insert(pointi);
2656 }
2657
2658 return pointi;
2659}
2660
2661
2663(
2664 const label pointi,
2665 const point& pt,
2666 const label zoneID,
2667 const bool inCell
2668)
2669{
2670 if (pointi < 0 || pointi >= points_.size())
2671 {
2673 << "illegal point label " << pointi << endl
2674 << "Valid point labels are 0 .. " << points_.size()-1
2675 << abort(FatalError);
2676 }
2677 if (pointRemoved(pointi) || pointMap_[pointi] == -1)
2678 {
2680 << "point " << pointi << " already marked for removal"
2681 << abort(FatalError);
2682 }
2683 points_[pointi] = pt;
2684
2685 if (zoneID >= 0)
2686 {
2687 pointZone_.set(pointi, zoneID);
2688 }
2689 else
2690 {
2691 pointZone_.erase(pointi);
2692 }
2693
2694 if (inCell)
2695 {
2696 retiredPoints_.erase(pointi);
2697 }
2698 else
2699 {
2700 retiredPoints_.insert(pointi);
2701 }
2702}
2703
2704
2706{
2707 if (newPoints.size() != points_.size())
2708 {
2710 << "illegal pointField size." << endl
2711 << "Size:" << newPoints.size() << endl
2712 << "Points in mesh:" << points_.size()
2713 << abort(FatalError);
2714 }
2715
2716 forAll(points_, pointi)
2717 {
2718 points_[pointi] = newPoints[pointi];
2719 }
2720}
2721
2722
2724(
2725 const label pointi,
2726 const label mergePointi
2727)
2728{
2729 if (pointi < 0 || pointi >= points_.size())
2730 {
2732 << "illegal point label " << pointi << endl
2733 << "Valid point labels are 0 .. " << points_.size()-1
2734 << abort(FatalError);
2735 }
2736
2737 if
2738 (
2739 strict_
2740 && (pointRemoved(pointi) || pointMap_[pointi] == -1)
2741 )
2742 {
2744 << "point " << pointi << " already marked for removal" << nl
2745 << "Point:" << points_[pointi] << " pointMap:" << pointMap_[pointi]
2746 << abort(FatalError);
2747 }
2748
2749 if (pointi == mergePointi)
2750 {
2752 << "Cannot remove/merge point " << pointi << " onto itself."
2753 << abort(FatalError);
2754 }
2755
2756 points_[pointi] = point::max;
2757 pointMap_[pointi] = -1;
2758 if (mergePointi >= 0)
2759 {
2760 reversePointMap_[pointi] = -mergePointi-2;
2761 }
2762 else
2763 {
2764 reversePointMap_[pointi] = -1;
2765 }
2766 pointZone_.erase(pointi);
2767 retiredPoints_.erase(pointi);
2768}
2769
2770
2772(
2773 const face& f,
2774 const label own,
2775 const label nei,
2776 const label masterPointID,
2777 const label masterEdgeID,
2778 const label masterFaceID,
2779 const bool flipFaceFlux,
2780 const label patchID,
2781 const label zoneID,
2782 const bool zoneFlip
2783)
2784{
2785 // Check validity
2786 if (debug)
2787 {
2788 checkFace(f, -1, own, nei, patchID, zoneID);
2789 }
2790
2791 label facei = faces_.size();
2792
2793 faces_.append(f);
2794 region_.append(patchID);
2795 faceOwner_.append(own);
2796 faceNeighbour_.append(nei);
2797
2798 if (masterPointID >= 0)
2799 {
2800 faceMap_.append(-1);
2801 faceFromPoint_.insert(facei, masterPointID);
2802 }
2803 else if (masterEdgeID >= 0)
2804 {
2805 faceMap_.append(-1);
2806 faceFromEdge_.insert(facei, masterEdgeID);
2807 }
2808 else if (masterFaceID >= 0)
2809 {
2810 faceMap_.append(masterFaceID);
2811 }
2812 else
2813 {
2814 // Allow inflate-from-nothing?
2815 //FatalErrorInFunction
2816 // << "Need to specify a master point, edge or face"
2817 // << "face:" << f << " own:" << own << " nei:" << nei
2818 // << abort(FatalError);
2819 faceMap_.append(-1);
2820 }
2821 reverseFaceMap_.append(facei);
2822
2823 flipFaceFlux_.set(facei, flipFaceFlux);
2824
2825 if (zoneID >= 0)
2826 {
2827 faceZone_.insert(facei, zoneID);
2828 }
2829 faceZoneFlip_.set(facei, zoneFlip);
2830
2831 return facei;
2832}
2833
2834
2836(
2837 const face& f,
2838 const label facei,
2839 const label own,
2840 const label nei,
2841 const bool flipFaceFlux,
2842 const label patchID,
2843 const label zoneID,
2844 const bool zoneFlip
2845)
2846{
2847 // Check validity
2848 if (debug)
2849 {
2850 checkFace(f, facei, own, nei, patchID, zoneID);
2851 }
2852
2853 faces_[facei] = f;
2854 faceOwner_[facei] = own;
2855 faceNeighbour_[facei] = nei;
2856 region_[facei] = patchID;
2857
2858 flipFaceFlux_.set(facei, flipFaceFlux);
2859 faceZoneFlip_.set(facei, zoneFlip);
2860
2861 if (zoneID >= 0)
2862 {
2863 faceZone_.set(facei, zoneID);
2864 }
2865 else
2866 {
2867 faceZone_.erase(facei);
2868 }
2869}
2870
2871
2873(
2874 const label facei,
2875 const label mergeFacei
2876)
2877{
2878 if (facei < 0 || facei >= faces_.size())
2879 {
2881 << "illegal face label " << facei << endl
2882 << "Valid face labels are 0 .. " << faces_.size()-1
2883 << abort(FatalError);
2884 }
2885
2886 if
2887 (
2888 strict_
2889 && (faceRemoved(facei) || faceMap_[facei] == -1)
2890 )
2891 {
2893 << "face " << facei
2894 << " already marked for removal"
2895 << abort(FatalError);
2896 }
2897
2898 faces_[facei].setSize(0);
2899 region_[facei] = -1;
2900 faceOwner_[facei] = -1;
2901 faceNeighbour_[facei] = -1;
2902 faceMap_[facei] = -1;
2903 if (mergeFacei >= 0)
2904 {
2905 reverseFaceMap_[facei] = -mergeFacei-2;
2906 }
2907 else
2908 {
2909 reverseFaceMap_[facei] = -1;
2910 }
2911 faceFromEdge_.erase(facei);
2912 faceFromPoint_.erase(facei);
2913 flipFaceFlux_.unset(facei);
2914 faceZoneFlip_.unset(facei);
2915 faceZone_.erase(facei);
2916}
2917
2918
2920(
2921 const label masterPointID,
2922 const label masterEdgeID,
2923 const label masterFaceID,
2924 const label masterCellID,
2925 const label zoneID
2926)
2927{
2928 label celli = cellMap_.size();
2929
2930 if (masterPointID >= 0)
2931 {
2932 cellMap_.append(-1);
2933 cellFromPoint_.insert(celli, masterPointID);
2934 }
2935 else if (masterEdgeID >= 0)
2936 {
2937 cellMap_.append(-1);
2938 cellFromEdge_.insert(celli, masterEdgeID);
2939 }
2940 else if (masterFaceID >= 0)
2941 {
2942 cellMap_.append(-1);
2943 cellFromFace_.insert(celli, masterFaceID);
2944 }
2945 else
2946 {
2947 cellMap_.append(masterCellID);
2948 }
2949 reverseCellMap_.append(celli);
2950 cellZone_.append(zoneID);
2951
2952 return celli;
2953}
2954
2955
2957(
2958 const label celli,
2959 const label zoneID
2960)
2961{
2962 cellZone_[celli] = zoneID;
2963}
2964
2965
2967(
2968 const label celli,
2969 const label mergeCelli
2970)
2971{
2972 if (celli < 0 || celli >= cellMap_.size())
2973 {
2975 << "illegal cell label " << celli << endl
2976 << "Valid cell labels are 0 .. " << cellMap_.size()-1
2977 << abort(FatalError);
2978 }
2979
2980 if (strict_ && cellMap_[celli] == -2)
2981 {
2983 << "cell " << celli
2984 << " already marked for removal"
2985 << abort(FatalError);
2986 }
2987
2988 cellMap_[celli] = -2;
2989 if (mergeCelli >= 0)
2990 {
2991 reverseCellMap_[celli] = -mergeCelli-2;
2992 }
2993 else
2994 {
2995 reverseCellMap_[celli] = -1;
2996 }
2997 cellFromPoint_.erase(celli);
2998 cellFromEdge_.erase(celli);
2999 cellFromFace_.erase(celli);
3000 cellZone_[celli] = -1;
3001}
3002
3003
3005(
3006 polyMesh& mesh,
3007 const labelUList& patchMap,
3008 const bool inflate,
3009 const bool syncParallel,
3010 const bool orderCells,
3011 const bool orderPoints
3012)
3013{
3014 if (debug)
3015 {
3016 Pout<< "polyTopoChange::changeMesh"
3017 << "(polyMesh&, const bool, const bool, const bool, const bool)"
3018 << endl;
3019 }
3020
3021 if (debug)
3022 {
3023 Pout<< "Old mesh:" << nl;
3024 writeMeshStats(mesh, Pout);
3025 }
3026
3027 // new mesh points
3028 pointField newPoints;
3029 // number of internal points
3030 label nInternalPoints;
3031 // patch slicing
3032 labelList patchSizes;
3033 labelList patchStarts;
3034 // inflate maps
3035 List<objectMap> pointsFromPoints;
3036 List<objectMap> facesFromPoints;
3037 List<objectMap> facesFromEdges;
3038 List<objectMap> facesFromFaces;
3039 List<objectMap> cellsFromPoints;
3040 List<objectMap> cellsFromEdges;
3041 List<objectMap> cellsFromFaces;
3042 List<objectMap> cellsFromCells;
3043 // old mesh info
3044 List<Map<label>> oldPatchMeshPointMaps;
3045 labelList oldPatchNMeshPoints;
3046 labelList oldPatchStarts;
3047 List<Map<label>> oldFaceZoneMeshPointMaps;
3048
3049 // Compact, reorder patch faces and calculate mesh/patch maps.
3050 compactAndReorder
3051 (
3052 mesh,
3053 patchMap,
3054 syncParallel,
3055 orderCells,
3056 orderPoints,
3057
3058 nInternalPoints,
3059 newPoints,
3060 patchSizes,
3061 patchStarts,
3062 pointsFromPoints,
3063 facesFromPoints,
3064 facesFromEdges,
3065 facesFromFaces,
3066 cellsFromPoints,
3067 cellsFromEdges,
3068 cellsFromFaces,
3069 cellsFromCells,
3070 oldPatchMeshPointMaps,
3071 oldPatchNMeshPoints,
3072 oldPatchStarts,
3073 oldFaceZoneMeshPointMaps
3074 );
3075
3076 const label nOldPoints(mesh.nPoints());
3077 const label nOldFaces(mesh.nFaces());
3078 const label nOldCells(mesh.nCells());
3079 autoPtr<scalarField> oldCellVolumes(new scalarField(mesh.cellVolumes()));
3080
3081
3082 // Change the mesh
3083 // ~~~~~~~~~~~~~~~
3084 // This will invalidate any addressing so better make sure you have
3085 // all the information you need!!!
3086
3087 if (inflate)
3088 {
3089 // Keep (renumbered) mesh points, store new points in map for inflation
3090 // (appended points (i.e. from nowhere) get value zero)
3091 pointField renumberedMeshPoints(newPoints.size());
3092
3093 forAll(pointMap_, newPointi)
3094 {
3095 const label oldPointi = pointMap_[newPointi];
3096
3097 if (oldPointi >= 0)
3098 {
3099 renumberedMeshPoints[newPointi] = mesh.points()[oldPointi];
3100 }
3101 else
3102 {
3103 renumberedMeshPoints[newPointi] = Zero;
3104 }
3105 }
3106
3108 (
3109 autoPtr<pointField>::New(std::move(renumberedMeshPoints)),
3110 autoPtr<faceList>::New(std::move(faces_)),
3111 autoPtr<labelList>::New(std::move(faceOwner_)),
3112 autoPtr<labelList>::New(std::move(faceNeighbour_)),
3113 patchSizes,
3114 patchStarts,
3115 syncParallel
3116 );
3117
3118 mesh.topoChanging(true);
3119 // Note: could already set moving flag as well
3120 // mesh.moving(true);
3121 }
3122 else
3123 {
3124 // Set new points.
3126 (
3127 autoPtr<pointField>::New(std::move(newPoints)),
3128 autoPtr<faceList>::New(std::move(faces_)),
3129 autoPtr<labelList>::New(std::move(faceOwner_)),
3130 autoPtr<labelList>::New(std::move(faceNeighbour_)),
3131 patchSizes,
3132 patchStarts,
3133 syncParallel
3134 );
3135 mesh.topoChanging(true);
3136 }
3137
3138 // Clear out primitives
3139 {
3140 retiredPoints_.clearStorage();
3141 region_.clearStorage();
3142 }
3143
3144
3145 if (debug)
3146 {
3147 // Some stats on changes
3148 label nAdd, nInflate, nMerge, nRemove;
3149 countMap(pointMap_, reversePointMap_, nAdd, nInflate, nMerge, nRemove);
3150 Pout<< "Points:"
3151 << " added(from point):" << nAdd
3152 << " added(from nothing):" << nInflate
3153 << " merged(into other point):" << nMerge
3154 << " removed:" << nRemove
3155 << nl;
3156
3157 countMap(faceMap_, reverseFaceMap_, nAdd, nInflate, nMerge, nRemove);
3158 Pout<< "Faces:"
3159 << " added(from face):" << nAdd
3160 << " added(inflated):" << nInflate
3161 << " merged(into other face):" << nMerge
3162 << " removed:" << nRemove
3163 << nl;
3164
3165 countMap(cellMap_, reverseCellMap_, nAdd, nInflate, nMerge, nRemove);
3166 Pout<< "Cells:"
3167 << " added(from cell):" << nAdd
3168 << " added(inflated):" << nInflate
3169 << " merged(into other cell):" << nMerge
3170 << " removed:" << nRemove
3171 << nl
3172 << endl;
3173 }
3174
3175 if (debug)
3176 {
3177 Pout<< "New mesh:" << nl;
3178 writeMeshStats(mesh, Pout);
3179 }
3180
3181
3182 // Zones
3183 // ~~~~~
3184
3185 // Inverse of point/face/cell zone addressing.
3186 // For every preserved point/face/cells in zone give the old position.
3187 // For added points, the index is set to -1
3188 labelListList pointZoneMap(mesh.pointZones().size());
3189 labelListList faceZoneFaceMap(mesh.faceZones().size());
3190 labelListList cellZoneMap(mesh.cellZones().size());
3191
3192 resetZones(mesh, mesh, pointZoneMap, faceZoneFaceMap, cellZoneMap);
3193
3194 // Clear zone info
3195 {
3196 pointZone_.clearStorage();
3197 faceZone_.clearStorage();
3198 faceZoneFlip_.clearStorage();
3199 cellZone_.clearStorage();
3200 }
3201
3202
3203 // Patch point renumbering
3204 // For every preserved point on a patch give the old position.
3205 // For added points, the index is set to -1
3206 labelListList patchPointMap(mesh.boundaryMesh().size());
3207 calcPatchPointMap
3208 (
3209 oldPatchMeshPointMaps,
3210 patchMap,
3212 patchPointMap
3213 );
3214
3215 // Create the face zone mesh point renumbering
3216 labelListList faceZonePointMap(mesh.faceZones().size());
3217 calcFaceZonePointMap(mesh, oldFaceZoneMeshPointMaps, faceZonePointMap);
3218
3219 labelHashSet flipFaceFluxSet(HashSetOps::used(flipFaceFlux_));
3220
3222 (
3223 mesh,
3224 nOldPoints,
3225 nOldFaces,
3226 nOldCells,
3227
3228 pointMap_,
3229 pointsFromPoints,
3230
3231 faceMap_,
3232 facesFromPoints,
3233 facesFromEdges,
3234 facesFromFaces,
3235
3236 cellMap_,
3237 cellsFromPoints,
3238 cellsFromEdges,
3239 cellsFromFaces,
3240 cellsFromCells,
3241
3242 reversePointMap_,
3243 reverseFaceMap_,
3244 reverseCellMap_,
3245
3246 flipFaceFluxSet,
3247
3248 patchPointMap,
3249
3250 pointZoneMap,
3251
3252 faceZonePointMap,
3253 faceZoneFaceMap,
3254 cellZoneMap,
3255
3256 newPoints, // if empty signals no inflation.
3257 oldPatchStarts,
3258 oldPatchNMeshPoints,
3259
3260 oldCellVolumes,
3261
3262 true // steal storage.
3263 );
3264
3265 // At this point all member DynamicList (pointMap_, cellMap_ etc.) will
3266 // be invalid.
3267}
3268
3269
3271(
3272 polyMesh& mesh,
3273 const bool inflate,
3274 const bool syncParallel,
3275 const bool orderCells,
3276 const bool orderPoints
3277)
3278{
3279 return changeMesh
3280 (
3281 mesh,
3283 inflate,
3284 syncParallel,
3285 orderCells,
3286 orderPoints
3287 );
3288}
3289
3290
3291// ************************************************************************* //
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
void setSize(const label n)
Alias for resize()
Definition: List.H:218
virtual const fileName & name() const
Get the name of the stream.
Definition: OSstream.H:107
void setSize(const label newLen)
Same as resize()
Definition: PtrList.H:151
void resize(const label newLen)
Adjust size of PtrList.
Definition: PtrList.C:103
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
bool found(const T &val, label pos=0) const
True if the value if found in the list.
Definition: UListI.H:265
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type test(const label i) const
Definition: UList.H:518
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
@ nonBlocking
"nonBlocking"
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
const DynamicList< point > & facePoints()
Returns the points of the cutting PLICface.
Definition: cutCellIso.C:173
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
void movePoints()
Update for new mesh geometry.
static const complex max
complex (VGREAT,VGREAT)
Definition: complex.H:293
const word & name() const noexcept
The patch name.
Class containing data for cell addition.
Definition: polyAddCell.H:51
label zoneID() const
Cell zone ID.
Definition: polyAddCell.H:182
label masterFaceID() const
Return master face ID.
Definition: polyAddCell.H:164
label masterPointID() const
Return master point ID.
Definition: polyAddCell.H:152
label masterCellID() const
Return master cell ID.
Definition: polyAddCell.H:170
label masterEdgeID() const
Return master edge ID.
Definition: polyAddCell.H:158
A face addition data class. A face can be inflated either from a point or from another face and can e...
Definition: polyAddFace.H:56
bool flipFaceFlux() const
Does the face flux need to be flipped.
Definition: polyAddFace.H:300
label owner() const
Return owner cell.
Definition: polyAddFace.H:246
label zoneID() const
Face zone ID.
Definition: polyAddFace.H:330
label masterFaceID() const
Return master face ID.
Definition: polyAddFace.H:294
label masterPointID() const
Return master point ID.
Definition: polyAddFace.H:282
const face & newFace() const
Return face.
Definition: polyAddFace.H:240
label patchID() const
Boundary patch ID.
Definition: polyAddFace.H:312
label zoneFlip() const
Face zone flip.
Definition: polyAddFace.H:336
label neighbour() const
Return neighbour cell.
Definition: polyAddFace.H:252
label masterEdgeID() const
Return master edge ID.
Definition: polyAddFace.H:288
Class containing data for point addition.
Definition: polyAddPoint.H:52
label zoneID() const
Point zone ID.
Definition: polyAddPoint.H:150
const point & newPoint() const
Point location.
Definition: polyAddPoint.H:126
bool inCell() const
Does the point support a cell.
Definition: polyAddPoint.H:156
label masterPointID() const
Master point label.
Definition: polyAddPoint.H:132
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
label start() const noexcept
The start label of boundary faces in the polyMesh face list.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:498
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1108
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1121
bool topoChanging() const noexcept
Is mesh topology changing.
Definition: polyMesh.H:548
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:504
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1127
const pointZoneMesh & pointZones() const noexcept
Return point zone mesh.
Definition: polyMesh.H:492
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1083
void resetPrimitives(autoPtr< pointField > &&points, autoPtr< faceList > &&faces, autoPtr< labelList > &&owner, autoPtr< labelList > &&neighbour, const labelUList &patchSizes, const labelUList &patchStarts, const bool validBoundary=true)
Reset mesh primitive data. Assumes all patch info correct.
Definition: polyMesh.C:718
Class describing modification of a cell.
bool removeFromZone() const
label zoneID() const
Cell zone ID.
label cellID() const
Cell ID.
Class describing modification of a face.
bool flipFaceFlux() const
Does the face flux need to be flipped.
label owner() const
Return owner cell ID.
label zoneID() const
Face zone ID.
const face & newFace() const
Return face.
label patchID() const
Boundary patch ID.
label zoneFlip() const
Face zone flip.
label neighbour() const
Return owner cell ID.
label faceID() const
Return master face ID.
Class describing modification of a point.
label zoneID() const
Point zone ID.
const point & newPoint() const
New point location.
bool inCell() const
Does the point support a cell.
label pointID() const
Point ID.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:364
Class containing data for cell removal.
label cellID() const
Return cell ID.
label mergeCellID() const
Return cell ID.
Class containing data for face removal.
label faceID() const
Return face ID.
label mergeFaceID() const
Return merge face ID.
Class containing data for point removal.
label mergePointID() const
label pointID() const
Return point ID.
Direct mesh changes based on v1.3 polyTopoChange syntax.
label addPoint(const point &pt, const label masterPointID, const label zoneID, const bool inCell)
Add point. Return new point label.
void removeCell(const label celli, const label mergeCelli)
Remove/merge cell.
void modifyPoint(const label pointi, const point &pt, const label zoneID, const bool inCell)
Modify coordinate.
void setCapacity(const label nPoints, const label nFaces, const label nCells)
label addCell(const label masterPointID, const label masterEdgeID, const label masterFaceID, const label masterCellID, const label zoneID)
Add cell. Return new cell label.
void addMesh(const polyMesh &mesh, const labelUList &patchMap, const labelUList &pointZoneMap, const labelUList &faceZoneMap, const labelUList &cellZoneMap)
void removePoint(const label pointi, const label mergePointi)
Remove/merge point.
void modifyCell(const label celli, const label zoneID)
Modify zone of cell.
void modifyFace(const face &f, const label facei, const label own, const label nei, const bool flipFaceFlux, const label patchID, const label zoneID, const bool zoneFlip)
Modify vertices or cell of face.
void removeFace(const label facei, const label mergeFacei)
Remove/merge face.
label addFace(const face &f, const label own, const label nei, const label masterPointID, const label masterEdgeID, const label masterFaceID, const bool flipFaceFlux, const label patchID, const label zoneID, const bool zoneFlip)
Add face to cells. Return new face label.
autoPtr< mapPolyMesh > changeMesh(polyMesh &mesh, const labelUList &patchMap, const bool inflate, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Inplace changes mesh without change of patches.
void clear()
Clear all storage.
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
const labelListList & pointCells() const
const scalarField & cellVolumes() 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 labelListList & pointFaces() const
const labelListList & edgeFaces() const
const labelListList & edgeCells() const
void compact()
Compact splitCells_. Removes all freeSplitCells_ elements.
A virtual base class for topological actions.
Definition: topoAction.H:52
setAction
Enumeration defining various actions.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
const polyBoundaryMesh & patches
faceListList boundary
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
const labelList nFaces(UPstream::listGatherValues< label >(aMesh.nFaces()))
const label nPatches
const pointField & points
label nPoints
const labelIOList & zoneID
#define WarningInFunction
Report a warning using Foam::Warning.
labelHashSet used(const bitSet &select)
Convert a bitset to a labelHashSet of the indices used.
Definition: HashOps.C:35
const wordList internal
Standard dimensioned field types (scalar, vector, tensor, etc)
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:38
List< label > labelList
A List of labels.
Definition: List.H:66
ZoneMesh< cellZone, polyMesh > cellZoneMesh
A ZoneMesh with the type cellZone.
ZoneMesh< pointZone, polyMesh > pointZoneMesh
A ZoneMesh with the type pointZone.
IntListType renumber(const labelUList &oldToNew, const IntListType &input)
Renumber the values (not the indices) of a list.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
constexpr label labelMax
Definition: label.H:61
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with the type faceZone.
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
List< List< bool > > boolListList
A List of boolList.
Definition: boolList.H:48
void inplaceReorder(const labelUList &oldToNew, ListType &input, const bool prune=false)
Inplace reorder the elements of a list.
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)
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
errorManip< error > abort(error &err)
Definition: errorManip.H:144
List< bool > boolList
A List of bools.
Definition: List.H:64
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
error FatalError
labelList invert(const label len, const labelUList &map)
Create an inverse one-to-one mapping.
Definition: ListOps.C:36
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
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
ListType reorder(const labelUList &oldToNew, const ListType &input, const bool prune=false)
Reorder the elements of a list.
void stableSort(UList< T > &list)
Stable sort the list.
Definition: UList.C:356
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
label newPointi
Definition: readKivaGrid.H:496
labelList f(nPoints)
labelList pointLabels(nPoints, -1)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278