polyMeshAdder.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "polyMeshAdder.H"
27 #include "mapAddedPolyMesh.H"
28 #include "IOobject.H"
29 #include "faceCoupleInfo.H"
30 #include "processorPolyPatch.H"
31 #include "SortableList.H"
32 #include "Time.H"
33 #include "globalMeshData.H"
34 #include "mergePoints.H"
35 #include "polyModifyFace.H"
36 #include "polyRemovePoint.H"
37 #include "polyTopoChange.H"
38 
39 #include "pointMesh.H"
40 #include "facePointPatch.H"
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 // Get index of patch in new set of patchnames/types
45 Foam::label Foam::polyMeshAdder::patchIndex
46 (
47  const polyPatch& p,
48  DynamicList<word>& allPatchNames,
49  DynamicList<word>& allPatchTypes
50 )
51 {
52  // Find the patch name on the list. If the patch is already there
53  // and patch types match, return index
54  const word& pType = p.type();
55  const word& pName = p.name();
56 
57  label patchi = findIndex(allPatchNames, pName);
58 
59  if (patchi == -1)
60  {
61  // Patch not found. Append to the list
62  allPatchNames.append(pName);
63  allPatchTypes.append(pType);
64 
65  return allPatchNames.size() - 1;
66  }
67  else if (allPatchTypes[patchi] == pType)
68  {
69  // Found name and types match
70  return patchi;
71  }
72  else
73  {
74  // Found the name, but type is different
75 
76  // Duplicate name is not allowed. Create a composite name from the
77  // patch name and case name
78  const word& caseName = p.boundaryMesh().mesh().time().caseName();
79 
80  allPatchNames.append(pName + "_" + caseName);
81  allPatchTypes.append(pType);
82 
83  Pout<< "label patchIndex(const polyPatch& p) : "
84  << "Patch " << p.index() << " named "
85  << pName << " in mesh " << caseName
86  << " already exists, but patch types"
87  << " do not match.\nCreating a composite name as "
88  << allPatchNames.last() << endl;
89 
90  return allPatchNames.size() - 1;
91  }
92 }
93 
94 
95 // Get index of zone in new set of zone names
96 Foam::label Foam::polyMeshAdder::zoneIndex
97 (
98  const word& curName,
99  DynamicList<word>& names
100 )
101 {
102  label zoneI = findIndex(names, curName);
103 
104  if (zoneI != -1)
105  {
106  return zoneI;
107  }
108  else
109  {
110  // Not found. Add new name to the list
111  names.append(curName);
112 
113  return names.size() - 1;
114  }
115 }
116 
117 
118 void Foam::polyMeshAdder::mergePatchNames
119 (
120  const polyBoundaryMesh& patches0,
121  const polyBoundaryMesh& patches1,
122 
123  DynamicList<word>& allPatchNames,
124  DynamicList<word>& allPatchTypes,
125 
126  labelList& from1ToAllPatches,
127  labelList& fromAllTo1Patches
128 )
129 {
130  // Insert the mesh0 patches and zones
131  allPatchNames.append(patches0.names());
132  allPatchTypes.append(patches0.types());
133 
134 
135  // Patches
136  // ~~~~~~~
137  // Patches from 0 are taken over as is; those from 1 get either merged
138  // (if they share name and type) or appended.
139  // Empty patches are filtered out much much later on.
140 
141  // Add mesh1 patches and build map both ways.
142  from1ToAllPatches.setSize(patches1.size());
143 
144  forAll(patches1, patchi)
145  {
146  from1ToAllPatches[patchi] = patchIndex
147  (
148  patches1[patchi],
149  allPatchNames,
150  allPatchTypes
151  );
152  }
153  allPatchTypes.shrink();
154  allPatchNames.shrink();
155 
156  // Invert 1 to all patch map
157  fromAllTo1Patches.setSize(allPatchNames.size());
158  fromAllTo1Patches = -1;
159 
160  forAll(from1ToAllPatches, i)
161  {
162  fromAllTo1Patches[from1ToAllPatches[i]] = i;
163  }
164 }
165 
166 
167 Foam::labelList Foam::polyMeshAdder::getPatchStarts
168 (
169  const polyBoundaryMesh& patches
170 )
171 {
172  labelList patchStarts(patches.size());
173  forAll(patches, patchi)
174  {
175  patchStarts[patchi] = patches[patchi].start();
176  }
177  return patchStarts;
178 }
179 
180 
181 Foam::labelList Foam::polyMeshAdder::getPatchSizes
182 (
183  const polyBoundaryMesh& patches
184 )
185 {
186  labelList patchSizes(patches.size());
187  forAll(patches, patchi)
188  {
189  patchSizes[patchi] = patches[patchi].size();
190  }
191  return patchSizes;
192 }
193 
194 
195 Foam::List<Foam::polyPatch*> Foam::polyMeshAdder::combinePatches
196 (
197  const polyMesh& mesh0,
198  const polyMesh& mesh1,
199  const polyBoundaryMesh& allBoundaryMesh,
200  const label nAllPatches,
201  const labelList& fromAllTo1Patches,
202 
203  const label nInternalFaces,
204  const labelList& nFaces,
205 
206  labelList& from0ToAllPatches,
207  labelList& from1ToAllPatches
208 )
209 {
210  const polyBoundaryMesh& patches0 = mesh0.boundaryMesh();
211  const polyBoundaryMesh& patches1 = mesh1.boundaryMesh();
212 
213 
214  // Compacted new patch list.
215  DynamicList<polyPatch*> allPatches(nAllPatches);
216 
217 
218  // Map from 0 to all patches (since gets compacted)
219  from0ToAllPatches.setSize(patches0.size());
220  from0ToAllPatches = -1;
221 
222  label startFacei = nInternalFaces;
223 
224  // Copy patches0 with new sizes. First patches always come from
225  // mesh0 and will always be present.
226  forAll(patches0, patchi)
227  {
228  // Originates from mesh0. Clone with new size & filter out empty
229  // patch.
230  label filteredPatchi;
231 
232  if (nFaces[patchi] == 0 && isA<processorPolyPatch>(patches0[patchi]))
233  {
234  // Pout<< "Removing zero sized mesh0 patch "
235  // << patches0[patchi].name() << endl;
236  filteredPatchi = -1;
237  }
238  else
239  {
240  filteredPatchi = allPatches.size();
241 
242  allPatches.append
243  (
244  patches0[patchi].clone
245  (
246  allBoundaryMesh,
247  filteredPatchi,
248  nFaces[patchi],
249  startFacei
250  ).ptr()
251  );
252  startFacei += nFaces[patchi];
253  }
254 
255  // Record new index in allPatches
256  from0ToAllPatches[patchi] = filteredPatchi;
257 
258  // Check if patch was also in mesh1 and update its addressing if so.
259  if (fromAllTo1Patches[patchi] != -1)
260  {
261  from1ToAllPatches[fromAllTo1Patches[patchi]] = filteredPatchi;
262  }
263  }
264 
265  // Copy unique patches of mesh1.
266  forAll(from1ToAllPatches, patchi)
267  {
268  label allPatchi = from1ToAllPatches[patchi];
269 
270  if (allPatchi >= patches0.size())
271  {
272  // Patch has not been merged with any mesh0 patch.
273 
274  label filteredPatchi;
275 
276  if
277  (
278  nFaces[allPatchi] == 0
279  && isA<processorPolyPatch>(patches1[patchi])
280  )
281  {
282  // Pout<< "Removing zero sized mesh1 patch "
283  // << patches1[patchi].name() << endl;
284  filteredPatchi = -1;
285  }
286  else
287  {
288  filteredPatchi = allPatches.size();
289 
290  allPatches.append
291  (
292  patches1[patchi].clone
293  (
294  allBoundaryMesh,
295  filteredPatchi,
296  nFaces[allPatchi],
297  startFacei
298  ).ptr()
299  );
300  startFacei += nFaces[allPatchi];
301  }
302 
303  from1ToAllPatches[patchi] = filteredPatchi;
304  }
305  }
306 
307  allPatches.shrink();
308 
309  return move(allPatches);
310 }
311 
312 
313 Foam::labelList Foam::polyMeshAdder::getFaceOrder
314 (
315  const cellList& cells,
316  const label nInternalFaces,
317  const labelList& owner,
318  const labelList& neighbour
319 )
320 {
321  labelList oldToNew(owner.size(), -1);
322 
323  // Leave boundary faces in order
324  for (label facei = nInternalFaces; facei < owner.size(); ++facei)
325  {
326  oldToNew[facei] = facei;
327  }
328 
329  // First unassigned face
330  label newFacei = 0;
331 
332  forAll(cells, celli)
333  {
334  const labelList& cFaces = cells[celli];
335 
336  SortableList<label> nbr(cFaces.size());
337 
338  forAll(cFaces, i)
339  {
340  label facei = cFaces[i];
341 
342  label nbrCelli = neighbour[facei];
343 
344  if (nbrCelli != -1)
345  {
346  // Internal face. Get cell on other side.
347  if (nbrCelli == celli)
348  {
349  nbrCelli = owner[facei];
350  }
351 
352  if (celli < nbrCelli)
353  {
354  // Celli is master
355  nbr[i] = nbrCelli;
356  }
357  else
358  {
359  // nbrCell is master. Let it handle this face.
360  nbr[i] = -1;
361  }
362  }
363  else
364  {
365  // External face. Do later.
366  nbr[i] = -1;
367  }
368  }
369 
370  nbr.sort();
371 
372  forAll(nbr, i)
373  {
374  if (nbr[i] != -1)
375  {
376  oldToNew[cFaces[nbr.indices()[i]]] = newFacei++;
377  }
378  }
379  }
380 
381 
382  // Check done all faces.
383  forAll(oldToNew, facei)
384  {
385  if (oldToNew[facei] == -1)
386  {
388  << "Did not determine new position"
389  << " for face " << facei
390  << abort(FatalError);
391  }
392  }
393 
394  return oldToNew;
395 }
396 
397 
398 // Extends face f with split points. cutEdgeToPoints gives for every
399 // edge the points introduced in between the endpoints.
400 void Foam::polyMeshAdder::insertVertices
401 (
402  const edgeLookup& cutEdgeToPoints,
403  const Map<label>& meshToMaster,
404  const labelList& masterToCutPoints,
405  const face& masterF,
406 
407  DynamicList<label>& workFace,
408  face& allF
409 )
410 {
411  workFace.clear();
412 
413  // Check any edge for being cut (check on the cut so takes account
414  // for any point merging on the cut)
415 
416  forAll(masterF, fp)
417  {
418  label v0 = masterF[fp];
419  label v1 = masterF.nextLabel(fp);
420 
421  // Copy existing face point
422  workFace.append(allF[fp]);
423 
424  // See if any edge between v0,v1
425 
426  Map<label>::const_iterator v0Fnd = meshToMaster.find(v0);
427  if (v0Fnd != meshToMaster.end())
428  {
429  Map<label>::const_iterator v1Fnd = meshToMaster.find(v1);
430  if (v1Fnd != meshToMaster.end())
431  {
432  // Get edge in cutPoint numbering
433  edge cutEdge
434  (
435  masterToCutPoints[v0Fnd()],
436  masterToCutPoints[v1Fnd()]
437  );
438 
439  edgeLookup::const_iterator iter = cutEdgeToPoints.find(cutEdge);
440 
441  if (iter != cutEdgeToPoints.end())
442  {
443  const edge& e = iter.key();
444  const labelList& addedPoints = iter();
445 
446  // cutPoints first in allPoints so no need for renumbering
447  if (e[0] == cutEdge[0])
448  {
449  forAll(addedPoints, i)
450  {
451  workFace.append(addedPoints[i]);
452  }
453  }
454  else
455  {
456  forAllReverse(addedPoints, i)
457  {
458  workFace.append(addedPoints[i]);
459  }
460  }
461  }
462  }
463  }
464  }
465 
466  if (workFace.size() != allF.size())
467  {
468  allF.transfer(workFace);
469  }
470 }
471 
472 
473 // Adds primitives (cells, faces, points)
474 // Cells:
475 // - all of mesh0
476 // - all of mesh1
477 // Faces:
478 // - all uncoupled of mesh0
479 // - all coupled faces
480 // - all uncoupled of mesh1
481 // Points:
482 // - all coupled
483 // - all uncoupled of mesh0
484 // - all uncoupled of mesh1
485 void Foam::polyMeshAdder::mergePrimitives
486 (
487  const polyMesh& mesh0,
488  const polyMesh& mesh1,
489  const faceCoupleInfo& coupleInfo,
490 
491  const label nAllPatches, // number of patches in the new mesh
492  const labelList& fromAllTo1Patches,
493  const labelList& from1ToAllPatches,
494 
495  pointField& allPoints,
496  labelList& from0ToAllPoints,
497  labelList& from1ToAllPoints,
498 
499  faceList& allFaces,
500  labelList& allOwner,
501  labelList& allNeighbour,
502  label& nInternalFaces,
503  labelList& nFacesPerPatch,
504  label& nCells,
505 
506  labelList& from0ToAllFaces,
507  labelList& from1ToAllFaces,
508  labelList& from1ToAllCells
509 )
510 {
511  const polyBoundaryMesh& patches0 = mesh0.boundaryMesh();
512  const polyBoundaryMesh& patches1 = mesh1.boundaryMesh();
513 
514  const primitiveFacePatch& cutFaces = coupleInfo.cutFaces();
515  const indirectPrimitivePatch& masterPatch = coupleInfo.masterPatch();
516  const indirectPrimitivePatch& slavePatch = coupleInfo.slavePatch();
517 
518 
519  // Points
520  // ~~~~~~
521 
522  // Storage for new points
523  allPoints.setSize(mesh0.nPoints() + mesh1.nPoints());
524  label allPointi = 0;
525 
526  from0ToAllPoints.setSize(mesh0.nPoints());
527  from0ToAllPoints = -1;
528  from1ToAllPoints.setSize(mesh1.nPoints());
529  from1ToAllPoints = -1;
530 
531  // Copy coupled points (on cut)
532  {
533  const pointField& cutPoints = coupleInfo.cutPoints();
534 
535  // const labelListList& cutToMasterPoints =
536  // coupleInfo.cutToMasterPoints();
537  labelListList cutToMasterPoints
538  (
540  (
541  cutPoints.size(),
542  coupleInfo.masterToCutPoints()
543  )
544  );
545 
546  // const labelListList& cutToSlavePoints =
547  // coupleInfo.cutToSlavePoints();
548  labelListList cutToSlavePoints
549  (
551  (
552  cutPoints.size(),
553  coupleInfo.slaveToCutPoints()
554  )
555  );
556 
557  forAll(cutPoints, i)
558  {
559  allPoints[allPointi] = cutPoints[i];
560 
561  // Mark all master and slave points referring to this point.
562 
563  const labelList& masterPoints = cutToMasterPoints[i];
564 
565  forAll(masterPoints, j)
566  {
567  label mesh0Pointi = masterPatch.meshPoints()[masterPoints[j]];
568  from0ToAllPoints[mesh0Pointi] = allPointi;
569  }
570 
571  const labelList& slavePoints = cutToSlavePoints[i];
572 
573  forAll(slavePoints, j)
574  {
575  label mesh1Pointi = slavePatch.meshPoints()[slavePoints[j]];
576  from1ToAllPoints[mesh1Pointi] = allPointi;
577  }
578  allPointi++;
579  }
580  }
581 
582  // Add uncoupled mesh0 points
583  forAll(mesh0.points(), pointi)
584  {
585  if (from0ToAllPoints[pointi] == -1)
586  {
587  allPoints[allPointi] = mesh0.points()[pointi];
588  from0ToAllPoints[pointi] = allPointi;
589  allPointi++;
590  }
591  }
592 
593  // Add uncoupled mesh1 points
594  forAll(mesh1.points(), pointi)
595  {
596  if (from1ToAllPoints[pointi] == -1)
597  {
598  allPoints[allPointi] = mesh1.points()[pointi];
599  from1ToAllPoints[pointi] = allPointi;
600  allPointi++;
601  }
602  }
603 
604  allPoints.setSize(allPointi);
605 
606 
607  // Faces
608  // ~~~~~
609 
610  // Sizes per patch
611  nFacesPerPatch.setSize(nAllPatches);
612  nFacesPerPatch = 0;
613 
614  // Storage for faces and owner/neighbour
615  allFaces.setSize(mesh0.nFaces() + mesh1.nFaces());
616  allOwner.setSize(allFaces.size());
617  allOwner = -1;
618  allNeighbour.setSize(allFaces.size());
619  allNeighbour = -1;
620  label allFacei = 0;
621 
622  from0ToAllFaces.setSize(mesh0.nFaces());
623  from0ToAllFaces = -1;
624  from1ToAllFaces.setSize(mesh1.nFaces());
625  from1ToAllFaces = -1;
626 
627  // Copy mesh0 internal faces (always uncoupled)
628  for (label facei = 0; facei < mesh0.nInternalFaces(); facei++)
629  {
630  allFaces[allFacei] = renumber(from0ToAllPoints, mesh0.faces()[facei]);
631  allOwner[allFacei] = mesh0.faceOwner()[facei];
632  allNeighbour[allFacei] = mesh0.faceNeighbour()[facei];
633  from0ToAllFaces[facei] = allFacei++;
634  }
635 
636  // Copy coupled faces. Every coupled face has an equivalent master and
637  // slave. Also uncount as boundary faces all the newly coupled faces.
638  const labelList& cutToMasterFaces = coupleInfo.cutToMasterFaces();
639  const labelList& cutToSlaveFaces = coupleInfo.cutToSlaveFaces();
640 
641  forAll(cutFaces, i)
642  {
643  label masterFacei = cutToMasterFaces[i];
644 
645  label mesh0Facei = masterPatch.addressing()[masterFacei];
646 
647  if (from0ToAllFaces[mesh0Facei] == -1)
648  {
649  // First occurrence of face
650  from0ToAllFaces[mesh0Facei] = allFacei;
651 
652  // External face becomes internal so uncount
653  label patch0 = patches0.whichPatch(mesh0Facei);
654  nFacesPerPatch[patch0]--;
655  }
656 
657  label slaveFacei = cutToSlaveFaces[i];
658 
659  label mesh1Facei = slavePatch.addressing()[slaveFacei];
660 
661  if (from1ToAllFaces[mesh1Facei] == -1)
662  {
663  from1ToAllFaces[mesh1Facei] = allFacei;
664 
665  label patch1 = patches1.whichPatch(mesh1Facei);
666  nFacesPerPatch[from1ToAllPatches[patch1]]--;
667  }
668 
669  // Copy cut face (since cutPoints are copied first no renumbering
670  // necessary)
671  allFaces[allFacei] = cutFaces[i];
672  allOwner[allFacei] = mesh0.faceOwner()[mesh0Facei];
673  allNeighbour[allFacei] = mesh1.faceOwner()[mesh1Facei] + mesh0.nCells();
674 
675  allFacei++;
676  }
677 
678  // Copy mesh1 internal faces (always uncoupled)
679  for (label facei = 0; facei < mesh1.nInternalFaces(); facei++)
680  {
681  allFaces[allFacei] = renumber(from1ToAllPoints, mesh1.faces()[facei]);
682  allOwner[allFacei] = mesh1.faceOwner()[facei] + mesh0.nCells();
683  allNeighbour[allFacei] = mesh1.faceNeighbour()[facei] + mesh0.nCells();
684  from1ToAllFaces[facei] = allFacei++;
685  }
686 
687  nInternalFaces = allFacei;
688 
689  // Copy (unmarked/uncoupled) external faces in new order.
690  for (label allPatchi = 0; allPatchi < nAllPatches; allPatchi++)
691  {
692  if (allPatchi < patches0.size())
693  {
694  // Patch is present in mesh0
695  const polyPatch& pp = patches0[allPatchi];
696 
697  nFacesPerPatch[allPatchi] += pp.size();
698 
699  label facei = pp.start();
700 
701  forAll(pp, i)
702  {
703  if (from0ToAllFaces[facei] == -1)
704  {
705  // Is uncoupled face since has not yet been dealt with
706  allFaces[allFacei] = renumber
707  (
708  from0ToAllPoints,
709  mesh0.faces()[facei]
710  );
711  allOwner[allFacei] = mesh0.faceOwner()[facei];
712  allNeighbour[allFacei] = -1;
713 
714  from0ToAllFaces[facei] = allFacei++;
715  }
716  facei++;
717  }
718  }
719  if (fromAllTo1Patches[allPatchi] != -1)
720  {
721  // Patch is present in mesh1
722  const polyPatch& pp = patches1[fromAllTo1Patches[allPatchi]];
723 
724  nFacesPerPatch[allPatchi] += pp.size();
725 
726  label facei = pp.start();
727 
728  forAll(pp, i)
729  {
730  if (from1ToAllFaces[facei] == -1)
731  {
732  // Is uncoupled face
733  allFaces[allFacei] = renumber
734  (
735  from1ToAllPoints,
736  mesh1.faces()[facei]
737  );
738  allOwner[allFacei] =
739  mesh1.faceOwner()[facei]
740  + mesh0.nCells();
741  allNeighbour[allFacei] = -1;
742 
743  from1ToAllFaces[facei] = allFacei++;
744  }
745  facei++;
746  }
747  }
748  }
749  allFaces.setSize(allFacei);
750  allOwner.setSize(allFacei);
751  allNeighbour.setSize(allFacei);
752 
753 
754  // So now we have all ok for one-to-one mapping.
755  // For split slace faces:
756  // - mesh consistent with slave side
757  // - mesh not consistent with owner side. It is not zipped up, the
758  // original faces need edges split.
759 
760  // Use brute force to prevent having to calculate addressing:
761  // - get map from master edge to split edges.
762  // - check all faces to find any edge that is split.
763  {
764  // From two cut-points to labels of cut-points in between.
765  // (in order: from e[0] to e[1]
766  const edgeLookup& cutEdgeToPoints = coupleInfo.cutEdgeToPoints();
767 
768  // Get map of master face (in mesh labels) that are in cut. These faces
769  // do not need to be renumbered.
770  labelHashSet masterCutFaces(cutToMasterFaces.size());
771  forAll(cutToMasterFaces, i)
772  {
773  label meshFacei = masterPatch.addressing()[cutToMasterFaces[i]];
774 
775  masterCutFaces.insert(meshFacei);
776  }
777 
778  DynamicList<label> workFace(100);
779 
780  forAll(from0ToAllFaces, face0)
781  {
782  if (!masterCutFaces.found(face0))
783  {
784  label allFacei = from0ToAllFaces[face0];
785 
786  insertVertices
787  (
788  cutEdgeToPoints,
789  masterPatch.meshPointMap(),
790  coupleInfo.masterToCutPoints(),
791  mesh0.faces()[face0],
792 
793  workFace,
794  allFaces[allFacei]
795  );
796  }
797  }
798 
799  // Same for slave face
800 
801  labelHashSet slaveCutFaces(cutToSlaveFaces.size());
802  forAll(cutToSlaveFaces, i)
803  {
804  label meshFacei = slavePatch.addressing()[cutToSlaveFaces[i]];
805 
806  slaveCutFaces.insert(meshFacei);
807  }
808 
809  forAll(from1ToAllFaces, face1)
810  {
811  if (!slaveCutFaces.found(face1))
812  {
813  label allFacei = from1ToAllFaces[face1];
814 
815  insertVertices
816  (
817  cutEdgeToPoints,
818  slavePatch.meshPointMap(),
819  coupleInfo.slaveToCutPoints(),
820  mesh1.faces()[face1],
821 
822  workFace,
823  allFaces[allFacei]
824  );
825  }
826  }
827  }
828 
829  // Now we have a full facelist and owner/neighbour addressing.
830 
831 
832  // Cells
833  // ~~~~~
834 
835  from1ToAllCells.setSize(mesh1.nCells());
836  from1ToAllCells = -1;
837 
838  forAll(mesh1.cells(), i)
839  {
840  from1ToAllCells[i] = i + mesh0.nCells();
841  }
842 
843  // Make cells (= cell-face addressing)
844  nCells = mesh0.nCells() + mesh1.nCells();
845  cellList allCells(nCells);
846  primitiveMesh::calcCells(allCells, allOwner, allNeighbour, nCells);
847 
848  // Reorder faces for upper-triangular order.
849  labelList oldToNew
850  (
851  getFaceOrder
852  (
853  allCells,
854  nInternalFaces,
855  allOwner,
856  allNeighbour
857  )
858  );
859 
860  inplaceReorder(oldToNew, allFaces);
861  inplaceReorder(oldToNew, allOwner);
862  inplaceReorder(oldToNew, allNeighbour);
863  inplaceRenumber(oldToNew, from0ToAllFaces);
864  inplaceRenumber(oldToNew, from1ToAllFaces);
865 }
866 
867 
868 void Foam::polyMeshAdder::mergePointZones
869 (
870  const label nAllPoints,
871  const pointZoneMesh& pz0,
872  const pointZoneMesh& pz1,
873  const labelList& from0ToAllPoints,
874  const labelList& from1ToAllPoints,
875 
876  DynamicList<word>& zoneNames,
877  labelList& from1ToAll,
878  List<DynamicList<label>>& pzPoints
879 )
880 {
881  zoneNames.setCapacity(pz0.size() + pz1.size());
882  zoneNames.append(pz0.names());
883 
884  from1ToAll.setSize(pz1.size());
885 
886  forAll(pz1, zoneI)
887  {
888  from1ToAll[zoneI] = zoneIndex(pz1[zoneI].name(), zoneNames);
889  }
890  zoneNames.shrink();
891 
892 
893 
894  // Zone(s) per point. Two levels: if only one zone
895  // stored in pointToZone. Any extra stored in additionalPointToZones.
896  // This is so we only allocate labelLists per point if absolutely
897  // necessary.
898  labelList pointToZone(nAllPoints, -1);
899  labelListList addPointToZones(nAllPoints);
900 
901  // mesh0 zones kept
902  forAll(pz0, zoneI)
903  {
904  const pointZone& pz = pz0[zoneI];
905 
906  forAll(pz, i)
907  {
908  label point0 = pz[i];
909  label allPointi = from0ToAllPoints[point0];
910 
911  if (pointToZone[allPointi] == -1)
912  {
913  pointToZone[allPointi] = zoneI;
914  }
915  else if (pointToZone[allPointi] != zoneI)
916  {
917  labelList& pZones = addPointToZones[allPointi];
918  if (findIndex(pZones, zoneI) == -1)
919  {
920  pZones.append(zoneI);
921  }
922  }
923  }
924  }
925 
926  // mesh1 zones renumbered
927  forAll(pz1, zoneI)
928  {
929  const pointZone& pz = pz1[zoneI];
930  const label allZoneI = from1ToAll[zoneI];
931 
932  forAll(pz, i)
933  {
934  label point1 = pz[i];
935  label allPointi = from1ToAllPoints[point1];
936 
937  if (pointToZone[allPointi] == -1)
938  {
939  pointToZone[allPointi] = allZoneI;
940  }
941  else if (pointToZone[allPointi] != allZoneI)
942  {
943  labelList& pZones = addPointToZones[allPointi];
944  if (findIndex(pZones, allZoneI) == -1)
945  {
946  pZones.append(allZoneI);
947  }
948  }
949  }
950  }
951 
952  // Extract back into zones
953 
954  // 1. Count
955  labelList nPoints(zoneNames.size(), 0);
956  forAll(pointToZone, allPointi)
957  {
958  label zoneI = pointToZone[allPointi];
959  if (zoneI != -1)
960  {
961  nPoints[zoneI]++;
962  }
963  }
964  forAll(addPointToZones, allPointi)
965  {
966  const labelList& pZones = addPointToZones[allPointi];
967  forAll(pZones, i)
968  {
969  nPoints[pZones[i]]++;
970  }
971  }
972 
973  // 2. Fill
974  pzPoints.setSize(zoneNames.size());
975  forAll(pzPoints, zoneI)
976  {
977  pzPoints[zoneI].setCapacity(nPoints[zoneI]);
978  }
979  forAll(pointToZone, allPointi)
980  {
981  label zoneI = pointToZone[allPointi];
982  if (zoneI != -1)
983  {
984  pzPoints[zoneI].append(allPointi);
985  }
986  }
987  forAll(addPointToZones, allPointi)
988  {
989  const labelList& pZones = addPointToZones[allPointi];
990  forAll(pZones, i)
991  {
992  pzPoints[pZones[i]].append(allPointi);
993  }
994  }
995  forAll(pzPoints, i)
996  {
997  pzPoints[i].shrink();
998  stableSort(pzPoints[i]);
999  }
1000 }
1001 
1002 
1003 void Foam::polyMeshAdder::mergeFaceZones
1004 (
1005  const labelList& allOwner,
1006 
1007  const polyMesh& mesh0,
1008  const polyMesh& mesh1,
1009 
1010  const labelList& from0ToAllFaces,
1011  const labelList& from1ToAllFaces,
1012 
1013  const labelList& from1ToAllCells,
1014 
1015  DynamicList<word>& zoneNames,
1016  labelList& from1ToAll,
1017  List<DynamicList<label>>& fzFaces,
1018  List<DynamicList<bool>>& fzFlips
1019 )
1020 {
1021  const faceZoneMesh& fz0 = mesh0.faceZones();
1022  const labelList& owner0 = mesh0.faceOwner();
1023  const faceZoneMesh& fz1 = mesh1.faceZones();
1024  const labelList& owner1 = mesh1.faceOwner();
1025 
1026 
1027  zoneNames.setCapacity(fz0.size() + fz1.size());
1028  zoneNames.append(fz0.names());
1029 
1030  from1ToAll.setSize(fz1.size());
1031 
1032  forAll(fz1, zoneI)
1033  {
1034  from1ToAll[zoneI] = zoneIndex(fz1[zoneI].name(), zoneNames);
1035  }
1036  zoneNames.shrink();
1037 
1038 
1039  // Zone(s) per face
1040  labelList faceToZone(allOwner.size(), -1);
1041  labelListList addFaceToZones(allOwner.size());
1042  boolList faceToFlip(allOwner.size(), false);
1043  boolListList addFaceToFlips(allOwner.size());
1044 
1045  // mesh0 zones kept
1046  forAll(fz0, zoneI)
1047  {
1048  const labelList& addressing = fz0[zoneI];
1049  const boolList& flipMap = fz0[zoneI].flipMap();
1050 
1051  forAll(addressing, i)
1052  {
1053  label face0 = addressing[i];
1054  bool flip0 = flipMap[i];
1055 
1056  label allFacei = from0ToAllFaces[face0];
1057  if (allFacei != -1)
1058  {
1059  // Check if orientation same
1060  label allCell0 = owner0[face0];
1061  if (allOwner[allFacei] != allCell0)
1062  {
1063  flip0 = !flip0;
1064  }
1065 
1066  if (faceToZone[allFacei] == -1)
1067  {
1068  faceToZone[allFacei] = zoneI;
1069  faceToFlip[allFacei] = flip0;
1070  }
1071  else if (faceToZone[allFacei] != zoneI)
1072  {
1073  labelList& fZones = addFaceToZones[allFacei];
1074  boolList& flipZones = addFaceToFlips[allFacei];
1075 
1076  if (findIndex(fZones, zoneI) == -1)
1077  {
1078  fZones.append(zoneI);
1079  flipZones.append(flip0);
1080  }
1081  }
1082  }
1083  }
1084  }
1085 
1086  // mesh1 zones renumbered
1087  forAll(fz1, zoneI)
1088  {
1089  const labelList& addressing = fz1[zoneI];
1090  const boolList& flipMap = fz1[zoneI].flipMap();
1091 
1092  const label allZoneI = from1ToAll[zoneI];
1093 
1094  forAll(addressing, i)
1095  {
1096  label face1 = addressing[i];
1097  bool flip1 = flipMap[i];
1098 
1099  label allFacei = from1ToAllFaces[face1];
1100  if (allFacei != -1)
1101  {
1102  // Check if orientation same
1103  label allCell1 = from1ToAllCells[owner1[face1]];
1104  if (allOwner[allFacei] != allCell1)
1105  {
1106  flip1 = !flip1;
1107  }
1108 
1109  if (faceToZone[allFacei] == -1)
1110  {
1111  faceToZone[allFacei] = allZoneI;
1112  faceToFlip[allFacei] = flip1;
1113  }
1114  else if (faceToZone[allFacei] != allZoneI)
1115  {
1116  labelList& fZones = addFaceToZones[allFacei];
1117  boolList& flipZones = addFaceToFlips[allFacei];
1118 
1119  if (findIndex(fZones, allZoneI) == -1)
1120  {
1121  fZones.append(allZoneI);
1122  flipZones.append(flip1);
1123  }
1124  }
1125  }
1126  }
1127  }
1128 
1129 
1130  // Extract back into zones
1131 
1132  // 1. Count
1133  labelList nFaces(zoneNames.size(), 0);
1134  forAll(faceToZone, allFacei)
1135  {
1136  label zoneI = faceToZone[allFacei];
1137  if (zoneI != -1)
1138  {
1139  nFaces[zoneI]++;
1140  }
1141  }
1142  forAll(addFaceToZones, allFacei)
1143  {
1144  const labelList& fZones = addFaceToZones[allFacei];
1145  forAll(fZones, i)
1146  {
1147  nFaces[fZones[i]]++;
1148  }
1149  }
1150 
1151  // 2. Fill
1152  fzFaces.setSize(zoneNames.size());
1153  fzFlips.setSize(zoneNames.size());
1154  forAll(fzFaces, zoneI)
1155  {
1156  fzFaces[zoneI].setCapacity(nFaces[zoneI]);
1157  fzFlips[zoneI].setCapacity(nFaces[zoneI]);
1158  }
1159  forAll(faceToZone, allFacei)
1160  {
1161  label zoneI = faceToZone[allFacei];
1162  bool flip = faceToFlip[allFacei];
1163  if (zoneI != -1)
1164  {
1165  fzFaces[zoneI].append(allFacei);
1166  fzFlips[zoneI].append(flip);
1167  }
1168  }
1169  forAll(addFaceToZones, allFacei)
1170  {
1171  const labelList& fZones = addFaceToZones[allFacei];
1172  const boolList& flipZones = addFaceToFlips[allFacei];
1173 
1174  forAll(fZones, i)
1175  {
1176  label zoneI = fZones[i];
1177  fzFaces[zoneI].append(allFacei);
1178  fzFlips[zoneI].append(flipZones[i]);
1179  }
1180  }
1181 
1182  forAll(fzFaces, i)
1183  {
1184  fzFaces[i].shrink();
1185  fzFlips[i].shrink();
1186 
1187  labelList order;
1188  sortedOrder(fzFaces[i], order);
1189  fzFaces[i] = UIndirectList<label>(fzFaces[i], order)();
1190  fzFlips[i] = UIndirectList<bool>(fzFlips[i], order)();
1191  }
1192 }
1193 
1194 
1195 void Foam::polyMeshAdder::mergeCellZones
1196 (
1197  const label nAllCells,
1198 
1199  const cellZoneMesh& cz0,
1200  const cellZoneMesh& cz1,
1201  const labelList& from1ToAllCells,
1202 
1203  DynamicList<word>& zoneNames,
1204  labelList& from1ToAll,
1205  List<DynamicList<label>>& czCells
1206 )
1207 {
1208  zoneNames.setCapacity(cz0.size() + cz1.size());
1209  zoneNames.append(cz0.names());
1210 
1211  from1ToAll.setSize(cz1.size());
1212  forAll(cz1, zoneI)
1213  {
1214  from1ToAll[zoneI] = zoneIndex(cz1[zoneI].name(), zoneNames);
1215  }
1216  zoneNames.shrink();
1217 
1218 
1219  // Zone(s) per cell. Two levels: if only one zone
1220  // stored in cellToZone. Any extra stored in additionalCellToZones.
1221  // This is so we only allocate labelLists per cell if absolutely
1222  // necessary.
1223  labelList cellToZone(nAllCells, -1);
1224  labelListList addCellToZones(nAllCells);
1225 
1226  // mesh0 zones kept
1227  forAll(cz0, zoneI)
1228  {
1229  const cellZone& cz = cz0[zoneI];
1230  forAll(cz, i)
1231  {
1232  label cell0 = cz[i];
1233 
1234  if (cellToZone[cell0] == -1)
1235  {
1236  cellToZone[cell0] = zoneI;
1237  }
1238  else if (cellToZone[cell0] != zoneI)
1239  {
1240  labelList& cZones = addCellToZones[cell0];
1241  if (findIndex(cZones, zoneI) == -1)
1242  {
1243  cZones.append(zoneI);
1244  }
1245  }
1246  }
1247  }
1248 
1249  // mesh1 zones renumbered
1250  forAll(cz1, zoneI)
1251  {
1252  const cellZone& cz = cz1[zoneI];
1253  const label allZoneI = from1ToAll[zoneI];
1254  forAll(cz, i)
1255  {
1256  label cell1 = cz[i];
1257  label allCelli = from1ToAllCells[cell1];
1258 
1259  if (cellToZone[allCelli] == -1)
1260  {
1261  cellToZone[allCelli] = allZoneI;
1262  }
1263  else if (cellToZone[allCelli] != allZoneI)
1264  {
1265  labelList& cZones = addCellToZones[allCelli];
1266  if (findIndex(cZones, allZoneI) == -1)
1267  {
1268  cZones.append(allZoneI);
1269  }
1270  }
1271  }
1272  }
1273 
1274  // Extract back into zones
1275 
1276  // 1. Count
1277  labelList nCells(zoneNames.size(), 0);
1278  forAll(cellToZone, allCelli)
1279  {
1280  label zoneI = cellToZone[allCelli];
1281  if (zoneI != -1)
1282  {
1283  nCells[zoneI]++;
1284  }
1285  }
1286  forAll(addCellToZones, allCelli)
1287  {
1288  const labelList& cZones = addCellToZones[allCelli];
1289  forAll(cZones, i)
1290  {
1291  nCells[cZones[i]]++;
1292  }
1293  }
1294 
1295  // 2. Fill
1296  czCells.setSize(zoneNames.size());
1297  forAll(czCells, zoneI)
1298  {
1299  czCells[zoneI].setCapacity(nCells[zoneI]);
1300  }
1301  forAll(cellToZone, allCelli)
1302  {
1303  label zoneI = cellToZone[allCelli];
1304  if (zoneI != -1)
1305  {
1306  czCells[zoneI].append(allCelli);
1307  }
1308  }
1309  forAll(addCellToZones, allCelli)
1310  {
1311  const labelList& cZones = addCellToZones[allCelli];
1312  forAll(cZones, i)
1313  {
1314  czCells[cZones[i]].append(allCelli);
1315  }
1316  }
1317  forAll(czCells, i)
1318  {
1319  czCells[i].shrink();
1320  stableSort(czCells[i]);
1321  }
1322 }
1323 
1324 
1325 void Foam::polyMeshAdder::mergeZones
1326 (
1327  const label nAllPoints,
1328  const labelList& allOwner,
1329  const label nAllCells,
1330 
1331  const polyMesh& mesh0,
1332  const polyMesh& mesh1,
1333  const labelList& from0ToAllPoints,
1334  const labelList& from0ToAllFaces,
1335 
1336  const labelList& from1ToAllPoints,
1337  const labelList& from1ToAllFaces,
1338  const labelList& from1ToAllCells,
1339 
1340  DynamicList<word>& pointZoneNames,
1341  List<DynamicList<label>>& pzPoints,
1342 
1343  DynamicList<word>& faceZoneNames,
1344  List<DynamicList<label>>& fzFaces,
1345  List<DynamicList<bool>>& fzFlips,
1346 
1347  DynamicList<word>& cellZoneNames,
1348  List<DynamicList<label>>& czCells
1349 )
1350 {
1351  labelList from1ToAllPZones;
1352  mergePointZones
1353  (
1354  nAllPoints,
1355  mesh0.pointZones(),
1356  mesh1.pointZones(),
1357  from0ToAllPoints,
1358  from1ToAllPoints,
1359 
1360  pointZoneNames,
1361  from1ToAllPZones,
1362  pzPoints
1363  );
1364 
1365  labelList from1ToAllFZones;
1366  mergeFaceZones
1367  (
1368  allOwner,
1369  mesh0,
1370  mesh1,
1371  from0ToAllFaces,
1372  from1ToAllFaces,
1373  from1ToAllCells,
1374 
1375  faceZoneNames,
1376  from1ToAllFZones,
1377  fzFaces,
1378  fzFlips
1379  );
1380 
1381  labelList from1ToAllCZones;
1382  mergeCellZones
1383  (
1384  nAllCells,
1385  mesh0.cellZones(),
1386  mesh1.cellZones(),
1387  from1ToAllCells,
1388 
1389  cellZoneNames,
1390  from1ToAllCZones,
1391  czCells
1392  );
1393 }
1394 
1395 
1396 void Foam::polyMeshAdder::addZones
1397 (
1398  const DynamicList<word>& pointZoneNames,
1399  const List<DynamicList<label>>& pzPoints,
1400 
1401  const DynamicList<word>& faceZoneNames,
1402  const List<DynamicList<label>>& fzFaces,
1403  const List<DynamicList<bool>>& fzFlips,
1404 
1405  const DynamicList<word>& cellZoneNames,
1406  const List<DynamicList<label>>& czCells,
1407 
1408  polyMesh& mesh
1409 )
1410 {
1411  List<pointZone*> pZones(pzPoints.size());
1412  forAll(pZones, i)
1413  {
1414  pZones[i] = new pointZone
1415  (
1416  pointZoneNames[i],
1417  pzPoints[i],
1418  i,
1419  mesh.pointZones()
1420  );
1421  }
1422 
1423  List<faceZone*> fZones(fzFaces.size());
1424  forAll(fZones, i)
1425  {
1426  fZones[i] = new faceZone
1427  (
1428  faceZoneNames[i],
1429  fzFaces[i],
1430  fzFlips[i],
1431  i,
1432  mesh.faceZones()
1433  );
1434  }
1435 
1436  List<cellZone*> cZones(czCells.size());
1437  forAll(cZones, i)
1438  {
1439  cZones[i] = new cellZone
1440  (
1441  cellZoneNames[i],
1442  czCells[i],
1443  i,
1444  mesh.cellZones()
1445  );
1446  }
1447 
1448  mesh.addZones
1449  (
1450  pZones,
1451  fZones,
1452  cZones
1453  );
1454 }
1455 
1456 
1457 
1458 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1459 
1460 // Returns new mesh and sets
1461 // - map from new cell/face/point/patch to either mesh0 or mesh1
1462 //
1463 // mesh0Faces/mesh1Faces: corresponding faces on both meshes.
1466  const IOobject& io,
1467  const polyMesh& mesh0,
1468  const polyMesh& mesh1,
1469  const faceCoupleInfo& coupleInfo,
1471 )
1472 {
1473  const polyBoundaryMesh& patches0 = mesh0.boundaryMesh();
1474  const polyBoundaryMesh& patches1 = mesh1.boundaryMesh();
1475 
1476 
1477  DynamicList<word> allPatchNames(patches0.size() + patches1.size());
1478  DynamicList<word> allPatchTypes(allPatchNames.size());
1479 
1480  // Patch maps
1481  labelList from1ToAllPatches(patches1.size());
1482  labelList fromAllTo1Patches(allPatchNames.size(), -1);
1483 
1484  mergePatchNames
1485  (
1486  patches0,
1487  patches1,
1488  allPatchNames,
1489  allPatchTypes,
1490  from1ToAllPatches,
1491  fromAllTo1Patches
1492  );
1493 
1494 
1495  // New points
1497 
1498  // Map from mesh0/1 points to allPoints.
1499  labelList from0ToAllPoints(mesh0.nPoints(), -1);
1500  labelList from1ToAllPoints(mesh1.nPoints(), -1);
1501 
1502  // New faces
1503  faceList allFaces;
1504  label nInternalFaces;
1505 
1506  // New cells
1507  labelList allOwner;
1508  labelList allNeighbour;
1509  label nCells;
1510 
1511  // Sizes per patch
1512  labelList nFaces(allPatchNames.size(), 0);
1513 
1514  // Maps
1515  labelList from0ToAllFaces(mesh0.nFaces(), -1);
1516  labelList from1ToAllFaces(mesh1.nFaces(), -1);
1517 
1518  // Map
1519  labelList from1ToAllCells(mesh1.nCells(), -1);
1520 
1521  mergePrimitives
1522  (
1523  mesh0,
1524  mesh1,
1525  coupleInfo,
1526 
1527  allPatchNames.size(),
1528  fromAllTo1Patches,
1529  from1ToAllPatches,
1530 
1531  allPoints,
1532  from0ToAllPoints,
1533  from1ToAllPoints,
1534 
1535  allFaces,
1536  allOwner,
1537  allNeighbour,
1538  nInternalFaces,
1539  nFaces,
1540  nCells,
1541 
1542  from0ToAllFaces,
1543  from1ToAllFaces,
1544  from1ToAllCells
1545  );
1546 
1547 
1548  // Zones
1549  // ~~~~~
1550 
1551  DynamicList<word> pointZoneNames;
1552  List<DynamicList<label>> pzPoints;
1553 
1554  DynamicList<word> faceZoneNames;
1555  List<DynamicList<label>> fzFaces;
1556  List<DynamicList<bool>> fzFlips;
1557 
1558  DynamicList<word> cellZoneNames;
1559  List<DynamicList<label>> czCells;
1560 
1561  mergeZones
1562  (
1563  allPoints.size(),
1564  allOwner,
1565  nCells,
1566 
1567  mesh0,
1568  mesh1,
1569 
1570  from0ToAllPoints,
1571  from0ToAllFaces,
1572 
1573  from1ToAllPoints,
1574  from1ToAllFaces,
1575  from1ToAllCells,
1576 
1577  pointZoneNames,
1578  pzPoints,
1579 
1580  faceZoneNames,
1581  fzFaces,
1582  fzFlips,
1583 
1584  cellZoneNames,
1585  czCells
1586  );
1587 
1588 
1589  // Patches
1590  // ~~~~~~~
1591 
1592  // Map from 0 to all patches (since gets compacted)
1593  labelList from0ToAllPatches(patches0.size(), -1);
1594 
1595  List<polyPatch*> allPatches
1596  (
1597  combinePatches
1598  (
1599  mesh0,
1600  mesh1,
1601  patches0, // Should be boundaryMesh() on new mesh.
1602  allPatchNames.size(),
1603  fromAllTo1Patches,
1604  mesh0.nInternalFaces()
1605  + mesh1.nInternalFaces()
1606  + coupleInfo.cutFaces().size(),
1607  nFaces,
1608 
1609  from0ToAllPatches,
1610  from1ToAllPatches
1611  )
1612  );
1613 
1614 
1615  // Map information
1616  // ~~~~~~~~~~~~~~~
1617 
1618  mapPtr.reset
1619  (
1620  new mapAddedPolyMesh
1621  (
1622  mesh0.nPoints(),
1623  mesh0.nFaces(),
1624  mesh0.nCells(),
1625 
1626  mesh1.nPoints(),
1627  mesh1.nFaces(),
1628  mesh1.nCells(),
1629 
1630  from0ToAllPoints,
1631  from0ToAllFaces,
1632  identity(mesh0.nCells()),
1633 
1634  from1ToAllPoints,
1635  from1ToAllFaces,
1636  from1ToAllCells,
1637 
1638  from0ToAllPatches,
1639  from1ToAllPatches,
1640  getPatchSizes(patches0),
1641  getPatchStarts(patches0)
1642  )
1643  );
1644 
1645 
1646 
1647  // Now we have extracted all information from all meshes.
1648  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1649 
1650  // Construct mesh
1651  autoPtr<polyMesh> tmesh
1652  (
1653  new polyMesh
1654  (
1655  io,
1656  move(allPoints),
1657  move(allFaces),
1658  move(allOwner),
1659  move(allNeighbour)
1660  )
1661  );
1662  polyMesh& mesh = tmesh();
1663 
1664  // Add zones to new mesh.
1665  addZones
1666  (
1667  pointZoneNames,
1668  pzPoints,
1669 
1670  faceZoneNames,
1671  fzFaces,
1672  fzFlips,
1673 
1674  cellZoneNames,
1675  czCells,
1676  mesh
1677  );
1678 
1679  // Add patches to new mesh
1680  mesh.addPatches(allPatches);
1681 
1682  return tmesh;
1683 }
1684 
1685 
1686 // Inplace add mesh1 to mesh0
1689  polyMesh& mesh0,
1690  const polyMesh& mesh1,
1691  const faceCoupleInfo& coupleInfo,
1692  const bool validBoundary
1693 )
1694 {
1695  const polyBoundaryMesh& patches0 = mesh0.boundaryMesh();
1696  const polyBoundaryMesh& patches1 = mesh1.boundaryMesh();
1697 
1698  DynamicList<word> allPatchNames(patches0.size() + patches1.size());
1699  DynamicList<word> allPatchTypes(allPatchNames.size());
1700 
1701  // Patch maps
1702  labelList from1ToAllPatches(patches1.size());
1703  labelList fromAllTo1Patches(allPatchNames.size(), -1);
1704 
1705  mergePatchNames
1706  (
1707  patches0,
1708  patches1,
1709  allPatchNames,
1710  allPatchTypes,
1711  from1ToAllPatches,
1712  fromAllTo1Patches
1713  );
1714 
1715 
1716  // New points
1718 
1719  // Map from mesh0/1 points to allPoints.
1720  labelList from0ToAllPoints(mesh0.nPoints(), -1);
1721  labelList from1ToAllPoints(mesh1.nPoints(), -1);
1722 
1723  // New faces
1724  faceList allFaces;
1725  labelList allOwner;
1726  labelList allNeighbour;
1727  label nInternalFaces;
1728  // Sizes per patch
1729  labelList nFaces(allPatchNames.size(), 0);
1730  label nCells;
1731 
1732  // Maps
1733  labelList from0ToAllFaces(mesh0.nFaces(), -1);
1734  labelList from1ToAllFaces(mesh1.nFaces(), -1);
1735  // Map
1736  labelList from1ToAllCells(mesh1.nCells(), -1);
1737 
1738  mergePrimitives
1739  (
1740  mesh0,
1741  mesh1,
1742  coupleInfo,
1743 
1744  allPatchNames.size(),
1745  fromAllTo1Patches,
1746  from1ToAllPatches,
1747 
1748  allPoints,
1749  from0ToAllPoints,
1750  from1ToAllPoints,
1751 
1752  allFaces,
1753  allOwner,
1754  allNeighbour,
1755  nInternalFaces,
1756  nFaces,
1757  nCells,
1758 
1759  from0ToAllFaces,
1760  from1ToAllFaces,
1761  from1ToAllCells
1762  );
1763 
1764 
1765  // Zones
1766  // ~~~~~
1767 
1768  DynamicList<word> pointZoneNames;
1769  List<DynamicList<label>> pzPoints;
1770 
1771  DynamicList<word> faceZoneNames;
1772  List<DynamicList<label>> fzFaces;
1773  List<DynamicList<bool>> fzFlips;
1774 
1775  DynamicList<word> cellZoneNames;
1776  List<DynamicList<label>> czCells;
1777 
1778  mergeZones
1779  (
1780  allPoints.size(),
1781  allOwner,
1782  nCells,
1783 
1784  mesh0,
1785  mesh1,
1786 
1787  from0ToAllPoints,
1788  from0ToAllFaces,
1789 
1790  from1ToAllPoints,
1791  from1ToAllFaces,
1792  from1ToAllCells,
1793 
1794  pointZoneNames,
1795  pzPoints,
1796 
1797  faceZoneNames,
1798  fzFaces,
1799  fzFlips,
1800 
1801  cellZoneNames,
1802  czCells
1803  );
1804 
1805 
1806  // Patches
1807  // ~~~~~~~
1808 
1809 
1810  // Store mesh0 patch info before modifying patches0.
1811  labelList mesh0PatchSizes(getPatchSizes(patches0));
1812  labelList mesh0PatchStarts(getPatchStarts(patches0));
1813 
1814  // Map from 0 to all patches (since gets compacted)
1815  labelList from0ToAllPatches(patches0.size(), -1);
1816 
1817  // Inplace extend mesh0 patches (note that patches0.size() now also
1818  // has changed)
1819  labelList patchSizes(allPatchNames.size());
1820  labelList patchStarts(allPatchNames.size());
1821 
1822  label startFacei = nInternalFaces;
1823 
1824  // Copy patches0 with new sizes. First patches always come from
1825  // mesh0 and will always be present.
1826  label allPatchi = 0;
1827 
1828  forAll(from0ToAllPatches, patch0)
1829  {
1830  // Originates from mesh0. Clone with new size & filter out empty
1831  // patch.
1832 
1833  if (nFaces[patch0] == 0 && isA<processorPolyPatch>(patches0[patch0]))
1834  {
1835  // Pout<< "Removing zero sized mesh0 patch "
1836  // << allPatchNames[patch0]
1837  // << endl;
1838  from0ToAllPatches[patch0] = -1;
1839 
1840  // Check if patch was also in mesh1 and update its addressing if so.
1841  if (fromAllTo1Patches[patch0] != -1)
1842  {
1843  from1ToAllPatches[fromAllTo1Patches[patch0]] = -1;
1844  }
1845  }
1846  else
1847  {
1848  patchSizes[allPatchi] = nFaces[patch0];
1849  patchStarts[allPatchi] = startFacei;
1850 
1851  // Record new index in allPatches
1852  from0ToAllPatches[patch0] = allPatchi;
1853 
1854  // Check if patch was also in mesh1 and update its addressing if so.
1855  if (fromAllTo1Patches[patch0] != -1)
1856  {
1857  from1ToAllPatches[fromAllTo1Patches[patch0]] = allPatchi;
1858  }
1859 
1860  startFacei += nFaces[patch0];
1861 
1862  allPatchi++;
1863  }
1864  }
1865 
1866  // Trim the existing patches
1867  {
1868  const label sz0 = from0ToAllPatches.size();
1869  labelList newToOld(sz0, sz0-1);
1870  label nNew = 0;
1871  forAll(from0ToAllPatches, patchi)
1872  {
1873  if (from0ToAllPatches[patchi] != -1)
1874  {
1875  newToOld[nNew++] = patchi;
1876  }
1877  }
1878  newToOld.setSize(nNew);
1879  mesh0.reorderPatches(newToOld, false);
1880  }
1881 
1882  // Copy unique patches of mesh1.
1883  forAll(from1ToAllPatches, patch1)
1884  {
1885  label uncompactAllPatchi = from1ToAllPatches[patch1];
1886 
1887  if (uncompactAllPatchi >= from0ToAllPatches.size())
1888  {
1889  // Patch has not been merged with any mesh0 patch.
1890 
1891  if
1892  (
1893  nFaces[uncompactAllPatchi] == 0
1894  && isA<processorPolyPatch>(patches1[patch1])
1895  )
1896  {
1897  // Pout<< "Removing zero sized mesh1 patch "
1898  // << allPatchNames[uncompactAllPatchi] << endl;
1899  from1ToAllPatches[patch1] = -1;
1900  }
1901  else
1902  {
1903  patchSizes[allPatchi] = nFaces[uncompactAllPatchi];
1904  patchStarts[allPatchi] = startFacei;
1905 
1906  // Clone. Note dummy size and start. Gets overwritten later in
1907  // resetPrimitives. This avoids getting temporarily illegal
1908  // SubList construction in polyPatch.
1909  mesh0.addPatch
1910  (
1911  allPatchi,
1912  patches1[patch1],
1913  dictionary(),
1914  "calculated",
1915  false
1916  );
1917 
1918  // Record new index in allPatches
1919  from1ToAllPatches[patch1] = allPatchi;
1920 
1921  startFacei += nFaces[uncompactAllPatchi];
1922 
1923  allPatchi++;
1924  }
1925  }
1926  }
1927  patchSizes.setSize(allPatchi);
1928  patchStarts.setSize(allPatchi);
1929 
1930 
1931  // Construct map information before changing mesh0 primitives
1932  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1933 
1935  (
1936  new mapAddedPolyMesh
1937  (
1938  mesh0.nPoints(),
1939  mesh0.nFaces(),
1940  mesh0.nCells(),
1941 
1942  mesh1.nPoints(),
1943  mesh1.nFaces(),
1944  mesh1.nCells(),
1945 
1946  from0ToAllPoints,
1947  from0ToAllFaces,
1948  identity(mesh0.nCells()),
1949 
1950  from1ToAllPoints,
1951  from1ToAllFaces,
1952  from1ToAllCells,
1953 
1954  from0ToAllPatches,
1955  from1ToAllPatches,
1956 
1957  mesh0PatchSizes,
1958  mesh0PatchStarts
1959  )
1960  );
1961 
1962 
1963 
1964  // Now we have extracted all information from all meshes.
1965  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1966 
1967  mesh0.resetMotion(); // delete any oldPoints.
1968  mesh0.resetPrimitives
1969  (
1970  move(allPoints),
1971  move(allFaces),
1972  move(allOwner),
1973  move(allNeighbour),
1974  patchSizes, // size
1975  patchStarts, // patchstarts
1976  validBoundary // boundary valid?
1977  );
1978 
1979  // Add zones to new mesh.
1980  mesh0.pointZones().clear();
1981  mesh0.faceZones().clear();
1982  mesh0.cellZones().clear();
1983  addZones
1984  (
1985  pointZoneNames,
1986  pzPoints,
1987 
1988  faceZoneNames,
1989  fzFaces,
1990  fzFlips,
1991 
1992  cellZoneNames,
1993  czCells,
1994  mesh0
1995  );
1996 
1997  return mapPtr;
1998 }
1999 
2000 
2003  const polyMesh& mesh,
2004  const scalar mergeDist
2005 )
2006 {
2007  const labelList& sharedPointLabels = mesh.globalData().sharedPointLabels();
2008  const labelList& sharedPointAddr = mesh.globalData().sharedPointAddr();
2009 
2010  // Because of adding the missing pieces e.g. when redistributing a mesh
2011  // it can be that there are multiple points on the same processor that
2012  // refer to the same shared point.
2013 
2014  // Invert point-to-shared addressing
2015  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2016 
2017  Map<labelList> sharedToMesh(sharedPointLabels.size());
2018 
2019  label nMultiple = 0;
2020 
2021  forAll(sharedPointLabels, i)
2022  {
2023  label pointi = sharedPointLabels[i];
2024 
2025  label sharedI = sharedPointAddr[i];
2026 
2027  Map<labelList>::iterator iter = sharedToMesh.find(sharedI);
2028 
2029  if (iter != sharedToMesh.end())
2030  {
2031  // sharedI already used by other point. Add this one.
2032 
2033  nMultiple++;
2034 
2035  labelList& connectedPointLabels = iter();
2036 
2037  label sz = connectedPointLabels.size();
2038 
2039  // Check just to make sure.
2040  if (findIndex(connectedPointLabels, pointi) != -1)
2041  {
2043  << "Duplicate point in sharedPoint addressing." << endl
2044  << "When trying to add point " << pointi << " on shared "
2045  << sharedI << " with connected points "
2046  << connectedPointLabels
2047  << abort(FatalError);
2048  }
2049 
2050  connectedPointLabels.setSize(sz+1);
2051  connectedPointLabels[sz] = pointi;
2052  }
2053  else
2054  {
2055  sharedToMesh.insert(sharedI, labelList(1, pointi));
2056  }
2057  }
2058 
2059 
2060  // Assign single master for every shared with multiple geometric points
2061  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2062 
2063  Map<label> pointToMaster(nMultiple);
2064 
2065  forAllConstIter(Map<labelList>, sharedToMesh, iter)
2066  {
2067  const labelList& connectedPointLabels = iter();
2068 
2069  // Pout<< "For shared:" << iter.key()
2070  // << " found points:" << connectedPointLabels
2071  // << " at coords:"
2072  // << pointField(mesh.points(), connectedPointLabels) << endl;
2073 
2074  if (connectedPointLabels.size() > 1)
2075  {
2076  const pointField connectedPoints
2077  (
2078  mesh.points(),
2079  connectedPointLabels
2080  );
2081 
2082  labelList toMergedPoints;
2083  label nUnique = Foam::mergePoints
2084  (
2085  connectedPoints,
2086  mergeDist,
2087  false,
2088  toMergedPoints
2089  );
2090 
2091  if (nUnique < connectedPoints.size())
2092  {
2093  // Invert toMergedPoints
2094  const labelListList mergeSets
2095  (
2097  (
2098  nUnique,
2099  toMergedPoints
2100  )
2101  );
2102 
2103  // Find master for valid merges
2104  forAll(mergeSets, setI)
2105  {
2106  const labelList& mergeSet = mergeSets[setI];
2107 
2108  if (mergeSet.size() > 1)
2109  {
2110  // Pick lowest numbered point
2111  label masterPointi = labelMax;
2112 
2113  forAll(mergeSet, i)
2114  {
2115  label pointi = connectedPointLabels[mergeSet[i]];
2116 
2117  masterPointi = min(masterPointi, pointi);
2118  }
2119 
2120  forAll(mergeSet, i)
2121  {
2122  label pointi = connectedPointLabels[mergeSet[i]];
2123 
2124  // Pout<< "Merging point " << pointi
2125  // << " at " << mesh.points()[pointi]
2126  // << " into master point "
2127  // << masterPointi
2128  // << " at " << mesh.points()[masterPointi]
2129  // << endl;
2130 
2131  pointToMaster.insert(pointi, masterPointi);
2132  }
2133  }
2134  }
2135  }
2136  }
2137  }
2138 
2139  //- Old: geometric merging. Causes problems for two close shared points.
2140  // labelList sharedToMerged;
2141  // label nUnique = Foam::mergePoints
2142  //(
2143  // pointField
2144  // (
2145  // mesh.points(),
2146  // sharedPointLabels
2147  // ),
2148  // mergeDist,
2149  // false,
2150  // sharedToMerged
2151  //);
2152  //
2155  //
2156  // Map<label> pointToMaster(10*sharedToMerged.size());
2157  //
2158  // if (nUnique < sharedPointLabels.size())
2159  //{
2160  // labelListList mergeSets
2161  // (
2162  // invertOneToMany
2163  // (
2164  // sharedToMerged.size(),
2165  // sharedToMerged
2166  // )
2167  // );
2168  //
2169  // label nMergeSets = 0;
2170  //
2171  // forAll(mergeSets, setI)
2172  // {
2173  // const labelList& mergeSet = mergeSets[setI];
2174  //
2175  // if (mergeSet.size() > 1)
2176  // {
2177  // // Take as master the shared point with the lowest mesh
2178  // // point label. (rather arbitrarily - could use max or
2179  // // any other one of the points)
2180  //
2181  // nMergeSets++;
2182  //
2183  // label masterI = labelMax;
2184  //
2185  // forAll(mergeSet, i)
2186  // {
2187  // label sharedI = mergeSet[i];
2188  //
2189  // masterI = min(masterI, sharedPointLabels[sharedI]);
2190  // }
2191  //
2192  // forAll(mergeSet, i)
2193  // {
2194  // label sharedI = mergeSet[i];
2195  //
2196  // pointToMaster.insert(sharedPointLabels[sharedI], masterI);
2197  // }
2198  // }
2199  // }
2200  //
2201  // // if (debug)
2202  // //{
2203  // // Pout<< "polyMeshAdder : merging:"
2204  // // << pointToMaster.size() << " into " << nMergeSets
2205  // // << " sets." << endl;
2206  // //}
2207  //}
2208 
2209  return pointToMaster;
2210 }
2211 
2212 
2215  const polyMesh& mesh,
2216  const Map<label>& pointToMaster,
2217  polyTopoChange& meshMod
2218 )
2219 {
2220  // Remove all non-master points.
2221  forAll(mesh.points(), pointi)
2222  {
2223  Map<label>::const_iterator iter = pointToMaster.find(pointi);
2224 
2225  if (iter != pointToMaster.end())
2226  {
2227  if (iter() != pointi)
2228  {
2229  meshMod.removePoint(pointi, iter());
2230  }
2231  }
2232  }
2233 
2234  // Modify faces for points. Note: could use pointFaces here but want to
2235  // avoid addressing calculation.
2236  const faceList& faces = mesh.faces();
2237 
2238  forAll(faces, facei)
2239  {
2240  const face& f = faces[facei];
2241 
2242  bool hasMerged = false;
2243 
2244  forAll(f, fp)
2245  {
2246  label pointi = f[fp];
2247 
2248  Map<label>::const_iterator iter = pointToMaster.find(pointi);
2249 
2250  if (iter != pointToMaster.end())
2251  {
2252  if (iter() != pointi)
2253  {
2254  hasMerged = true;
2255  break;
2256  }
2257  }
2258  }
2259 
2260  if (hasMerged)
2261  {
2262  face newF(f);
2263 
2264  forAll(f, fp)
2265  {
2266  label pointi = f[fp];
2267 
2268  Map<label>::const_iterator iter = pointToMaster.find(pointi);
2269 
2270  if (iter != pointToMaster.end())
2271  {
2272  newF[fp] = iter();
2273  }
2274  }
2275 
2276  label patchID = mesh.boundaryMesh().whichPatch(facei);
2277  label nei = (patchID == -1 ? mesh.faceNeighbour()[facei] : -1);
2278  label zoneID = mesh.faceZones().whichZone(facei);
2279  bool zoneFlip = false;
2280 
2281  if (zoneID >= 0)
2282  {
2283  const faceZone& fZone = mesh.faceZones()[zoneID];
2284  zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
2285  }
2286 
2287  meshMod.setAction
2288  (
2290  (
2291  newF, // modified face
2292  facei, // label of face
2293  mesh.faceOwner()[facei], // owner
2294  nei, // neighbour
2295  false, // face flip
2296  patchID, // patch for face
2297  false, // remove from zone
2298  zoneID, // zone for face
2299  zoneFlip // face flip in zone
2300  )
2301  );
2302  }
2303  }
2304 }
2305 
2306 
2307 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
ListType renumber(const labelUList &oldToNew, const ListType &)
Renumber the values (not the indices) of a list.
const labelList & sharedPointLabels() const
Return indices of local points that are globally shared.
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
void inplaceReorder(const labelUList &oldToNew, ListType &)
Inplace reorder the elements of a list.
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with the type faceZone.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
List< List< bool > > boolListList
Definition: boolList.H:51
void sortedOrder(const UList< T > &, labelList &order)
Generate the (stable) sort order for the list.
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:476
Class describing modification of a face.
void reset(T *=nullptr)
If object pointer already set, delete object and set to given.
Definition: autoPtrI.H:114
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1175
label nFaces() const
labelListList invertOneToMany(const label len, const labelUList &)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:67
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:248
List< face > faceList
Definition: faceListFwd.H:43
label nCells() const
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: UList.H:446
Class containing mesh-to-mesh mapping information after a mesh addition where we add a mesh (&#39;added m...
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
HashTable< label, label, Hash< label > >::const_iterator const_iterator
Definition: Map.H:59
label whichFace(const label globalCellID) const
Helper function to re-direct to zone::localID(...)
Definition: faceZone.C:306
const primitiveFacePatch & cutFaces() const
Addressing engine for combined set of faces.
void removePoint(const label, const label)
Remove/merge point.
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
Foam::indirectPrimitivePatch.
Container for information needed to couple to meshes. When constructed from two meshes and a geometri...
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1131
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:142
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:211
An STL-conforming iterator.
Definition: HashTable.H:426
void clear()
Delete object (if the pointer is valid) and set pointer to.
Definition: autoPtrI.H:126
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
T clone(const T &t)
Definition: List.H:54
static void mergePoints(const polyMesh &, const Map< label > &pointToMaster, polyTopoChange &meshMod)
Helper: Merge points.
label nPoints
friend class const_iterator
Declare friendship with the const_iterator.
Definition: HashTable.H:197
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1169
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1394
static const label labelMax
Definition: label.H:62
List< label > labelList
A List of labels.
Definition: labelList.H:56
ZoneMesh< pointZone, polyMesh > pointZoneMesh
A ZoneMesh with the type pointZone.
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1156
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
errorManip< error > abort(error &err)
Definition: errorManip.H:131
void addPatches(const List< polyPatch *> &, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:909
Foam::polyBoundaryMesh.
Merge points. See below.
labelList f(nPoints)
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
PrimitivePatch< List< face >, const pointField & > primitiveFacePatch
Foam::primitiveFacePatch.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: ZoneMesh.C:221
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
void setSize(const label)
Reset size of List.
Definition: List.C:281
label patchi
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Direct mesh changes based on v1.3 polyTopoChange syntax.
label mergePoints(const UList< Type > &points, const scalar mergeTol, const bool verbose, labelList &pointMap, const Type &origin=Type::zero)
Sorts and merges points. All points closer than/equal mergeTol get merged.
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
label nPoints() const
HashTable< labelList, edge, Hash< edge > > edgeLookup
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
void stableSort(UList< T > &)
Definition: UList.C:129
static autoPtr< polyMesh > add(const IOobject &io, const polyMesh &mesh0, const polyMesh &mesh1, const faceCoupleInfo &coupleInfo, autoPtr< mapAddedPolyMesh > &mapPtr)
Add two polyMeshes. Returns new polyMesh and map construct.
const labelList & sharedPointAddr() const
Return addressing into the complete globally shared points.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
List< cell > cellList
list of cells
Definition: cellList.H:42
ZoneMesh< cellZone, polyMesh > cellZoneMesh
A ZoneMesh with the type cellZone.
static Map< label > findSharedPoints(const polyMesh &, const scalar mergeTol)
Find topologically and geometrically shared points.
IOporosityModelList pZones(mesh)
A HashTable to objects of type <T> with a label key.
Definition: Map.H:49
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.