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