Loading...
Searching...
No Matches
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,2024 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"
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{
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 const auto& czs = mesh.cellZones();
206 labelList cellZoneSizes(czs.size(), 0);
207 for (const auto& cz : czs)
208 {
209 cellZoneSizes[cz.index()] = cz.size();
210 }
211 const auto& fzs = mesh.faceZones();
212 labelList faceZoneSizes(fzs.size(), 0);
213 for (const auto& fz : fzs)
214 {
215 faceZoneSizes[fz.index()] = fz.size();
216 }
217 const auto& pzs = mesh.pointZones();
218 labelList pointZoneSizes(pzs.size(), 0);
219 for (const auto& pz : pzs)
220 {
221 pointZoneSizes[pz.index()] = pz.size();
222 }
223
224 os << " Points : " << mesh.nPoints() << nl
225 << " Faces : " << mesh.nFaces() << nl
226 << " Cells : " << mesh.nCells() << nl
227 << " PatchSizes : " << flatOutput(patchSizes) << nl
228 << " PatchStarts : " << flatOutput(patchStarts) << nl
229 << " cZoneSizes : " << flatOutput(cellZoneSizes) << nl
230 << " fZoneSizes : " << flatOutput(faceZoneSizes) << nl
231 << " pZoneSizes : " << flatOutput(pointZoneSizes) << nl
232 << endl;
233}
234
235
236void Foam::polyTopoChange::getMergeSets
237(
238 const labelUList& reverseCellMap,
239 const labelUList& cellMap,
240 List<objectMap>& cellsFromCells
241)
242{
243 // Per new cell the number of old cells that have been merged into it
244 labelList nMerged(cellMap.size(), 1);
245
246 forAll(reverseCellMap, oldCelli)
247 {
248 const label newCelli = reverseCellMap[oldCelli];
249
250 if (newCelli < -1)
251 {
252 label mergeCelli = -newCelli-2;
253
254 nMerged[mergeCelli]++;
255 }
256 }
257
258 // From merged cell to set index
259 labelList cellToMergeSet(cellMap.size(), -1);
260
261 label nSets = 0;
262
263 forAll(nMerged, celli)
264 {
265 if (nMerged[celli] > 1)
266 {
267 cellToMergeSet[celli] = nSets++;
268 }
269 }
270
271 // Collect cell labels.
272 // Each objectMap will have
273 // - index : new mesh cell label
274 // - masterObjects : list of old cells that have been merged. Element 0
275 // will be the original destination cell label.
276
277 cellsFromCells.setSize(nSets);
278
279 forAll(reverseCellMap, oldCelli)
280 {
281 const label newCelli = reverseCellMap[oldCelli];
282
283 if (newCelli < -1)
284 {
285 const label mergeCelli = -newCelli-2;
286
287 // oldCelli was merged into mergeCelli
288
289 const label setI = cellToMergeSet[mergeCelli];
290
291 objectMap& mergeSet = cellsFromCells[setI];
292
293 if (mergeSet.masterObjects().empty())
294 {
295 // First occurrence of master cell mergeCelli
296
297 mergeSet.index() = mergeCelli;
298 mergeSet.masterObjects().setSize(nMerged[mergeCelli]);
299
300 // old master label
301 mergeSet.masterObjects()[0] = cellMap[mergeCelli];
302
303 // old slave label
304 mergeSet.masterObjects()[1] = oldCelli;
305
306 nMerged[mergeCelli] = 2;
307 }
308 else
309 {
310 mergeSet.masterObjects()[nMerged[mergeCelli]++] = oldCelli;
311 }
312 }
313 }
314}
315
316
317bool Foam::polyTopoChange::hasValidPoints(const face& f) const
318{
319 for (const label fp : f)
320 {
321 if (fp < 0 || fp >= points_.size())
322 {
323 return false;
324 }
325 }
326 return true;
327}
328
329
330Foam::pointField Foam::polyTopoChange::facePoints(const face& f) const
331{
333 forAll(f, fp)
334 {
335 if (f[fp] < 0 && f[fp] >= points_.size())
336 {
338 << "Problem." << abort(FatalError);
339 }
340 points[fp] = points_[f[fp]];
341 }
342 return points;
343}
344
345
346void Foam::polyTopoChange::checkFace
347(
348 const face& f,
349 const label facei,
350 const label own,
351 const label nei,
352 const label patchi,
353 const label zoneI
354) const
355{
356 if (nei == -1)
357 {
358 if (own == -1 && zoneI != -1)
359 {
360 // retired face
361 }
362 else if (patchi == -1 || patchi >= nPatches_)
363 {
365 << "Face has no neighbour (so external) but does not have"
366 << " a valid patch" << nl
367 << "f:" << f
368 << " facei(-1 if added face):" << facei
369 << " own:" << own << " nei:" << nei
370 << " patchi:" << patchi << nl;
371 if (hasValidPoints(f))
372 {
374 << "points (removed points marked with "
375 << vector::max << ") " << facePoints(f);
376 }
378 }
379 }
380 else
381 {
382 if (patchi != -1)
383 {
385 << "Cannot both have valid patchi and neighbour" << nl
386 << "f:" << f
387 << " facei(-1 if added face):" << facei
388 << " own:" << own << " nei:" << nei
389 << " patchi:" << patchi << nl;
390 if (hasValidPoints(f))
391 {
393 << "points (removed points marked with "
394 << vector::max << ") : " << facePoints(f);
395 }
397 }
398
399 if (nei <= own)
400 {
402 << "Owner cell label should be less than neighbour cell label"
403 << nl
404 << "f:" << f
405 << " facei(-1 if added face):" << facei
406 << " own:" << own << " nei:" << nei
407 << " patchi:" << patchi << nl;
408 if (hasValidPoints(f))
409 {
411 << "points (removed points marked with "
412 << vector::max << ") : " << facePoints(f);
413 }
415 }
416 }
417
418 if (f.size() < 3 || f.found(-1))
419 {
421 << "Illegal vertices in face"
422 << nl
423 << "f:" << f
424 << " facei(-1 if added face):" << facei
425 << " own:" << own << " nei:" << nei
426 << " patchi:" << patchi << nl;
427 if (hasValidPoints(f))
428 {
430 << "points (removed points marked with "
431 << vector::max << ") : " << facePoints(f);
432 }
434 }
435 if (facei >= 0 && facei < faces_.size() && faceRemoved(facei))
436 {
438 << "Face already marked for removal"
439 << nl
440 << "f:" << f
441 << " facei(-1 if added face):" << facei
442 << " own:" << own << " nei:" << nei
443 << " patchi:" << patchi << nl;
444 if (hasValidPoints(f))
445 {
447 << "points (removed points marked with "
448 << vector::max << ") : " << facePoints(f);
449 }
451 }
452 forAll(f, fp)
453 {
454 if (f[fp] < points_.size() && pointRemoved(f[fp]))
455 {
457 << "Face uses removed vertices"
458 << nl
459 << "f:" << f
460 << " facei(-1 if added face):" << facei
461 << " own:" << own << " nei:" << nei
462 << " patchi:" << patchi << nl;
463 if (hasValidPoints(f))
464 {
466 << "points (removed points marked with "
467 << vector::max << ") : " << facePoints(f);
468 }
470 }
471 }
472}
473
474
475void Foam::polyTopoChange::makeCells
476(
477 const label nActiveFaces,
478 labelList& cellFaces,
479 labelList& cellFaceOffsets
480) const
481{
482 cellFaces.setSize(2*nActiveFaces);
483 cellFaceOffsets.setSize(cellMap_.size() + 1);
484
485 // Faces per cell
486 labelList nNbrs(cellMap_.size(), Zero);
487
488 // 1. Count faces per cell
489
490 for (label facei = 0; facei < nActiveFaces; facei++)
491 {
492 if (faceOwner_[facei] < 0)
493 {
494 pointField newPoints;
495 if (facei < faces_.size())
496 {
497 const face& f = faces_[facei];
498 newPoints.setSize(f.size(), vector::max);
499 forAll(f, fp)
500 {
501 if (f[fp] < points_.size())
502 {
503 newPoints[fp] = points_[f[fp]];
504 }
505 }
506 }
507
508
510 << "Face " << facei << " is active but its owner has"
511 << " been deleted. This is usually due to deleting cells"
512 << " without modifying exposed faces to be boundary faces."
513 << exit(FatalError);
514 }
515 nNbrs[faceOwner_[facei]]++;
516 }
517 for (label facei = 0; facei < nActiveFaces; facei++)
518 {
519 if (faceNeighbour_[facei] >= 0)
520 {
521 nNbrs[faceNeighbour_[facei]]++;
522 }
523 }
524
525 // 2. Calculate offsets
526
527 cellFaceOffsets[0] = 0;
528 forAll(nNbrs, celli)
529 {
530 cellFaceOffsets[celli+1] = cellFaceOffsets[celli] + nNbrs[celli];
531 }
532
533 // 3. Fill faces per cell
534
535 // reset the whole list to use as counter
536 nNbrs = 0;
537
538 for (label facei = 0; facei < nActiveFaces; facei++)
539 {
540 label celli = faceOwner_[facei];
541
542 cellFaces[cellFaceOffsets[celli] + nNbrs[celli]++] = facei;
543 }
544
545 for (label facei = 0; facei < nActiveFaces; facei++)
546 {
547 label celli = faceNeighbour_[facei];
548
549 if (celli >= 0)
550 {
551 cellFaces[cellFaceOffsets[celli] + nNbrs[celli]++] = facei;
552 }
553 }
554
555 // Last offset points to beyond end of cellFaces.
556 cellFaces.setSize(cellFaceOffsets[cellMap_.size()]);
557}
558
559
560// Create cell-cell addressing. Called after compaction (but before ordering)
561// of faces
562void Foam::polyTopoChange::makeCellCells
563(
564 const label nActiveFaces,
565 CompactListList<label>& cellCells
566) const
567{
568 // Neighbours per cell
569 labelList nNbrs(cellMap_.size(), Zero);
570
571 // 1. Count neighbours (through internal faces) per cell
572
573 for (label facei = 0; facei < nActiveFaces; facei++)
574 {
575 if (faceNeighbour_[facei] >= 0)
576 {
577 nNbrs[faceOwner_[facei]]++;
578 nNbrs[faceNeighbour_[facei]]++;
579 }
580 }
581
582 // 2. Construct csr
583 cellCells.setSize(nNbrs);
584
585
586 // 3. Fill faces per cell
587
588 // reset the whole list to use as counter
589 nNbrs = 0;
590
591 for (label facei = 0; facei < nActiveFaces; facei++)
592 {
593 label nei = faceNeighbour_[facei];
594
595 if (nei >= 0)
596 {
597 label own = faceOwner_[facei];
598 cellCells(own, nNbrs[own]++) = nei;
599 cellCells(nei, nNbrs[nei]++) = own;
600 }
601 }
602}
603
604
605// Cell ordering (based on bandCompression).
606// Handles removed cells. Returns number of remaining cells.
607Foam::label Foam::polyTopoChange::getCellOrder
608(
609 const CompactListList<label>& cellCellAddressing,
610 labelList& oldToNew
611) const
612{
613 const label nOldCells(cellCellAddressing.size());
614
615 // Which cells are visited/unvisited
616 bitSet unvisited(nOldCells, true);
617
618 // Exclude removed cells
619 for (label celli = 0; celli < nOldCells; ++celli)
620 {
621 if (cellRemoved(celli))
622 {
623 unvisited.unset(celli);
624 }
625 }
626
627 // The new output order
628 labelList newOrder(unvisited.count());
629
630
631 // Various work arrays
632 // ~~~~~~~~~~~~~~~~~~~
633
634 //- Neighbour cells
635 DynamicList<label> nbrCells;
636
637 //- Neighbour ordering
638 DynamicList<label> nbrOrder;
639
640 //- Corresponding weights for neighbour cells
641 DynamicList<label> weights;
642
643 // FIFO buffer for renumbering.
644 CircularBuffer<label> queuedCells(1024);
645
646
647 label cellInOrder = 0;
648
649 while (true)
650 {
651 // For a disconnected region find the lowest connected cell.
652 label currCelli = -1;
653 label minCount = labelMax;
654
655 for (const label celli : unvisited)
656 {
657 const label nbrCount = cellCellAddressing[celli].size();
658
659 if (minCount > nbrCount)
660 {
661 minCount = nbrCount;
662 currCelli = celli;
663 }
664 }
665
666 if (currCelli == -1)
667 {
668 break;
669 }
670
671
672 // Starting from currCelli - walk breadth-first
673
674 queuedCells.push_back(currCelli);
675
676 // Loop through queuedCells list. Add the first cell into the
677 // cell order if it has not already been visited and ask for its
678 // neighbours. If the neighbour in question has not been visited,
679 // add it to the end of the queuedCells list
680
681 while (!queuedCells.empty())
682 {
683 // Process as FIFO
684 currCelli = queuedCells.front();
685 queuedCells.pop_front();
686
687 if (unvisited.test(currCelli))
688 {
689 // First visit...
690 unvisited.unset(currCelli);
691
692 // Add into cellOrder
693 newOrder[cellInOrder] = currCelli;
694 ++cellInOrder;
695
696 // Find if the neighbours have been visited
697 const auto& neighbours = cellCellAddressing[currCelli];
698
699 // Add in increasing order of connectivity
700
701 // 1. Count neighbours of unvisited neighbours
702 nbrCells.clear();
703 weights.clear();
704
705 for (const label nbr : neighbours)
706 {
707 const label nbrCount = cellCellAddressing[nbr].size();
708
709 if (unvisited.test(nbr))
710 {
711 // Not visited (or removed), add to the list
712 nbrCells.push_back(nbr);
713 weights.push_back(nbrCount);
714 }
715 }
716
717 // Resize DynamicList prior to sortedOrder
718 nbrOrder.resize_nocopy(weights.size());
719
720 // 2. Ascending order
721 Foam::sortedOrder(weights, nbrOrder);
722
723 // 3. Add to FIFO in sorted order
724 for (const label nbrIdx : nbrOrder)
725 {
726 queuedCells.push_back(nbrCells[nbrIdx]);
727 }
728 }
729 }
730 }
731
732 // Now we have new-to-old in newOrder.
733 newOrder.resize(cellInOrder); // Extra safety, but should be a no-op
734
735 // Invert to get old-to-new. Make sure removed (i.e. unmapped) cells are -1.
736 oldToNew = invert(nOldCells, newOrder);
737
738 return cellInOrder;
739}
740
741
742// Determine order for faces:
743// - upper-triangular order for internal faces
744// - external faces after internal faces and in patch order.
745void Foam::polyTopoChange::getFaceOrder
746(
747 const label nActiveFaces,
748 const labelUList& cellFaces,
749 const labelUList& cellFaceOffsets,
750
751 labelList& oldToNew,
752 labelList& patchSizes,
753 labelList& patchStarts
754) const
755{
756 oldToNew.setSize(faceOwner_.size());
757 oldToNew = -1;
758
759 // First unassigned face
760 label newFacei = 0;
761
762 labelList nbr;
763 labelList order;
764
765 forAll(cellMap_, celli)
766 {
767 label startOfCell = cellFaceOffsets[celli];
768 label nFaces = cellFaceOffsets[celli+1] - startOfCell;
769
770 // Neighbouring cells
771 nbr.setSize(nFaces);
772
773 for (label i = 0; i < nFaces; i++)
774 {
775 label facei = cellFaces[startOfCell + i];
776
777 label nbrCelli = faceNeighbour_[facei];
778
779 if (facei >= nActiveFaces)
780 {
781 // Retired face.
782 nbr[i] = -1;
783 }
784 else if (nbrCelli != -1)
785 {
786 // Internal face. Get cell on other side.
787 if (nbrCelli == celli)
788 {
789 nbrCelli = faceOwner_[facei];
790 }
791
792 if (celli < nbrCelli)
793 {
794 // Celli is master
795 nbr[i] = nbrCelli;
796 }
797 else
798 {
799 // nbrCell is master. Let it handle this face.
800 nbr[i] = -1;
801 }
802 }
803 else
804 {
805 // External face. Do later.
806 nbr[i] = -1;
807 }
808 }
809
810 sortedOrder(nbr, order);
811
812 for (const label index : order)
813 {
814 if (nbr[index] != -1)
815 {
816 oldToNew[cellFaces[startOfCell + index]] = newFacei++;
817 }
818 }
819 }
820
821
822 // Pick up all patch faces in patch face order.
823 patchStarts.setSize(nPatches_);
824 patchStarts = 0;
825 patchSizes.setSize(nPatches_);
826 patchSizes = 0;
827
828 if (nPatches_ > 0)
829 {
830 patchStarts[0] = newFacei;
831
832 for (label facei = 0; facei < nActiveFaces; facei++)
833 {
834 if (region_[facei] >= 0)
835 {
836 patchSizes[region_[facei]]++;
837 }
838 }
839
840 label facei = patchStarts[0];
841
842 forAll(patchStarts, patchi)
843 {
844 patchStarts[patchi] = facei;
845 facei += patchSizes[patchi];
846 }
847 }
848
849 //if (debug)
850 //{
851 // Pout<< "patchSizes:" << patchSizes << nl
852 // << "patchStarts:" << patchStarts << endl;
853 //}
854
855 labelList workPatchStarts(patchStarts);
856
857 for (label facei = 0; facei < nActiveFaces; facei++)
858 {
859 if (region_[facei] >= 0)
860 {
861 oldToNew[facei] = workPatchStarts[region_[facei]]++;
862 }
863 }
864
865 // Retired faces.
866 for (label facei = nActiveFaces; facei < oldToNew.size(); facei++)
867 {
868 oldToNew[facei] = facei;
869 }
870
871 // Check done all faces.
872 forAll(oldToNew, facei)
873 {
874 if (oldToNew[facei] == -1)
875 {
877 << "Did not determine new position"
878 << " for face " << facei
879 << " owner " << faceOwner_[facei]
880 << " neighbour " << faceNeighbour_[facei]
881 << " region " << region_[facei] << endl
882 << "This is usually caused by not specifying a patch for"
883 << " a boundary face." << nl
884 << "Switch on the polyTopoChange::debug flag to catch"
885 << " this error earlier." << nl;
886 if (hasValidPoints(faces_[facei]))
887 {
889 << "points (removed points marked with "
890 << vector::max << ") " << facePoints(faces_[facei]);
891 }
893 }
894 }
895}
896
897
898// Reorder and compact faces according to map.
899void Foam::polyTopoChange::reorderCompactFaces
900(
901 const label newSize,
902 const labelUList& oldToNew
903)
904{
905 reorder(oldToNew, faces_);
906 faces_.setCapacity(newSize);
907
908 reorder(oldToNew, region_);
909 region_.setCapacity(newSize);
910
911 reorder(oldToNew, faceOwner_);
912 faceOwner_.setCapacity(newSize);
913
914 reorder(oldToNew, faceNeighbour_);
915 faceNeighbour_.setCapacity(newSize);
916
917 // Update faceMaps.
918 reorder(oldToNew, faceMap_);
919 faceMap_.setCapacity(newSize);
920
921 renumberReverseMap(oldToNew, reverseFaceMap_);
922
923 renumberKey(oldToNew, faceFromPoint_);
924 renumberKey(oldToNew, faceFromEdge_);
925 inplaceReorder(oldToNew, flipFaceFlux_);
926 flipFaceFlux_.setCapacity(newSize);
927 renumberKey(oldToNew, faceZone_);
928 inplaceReorder(oldToNew, faceZoneFlip_);
929 faceZoneFlip_.setCapacity(newSize);
930 if (faceAdditionalZones_.size())
931 {
932 // Extend to number of faces so oldToNew can be used
933 faceAdditionalZones_.setSize(newSize);
934 inplaceReorder(oldToNew, faceAdditionalZones_);
935 //faceAdditionalZones_.setCapacity(newSize);
936 }
937}
938
941{
942 points_.shrink();
943 pointMap_.shrink();
944 reversePointMap_.shrink();
945 pointAdditionalZones_.shrink();
946
947 faces_.shrink();
948 region_.shrink();
949 faceOwner_.shrink();
950 faceNeighbour_.shrink();
951 faceMap_.shrink();
952 reverseFaceMap_.shrink();
953 faceAdditionalZones_.shrink();
954
955 cellMap_.shrink();
956 reverseCellMap_.shrink();
957 cellZone_.shrink();
958 cellAdditionalZones_.shrink();
959}
960
961
962// Compact all and orders points and faces:
963// - points into internal followed by external points
964// - internalfaces upper-triangular
965// - externalfaces after internal ones.
966void Foam::polyTopoChange::compact
967(
968 const bool orderCells,
969 const bool orderPoints,
970 label& nInternalPoints,
971 labelList& patchSizes,
972 labelList& patchStarts
973)
974{
975 shrink();
976
977 // Compact points
978 label nActivePoints = 0;
979 {
980 labelList localPointMap(points_.size(), -1);
981 label newPointi = 0;
982
983 if (!orderPoints)
984 {
985 nInternalPoints = -1;
986
987 forAll(points_, pointi)
988 {
989 if (!pointRemoved(pointi) && !retiredPoints_.found(pointi))
990 {
991 localPointMap[pointi] = newPointi++;
992 }
993 }
994 nActivePoints = newPointi;
995 }
996 else
997 {
998 forAll(points_, pointi)
999 {
1000 if (!pointRemoved(pointi) && !retiredPoints_.found(pointi))
1001 {
1002 nActivePoints++;
1003 }
1004 }
1005
1006 // Mark boundary points
1007 forAll(faceOwner_, facei)
1008 {
1009 if
1010 (
1011 !faceRemoved(facei)
1012 && faceOwner_[facei] >= 0
1013 && faceNeighbour_[facei] < 0
1014 )
1015 {
1016 // Valid boundary face
1017 const face& f = faces_[facei];
1018
1019 forAll(f, fp)
1020 {
1021 label pointi = f[fp];
1022
1023 if (localPointMap[pointi] == -1)
1024 {
1025 if
1026 (
1027 pointRemoved(pointi)
1028 || retiredPoints_.found(pointi)
1029 )
1030 {
1032 << "Removed or retired point " << pointi
1033 << " in face " << f
1034 << " at position " << facei << endl
1035 << "Probably face has not been adapted for"
1036 << " removed points." << abort(FatalError);
1037 }
1038 localPointMap[pointi] = newPointi++;
1039 }
1040 }
1041 }
1042 }
1043
1044 label nBoundaryPoints = newPointi;
1045 nInternalPoints = nActivePoints - nBoundaryPoints;
1046
1047 // Move the boundary addressing up
1048 forAll(localPointMap, pointi)
1049 {
1050 if (localPointMap[pointi] != -1)
1051 {
1052 localPointMap[pointi] += nInternalPoints;
1053 }
1054 }
1055
1056 newPointi = 0;
1057
1058 // Mark internal points
1059 forAll(faceOwner_, facei)
1060 {
1061 if
1062 (
1063 !faceRemoved(facei)
1064 && faceOwner_[facei] >= 0
1065 && faceNeighbour_[facei] >= 0
1066 )
1067 {
1068 // Valid internal face
1069 const face& f = faces_[facei];
1070
1071 forAll(f, fp)
1072 {
1073 label pointi = f[fp];
1074
1075 if (localPointMap[pointi] == -1)
1076 {
1077 if
1078 (
1079 pointRemoved(pointi)
1080 || retiredPoints_.found(pointi)
1081 )
1082 {
1084 << "Removed or retired point " << pointi
1085 << " in face " << f
1086 << " at position " << facei << endl
1087 << "Probably face has not been adapted for"
1088 << " removed points." << abort(FatalError);
1089 }
1090 localPointMap[pointi] = newPointi++;
1091 }
1092 }
1093 }
1094 }
1095
1096 if (newPointi != nInternalPoints)
1097 {
1099 << "Problem." << abort(FatalError);
1100 }
1101 newPointi = nActivePoints;
1102 }
1103
1104 for (const label pointi : retiredPoints_)
1105 {
1106 localPointMap[pointi] = newPointi++;
1107 }
1108
1109
1110 if (debug)
1111 {
1112 Pout<< "Points : active:" << nActivePoints
1113 << " removed:" << points_.size()-newPointi << endl;
1114 }
1115
1116 reorder(localPointMap, points_);
1117 points_.setCapacity(newPointi);
1118
1119 // Update pointMaps
1120 reorder(localPointMap, pointMap_);
1121 pointMap_.setCapacity(newPointi);
1122 renumberReverseMap(localPointMap, reversePointMap_);
1123
1124 renumberKey(localPointMap, pointZone_);
1125 renumber(localPointMap, retiredPoints_);
1126 if (pointAdditionalZones_.size())
1127 {
1128 reorder(localPointMap, pointAdditionalZones_);
1129 }
1130
1131 // Use map to relabel face vertices
1132 forAll(faces_, facei)
1133 {
1134 face& f = faces_[facei];
1135
1136 //labelList oldF(f);
1137 renumberCompact(localPointMap, f);
1138
1139 if (!faceRemoved(facei) && f.size() < 3)
1140 {
1142 << "Created illegal face " << f
1143 //<< " from face " << oldF
1144 << " at position:" << facei
1145 << " when filtering removed points"
1146 << abort(FatalError);
1147 }
1148 }
1149 }
1150
1151
1152 // Compact faces.
1153 {
1154 labelList localFaceMap(faces_.size(), -1);
1155 label newFacei = 0;
1156
1157 forAll(faces_, facei)
1158 {
1159 if (!faceRemoved(facei) && faceOwner_[facei] >= 0)
1160 {
1161 localFaceMap[facei] = newFacei++;
1162 }
1163 }
1164 nActiveFaces_ = newFacei;
1165
1166 forAll(faces_, facei)
1167 {
1168 if (!faceRemoved(facei) && faceOwner_[facei] < 0)
1169 {
1170 // Retired face
1171 localFaceMap[facei] = newFacei++;
1172 }
1173 }
1174
1175 if (debug)
1176 {
1177 Pout<< "Faces : active:" << nActiveFaces_
1178 << " removed:" << faces_.size()-newFacei << endl;
1179 }
1180
1181 // Reorder faces.
1182 reorderCompactFaces(newFacei, localFaceMap);
1183 }
1184
1185 // Compact cells.
1186 {
1187 labelList localCellMap;
1188 label newCelli;
1189
1190 if (orderCells)
1191 {
1192 // Construct cellCell addressing
1193 CompactListList<label> cellCells;
1194 makeCellCells(nActiveFaces_, cellCells);
1195
1196 // Cell ordering (based on bandCompression). Handles removed cells.
1197 newCelli = getCellOrder(cellCells, localCellMap);
1198 }
1199 else
1200 {
1201 // Compact out removed cells
1202 localCellMap.setSize(cellMap_.size());
1203 localCellMap = -1;
1204
1205 newCelli = 0;
1206 forAll(cellMap_, celli)
1207 {
1208 if (!cellRemoved(celli))
1209 {
1210 localCellMap[celli] = newCelli++;
1211 }
1212 }
1213 }
1214
1215 if (debug)
1216 {
1217 Pout<< "Cells : active:" << newCelli
1218 << " removed:" << cellMap_.size()-newCelli << endl;
1219 }
1220
1221 // Renumber -if cells reordered or -if cells removed
1222 if (orderCells || (newCelli != cellMap_.size()))
1223 {
1224 reorder(localCellMap, cellMap_);
1225 cellMap_.setCapacity(newCelli);
1226 renumberReverseMap(localCellMap, reverseCellMap_);
1227
1228 reorder(localCellMap, cellZone_);
1229 cellZone_.setCapacity(newCelli);
1230 if (cellAdditionalZones_.size())
1231 {
1232 reorder(localCellMap, cellAdditionalZones_);
1233 cellAdditionalZones_.setCapacity(newCelli);
1234 }
1235
1236 renumberKey(localCellMap, cellFromPoint_);
1237 renumberKey(localCellMap, cellFromEdge_);
1238 renumberKey(localCellMap, cellFromFace_);
1239
1240 // Renumber owner/neighbour. Take into account if neighbour suddenly
1241 // gets lower cell than owner.
1242 forAll(faceOwner_, facei)
1243 {
1244 label own = faceOwner_[facei];
1245 label nei = faceNeighbour_[facei];
1246
1247 if (own >= 0)
1248 {
1249 // Update owner
1250 faceOwner_[facei] = localCellMap[own];
1251
1252 if (nei >= 0)
1253 {
1254 // Update neighbour.
1255 faceNeighbour_[facei] = localCellMap[nei];
1256
1257 // Check if face needs reversing.
1258 if
1259 (
1260 faceNeighbour_[facei] >= 0
1261 && faceNeighbour_[facei] < faceOwner_[facei]
1262 )
1263 {
1264 faces_[facei].flip();
1265 std::swap(faceOwner_[facei], faceNeighbour_[facei]);
1266 flipFaceFlux_.flip(facei);
1267 faceZoneFlip_.flip(facei);
1268 if (facei < faceAdditionalZones_.size())
1269 {
1270 for (auto& zas : faceAdditionalZones_[facei])
1271 {
1272 // Flip sign
1273 zas = -zas;
1274 }
1275 }
1276 }
1277 }
1278 }
1279 else if (nei >= 0)
1280 {
1281 // Update neighbour.
1282 faceNeighbour_[facei] = localCellMap[nei];
1283 }
1284 }
1285 }
1286 }
1287
1288 // Reorder faces into upper-triangular and patch ordering
1289 {
1290 // Create cells (packed storage)
1291 labelList cellFaces;
1292 labelList cellFaceOffsets;
1293 makeCells(nActiveFaces_, cellFaces, cellFaceOffsets);
1294
1295 // Do upper triangular order and patch sorting
1296 labelList localFaceMap;
1297 getFaceOrder
1298 (
1299 nActiveFaces_,
1300 cellFaces,
1301 cellFaceOffsets,
1302
1303 localFaceMap,
1304 patchSizes,
1305 patchStarts
1306 );
1307
1308 // Reorder faces.
1309 reorderCompactFaces(localFaceMap.size(), localFaceMap);
1310 }
1311}
1312
1313
1314// Find faces to interpolate to create value for new face. Only used if
1315// face was inflated from edge or point. Internal faces should only be
1316// created from internal faces, external faces only from external faces
1317// (and ideally the same patch)
1318// Is bit problematic if there are no faces to select, i.e. in polyDualMesh
1319// an internal face can be created from a boundary edge with no internal
1320// faces connected to it.
1321Foam::labelList Foam::polyTopoChange::selectFaces
1322(
1323 const primitiveMesh& mesh,
1324 const labelUList& faceLabels,
1325 const bool internalFacesOnly
1326)
1327{
1328 label nFaces = 0;
1329
1330 forAll(faceLabels, i)
1331 {
1332 label facei = faceLabels[i];
1333
1334 if (internalFacesOnly == mesh.isInternalFace(facei))
1335 {
1336 nFaces++;
1337 }
1338 }
1339
1340 labelList collectedFaces;
1341
1342 if (nFaces == 0)
1343 {
1344 // Did not find any faces of the correct type so just use any old
1345 // face.
1346 collectedFaces = faceLabels;
1347 }
1348 else
1349 {
1350 collectedFaces.setSize(nFaces);
1351
1352 nFaces = 0;
1353
1354 forAll(faceLabels, i)
1355 {
1356 label facei = faceLabels[i];
1357
1358 if (internalFacesOnly == mesh.isInternalFace(facei))
1359 {
1360 collectedFaces[nFaces++] = facei;
1361 }
1362 }
1363 }
1364
1365 return collectedFaces;
1366}
1367
1368
1369// Calculate pointMap per patch (so from patch point label to old patch point
1370// label)
1371void Foam::polyTopoChange::calcPatchPointMap
1372(
1373 const UList<Map<label>>& oldPatchMeshPointMaps,
1374 const labelUList& patchMap,
1376 labelListList& patchPointMap
1377) const
1378{
1379 patchPointMap.setSize(patchMap.size());
1380
1381 forAll(patchMap, patchi)
1382 {
1383 const label oldPatchi = patchMap[patchi];
1384
1385 if (oldPatchi != -1)
1386 {
1387 const labelList& meshPoints = boundary[patchi].meshPoints();
1388
1389 const Map<label>& oldMeshPointMap =
1390 oldPatchMeshPointMaps[oldPatchi];
1391
1392 labelList& curPatchPointRnb = patchPointMap[patchi];
1393
1394 curPatchPointRnb.setSize(meshPoints.size());
1395
1396 forAll(meshPoints, i)
1397 {
1398 if (meshPoints[i] < pointMap_.size())
1399 {
1400 // Check if old point was part of same patch
1401 curPatchPointRnb[i] = oldMeshPointMap.lookup
1402 (
1403 pointMap_[meshPoints[i]],
1404 -1
1405 );
1406 }
1407 else
1408 {
1409 curPatchPointRnb[i] = -1;
1410 }
1411 }
1412 }
1413 }
1414}
1415
1416
1417void Foam::polyTopoChange::calcFaceInflationMaps
1418(
1419 const polyMesh& mesh,
1420 List<objectMap>& facesFromPoints,
1421 List<objectMap>& facesFromEdges,
1422 List<objectMap>& facesFromFaces
1423) const
1424{
1425 // Faces inflated from points
1426 // ~~~~~~~~~~~~~~~~~~~~~~~~~~
1427
1428 facesFromPoints.setSize(faceFromPoint_.size());
1429
1430 if (faceFromPoint_.size())
1431 {
1432 label nFacesFromPoints = 0;
1433
1434 // Collect all still existing faces connected to this point.
1435 forAllConstIters(faceFromPoint_, iter)
1436 {
1437 const label facei = iter.key();
1438 const label pointi = iter.val();
1439
1440 // Get internal or patch faces using point on old mesh
1441 const bool internal = (region_[facei] == -1);
1442
1443 facesFromPoints[nFacesFromPoints++] = objectMap
1444 (
1445 facei,
1446 selectFaces
1447 (
1448 mesh,
1449 mesh.pointFaces()[pointi],
1450 internal
1451 )
1452 );
1453 }
1454 }
1455
1456
1457 // Faces inflated from edges
1458 // ~~~~~~~~~~~~~~~~~~~~~~~~~
1459
1460 facesFromEdges.setSize(faceFromEdge_.size());
1461
1462 if (faceFromEdge_.size())
1463 {
1464 label nFacesFromEdges = 0;
1465
1466 // Collect all still existing faces connected to this edge.
1467 forAllConstIters(faceFromEdge_, iter)
1468 {
1469 const label facei = iter.key();
1470 const label edgei = iter.val();
1471
1472 // Get internal or patch faces using edge on old mesh
1473 const bool internal = (region_[facei] == -1);
1474
1475 facesFromEdges[nFacesFromEdges++] = objectMap
1476 (
1477 facei,
1478 selectFaces
1479 (
1480 mesh,
1481 mesh.edgeFaces(edgei),
1482 internal
1483 )
1484 );
1485 }
1486 }
1487
1488
1489 // Faces from face merging
1490 // ~~~~~~~~~~~~~~~~~~~~~~~
1491
1492 getMergeSets
1493 (
1494 reverseFaceMap_,
1495 faceMap_,
1496 facesFromFaces
1497 );
1498}
1499
1500
1501void Foam::polyTopoChange::calcCellInflationMaps
1502(
1503 const polyMesh& mesh,
1504 List<objectMap>& cellsFromPoints,
1505 List<objectMap>& cellsFromEdges,
1506 List<objectMap>& cellsFromFaces,
1507 List<objectMap>& cellsFromCells
1508) const
1509{
1510 cellsFromPoints.setSize(cellFromPoint_.size());
1511
1512 if (cellFromPoint_.size())
1513 {
1514 label nCellsFromPoints = 0;
1515
1516 // Collect all still existing cells connected to this point.
1517 forAllConstIters(cellFromPoint_, iter)
1518 {
1519 const label celli = iter.key();
1520 const label pointi = iter.val();
1521
1522 cellsFromPoints[nCellsFromPoints++] = objectMap
1523 (
1524 celli,
1525 mesh.pointCells()[pointi]
1526 );
1527 }
1528 }
1529
1530
1531 cellsFromEdges.setSize(cellFromEdge_.size());
1532
1533 if (cellFromEdge_.size())
1534 {
1535 label nCellsFromEdges = 0;
1536
1537 // Collect all still existing cells connected to this edge.
1538 forAllConstIters(cellFromEdge_, iter)
1539 {
1540 const label celli = iter.key();
1541 const label edgei = iter.val();
1542
1543 cellsFromEdges[nCellsFromEdges++] = objectMap
1544 (
1545 celli,
1546 mesh.edgeCells()[edgei]
1547 );
1548 }
1549 }
1550
1551
1552 cellsFromFaces.setSize(cellFromFace_.size());
1553
1554 if (cellFromFace_.size())
1555 {
1556 label nCellsFromFaces = 0;
1557
1558 labelList twoCells(2);
1559
1560 // Collect all still existing faces connected to this point.
1561 forAllConstIters(cellFromFace_, iter)
1562 {
1563 const label celli = iter.key();
1564 const label oldFacei = iter.val();
1565
1566 if (mesh.isInternalFace(oldFacei))
1567 {
1568 twoCells[0] = mesh.faceOwner()[oldFacei];
1569 twoCells[1] = mesh.faceNeighbour()[oldFacei];
1570 cellsFromFaces[nCellsFromFaces++] = objectMap
1571 (
1572 celli,
1573 twoCells
1574 );
1575 }
1576 else
1577 {
1578 cellsFromFaces[nCellsFromFaces++] = objectMap
1579 (
1580 celli,
1581 labelList(1, mesh.faceOwner()[oldFacei])
1582 );
1583 }
1584 }
1585 }
1586
1587
1588 // Cells from cell merging
1589 // ~~~~~~~~~~~~~~~~~~~~~~~
1590
1591 getMergeSets
1592 (
1593 reverseCellMap_,
1594 cellMap_,
1595 cellsFromCells
1596 );
1597}
1598
1599
1600void Foam::polyTopoChange::resetZones
1601(
1602 const polyMesh& mesh,
1603 polyMesh& newMesh,
1604 labelListList& pointZoneMap,
1605 labelListList& faceZoneFaceMap,
1606 labelListList& cellZoneMap
1607) const
1608{
1609 // pointZones
1610 // ~~~~~~~~~~
1611
1612 pointZoneMap.setSize(mesh.pointZones().size());
1613 {
1614 const pointZoneMesh& pointZones = mesh.pointZones();
1615
1616 // Count points per zone
1617
1618 labelList nPoints(pointZones.size(), Zero);
1619
1620 forAllConstIters(pointZone_, iter)
1621 {
1622 const label pointi = iter.key();
1623 const label zonei = iter.val();
1624
1625 if (zonei < 0 || zonei >= pointZones.size())
1626 {
1628 << "Illegal zoneID " << zonei << " for point "
1629 << pointi << " coord " << mesh.points()[pointi]
1630 << abort(FatalError);
1631 }
1632 nPoints[zonei]++;
1633 }
1634 forAll(pointAdditionalZones_, pointi)
1635 {
1636 for (const label zonei : pointAdditionalZones_[pointi])
1637 {
1638 if (zonei < 0 || zonei >= pointZones.size())
1639 {
1641 << "Illegal zoneID " << zonei << " for point "
1642 << pointi << abort(FatalError);
1643 }
1644 nPoints[zonei]++;
1645 }
1646 }
1647
1648 // Distribute points per zone
1649
1650 labelListList addressing(pointZones.size());
1651 forAll(addressing, zonei)
1652 {
1653 addressing[zonei].setSize(nPoints[zonei]);
1654 }
1655 nPoints = 0;
1656
1657 forAllConstIters(pointZone_, iter)
1658 {
1659 const label pointi = iter.key();
1660 const label zonei = iter.val();
1661
1662 addressing[zonei][nPoints[zonei]++] = pointi;
1663 }
1664 forAll(pointAdditionalZones_, pointi)
1665 {
1666 for (const label zonei : pointAdditionalZones_[pointi])
1667 {
1668 addressing[zonei][nPoints[zonei]++] = pointi;
1669 }
1670 }
1671 // Sort the addressing
1672 forAll(addressing, zonei)
1673 {
1674 stableSort(addressing[zonei]);
1675 }
1676
1677 // So now we both have old zones and the new addressing.
1678 // Invert the addressing to get pointZoneMap.
1679 forAll(addressing, zonei)
1680 {
1681 const pointZone& oldZone = pointZones[zonei];
1682 const labelList& newZoneAddr = addressing[zonei];
1683
1684 labelList& curPzRnb = pointZoneMap[zonei];
1685 curPzRnb.setSize(newZoneAddr.size());
1686
1687 forAll(newZoneAddr, i)
1688 {
1689 if (newZoneAddr[i] < pointMap_.size())
1690 {
1691 curPzRnb[i] = oldZone.whichPoint(pointMap_[newZoneAddr[i]]);
1692 }
1693 else
1694 {
1695 curPzRnb[i] = -1;
1696 }
1697 }
1698 }
1699
1700 // Reset the addressing on the zone
1701 newMesh.pointZones().clearAddressing();
1702 forAll(newMesh.pointZones(), zonei)
1703 {
1704 if (debug)
1705 {
1706 Pout<< "pointZone:" << zonei
1707 << " name:" << newMesh.pointZones()[zonei].name()
1708 << " size:" << addressing[zonei].size()
1709 << endl;
1710 }
1711
1712 newMesh.pointZones()[zonei] = addressing[zonei];
1713 }
1714 }
1715
1716
1717 // faceZones
1718 // ~~~~~~~~~
1719
1720 faceZoneFaceMap.setSize(mesh.faceZones().size());
1721 {
1722 const faceZoneMesh& faceZones = mesh.faceZones();
1723
1724 labelList nFaces(faceZones.size(), Zero);
1725
1726 forAllConstIters(faceZone_, iter)
1727 {
1728 const label facei = iter.key();
1729 const label zonei = iter.val();
1730
1731 if (zonei < 0 || zonei >= faceZones.size())
1732 {
1734 << "Illegal zoneID " << zonei << " for face "
1735 << facei
1736 << abort(FatalError);
1737 }
1738 nFaces[zonei]++;
1739 }
1740 forAll(faceAdditionalZones_, facei)
1741 {
1742 for (const label zoneAndSign : faceAdditionalZones_[facei])
1743 {
1744 const label zonei = mag(zoneAndSign)-1;
1745 if (zonei < 0 || zonei >= faceZones.size())
1746 {
1748 << "Illegal zoneID " << zonei << " for face "
1749 << facei << abort(FatalError);
1750 }
1751 nFaces[zonei]++;
1752 }
1753 }
1754
1755 labelListList addressing(faceZones.size());
1756 boolListList flipMode(faceZones.size());
1757
1758 forAll(addressing, zonei)
1759 {
1760 addressing[zonei].setSize(nFaces[zonei]);
1761 flipMode[zonei].setSize(nFaces[zonei]);
1762 }
1763 nFaces = 0;
1764
1765 forAllConstIters(faceZone_, iter)
1766 {
1767 const label facei = iter.key();
1768 const label zonei = iter.val();
1769
1770 const label index = nFaces[zonei]++;
1771 addressing[zonei][index] = facei;
1772 flipMode[zonei][index] = faceZoneFlip_.test(facei);
1773
1774 if (facei < faceAdditionalZones_.size())
1775 {
1776 for (const label zoneAndSign : faceAdditionalZones_[facei])
1777 {
1778 const label zonei = mag(zoneAndSign)-1;
1779 const bool flip = (zoneAndSign > 0);
1780
1781 const label index = nFaces[zonei]++;
1782 addressing[zonei][index] = facei;
1783 flipMode[zonei][index] = flip;
1784 }
1785 }
1786 }
1787
1788 // Sort the addressing
1789 forAll(addressing, zonei)
1790 {
1791 labelList newToOld(sortedOrder(addressing[zonei]));
1792 {
1793 labelList newAddressing(addressing[zonei].size());
1794 forAll(newAddressing, i)
1795 {
1796 newAddressing[i] = addressing[zonei][newToOld[i]];
1797 }
1798 addressing[zonei].transfer(newAddressing);
1799 }
1800 {
1801 boolList newFlipMode(flipMode[zonei].size());
1802 forAll(newFlipMode, i)
1803 {
1804 newFlipMode[i] = flipMode[zonei][newToOld[i]];
1805 }
1806 flipMode[zonei].transfer(newFlipMode);
1807 }
1808 }
1809
1810 // So now we both have old zones and the new addressing.
1811 // Invert the addressing to get faceZoneFaceMap.
1812 forAll(addressing, zonei)
1813 {
1814 const faceZone& oldZone = faceZones[zonei];
1815 const labelList& newZoneAddr = addressing[zonei];
1816
1817 labelList& curFzFaceRnb = faceZoneFaceMap[zonei];
1818
1819 curFzFaceRnb.setSize(newZoneAddr.size());
1820
1821 forAll(newZoneAddr, i)
1822 {
1823 if (newZoneAddr[i] < faceMap_.size())
1824 {
1825 curFzFaceRnb[i] =
1826 oldZone.whichFace(faceMap_[newZoneAddr[i]]);
1827 }
1828 else
1829 {
1830 curFzFaceRnb[i] = -1;
1831 }
1832 }
1833 }
1834
1835
1836 // Reset the addressing on the zone
1837 newMesh.faceZones().clearAddressing();
1838 forAll(newMesh.faceZones(), zonei)
1839 {
1840 if (debug)
1841 {
1842 Pout<< "faceZone:" << zonei
1843 << " name:" << newMesh.faceZones()[zonei].name()
1844 << " size:" << addressing[zonei].size()
1845 << endl;
1846 }
1847
1848 newMesh.faceZones()[zonei].resetAddressing
1849 (
1850 addressing[zonei],
1851 flipMode[zonei]
1852 );
1853 }
1854 }
1855
1856
1857 // cellZones
1858 // ~~~~~~~~~
1859
1860 cellZoneMap.setSize(mesh.cellZones().size());
1861 {
1862 const cellZoneMesh& cellZones = mesh.cellZones();
1863
1864 labelList nCells(cellZones.size(), Zero);
1865
1866 forAll(cellZone_, celli)
1867 {
1868 const label zonei = cellZone_[celli];
1869
1870 if (zonei >= cellZones.size())
1871 {
1873 << "Illegal zoneID " << zonei << " for cell "
1874 << celli << abort(FatalError);
1875 }
1876
1877 if (zonei >= 0)
1878 {
1879 nCells[zonei]++;
1880 }
1881 }
1882 forAll(cellAdditionalZones_, celli)
1883 {
1884 for (const label zonei : cellAdditionalZones_[celli])
1885 {
1886 if (zonei < 0 || zonei >= cellZones.size())
1887 {
1889 << "Illegal zoneID " << zonei << " for cell "
1890 << celli << abort(FatalError);
1891 }
1892 nCells[zonei]++;
1893 }
1894 }
1895
1896
1897 labelListList addressing(cellZones.size());
1898 forAll(addressing, zonei)
1899 {
1900 addressing[zonei].setSize(nCells[zonei]);
1901 }
1902 nCells = 0;
1903
1904 forAll(cellZone_, celli)
1905 {
1906 const label zonei = cellZone_[celli];
1907
1908 if (zonei >= 0)
1909 {
1910 addressing[zonei][nCells[zonei]++] = celli;
1911 }
1912 }
1913 forAll(cellAdditionalZones_, celli)
1914 {
1915 for (const label zonei : cellAdditionalZones_[celli])
1916 {
1917 addressing[zonei][nCells[zonei]++] = celli;
1918 }
1919 }
1920 // Sort the addressing
1921 forAll(addressing, zonei)
1922 {
1923 stableSort(addressing[zonei]);
1924 }
1925
1926 // So now we both have old zones and the new addressing.
1927 // Invert the addressing to get cellZoneMap.
1928 forAll(addressing, zonei)
1929 {
1930 const cellZone& oldZone = cellZones[zonei];
1931 const labelList& newZoneAddr = addressing[zonei];
1932
1933 labelList& curCellRnb = cellZoneMap[zonei];
1934
1935 curCellRnb.setSize(newZoneAddr.size());
1936
1937 forAll(newZoneAddr, i)
1938 {
1939 if (newZoneAddr[i] < cellMap_.size())
1940 {
1941 curCellRnb[i] =
1942 oldZone.whichCell(cellMap_[newZoneAddr[i]]);
1943 }
1944 else
1945 {
1946 curCellRnb[i] = -1;
1947 }
1948 }
1949 }
1950
1951 // Reset the addressing on the zone
1952 newMesh.cellZones().clearAddressing();
1953 forAll(newMesh.cellZones(), zonei)
1954 {
1955 if (debug)
1956 {
1957 Pout<< "cellZone:" << zonei
1958 << " name:" << newMesh.cellZones()[zonei].name()
1959 << " size:" << addressing[zonei].size()
1960 << endl;
1961 }
1962
1963 newMesh.cellZones()[zonei] = addressing[zonei];
1964 }
1965 }
1966}
1967
1968
1969void Foam::polyTopoChange::calcFaceZonePointMap
1970(
1971 const polyMesh& mesh,
1972 const UList<Map<label>>& oldFaceZoneMeshPointMaps,
1973 labelListList& faceZonePointMap
1974) const
1975{
1976 const faceZoneMesh& faceZones = mesh.faceZones();
1977
1978 faceZonePointMap.setSize(faceZones.size());
1979
1980 forAll(faceZones, zonei)
1981 {
1982 const faceZone& newZone = faceZones[zonei];
1983
1984 const labelList& newZoneMeshPoints = newZone.patch().meshPoints();
1985
1986 const Map<label>& oldZoneMeshPointMap = oldFaceZoneMeshPointMaps[zonei];
1987
1988 labelList& curFzPointRnb = faceZonePointMap[zonei];
1989
1990 curFzPointRnb.setSize(newZoneMeshPoints.size());
1991
1992 forAll(newZoneMeshPoints, pointi)
1993 {
1994 if (newZoneMeshPoints[pointi] < pointMap_.size())
1995 {
1996 curFzPointRnb[pointi] =
1997 oldZoneMeshPointMap.lookup
1998 (
1999 pointMap_[newZoneMeshPoints[pointi]],
2000 -1
2001 );
2002 }
2003 else
2004 {
2005 curFzPointRnb[pointi] = -1;
2006 }
2007 }
2008 }
2009}
2010
2011
2012void Foam::polyTopoChange::reorderCoupledFaces
2013(
2014 const bool syncParallel,
2015 const polyBoundaryMesh& oldBoundary,
2016 const labelUList& patchMap, // new to old patches
2017 const labelUList& patchStarts,
2018 const labelUList& patchSizes,
2019 const pointField& points
2020)
2021{
2022 // Mapping for faces (old to new). Extends over all mesh faces for
2023 // convenience (could be just the external faces)
2024 labelList oldToNew(identity(faces_.size()));
2025
2026 // Rotation on new faces.
2027 labelList rotation(faces_.size(), Zero);
2028
2029 // Allocate unique tag for all comms
2030 const int oldTag = UPstream::incrMsgType();
2031
2032 PstreamBuffers pBufs;
2033
2034 // Send ordering
2035 forAll(patchMap, patchi)
2036 {
2037 const label oldPatchi = patchMap[patchi];
2038
2039 if
2040 (
2041 oldPatchi != -1
2042 && (syncParallel || !isA<processorPolyPatch>(oldBoundary[oldPatchi]))
2043 )
2044 {
2045 oldBoundary[oldPatchi].initOrder
2046 (
2047 pBufs,
2049 (
2051 (
2052 faces_,
2053 patchSizes[patchi],
2054 patchStarts[patchi]
2055 ),
2056 points
2057 )
2058 );
2059 }
2060 }
2061
2062 if (syncParallel)
2063 {
2064 pBufs.finishedSends();
2065 }
2066
2067 // Receive and calculate ordering
2068
2069 bool anyChanged = false;
2070
2071 forAll(patchMap, patchi)
2072 {
2073 const label oldPatchi = patchMap[patchi];
2074
2075 if
2076 (
2077 oldPatchi != -1
2078 && (syncParallel || !isA<processorPolyPatch>(oldBoundary[oldPatchi]))
2079 )
2080 {
2081 labelList patchFaceMap(patchSizes[patchi], -1);
2082 labelList patchFaceRotation(patchSizes[patchi], Zero);
2083
2084 const bool changed = oldBoundary[oldPatchi].order
2085 (
2086 pBufs,
2088 (
2090 (
2091 faces_,
2092 patchSizes[patchi],
2093 patchStarts[patchi]
2094 ),
2095 points
2096 ),
2097 patchFaceMap,
2098 patchFaceRotation
2099 );
2100
2101 if (changed)
2102 {
2103 // Merge patch face reordering into mesh face reordering table
2104 const label start = patchStarts[patchi];
2105
2106 forAll(patchFaceMap, patchFacei)
2107 {
2108 oldToNew[patchFacei + start] =
2109 start + patchFaceMap[patchFacei];
2110 }
2111
2112 forAll(patchFaceRotation, patchFacei)
2113 {
2114 rotation[patchFacei + start] =
2115 patchFaceRotation[patchFacei];
2116 }
2117
2118 anyChanged = true;
2119 }
2120 }
2121 }
2122
2123 if (syncParallel ? returnReduceOr(anyChanged) : anyChanged)
2124 {
2125 // Reorder faces according to oldToNew.
2126 reorderCompactFaces(oldToNew.size(), oldToNew);
2127
2128 // Rotate faces (rotation is already in new face indices).
2129 forAll(rotation, facei)
2130 {
2131 if (rotation[facei] != 0)
2132 {
2133 inplaceRotateList<List, label>(faces_[facei], rotation[facei]);
2134 }
2135 }
2136 }
2137
2138 // Reset tag
2139 UPstream::msgType(oldTag);
2140}
2141
2142
2143void Foam::polyTopoChange::compactAndReorder
2144(
2145 const polyMesh& mesh,
2146 const labelUList& patchMap, // from new to old patch
2147 const bool syncParallel,
2148 const bool orderCells,
2149 const bool orderPoints,
2150
2151 label& nInternalPoints,
2152 pointField& newPoints,
2153 labelList& patchSizes,
2154 labelList& patchStarts,
2155 List<objectMap>& pointsFromPoints,
2156 List<objectMap>& facesFromPoints,
2157 List<objectMap>& facesFromEdges,
2158 List<objectMap>& facesFromFaces,
2159 List<objectMap>& cellsFromPoints,
2160 List<objectMap>& cellsFromEdges,
2161 List<objectMap>& cellsFromFaces,
2162 List<objectMap>& cellsFromCells,
2163 List<Map<label>>& oldPatchMeshPointMaps,
2164 labelList& oldPatchNMeshPoints,
2165 labelList& oldPatchStarts,
2166 List<Map<label>>& oldFaceZoneMeshPointMaps
2167)
2168{
2169 if (patchMap.size() != nPatches_)
2170 {
2172 << "polyTopoChange was constructed with a mesh with "
2173 << nPatches_ << " patches." << endl
2174 << "The mesh now provided has a different number of patches "
2175 << patchMap.size()
2176 << " which is illegal" << endl
2177 << abort(FatalError);
2178 }
2179
2180 // Remove any holes from points/faces/cells and sort faces.
2181 // Sets nActiveFaces_.
2182 compact(orderCells, orderPoints, nInternalPoints, patchSizes, patchStarts);
2183
2184 // Transfer points to pointField. points_ are now cleared!
2185 // Only done since e.g. reorderCoupledFaces requires pointField.
2186 newPoints.transfer(points_);
2187
2188 // Reorder any coupled faces
2189 reorderCoupledFaces
2190 (
2191 syncParallel,
2193 patchMap,
2194 patchStarts,
2195 patchSizes,
2196 newPoints
2197 );
2198
2199
2200 // Calculate inflation/merging maps
2201 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2202 // These are for the new face(/point/cell) the old faces whose value
2203 // needs to be
2204 // averaged/summed to get the new value. These old faces come either from
2205 // merged old faces (face remove into other face),
2206 // the old edgeFaces (inflate from edge) or the old pointFaces (inflate
2207 // from point). As an additional complexity will use only internal faces
2208 // to create new value for internal face and vice versa only patch
2209 // faces to to create patch face value.
2210
2211 // For point only point merging
2212 getMergeSets
2213 (
2214 reversePointMap_,
2215 pointMap_,
2216 pointsFromPoints
2217 );
2218
2219 calcFaceInflationMaps
2220 (
2221 mesh,
2222 facesFromPoints,
2223 facesFromEdges,
2224 facesFromFaces
2225 );
2226
2227 calcCellInflationMaps
2228 (
2229 mesh,
2230 cellsFromPoints,
2231 cellsFromEdges,
2232 cellsFromFaces,
2233 cellsFromCells
2234 );
2235
2236 // Clear inflation info
2237 {
2238 faceFromPoint_.clearStorage();
2239 faceFromEdge_.clearStorage();
2240
2241 cellFromPoint_.clearStorage();
2242 cellFromEdge_.clearStorage();
2243 cellFromFace_.clearStorage();
2244 }
2245
2246
2248
2249 // Grab patch mesh point maps
2250 oldPatchMeshPointMaps.setSize(boundary.size());
2251 oldPatchNMeshPoints.setSize(boundary.size());
2252 oldPatchStarts.setSize(boundary.size());
2253
2254 forAll(boundary, patchi)
2255 {
2256 // Copy old face zone mesh point maps
2257 oldPatchMeshPointMaps[patchi] = boundary[patchi].meshPointMap();
2258 oldPatchNMeshPoints[patchi] = boundary[patchi].meshPoints().size();
2259 oldPatchStarts[patchi] = boundary[patchi].start();
2260 }
2261
2262 // Grab old face zone mesh point maps.
2263 // These need to be saved before resetting the mesh and are used
2264 // later on to calculate the faceZone pointMaps.
2265 oldFaceZoneMeshPointMaps.setSize(mesh.faceZones().size());
2266
2267 forAll(mesh.faceZones(), zonei)
2268 {
2269 const faceZone& oldZone = mesh.faceZones()[zonei];
2270
2271 oldFaceZoneMeshPointMaps[zonei] = oldZone.patch().meshPointMap();
2272 }
2273}
2274
2275
2276// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
2277
2278// Construct from components
2279Foam::polyTopoChange::polyTopoChange(const label nPatches, const bool strict)
2280:
2281 strict_(strict),
2282 nPatches_(nPatches),
2283 points_(0),
2284 pointMap_(0),
2285 reversePointMap_(0),
2286 pointZone_(0),
2287 retiredPoints_(0),
2288 faces_(0),
2289 region_(0),
2290 faceOwner_(0),
2291 faceNeighbour_(0),
2292 faceMap_(0),
2293 reverseFaceMap_(0),
2294 faceFromPoint_(0),
2295 faceFromEdge_(0),
2296 flipFaceFlux_(0),
2297 faceZone_(0),
2298 faceZoneFlip_(0),
2299 nActiveFaces_(0),
2300 cellMap_(0),
2301 reverseCellMap_(0),
2302 cellFromPoint_(0),
2303 cellFromEdge_(0),
2304 cellFromFace_(0),
2305 cellZone_(0)
2306{}
2307
2308
2309// Construct from components
2311(
2312 const polyMesh& mesh,
2313 const bool strict
2314)
2315:
2316 strict_(strict),
2317 nPatches_(0),
2318 points_(0),
2319 pointMap_(0),
2320 reversePointMap_(0),
2321 pointZone_(0),
2322 retiredPoints_(0),
2323 faces_(0),
2324 region_(0),
2325 faceOwner_(0),
2326 faceNeighbour_(0),
2327 faceMap_(0),
2328 reverseFaceMap_(0),
2329 faceFromPoint_(0),
2330 faceFromEdge_(0),
2331 flipFaceFlux_(0),
2332 faceZone_(0),
2333 faceZoneFlip_(0),
2334 nActiveFaces_(0),
2335 cellMap_(0),
2336 reverseCellMap_(0),
2337 cellFromPoint_(0),
2338 cellFromEdge_(0),
2339 cellFromFace_(0),
2340 cellZone_(0)
2341{
2342 addMesh
2343 (
2344 mesh,
2345 identity(mesh.boundaryMesh().size()),
2346 identity(mesh.pointZones().size()),
2347 identity(mesh.faceZones().size()),
2348 identity(mesh.cellZones().size())
2349 );
2350}
2351
2352
2353// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
2356{
2357 points_.clearStorage();
2358 pointMap_.clearStorage();
2359 reversePointMap_.clearStorage();
2360 pointZone_.clearStorage();
2361 pointAdditionalZones_.clearStorage();
2362 retiredPoints_.clearStorage();
2363
2364 faces_.clearStorage();
2365 region_.clearStorage();
2366 faceOwner_.clearStorage();
2367 faceNeighbour_.clearStorage();
2368 faceMap_.clearStorage();
2369 reverseFaceMap_.clearStorage();
2370 faceFromPoint_.clearStorage();
2371 faceFromEdge_.clearStorage();
2372 flipFaceFlux_.clearStorage();
2373 faceZone_.clearStorage();
2374 faceZoneFlip_.clearStorage();
2375 faceAdditionalZones_.clearStorage();
2376 nActiveFaces_ = 0;
2377
2378 cellMap_.clearStorage();
2379 reverseCellMap_.clearStorage();
2380 cellZone_.clearStorage();
2381 cellAdditionalZones_.clearStorage();
2382 cellFromPoint_.clearStorage();
2383 cellFromEdge_.clearStorage();
2384 cellFromFace_.clearStorage();
2385}
2386
2389(
2390 const polyMesh& mesh,
2391 const labelUList& patchMap,
2392 const labelUList& pointZoneMap,
2393 const labelUList& faceZoneMap,
2394 const labelUList& cellZoneMap
2395)
2396{
2397 label maxRegion = nPatches_ - 1;
2398 for (const label regioni : patchMap)
2399 {
2400 maxRegion = max(maxRegion, regioni);
2401 }
2402 nPatches_ = maxRegion + 1;
2403
2404
2405 // Add points
2406 {
2407 const pointField& points = mesh.points();
2408 const pointZoneMesh& pointZones = mesh.pointZones();
2409
2410 // Extend
2411 points_.setCapacity(points_.size() + points.size());
2412 pointMap_.setCapacity(pointMap_.size() + points.size());
2413 reversePointMap_.setCapacity(reversePointMap_.size() + points.size());
2414 pointZone_.resize(pointZone_.size() + points.size()/128);
2415
2416 // Add points without pointZone
2417 labelList newPointID(points.size());
2418 forAll(newPointID, pointi)
2419 {
2420 // Add point from point
2421 newPointID[pointi] = addPoint
2422 (
2423 points[pointi],
2424 pointi,
2425 -1,
2426 true
2427 );
2428 }
2429
2430 // Modify for pointZones
2431 DynamicList<label> zones;
2432 forAll(pointZones, zonei)
2433 {
2434 for (const label pointi : pointZones[zonei])
2435 {
2436 // Get previous zones
2437 this->pointZones(newPointID[pointi], zones);
2438 for (auto& zonei : zones)
2439 {
2440 zonei = pointZoneMap[zonei];
2441 }
2442 // Add new zone
2443 zones.append(pointZoneMap[zonei]);
2445 (
2446 newPointID[pointi],
2447 points_[newPointID[pointi]],
2448 zones,
2449 true
2450 );
2451 }
2452 }
2453 }
2454
2455 // Add cells
2456 {
2457 const cellZoneMesh& cellZones = mesh.cellZones();
2458
2459 // Resize
2460
2461 // Note: polyMesh does not allow retired cells anymore. So allCells
2462 // always equals nCells
2463 label nAllCells = mesh.nCells();
2464
2465 cellMap_.setCapacity(cellMap_.size() + nAllCells);
2466 reverseCellMap_.setCapacity(reverseCellMap_.size() + nAllCells);
2467 cellFromPoint_.resize(cellFromPoint_.size() + nAllCells/128);
2468 cellFromEdge_.resize(cellFromEdge_.size() + nAllCells/128);
2469 cellFromFace_.resize(cellFromFace_.size() + nAllCells/128);
2470 cellZone_.setCapacity(cellZone_.size() + nAllCells);
2471
2472
2473 // Add cells without cellZone
2474 labelList newCellID(nAllCells);
2475 for (label celli = 0; celli < nAllCells; celli++)
2476 {
2477 // Add cell from cell
2478 newCellID[celli] = addCell(-1, -1, -1, celli, -1);
2479 }
2480
2481 // Modify for cellZones
2482 DynamicList<label> zones;
2483 forAll(cellZones, zonei)
2484 {
2485 for (const label celli : cellZones[zonei])
2486 {
2487 // Get previous zones
2488 this->cellZones(newCellID[celli], zones);
2489 for (auto& zonei : zones)
2490 {
2491 zonei = cellZoneMap[zonei];
2492 }
2493 // Add new zone
2494 zones.append(cellZoneMap[zonei]);
2495 modifyCell(newCellID[celli], zones);
2496 }
2497 }
2498 }
2499
2500 // Add faces
2501 {
2502 const polyBoundaryMesh& patches = mesh.boundaryMesh();
2503 const faceList& faces = mesh.faces();
2504 const labelList& faceOwner = mesh.faceOwner();
2505 const labelList& faceNeighbour = mesh.faceNeighbour();
2506 const faceZoneMesh& faceZones = mesh.faceZones();
2507
2508 // Resize
2509 label nAllFaces = mesh.faces().size();
2510
2511 faces_.setCapacity(faces_.size() + nAllFaces);
2512 region_.setCapacity(region_.size() + nAllFaces);
2513 faceOwner_.setCapacity(faceOwner_.size() + nAllFaces);
2514 faceNeighbour_.setCapacity(faceNeighbour_.size() + nAllFaces);
2515 faceMap_.setCapacity(faceMap_.size() + nAllFaces);
2516 reverseFaceMap_.setCapacity(reverseFaceMap_.size() + nAllFaces);
2517 faceFromPoint_.resize(faceFromPoint_.size() + nAllFaces/128);
2518 faceFromEdge_.resize(faceFromEdge_.size() + nAllFaces/128);
2519 flipFaceFlux_.setCapacity(faces_.size() + nAllFaces);
2520 faceZone_.resize(faceZone_.size() + nAllFaces/128);
2521 faceZoneFlip_.setCapacity(faceMap_.capacity());
2522
2523
2524 // Add faces (in mesh order) without faceZone
2525 labelList newFaceID(nAllFaces);
2526
2527 // 1. Internal faces
2528 for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
2529 {
2530 newFaceID[facei] = addFace
2531 (
2532 faces[facei],
2533 faceOwner[facei],
2534 faceNeighbour[facei],
2535 -1, // masterPointID
2536 -1, // masterEdgeID
2537 facei, // masterFaceID
2538 false, // flipFaceFlux
2539 -1, // patchID
2540 -1, // zoneID
2541 false // zoneFlip
2542 );
2543 }
2544
2545 // 2. Patch faces
2546 forAll(patches, patchi)
2547 {
2548 const polyPatch& pp = patches[patchi];
2549
2550 if (pp.start() != faces_.size())
2551 {
2553 << "Problem : "
2554 << "Patch " << pp.name() << " starts at " << pp.start()
2555 << endl
2556 << "Current face counter at " << faces_.size() << endl
2557 << "Are patches in incremental order?"
2558 << abort(FatalError);
2559 }
2560 forAll(pp, patchFacei)
2561 {
2562 const label facei = pp.start() + patchFacei;
2563
2564 newFaceID[facei] = addFace
2565 (
2566 faces[facei],
2567 faceOwner[facei],
2568 -1, // neighbour
2569 -1, // masterPointID
2570 -1, // masterEdgeID
2571 facei, // masterFaceID
2572 false, // flipFaceFlux
2573 patchMap[patchi], // patchID
2574 -1, // zoneID
2575 false // zoneFlip
2576 );
2577 }
2578 }
2579
2580
2581 // Modify for faceZones
2582
2583 DynamicList<label> zones;
2584 DynamicList<bool> flips;
2585 forAll(faceZones, zonei)
2586 {
2587 const auto& fz = faceZones[zonei];
2588
2589 forAll(fz, i)
2590 {
2591 const label facei = fz[i];
2592 const bool flip = fz.flipMap()[i];
2593
2594 // Modify face index for combined zones
2595 const label newFacei = newFaceID[facei];
2596
2597 // Get previous zones
2598 this->faceZones(newFacei, zones, flips);
2599 for (auto& zonei : zones)
2600 {
2601 zonei = faceZoneMap[zonei];
2602 }
2603 // Add new zone
2604 zones.append(faceZoneMap[zonei]);
2605 flips.append(flip);
2606
2608 (
2609 faces_[newFacei],
2610 newFacei,
2611 faceOwner_[newFacei],
2612 faceNeighbour_[newFacei],
2613 flipFaceFlux_(newFacei),
2614 region_[newFacei],
2615 zones,
2616 flips
2617 );
2618 }
2619 }
2620 }
2621}
2622
2625(
2626 const label nPoints,
2627 const label nFaces,
2628 const label nCells
2629)
2630{
2631 points_.setCapacity(nPoints);
2632 pointMap_.setCapacity(nPoints);
2633 reversePointMap_.setCapacity(nPoints);
2634 pointZone_.resize(pointZone_.size() + nPoints/128);
2635
2636 faces_.setCapacity(nFaces);
2637 region_.setCapacity(nFaces);
2638 faceOwner_.setCapacity(nFaces);
2639 faceNeighbour_.setCapacity(nFaces);
2640 faceMap_.setCapacity(nFaces);
2641 reverseFaceMap_.setCapacity(nFaces);
2642 faceFromPoint_.resize(faceFromPoint_.size() + nFaces/128);
2643 faceFromEdge_.resize(faceFromEdge_.size() + nFaces/128);
2644 flipFaceFlux_.setCapacity(nFaces);
2645 faceZone_.resize(faceZone_.size() + nFaces/128);
2646 faceZoneFlip_.setCapacity(nFaces);
2647
2648 cellMap_.setCapacity(nCells);
2649 reverseCellMap_.setCapacity(nCells);
2650 cellFromPoint_.resize(cellFromPoint_.size() + nCells/128);
2651 cellFromEdge_.resize(cellFromEdge_.size() + nCells/128);
2652 cellFromFace_.resize(cellFromFace_.size() + nCells/128);
2653 cellZone_.setCapacity(nCells);
2654}
2655
2657Foam::label Foam::polyTopoChange::setAction(const topoAction& action)
2658{
2659 if (isType<polyAddPoint>(action))
2660 {
2661 const polyAddPoint& pap = refCast<const polyAddPoint>(action);
2662
2663 return addPoint
2664 (
2665 pap.newPoint(),
2666 pap.masterPointID(),
2667 pap.zoneID(),
2668 pap.inCell()
2669 );
2670 }
2671 else if (isType<polyModifyPoint>(action))
2672 {
2674
2676 (
2677 pmp.pointID(),
2678 pmp.newPoint(),
2679 pmp.zoneID(),
2680 pmp.inCell()
2681 );
2682
2683 return -1;
2684 }
2685 else if (isType<polyRemovePoint>(action))
2686 {
2688
2689 removePoint(prp.pointID(), prp.mergePointID());
2690
2691 return -1;
2692 }
2693 else if (isType<polyAddFace>(action))
2694 {
2695 const polyAddFace& paf = refCast<const polyAddFace>(action);
2696
2697 return addFace
2698 (
2699 paf.newFace(),
2700 paf.owner(),
2701 paf.neighbour(),
2702 paf.masterPointID(),
2703 paf.masterEdgeID(),
2704 paf.masterFaceID(),
2705 paf.flipFaceFlux(),
2706 paf.patchID(),
2707 paf.zoneID(),
2708 paf.zoneFlip()
2709 );
2710 }
2711 else if (isType<polyModifyFace>(action))
2712 {
2713 const polyModifyFace& pmf = refCast<const polyModifyFace>(action);
2714
2716 (
2717 pmf.newFace(),
2718 pmf.faceID(),
2719 pmf.owner(),
2720 pmf.neighbour(),
2721 pmf.flipFaceFlux(),
2722 pmf.patchID(),
2723 pmf.zoneID(),
2724 pmf.zoneFlip()
2725 );
2726
2727 return -1;
2728 }
2729 else if (isType<polyRemoveFace>(action))
2730 {
2731 const polyRemoveFace& prf = refCast<const polyRemoveFace>(action);
2732
2733 removeFace(prf.faceID(), prf.mergeFaceID());
2734
2735 return -1;
2736 }
2737 else if (isType<polyAddCell>(action))
2738 {
2739 const polyAddCell& pac = refCast<const polyAddCell>(action);
2740
2741 return addCell
2742 (
2743 pac.masterPointID(),
2744 pac.masterEdgeID(),
2745 pac.masterFaceID(),
2746 pac.masterCellID(),
2747 pac.zoneID()
2748 );
2749 }
2750 else if (isType<polyModifyCell>(action))
2751 {
2752 const polyModifyCell& pmc = refCast<const polyModifyCell>(action);
2753
2754 if (pmc.removeFromZone())
2755 {
2756 modifyCell(pmc.cellID(), -1);
2757 }
2758 else
2759 {
2760 modifyCell(pmc.cellID(), pmc.zoneID());
2761 }
2762
2763 return -1;
2764 }
2765 else if (isType<polyRemoveCell>(action))
2766 {
2767 const polyRemoveCell& prc = refCast<const polyRemoveCell>(action);
2768
2769 removeCell(prc.cellID(), prc.mergeCellID());
2770
2771 return -1;
2772 }
2773 else
2774 {
2776 << "Unknown type of topoChange: " << action.type()
2777 << abort(FatalError);
2778
2779 // Dummy return to keep compiler happy
2780 return -1;
2781 }
2782}
2783
2786(
2787 const point& pt,
2788 const label masterPointID,
2789 const label zoneID,
2790 const bool inCell
2791)
2792{
2793 const label pointi = points_.size();
2794
2795 points_.append(pt);
2796 pointMap_.append(masterPointID);
2797 reversePointMap_.append(pointi);
2798
2799 if (zoneID >= 0)
2800 {
2801 pointZone_.insert(pointi, zoneID);
2802 }
2803 // Delay extending the pointAdditionalZones
2804
2805 if (!inCell)
2806 {
2807 retiredPoints_.insert(pointi);
2808 }
2809
2810 return pointi;
2811}
2812
2815(
2816 const point& pt,
2817 const label masterPointID,
2818 const labelUList& zoneIDs,
2819 const bool inCell
2820)
2821{
2822 const label pointi = points_.size();
2823
2824 points_.append(pt);
2825 pointMap_.append(masterPointID);
2826 reversePointMap_.append(pointi);
2827
2828 if (zoneIDs.size())
2829 {
2830 const label minIndex = findMin(zoneIDs);
2831
2832 pointZone_.set(pointi, zoneIDs[minIndex]);
2833
2834 // Clear any additional storage
2835 if (pointi < pointAdditionalZones_.size())
2836 {
2837 pointAdditionalZones_[pointi].clear();
2838 }
2839 // (auto-vivify and) store additional zones
2840 forAll(zoneIDs, i)
2841 {
2842 if (i != minIndex)
2843 {
2844 if (zoneIDs[i] == zoneIDs[minIndex])
2845 {
2846 FatalErrorInFunction << "Duplicates in zones "
2847 << flatOutput(zoneIDs) << " for point " << pointi
2848 << exit(FatalError);
2849 }
2850 pointAdditionalZones_(pointi).push_uniq(zoneIDs[i]);
2851 }
2852 }
2853 }
2854
2855 if (!inCell)
2856 {
2857 retiredPoints_.insert(pointi);
2858 }
2859
2860 return pointi;
2861}
2862
2865(
2866 const label pointi,
2867 const point& pt,
2868 const label zoneID,
2869 const bool inCell,
2870 const bool multiZone
2871)
2872{
2873 if (pointi < 0 || pointi >= points_.size())
2874 {
2876 << "illegal point label " << pointi << endl
2877 << "Valid point labels are 0 .. " << points_.size()-1
2878 << abort(FatalError);
2879 }
2880 if (pointRemoved(pointi) || pointMap_[pointi] == -1)
2881 {
2883 << "point " << pointi << " already marked for removal"
2884 << abort(FatalError);
2885 }
2886 points_[pointi] = pt;
2887
2888 if (!multiZone)
2889 {
2890 if (zoneID >= 0)
2891 {
2892 pointZone_.set(pointi, zoneID);
2893 }
2894 else
2895 {
2896 pointZone_.erase(pointi);
2897 }
2898 if (pointi < pointAdditionalZones_.size())
2899 {
2900 pointAdditionalZones_[pointi].clear();
2901 }
2902 }
2903 else
2904 {
2905 auto fnd = pointZone_.find(pointi);
2906 if (!fnd)
2907 {
2908 pointZone_.set(pointi, zoneID);
2909
2910 // Check that no leftover additional zones
2911 if (pointAdditionalZones_(pointi).size())
2912 {
2914 << "Additional zones for point:"
2915 << pointAdditionalZones_(pointi)
2916 << abort(FatalError);
2917 }
2918 }
2919 else if (zoneID == -1)
2920 {
2921 // Clear
2922 pointZone_.erase(fnd);
2923 if (pointi < pointAdditionalZones_.size())
2924 {
2925 pointAdditionalZones_[pointi].clear();
2926 }
2927 }
2928 else if (zoneID != fnd())
2929 {
2930 // New zone. Store in additional storage.
2931 pointAdditionalZones_(pointi).push_uniq(zoneID);
2932 }
2933 }
2934
2935 if (inCell)
2936 {
2937 retiredPoints_.erase(pointi);
2938 }
2939 else
2940 {
2941 retiredPoints_.insert(pointi);
2942 }
2943}
2944
2947(
2948 const label pointi,
2949 const point& pt,
2950 const labelUList& zoneIDs,
2951 const bool inCell
2952)
2953{
2954 if (pointi < 0 || pointi >= points_.size())
2955 {
2957 << "illegal point label " << pointi << endl
2958 << "Valid point labels are 0 .. " << points_.size()-1
2959 << abort(FatalError);
2960 }
2961 if (pointRemoved(pointi) || pointMap_[pointi] == -1)
2962 {
2964 << "point " << pointi << " already marked for removal"
2965 << abort(FatalError);
2966 }
2967 points_[pointi] = pt;
2968
2969 if (zoneIDs.size())
2970 {
2971 if (zoneIDs.found(-1))
2972 {
2974 << "zones cannot contain -1 : " << flatOutput(zoneIDs)
2975 << " for point:" << pointi
2976 << abort(FatalError);
2977 }
2978
2979 pointZone_.set(pointi, zoneIDs[0]);
2980 if (pointi < pointAdditionalZones_.size())
2981 {
2982 pointAdditionalZones_[pointi].clear();
2983 }
2984 for (label i = 1; i < zoneIDs.size(); ++i)
2985 {
2986 pointAdditionalZones_(pointi).push_uniq(zoneIDs[i]);
2987 }
2988 }
2989 else
2990 {
2991 pointZone_.erase(pointi);
2992 if (pointi < pointAdditionalZones_.size())
2993 {
2994 pointAdditionalZones_[pointi].clear();
2995 }
2996 }
2997
2998 if (inCell)
2999 {
3000 retiredPoints_.erase(pointi);
3001 }
3002 else
3003 {
3004 retiredPoints_.insert(pointi);
3005 }
3006}
3008void Foam::polyTopoChange::movePoints(const pointField& newPoints)
3009{
3010 if (newPoints.size() != points_.size())
3011 {
3013 << "illegal pointField size." << endl
3014 << "Size:" << newPoints.size() << endl
3015 << "Points in mesh:" << points_.size()
3016 << abort(FatalError);
3017 }
3018
3019 forAll(points_, pointi)
3020 {
3021 points_[pointi] = newPoints[pointi];
3022 }
3023}
3024
3027(
3028 const label pointi,
3029 const label mergePointi
3030)
3031{
3032 if (pointi < 0 || pointi >= points_.size())
3033 {
3035 << "illegal point label " << pointi << endl
3036 << "Valid point labels are 0 .. " << points_.size()-1
3037 << abort(FatalError);
3038 }
3039
3040 if
3041 (
3042 strict_
3043 && (pointRemoved(pointi) || pointMap_[pointi] == -1)
3044 )
3045 {
3047 << "point " << pointi << " already marked for removal" << nl
3048 << "Point:" << points_[pointi] << " pointMap:" << pointMap_[pointi]
3049 << abort(FatalError);
3050 }
3051
3052 if (pointi == mergePointi)
3053 {
3055 << "Cannot remove/merge point " << pointi << " onto itself."
3056 << abort(FatalError);
3057 }
3058
3059 points_[pointi] = point::max;
3060 pointMap_[pointi] = -1;
3061 if (mergePointi >= 0)
3062 {
3063 reversePointMap_[pointi] = -mergePointi-2;
3064 }
3065 else
3066 {
3067 reversePointMap_[pointi] = -1;
3068 }
3069 pointZone_.erase(pointi);
3070 if (pointi < pointAdditionalZones_.size())
3071 {
3072 pointAdditionalZones_[pointi].clear();
3073 }
3074 retiredPoints_.erase(pointi);
3075}
3076
3079(
3080 const face& f,
3081 const label own,
3082 const label nei,
3083 const label masterPointID,
3084 const label masterEdgeID,
3085 const label masterFaceID,
3086 const bool flipFaceFlux,
3087 const label patchID,
3088 const label zoneID,
3089 const bool zoneFlip
3090)
3091{
3092 // Check validity
3093 if (debug)
3094 {
3095 checkFace(f, -1, own, nei, patchID, zoneID);
3096 }
3097
3098 label facei = faces_.size();
3099
3100 faces_.append(f);
3101 region_.append(patchID);
3102 faceOwner_.append(own);
3103 faceNeighbour_.append(nei);
3104
3105 if (masterPointID >= 0)
3106 {
3107 faceMap_.append(-1);
3108 faceFromPoint_.insert(facei, masterPointID);
3109 }
3110 else if (masterEdgeID >= 0)
3111 {
3112 faceMap_.append(-1);
3113 faceFromEdge_.insert(facei, masterEdgeID);
3114 }
3115 else if (masterFaceID >= 0)
3116 {
3117 faceMap_.append(masterFaceID);
3118 }
3119 else
3120 {
3121 // Allow inflate-from-nothing?
3122 //FatalErrorInFunction
3123 // << "Need to specify a master point, edge or face"
3124 // << "face:" << f << " own:" << own << " nei:" << nei
3125 // << abort(FatalError);
3126 faceMap_.append(-1);
3127 }
3128 reverseFaceMap_.append(facei);
3129
3130 flipFaceFlux_.set(facei, flipFaceFlux);
3131
3132 if (zoneID >= 0)
3133 {
3134 faceZone_.insert(facei, zoneID);
3135 }
3136 faceZoneFlip_.set(facei, zoneFlip);
3137 // Delay extending the faceAdditionalZones
3138
3139 return facei;
3140}
3141
3144(
3145 const face& f,
3146 const label own,
3147 const label nei,
3148 const label masterPointID,
3149 const label masterEdgeID,
3150 const label masterFaceID,
3151 const bool flipFaceFlux,
3152 const label patchID,
3153 const labelUList& zoneIDs,
3154 const UList<bool>& zoneFlips
3155)
3156{
3157 // Check validity
3158 if (debug)
3159 {
3160 // Note: no check on zones
3161 checkFace(f, -1, own, nei, patchID, -1);
3162 }
3163
3164 label facei = faces_.size();
3165
3166 faces_.append(f);
3167 region_.append(patchID);
3168 faceOwner_.append(own);
3169 faceNeighbour_.append(nei);
3170
3171 if (masterPointID >= 0)
3172 {
3173 faceMap_.append(-1);
3174 faceFromPoint_.insert(facei, masterPointID);
3175 }
3176 else if (masterEdgeID >= 0)
3177 {
3178 faceMap_.append(-1);
3179 faceFromEdge_.insert(facei, masterEdgeID);
3180 }
3181 else if (masterFaceID >= 0)
3182 {
3183 faceMap_.append(masterFaceID);
3184 }
3185 else
3186 {
3187 // Allow inflate-from-nothing?
3188 //FatalErrorInFunction
3189 // << "Need to specify a master point, edge or face"
3190 // << "face:" << f << " own:" << own << " nei:" << nei
3191 // << abort(FatalError);
3192 faceMap_.append(-1);
3193 }
3194 reverseFaceMap_.append(facei);
3195
3196 flipFaceFlux_.set(facei, flipFaceFlux);
3197
3198 if (zoneIDs.size())
3199 {
3200 const label minIndex = findMin(zoneIDs);
3201
3202 faceZone_.set(facei, zoneIDs[minIndex]);
3203 faceZoneFlip_.set(facei, zoneFlips[minIndex]);
3204
3205 // Clear any additional storage
3206 if (facei < faceAdditionalZones_.size())
3207 {
3208 faceAdditionalZones_[facei].clear();
3209 }
3210 // (auto-vivify and) store additional zones
3211 forAll(zoneIDs, i)
3212 {
3213 if (i != minIndex)
3214 {
3215 if (zoneIDs[i] == zoneIDs[minIndex])
3216 {
3217 FatalErrorInFunction << "Duplicates in zones "
3218 << flatOutput(zoneIDs) << " for face " << facei
3219 << exit(FatalError);
3220 }
3221 const label zoneAndSign
3222 (
3223 zoneFlips[i]
3224 ? zoneIDs[i]+1
3225 : -(zoneIDs[i]+1)
3226 );
3227 faceAdditionalZones_(facei).push_uniq(zoneAndSign);
3228 }
3229 }
3230 }
3231
3232 return facei;
3233}
3234
3237(
3238 const face& f,
3239 const label facei,
3240 const label own,
3241 const label nei,
3242 const bool flipFaceFlux,
3243 const label patchID,
3244 const label zoneID,
3245 const bool zoneFlip,
3246 const bool multiZone
3247)
3248{
3249 // Check validity
3250 if (debug)
3251 {
3252 checkFace(f, facei, own, nei, patchID, zoneID);
3253 }
3254
3255 faces_[facei] = f;
3256 faceOwner_[facei] = own;
3257 faceNeighbour_[facei] = nei;
3258 region_[facei] = patchID;
3259
3260 flipFaceFlux_.set(facei, flipFaceFlux);
3261
3262 // Note: same logic as modifyCell except storage is now sparse instead
3263 // of marked with -1
3264
3265 if (!multiZone)
3266 {
3267 if (zoneID == -1)
3268 {
3269 faceZone_.erase(facei);
3270 faceZoneFlip_.set(facei, zoneFlip);
3271 }
3272 else
3273 {
3274 faceZone_.set(facei, zoneID);
3275 faceZoneFlip_.set(facei, zoneFlip);
3276 }
3277 if (facei < faceAdditionalZones_.size())
3278 {
3279 faceAdditionalZones_[facei].clear();
3280 }
3281 }
3282 else
3283 {
3284 auto fnd = faceZone_.find(facei);
3285 if (!fnd)
3286 {
3287 faceZone_.set(facei, zoneID);
3288 faceZoneFlip_.set(facei, zoneFlip);
3289
3290 // Check that no leftover additional zones
3291 if (faceAdditionalZones_(facei).size())
3292 {
3294 << "Additional zones for face:"
3295 << faceAdditionalZones_(facei)
3296 << abort(FatalError);
3297 }
3298 }
3299 else if (zoneID == -1)
3300 {
3301 // Clear
3302 faceZone_.erase(fnd);
3303 faceZoneFlip_.set(facei, false);
3304 if (facei < faceAdditionalZones_.size())
3305 {
3306 faceAdditionalZones_[facei].clear();
3307 }
3308 }
3309 else if (zoneID != fnd())
3310 {
3311 // New zone. Store in additional storage.
3312 const label zoneAndSign
3313 (
3314 zoneFlip
3315 ? zoneID+1
3316 : -(zoneID+1)
3317 );
3318 faceAdditionalZones_(facei).push_uniq(zoneAndSign);
3319 }
3320 }
3321}
3322
3325(
3326 const face& f,
3327 const label facei,
3328 const label own,
3329 const label nei,
3330 const bool flipFaceFlux,
3331 const label patchID,
3332 const labelUList& zoneIDs,
3333 const UList<bool>& zoneFlips
3334)
3335{
3336 // Check validity
3337 if (debug)
3338 {
3339 // Note: no check on zones
3340 checkFace(f, facei, own, nei, patchID, -1);
3341 }
3342
3343 faces_[facei] = f;
3344 faceOwner_[facei] = own;
3345 faceNeighbour_[facei] = nei;
3346 region_[facei] = patchID;
3347
3348 flipFaceFlux_.set(facei, flipFaceFlux);
3349
3350 // Note: same logic as modifyCell except storage is now sparse instead
3351 // of marked with -1
3352
3353 if (zoneIDs.size())
3354 {
3355 if (zoneIDs.found(-1))
3356 {
3358 << "zones cannot contain -1 : " << flatOutput(zoneIDs)
3359 << " for face:" << facei
3360 << abort(FatalError);
3361 }
3362
3363 faceZone_.set(facei, zoneIDs[0]);
3364 faceZoneFlip_.set(facei, zoneFlips[0]);
3365 if (facei < faceAdditionalZones_.size())
3366 {
3367 faceAdditionalZones_[facei].clear();
3368 }
3369 for (label i = 1; i < zoneIDs.size(); ++i)
3370 {
3371 const label zoneAndSign
3372 (
3373 zoneFlips[i]
3374 ? zoneIDs[i]+1
3375 : -(zoneIDs[i]+1)
3376 );
3377 faceAdditionalZones_(facei).push_uniq(zoneAndSign);
3378 }
3379 }
3380 else
3381 {
3382 faceZone_.erase(facei);
3383 faceZoneFlip_.set(facei, false);
3384 if (facei < faceAdditionalZones_.size())
3385 {
3386 faceAdditionalZones_[facei].clear();
3387 }
3388 }
3389}
3390
3393(
3394 const label facei,
3395 const label mergeFacei
3396)
3397{
3398 if (facei < 0 || facei >= faces_.size())
3399 {
3401 << "illegal face label " << facei << endl
3402 << "Valid face labels are 0 .. " << faces_.size()-1
3403 << abort(FatalError);
3404 }
3405
3406 if
3407 (
3408 strict_
3409 && (faceRemoved(facei) || faceMap_[facei] == -1)
3410 )
3411 {
3413 << "face " << facei
3414 << " already marked for removal"
3415 << abort(FatalError);
3416 }
3417
3418 faces_[facei].clear();
3419 region_[facei] = -1;
3420 faceOwner_[facei] = -1;
3421 faceNeighbour_[facei] = -1;
3422 faceMap_[facei] = -1;
3423 if (mergeFacei >= 0)
3424 {
3425 reverseFaceMap_[facei] = -mergeFacei-2;
3426 }
3427 else
3428 {
3429 reverseFaceMap_[facei] = -1;
3430 }
3431 faceFromEdge_.erase(facei);
3432 faceFromPoint_.erase(facei);
3433 flipFaceFlux_.unset(facei);
3434 faceZoneFlip_.unset(facei);
3435 faceZone_.erase(facei);
3436 if (facei < faceAdditionalZones_.size())
3437 {
3438 faceAdditionalZones_[facei].clear();
3439 }
3440}
3441
3444(
3445 const label masterPointID,
3446 const label masterEdgeID,
3447 const label masterFaceID,
3448 const label masterCellID,
3449 const label zoneID
3450)
3451{
3452 label celli = cellMap_.size();
3453
3454 if (masterPointID >= 0)
3455 {
3456 cellMap_.append(-1);
3457 cellFromPoint_.insert(celli, masterPointID);
3458 }
3459 else if (masterEdgeID >= 0)
3460 {
3461 cellMap_.append(-1);
3462 cellFromEdge_.insert(celli, masterEdgeID);
3463 }
3464 else if (masterFaceID >= 0)
3465 {
3466 cellMap_.append(-1);
3467 cellFromFace_.insert(celli, masterFaceID);
3468 }
3469 else
3470 {
3471 cellMap_.append(masterCellID);
3472 }
3473 reverseCellMap_.append(celli);
3474
3475 cellZone_.append(zoneID);
3476 // Delay extending the cellAdditionalZones
3477
3478 return celli;
3479}
3480
3483(
3484 const label masterPointID,
3485 const label masterEdgeID,
3486 const label masterFaceID,
3487 const label masterCellID,
3488 const labelUList& zoneIDs
3489)
3490{
3491 label celli = cellMap_.size();
3492
3493 if (masterPointID >= 0)
3494 {
3495 cellMap_.append(-1);
3496 cellFromPoint_.insert(celli, masterPointID);
3497 }
3498 else if (masterEdgeID >= 0)
3499 {
3500 cellMap_.append(-1);
3501 cellFromEdge_.insert(celli, masterEdgeID);
3502 }
3503 else if (masterFaceID >= 0)
3504 {
3505 cellMap_.append(-1);
3506 cellFromFace_.insert(celli, masterFaceID);
3507 }
3508 else
3509 {
3510 cellMap_.append(masterCellID);
3511 }
3512 reverseCellMap_.append(celli);
3513
3514 if (zoneIDs.size())
3515 {
3516 cellZone_.append(zoneIDs[0]);
3517
3518 // Clear any additional storage
3519 if (celli < cellAdditionalZones_.size())
3520 {
3521 cellAdditionalZones_[celli].clear();
3522 }
3523 // (auto-vivify and) store additional zones
3524 for (label i = 1; i < zoneIDs.size(); ++i)
3525 {
3526 cellAdditionalZones_(celli).push_uniq(zoneIDs[i]);
3527 }
3528 }
3529 else
3530 {
3531 cellZone_.append(-1);
3532 }
3533
3534 return celli;
3535}
3536
3539(
3540 const label celli,
3541 const label zoneID,
3542 const bool multiZone
3543)
3544{
3545 if (!multiZone)
3546 {
3547 cellZone_[celli] = zoneID;
3548 if (celli < cellAdditionalZones_.size())
3549 {
3550 cellAdditionalZones_[celli].clear();
3551 }
3552 }
3553 else
3554 {
3555 if (cellZone_[celli] == -1)
3556 {
3557 cellZone_[celli] = zoneID;
3558
3559 // Check that no leftover additional zones
3560 if (cellAdditionalZones_(celli).size())
3561 {
3563 << "Additional zones for cell:"
3564 << cellAdditionalZones_(celli)
3565 << abort(FatalError);
3566 }
3567 }
3568 else
3569 {
3570 if (zoneID == -1)
3571 {
3572 // Clear
3573 cellZone_[celli] = zoneID;
3574 if (celli < cellAdditionalZones_.size())
3575 {
3576 cellAdditionalZones_[celli].clear();
3577 }
3578 }
3579 else if (zoneID != cellZone_[celli])
3580 {
3581 // New zone. Store in additional storage.
3582 cellAdditionalZones_(celli).push_uniq(zoneID);
3583 }
3584 }
3585 }
3586}
3587
3590(
3591 const label celli,
3592 const labelUList& zoneIDs
3593)
3594{
3595 if (celli < 0 || celli >= cellMap_.size())
3596 {
3598 << "illegal cell label " << celli << endl
3599 << "Valid cell labels are 0 .. " << cellMap_.size()-1
3600 << abort(FatalError);
3601 }
3602
3603 if (zoneIDs.size())
3604 {
3605 if (zoneIDs.found(-1))
3606 {
3608 << "zones cannot contain -1 : " << flatOutput(zoneIDs)
3609 << " for cell:" << celli
3610 << abort(FatalError);
3611 }
3612
3613 cellZone_[celli] = zoneIDs[0];
3614 if (celli < cellAdditionalZones_.size())
3615 {
3616 cellAdditionalZones_[celli].clear();
3617 }
3618 for (label i = 1; i < zoneIDs.size(); ++i)
3619 {
3620 cellAdditionalZones_(celli).push_uniq(zoneIDs[i]);
3621 }
3622 }
3623 else
3624 {
3625 cellZone_[celli] = -1;
3626 if (celli < cellAdditionalZones_.size())
3627 {
3628 cellAdditionalZones_[celli].clear();
3629 }
3630 }
3631}
3632
3635(
3636 const label celli,
3637 const label mergeCelli
3638)
3639{
3640 if (celli < 0 || celli >= cellMap_.size())
3641 {
3643 << "illegal cell label " << celli << endl
3644 << "Valid cell labels are 0 .. " << cellMap_.size()-1
3645 << abort(FatalError);
3646 }
3647
3648 if (strict_ && cellMap_[celli] == -2)
3649 {
3651 << "cell " << celli
3652 << " already marked for removal"
3653 << abort(FatalError);
3654 }
3655
3656 cellMap_[celli] = -2;
3657 if (mergeCelli >= 0)
3658 {
3659 reverseCellMap_[celli] = -mergeCelli-2;
3660 }
3661 else
3662 {
3663 reverseCellMap_[celli] = -1;
3664 }
3665 cellFromPoint_.erase(celli);
3666 cellFromEdge_.erase(celli);
3667 cellFromFace_.erase(celli);
3668 cellZone_[celli] = -1;
3669 if (celli < cellAdditionalZones_.size())
3670 {
3671 cellAdditionalZones_[celli].clear();
3672 }
3673}
3674
3677(
3678 const label pointi,
3679 DynamicList<label>& zones
3680) const
3681{
3682 if (pointi < 0 || pointi >= pointMap_.size())
3683 {
3685 << "illegal point label " << pointi << endl
3686 << "Valid point labels are 0 .. " << pointMap_.size()-1
3687 << abort(FatalError);
3688 }
3689 zones.clear();
3690
3691 const auto fnd = pointZone_.find(pointi);
3692 if (fnd)
3693 {
3694 zones.push_back(fnd());
3695 if (pointi < pointAdditionalZones_.size())
3696 {
3697 for (const label zonei : pointAdditionalZones_[pointi])
3698 {
3699 zones.push_uniq(zonei);
3700 }
3701 }
3702 }
3703 return zones.size();
3704}
3705
3708(
3709 const label facei,
3710 DynamicList<label>& zones,
3711 DynamicList<bool>& flips
3712) const
3713{
3714 if (facei < 0 || facei >= faceMap_.size())
3715 {
3717 << "illegal face label " << facei << endl
3718 << "Valid face labels are 0 .. " << faceMap_.size()-1
3719 << abort(FatalError);
3720 }
3721 zones.clear();
3722 flips.clear();
3723
3724 const auto fnd = faceZone_.find(facei);
3725 if (fnd)
3726 {
3727 zones.push_back(fnd());
3728 flips.push_back(flipFaceFlux_[facei]);
3729 if (facei < faceAdditionalZones_.size())
3730 {
3731 for (const label zoneAndSign : faceAdditionalZones_[facei])
3732 {
3733 const label zonei = mag(zoneAndSign)-1;
3734 if (!zones.found(zonei))
3735 {
3736 zones.push_back(zonei);
3737 flips.push_back((zoneAndSign > 0));
3738 }
3739 }
3740 }
3741 }
3742 return zones.size();
3743}
3744
3747(
3748 const label celli,
3749 DynamicList<label>& zones
3750) const
3751{
3752 if (celli < 0 || celli >= cellMap_.size())
3753 {
3755 << "illegal cell label " << celli << endl
3756 << "Valid cell labels are 0 .. " << cellMap_.size()-1
3757 << abort(FatalError);
3758 }
3759 zones.clear();
3760 if (cellZone_[celli] != -1)
3761 {
3762 zones.push_back(cellZone_[celli]);
3763 }
3764 if (celli < cellAdditionalZones_.size())
3765 {
3766 for (const label zonei : cellAdditionalZones_[celli])
3767 {
3768 zones.push_uniq(zonei);
3769 }
3770 }
3771 return zones.size();
3772}
3773
3776(
3777 polyMesh& mesh,
3778 const labelUList& patchMap,
3779 const bool inflate,
3780 const bool syncParallel,
3781 const bool orderCells,
3782 const bool orderPoints
3783)
3784{
3785 if (debug)
3786 {
3787 Pout<< "polyTopoChange::changeMesh"
3788 << "(polyMesh&, const bool, const bool, const bool, const bool)"
3789 << endl;
3790 }
3791
3792 if (debug)
3793 {
3794 Pout<< "Old mesh:" << nl;
3795 writeMeshStats(mesh, Pout);
3796 }
3797
3798 // new mesh points
3799 pointField newPoints;
3800 // number of internal points
3801 label nInternalPoints;
3802 // patch slicing
3803 labelList patchSizes;
3804 labelList patchStarts;
3805 // inflate maps
3806 List<objectMap> pointsFromPoints;
3807 List<objectMap> facesFromPoints;
3808 List<objectMap> facesFromEdges;
3809 List<objectMap> facesFromFaces;
3810 List<objectMap> cellsFromPoints;
3811 List<objectMap> cellsFromEdges;
3812 List<objectMap> cellsFromFaces;
3813 List<objectMap> cellsFromCells;
3814 // old mesh info
3815 List<Map<label>> oldPatchMeshPointMaps;
3816 labelList oldPatchNMeshPoints;
3817 labelList oldPatchStarts;
3818 List<Map<label>> oldFaceZoneMeshPointMaps;
3819
3820 // Compact, reorder patch faces and calculate mesh/patch maps.
3821 compactAndReorder
3822 (
3823 mesh,
3824 patchMap,
3825 syncParallel,
3826 orderCells,
3827 orderPoints,
3828
3829 nInternalPoints,
3830 newPoints,
3831 patchSizes,
3832 patchStarts,
3833 pointsFromPoints,
3834 facesFromPoints,
3835 facesFromEdges,
3836 facesFromFaces,
3837 cellsFromPoints,
3838 cellsFromEdges,
3839 cellsFromFaces,
3840 cellsFromCells,
3841 oldPatchMeshPointMaps,
3842 oldPatchNMeshPoints,
3843 oldPatchStarts,
3844 oldFaceZoneMeshPointMaps
3845 );
3846
3847 const label nOldPoints(mesh.nPoints());
3848 const label nOldFaces(mesh.nFaces());
3849 const label nOldCells(mesh.nCells());
3850 autoPtr<scalarField> oldCellVolumes(new scalarField(mesh.cellVolumes()));
3851
3852
3853 // Change the mesh
3854 // ~~~~~~~~~~~~~~~
3855 // This will invalidate any addressing so better make sure you have
3856 // all the information you need!!!
3857
3858 if (inflate)
3859 {
3860 // Keep (renumbered) mesh points, store new points in map for inflation
3861 // (appended points (i.e. from nowhere) get value zero)
3862 pointField renumberedMeshPoints(newPoints.size());
3863
3864 forAll(pointMap_, newPointi)
3865 {
3866 const label oldPointi = pointMap_[newPointi];
3867
3868 if (oldPointi >= 0)
3869 {
3870 renumberedMeshPoints[newPointi] = mesh.points()[oldPointi];
3871 }
3872 else
3873 {
3874 renumberedMeshPoints[newPointi] = Zero;
3875 }
3876 }
3877
3878 mesh.resetPrimitives
3879 (
3880 autoPtr<pointField>::New(std::move(renumberedMeshPoints)),
3881 autoPtr<faceList>::New(std::move(faces_)),
3882 autoPtr<labelList>::New(std::move(faceOwner_)),
3883 autoPtr<labelList>::New(std::move(faceNeighbour_)),
3884 patchSizes,
3885 patchStarts,
3886 syncParallel
3887 );
3888
3889 mesh.topoChanging(true);
3890 // Note: could already set moving flag as well
3891 // mesh.moving(true);
3892 }
3893 else
3894 {
3895 // Set new points.
3896 mesh.resetPrimitives
3897 (
3898 autoPtr<pointField>::New(std::move(newPoints)),
3899 autoPtr<faceList>::New(std::move(faces_)),
3900 autoPtr<labelList>::New(std::move(faceOwner_)),
3901 autoPtr<labelList>::New(std::move(faceNeighbour_)),
3902 patchSizes,
3903 patchStarts,
3904 syncParallel
3905 );
3906 mesh.topoChanging(true);
3907 }
3908
3909 // Clear out primitives
3910 {
3911 retiredPoints_.clearStorage();
3912 region_.clearStorage();
3913 }
3914
3915
3916 if (debug)
3917 {
3918 // Some stats on changes
3919 label nAdd, nInflate, nMerge, nRemove;
3920 countMap(pointMap_, reversePointMap_, nAdd, nInflate, nMerge, nRemove);
3921 Pout<< "Points:"
3922 << " added(from point):" << nAdd
3923 << " added(from nothing):" << nInflate
3924 << " merged(into other point):" << nMerge
3925 << " removed:" << nRemove
3926 << nl;
3927
3928 countMap(faceMap_, reverseFaceMap_, nAdd, nInflate, nMerge, nRemove);
3929 Pout<< "Faces:"
3930 << " added(from face):" << nAdd
3931 << " added(inflated):" << nInflate
3932 << " merged(into other face):" << nMerge
3933 << " removed:" << nRemove
3934 << nl;
3935
3936 countMap(cellMap_, reverseCellMap_, nAdd, nInflate, nMerge, nRemove);
3937 Pout<< "Cells:"
3938 << " added(from cell):" << nAdd
3939 << " added(inflated):" << nInflate
3940 << " merged(into other cell):" << nMerge
3941 << " removed:" << nRemove
3942 << nl
3943 << endl;
3944 }
3945
3946 // Zones
3947 // ~~~~~
3948
3949 // Inverse of point/face/cell zone addressing.
3950 // For every preserved point/face/cells in zone give the old position.
3951 // For added points, the index is set to -1
3952 labelListList pointZoneMap(mesh.pointZones().size());
3953 labelListList faceZoneFaceMap(mesh.faceZones().size());
3954 labelListList cellZoneMap(mesh.cellZones().size());
3955
3956 resetZones(mesh, mesh, pointZoneMap, faceZoneFaceMap, cellZoneMap);
3957
3958 if (debug)
3959 {
3960 Pout<< "New mesh:" << nl;
3961 writeMeshStats(mesh, Pout);
3962 }
3963
3964 // Clear zone info
3965 {
3966 pointZone_.clearStorage();
3967 faceZone_.clearStorage();
3968 faceZoneFlip_.clearStorage();
3969 cellZone_.clearStorage();
3970 }
3971
3972
3973 // Patch point renumbering
3974 // For every preserved point on a patch give the old position.
3975 // For added points, the index is set to -1
3976 labelListList patchPointMap(mesh.boundaryMesh().size());
3977 calcPatchPointMap
3978 (
3979 oldPatchMeshPointMaps,
3980 patchMap,
3981 mesh.boundaryMesh(),
3982 patchPointMap
3983 );
3984
3985 // Create the face zone mesh point renumbering
3986 labelListList faceZonePointMap(mesh.faceZones().size());
3987 calcFaceZonePointMap(mesh, oldFaceZoneMeshPointMaps, faceZonePointMap);
3988
3989 labelHashSet flipFaceFluxSet(HashSetOps::used(flipFaceFlux_));
3990
3992 (
3993 mesh,
3994 nOldPoints,
3995 nOldFaces,
3996 nOldCells,
3997
3998 pointMap_,
3999 pointsFromPoints,
4000
4001 faceMap_,
4002 facesFromPoints,
4003 facesFromEdges,
4004 facesFromFaces,
4005
4006 cellMap_,
4007 cellsFromPoints,
4008 cellsFromEdges,
4009 cellsFromFaces,
4010 cellsFromCells,
4011
4012 reversePointMap_,
4013 reverseFaceMap_,
4014 reverseCellMap_,
4015
4016 flipFaceFluxSet,
4017
4018 patchPointMap,
4019
4020 pointZoneMap,
4021
4022 faceZonePointMap,
4023 faceZoneFaceMap,
4024 cellZoneMap,
4025
4026 newPoints, // if empty signals no inflation.
4027 oldPatchStarts,
4028 oldPatchNMeshPoints,
4029
4030 oldCellVolumes,
4031
4032 true // steal storage.
4033 );
4034
4035 // At this point all member DynamicList (pointMap_, cellMap_ etc.) will
4036 // be invalid.
4037}
4038
4041(
4042 polyMesh& mesh,
4043 const bool inflate,
4044 const bool syncParallel,
4045 const bool orderCells,
4046 const bool orderPoints
4047)
4048{
4049 return changeMesh
4050 (
4051 mesh,
4052 identity(mesh.boundaryMesh().size()),
4053 inflate,
4054 syncParallel,
4055 orderCells,
4056 orderPoints
4057 );
4058}
4059
4060
4061// ************************************************************************* //
labelList faceLabels(nFaceLabels)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
A simple list of objects of type <T> that is intended to be used as a circular buffer (eg,...
A packed storage of objects of type <T> using an offset table for access.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
void append(const T &val)
Copy append an element to the end of this list.
label push_uniq(const T &val)
Append an element if not already in the list.
void push_back(const T &val)
Copy append an element to the end of this list.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void setSize(label n)
Alias for resize().
Definition List.H:536
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Buffers for inter-processor communications streams (UOPstream, UIPstream).
void setSize(const label n)
Same as resize().
Definition PtrList.H:357
A non-owning sub-view of a List (allocated or unallocated storage).
Definition SubList.H:61
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
bool found(const T &val, label pos=0) const
Same as contains().
Definition UList.H:983
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static int & msgType() noexcept
Message tag of standard messages.
Definition UPstream.H:1926
static int incrMsgType(int val=1) noexcept
Increment the message tag for standard messages.
Definition UPstream.H:1948
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
static const Form max
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition autoPtr.H:178
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
A subset of mesh cells.
Definition cellZone.H:61
A subset of mesh faces organised as a primitive patch.
Definition faceZone.H:63
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
An objectMap is a pair of labels defining the mapping of an object from another object,...
Definition objectMap.H:57
A subset of mesh points.
Definition pointZone.H:64
Class containing data for cell addition.
Definition polyAddCell.H:45
label zoneID() const
Cell zone ID.
label masterFaceID() const
Return master face ID.
label masterPointID() const
Return master point ID.
label masterCellID() const
Return master cell ID.
label masterEdgeID() const
Return master edge ID.
A face addition data class. A face can be inflated either from a point or from another face and can e...
Definition polyAddFace.H:50
bool flipFaceFlux() const
Does the face flux need to be flipped.
label owner() const
Return owner cell.
label zoneID() const
Face zone ID.
label masterFaceID() const
Return master face ID.
label masterPointID() const
Return master point ID.
const face & newFace() const
Return face.
label patchID() const
Boundary patch ID.
label zoneFlip() const
Face zone flip.
label neighbour() const
Return neighbour cell.
label masterEdgeID() const
Return master edge ID.
Class containing data for point addition.
label zoneID() const
Point zone ID.
const point & newPoint() const
Point location.
bool inCell() const
Does the point support a cell.
label masterPointID() const
Master point label.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
label start() const noexcept
The start label of boundary faces in the polyMesh face list.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition polyMesh.H:671
virtual const labelList & faceOwner() const
Return face owner.
Definition polyMesh.C:1101
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition polyMesh.H:679
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition polyMesh.C:1107
const pointZoneMesh & pointZones() const noexcept
Return point zone mesh.
Definition polyMesh.H:663
virtual const pointField & points() const
Return raw points.
Definition polyMesh.C:1063
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:73
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 pointID() const
Return point ID.
Direct mesh changes based on v1.3 polyTopoChange syntax.
void movePoints(const pointField &newPoints)
Move all points. Incompatible with other topology changes.
void shrink()
Shrink storage (does not remove any elements; just compacts dynamic lists.
label faceZones(const label facei, DynamicList< label > &zones, DynamicList< bool > &flips) const
Get current faceZone(s). Return number of zones.
const DynamicList< face > & faces() const
polyTopoChange(const label nPatches, const bool strict=true)
Construct without mesh. Either specify nPatches or use setNumPatches before trying to make a mesh (ma...
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 modifyCell(const label celli, const label zoneID, const bool multiZone=false)
Modify zone of cell. Optionally add to zone.
bool pointRemoved(const label pointi) const
Is point removed? Considered removed if point is GREAT.
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.
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, const bool multiZone=false)
Modify vertices or cell of face.
label pointZones(const label pointi, DynamicList< label > &zones) const
Get current cellZone(s). Return number of zones.
void setCapacity(const label nPoints, const label nFaces, const label nCells)
Explicitly pre-size the dynamic storage for expected mesh size for if construct-without-mesh.
label addCell(const label masterPointID, const label masterEdgeID, const label masterFaceID, const label masterCellID, const label zoneID)
Add cell. Return new cell label.
bool faceRemoved(const label facei) const
Is face removed? Considered removed if face is empty.
void addMesh(const polyMesh &mesh, const labelUList &patchMap, const labelUList &pointZoneMap, const labelUList &faceZoneMap, const labelUList &cellZoneMap)
Add all points/faces/cells of mesh. Additional offset for patch or zone ids.
void removePoint(const label pointi, const label mergePointi)
Remove/merge point.
void removeFace(const label facei, const label mergeFacei)
Remove/merge face.
void modifyPoint(const label pointi, const point &pt, const label zoneID, const bool inCell, const bool multiZone=false)
Modify coordinate.
const DynamicList< label > & faceOwner() const
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.
const DynamicList< point > & points() const
Points. Shrunk after constructing mesh (or calling of compact()).
void clear()
Clear all storage.
const DynamicList< label > & faceNeighbour() const
label cellZones(const label celli, DynamicList< label > &zones) const
Get current cellZone(s). Return number of zones.
Cell-face mesh analysis engine.
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
const labelListList & pointCells() 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
A virtual base class for topological actions.
Definition topoAction.H:48
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
const polyBoundaryMesh & patches
faceListList boundary
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
label nPoints
const labelIOList & zoneIDs
Definition correctPhi.H:59
labelHashSet used(const bitSet &select)
Convert a bitset to a labelHashSet of the indices used.
Definition HashOps.C:26
Namespace for handling debugging switches.
Definition debug.C:45
const wordList internal
Standard volume internal field types (scalar, vector, tensor, etc).
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
ZoneMesh< pointZone, polyMesh > pointZoneMesh
A ZoneMesh with pointZone content on a polyMesh.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
label findMin(const ListType &input, label start=0)
Linear search for the index of the min element, similar to std::min_element but for lists and returns...
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field.
IntListType renumber(const labelUList &oldToNew, const IntListType &input)
Renumber the values within a list.
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with faceZone content on a polyMesh.
List< face > faceList
List of faces.
Definition faceListFwd.H:41
void inplaceRotateList(ListType< DataType > &list, label n)
Inplace reversal of a list using the Reversal Block Swapping algorithm.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
constexpr label labelMax
Definition label.H:55
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
void inplaceReorder(const labelUList &oldToNew, ListType &input, const bool prune=false)
Inplace reorder the elements of a list.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
ZoneMesh< cellZone, polyMesh > cellZoneMesh
A ZoneMesh with cellZone content on a polyMesh.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition FlatOutput.H:217
List< List< bool > > boolListList
List of boolList.
Definition boolList.H:36
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
errorManip< error > abort(error &err)
Definition errorManip.H:139
vector point
Point is a vector.
Definition point.H:37
List< bool > boolList
A List of bools.
Definition List.H:60
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
bool isType(const U &obj)
Check if typeid of the object and Type are identical.
Definition typeInfo.H:112
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
labelList invert(const label len, const labelUList &map)
Create an inverse one-to-one mapping.
Definition ListOps.C:28
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
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:299
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
label newPointi
labelList f(nPoints)
label nPatches
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235