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-2021 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 "polyModifyFace.H"
35 #include "polyRemovePoint.H"
36 #include "polyTopoChange.H"
37 
38 #include "pointMesh.H"
39 #include "facePointPatch.H"
40 
41 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
42 
43 // Get index of patch in new set of patchnames/types
44 Foam::label Foam::polyMeshAdder::patchIndex
45 (
46  const polyPatch& p,
47  DynamicList<word>& allPatchNames,
48  DynamicList<word>& allPatchTypes
49 )
50 {
51  // Find the patch name on the list. If the patch is already there
52  // and patch types match, return index
53  const word& pType = p.type();
54  const word& pName = p.name();
55 
56  label patchi = findIndex(allPatchNames, pName);
57 
58  if (patchi == -1)
59  {
60  // Patch not found. Append to the list
61  allPatchNames.append(pName);
62  allPatchTypes.append(pType);
63 
64  return allPatchNames.size() - 1;
65  }
66  else if (allPatchTypes[patchi] == pType)
67  {
68  // Found name and types match
69  return patchi;
70  }
71  else
72  {
73  // Found the name, but type is different
74 
75  // Duplicate name is not allowed. Create a composite name from the
76  // patch name and case name
77  const word& caseName = p.boundaryMesh().mesh().time().caseName();
78 
79  allPatchNames.append(pName + "_" + caseName);
80  allPatchTypes.append(pType);
81 
82  Pout<< "label patchIndex(const polyPatch& p) : "
83  << "Patch " << p.index() << " named "
84  << pName << " in mesh " << caseName
85  << " already exists, but patch types"
86  << " do not match.\nCreating a composite name as "
87  << allPatchNames.last() << endl;
88 
89  return allPatchNames.size() - 1;
90  }
91 }
92 
93 
94 // Get index of zone in new set of zone names
95 Foam::label Foam::polyMeshAdder::zoneIndex
96 (
97  const word& curName,
98  DynamicList<word>& names
99 )
100 {
101  label zoneI = findIndex(names, curName);
102 
103  if (zoneI != -1)
104  {
105  return zoneI;
106  }
107  else
108  {
109  // Not found. Add new name to the list
110  names.append(curName);
111 
112  return names.size() - 1;
113  }
114 }
115 
116 
117 void Foam::polyMeshAdder::mergePatchNames
118 (
119  const polyBoundaryMesh& patches0,
120  const polyBoundaryMesh& patches1,
121 
122  DynamicList<word>& allPatchNames,
123  DynamicList<word>& allPatchTypes,
124 
125  labelList& from1ToAllPatches,
126  labelList& fromAllTo1Patches
127 )
128 {
129  // Insert the mesh0 patches and zones
130  allPatchNames.append(patches0.names());
131  allPatchTypes.append(patches0.types());
132 
133 
134  // Patches
135  // ~~~~~~~
136  // Patches from 0 are taken over as is; those from 1 get either merged
137  // (if they share name and type) or appended.
138  // Empty patches are filtered out much much later on.
139 
140  // Add mesh1 patches and build map both ways.
141  from1ToAllPatches.setSize(patches1.size());
142 
143  forAll(patches1, patchi)
144  {
145  from1ToAllPatches[patchi] = patchIndex
146  (
147  patches1[patchi],
148  allPatchNames,
149  allPatchTypes
150  );
151  }
152  allPatchTypes.shrink();
153  allPatchNames.shrink();
154 
155  // Invert 1 to all patch map
156  fromAllTo1Patches.setSize(allPatchNames.size());
157  fromAllTo1Patches = -1;
158 
159  forAll(from1ToAllPatches, i)
160  {
161  fromAllTo1Patches[from1ToAllPatches[i]] = i;
162  }
163 }
164 
165 
166 Foam::labelList Foam::polyMeshAdder::getPatchStarts
167 (
168  const polyBoundaryMesh& patches
169 )
170 {
171  labelList patchStarts(patches.size());
172  forAll(patches, patchi)
173  {
174  patchStarts[patchi] = patches[patchi].start();
175  }
176  return patchStarts;
177 }
178 
179 
180 Foam::labelList Foam::polyMeshAdder::getPatchSizes
181 (
182  const polyBoundaryMesh& patches
183 )
184 {
185  labelList patchSizes(patches.size());
186  forAll(patches, patchi)
187  {
188  patchSizes[patchi] = patches[patchi].size();
189  }
190  return patchSizes;
191 }
192 
193 
194 Foam::List<Foam::polyPatch*> Foam::polyMeshAdder::combinePatches
195 (
196  const polyMesh& mesh0,
197  const polyMesh& mesh1,
198  const polyBoundaryMesh& allBoundaryMesh,
199  const label nAllPatches,
200  const labelList& fromAllTo1Patches,
201 
202  const label nInternalFaces,
203  const labelList& nFaces,
204 
205  labelList& from0ToAllPatches,
206  labelList& from1ToAllPatches
207 )
208 {
209  const polyBoundaryMesh& patches0 = mesh0.boundaryMesh();
210  const polyBoundaryMesh& patches1 = mesh1.boundaryMesh();
211 
212 
213  // Compacted new patch list.
214  DynamicList<polyPatch*> allPatches(nAllPatches);
215 
216 
217  // Map from 0 to all patches (since gets compacted)
218  from0ToAllPatches.setSize(patches0.size());
219  from0ToAllPatches = -1;
220 
221  label startFacei = nInternalFaces;
222 
223  // Copy patches0 with new sizes. First patches always come from
224  // mesh0 and will always be present.
225  forAll(patches0, patchi)
226  {
227  // Originates from mesh0. Clone with new size & filter out empty
228  // patch.
229  label filteredPatchi;
230 
231  if (nFaces[patchi] == 0 && isA<processorPolyPatch>(patches0[patchi]))
232  {
233  // Pout<< "Removing zero sized mesh0 patch "
234  // << patches0[patchi].name() << endl;
235  filteredPatchi = -1;
236  }
237  else
238  {
239  filteredPatchi = allPatches.size();
240 
241  allPatches.append
242  (
243  patches0[patchi].clone
244  (
245  allBoundaryMesh,
246  filteredPatchi,
247  nFaces[patchi],
248  startFacei
249  ).ptr()
250  );
251  startFacei += nFaces[patchi];
252  }
253 
254  // Record new index in allPatches
255  from0ToAllPatches[patchi] = filteredPatchi;
256 
257  // Check if patch was also in mesh1 and update its addressing if so.
258  if (fromAllTo1Patches[patchi] != -1)
259  {
260  from1ToAllPatches[fromAllTo1Patches[patchi]] = filteredPatchi;
261  }
262  }
263 
264  // Copy unique patches of mesh1.
265  forAll(from1ToAllPatches, patchi)
266  {
267  label allPatchi = from1ToAllPatches[patchi];
268 
269  if (allPatchi >= patches0.size())
270  {
271  // Patch has not been merged with any mesh0 patch.
272 
273  label filteredPatchi;
274 
275  if
276  (
277  nFaces[allPatchi] == 0
278  && isA<processorPolyPatch>(patches1[patchi])
279  )
280  {
281  // Pout<< "Removing zero sized mesh1 patch "
282  // << patches1[patchi].name() << endl;
283  filteredPatchi = -1;
284  }
285  else
286  {
287  filteredPatchi = allPatches.size();
288 
289  allPatches.append
290  (
291  patches1[patchi].clone
292  (
293  allBoundaryMesh,
294  filteredPatchi,
295  nFaces[allPatchi],
296  startFacei
297  ).ptr()
298  );
299  startFacei += nFaces[allPatchi];
300  }
301 
302  from1ToAllPatches[patchi] = filteredPatchi;
303  }
304  }
305 
306  allPatches.shrink();
307 
308  return move(allPatches);
309 }
310 
311 
312 Foam::labelList Foam::polyMeshAdder::getFaceOrder
313 (
314  const cellList& cells,
315  const label nInternalFaces,
316  const labelList& owner,
317  const labelList& neighbour
318 )
319 {
320  labelList oldToNew(owner.size(), -1);
321 
322  // Leave boundary faces in order
323  for (label facei = nInternalFaces; facei < owner.size(); ++facei)
324  {
325  oldToNew[facei] = facei;
326  }
327 
328  // First unassigned face
329  label newFacei = 0;
330 
331  forAll(cells, celli)
332  {
333  const labelList& cFaces = cells[celli];
334 
335  SortableList<label> nbr(cFaces.size());
336 
337  forAll(cFaces, i)
338  {
339  label facei = cFaces[i];
340 
341  label nbrCelli = neighbour[facei];
342 
343  if (nbrCelli != -1)
344  {
345  // Internal face. Get cell on other side.
346  if (nbrCelli == celli)
347  {
348  nbrCelli = owner[facei];
349  }
350 
351  if (celli < nbrCelli)
352  {
353  // Celli is master
354  nbr[i] = nbrCelli;
355  }
356  else
357  {
358  // nbrCell is master. Let it handle this face.
359  nbr[i] = -1;
360  }
361  }
362  else
363  {
364  // External face. Do later.
365  nbr[i] = -1;
366  }
367  }
368 
369  nbr.sort();
370 
371  forAll(nbr, i)
372  {
373  if (nbr[i] != -1)
374  {
375  oldToNew[cFaces[nbr.indices()[i]]] = newFacei++;
376  }
377  }
378  }
379 
380 
381  // Check done all faces.
382  forAll(oldToNew, facei)
383  {
384  if (oldToNew[facei] == -1)
385  {
387  << "Did not determine new position"
388  << " for face " << facei
389  << abort(FatalError);
390  }
391  }
392 
393  return oldToNew;
394 }
395 
396 
397 // Adds primitives (cells, faces, points)
398 // Cells:
399 // - all of mesh0
400 // - all of mesh1
401 // Faces:
402 // - all uncoupled of mesh0
403 // - all coupled faces
404 // - all uncoupled of mesh1
405 // Points:
406 // - all coupled
407 // - all uncoupled of mesh0
408 // - all uncoupled of mesh1
409 void Foam::polyMeshAdder::mergePrimitives
410 (
411  const polyMesh& mesh0,
412  const polyMesh& mesh1,
413  const faceCoupleInfo& coupleInfo,
414 
415  const label nAllPatches, // number of patches in the new mesh
416  const labelList& fromAllTo1Patches,
417  const labelList& from1ToAllPatches,
418 
419  pointField& allPoints,
420  labelList& from0ToAllPoints,
421  labelList& from1ToAllPoints,
422 
423  faceList& allFaces,
424  labelList& allOwner,
425  labelList& allNeighbour,
426  label& nInternalFaces,
427  labelList& nFacesPerPatch,
428  label& nCells,
429 
430  labelList& from0ToAllFaces,
431  labelList& from1ToAllFaces,
432  labelList& from1ToAllCells
433 )
434 {
435  const polyBoundaryMesh& patches0 = mesh0.boundaryMesh();
436  const polyBoundaryMesh& patches1 = mesh1.boundaryMesh();
437 
438  const indirectPrimitivePatch& masterPatch = coupleInfo.masterPatch();
439  const indirectPrimitivePatch& slavePatch = coupleInfo.slavePatch();
440 
441 
442  // Points
443  // ~~~~~~
444 
445  // Storage for new points
446  allPoints.setSize(mesh0.nPoints() + mesh1.nPoints());
447  label allPointi = 0;
448 
449  from0ToAllPoints.setSize(mesh0.nPoints());
450  from0ToAllPoints = -1;
451  from1ToAllPoints.setSize(mesh1.nPoints());
452  from1ToAllPoints = -1;
453 
454  // Copy coupled points
455  {
456  const labelListList coupleToMasterPoints
457  (
458  coupleInfo.coupleToMasterPoints()
459  );
460  const labelListList coupleToSlavePoints
461  (
462  coupleInfo.coupleToSlavePoints()
463  );
464 
465  forAll(coupleToMasterPoints, couplePointi)
466  {
467  const labelList& masterPoints = coupleToMasterPoints[couplePointi];
468 
469  forAll(masterPoints, j)
470  {
471  label mesh0Pointi = masterPatch.meshPoints()[masterPoints[j]];
472  from0ToAllPoints[mesh0Pointi] = allPointi;
473  allPoints[allPointi] = mesh0.points()[mesh0Pointi];
474  }
475 
476  const labelList& slavePoints = coupleToSlavePoints[couplePointi];
477 
478  forAll(slavePoints, j)
479  {
480  label mesh1Pointi = slavePatch.meshPoints()[slavePoints[j]];
481  from1ToAllPoints[mesh1Pointi] = allPointi;
482  allPoints[allPointi] = mesh1.points()[mesh1Pointi];
483  }
484 
485  allPointi++;
486  }
487  }
488 
489  // Add uncoupled mesh0 points
490  forAll(mesh0.points(), pointi)
491  {
492  if (from0ToAllPoints[pointi] == -1)
493  {
494  allPoints[allPointi] = mesh0.points()[pointi];
495  from0ToAllPoints[pointi] = allPointi;
496  allPointi++;
497  }
498  }
499 
500  // Add uncoupled mesh1 points
501  forAll(mesh1.points(), pointi)
502  {
503  if (from1ToAllPoints[pointi] == -1)
504  {
505  allPoints[allPointi] = mesh1.points()[pointi];
506  from1ToAllPoints[pointi] = allPointi;
507  allPointi++;
508  }
509  }
510 
511  allPoints.setSize(allPointi);
512 
513 
514  // Faces
515  // ~~~~~
516 
517  // Sizes per patch
518  nFacesPerPatch.setSize(nAllPatches);
519  nFacesPerPatch = 0;
520 
521  // Storage for faces and owner/neighbour
522  allFaces.setSize(mesh0.nFaces() + mesh1.nFaces());
523  allOwner.setSize(allFaces.size());
524  allOwner = -1;
525  allNeighbour.setSize(allFaces.size());
526  allNeighbour = -1;
527  label allFacei = 0;
528 
529  from0ToAllFaces.setSize(mesh0.nFaces());
530  from0ToAllFaces = -1;
531  from1ToAllFaces.setSize(mesh1.nFaces());
532  from1ToAllFaces = -1;
533 
534  // Copy mesh0 internal faces (always uncoupled)
535  for (label facei = 0; facei < mesh0.nInternalFaces(); facei++)
536  {
537  allFaces[allFacei] = renumber(from0ToAllPoints, mesh0.faces()[facei]);
538  allOwner[allFacei] = mesh0.faceOwner()[facei];
539  allNeighbour[allFacei] = mesh0.faceNeighbour()[facei];
540  from0ToAllFaces[facei] = allFacei++;
541  }
542 
543  // Copy coupled faces. Every coupled face has an equivalent master and
544  // slave. Also uncount as boundary faces all the newly coupled faces.
545  forAll(masterPatch, coupleFacei)
546  {
547  const label mesh0Facei = masterPatch.addressing()[coupleFacei];
548 
549  if (from0ToAllFaces[mesh0Facei] == -1)
550  {
551  // First occurrence of face
552  from0ToAllFaces[mesh0Facei] = allFacei;
553 
554  // External face becomes internal so uncount
555  label patch0 = patches0.whichPatch(mesh0Facei);
556  nFacesPerPatch[patch0]--;
557  }
558 
559  const label mesh1Facei = slavePatch.addressing()[coupleFacei];
560 
561  if (from1ToAllFaces[mesh1Facei] == -1)
562  {
563  from1ToAllFaces[mesh1Facei] = allFacei;
564 
565  label patch1 = patches1.whichPatch(mesh1Facei);
566  nFacesPerPatch[from1ToAllPatches[patch1]]--;
567  }
568 
569  // Copy cut face (since cutPoints are copied first no renumbering
570  // necessary)
571  allFaces[allFacei] = coupleInfo.coupleFace(coupleFacei);
572  allOwner[allFacei] = mesh0.faceOwner()[mesh0Facei];
573  allNeighbour[allFacei] = mesh1.faceOwner()[mesh1Facei] + mesh0.nCells();
574 
575  allFacei++;
576  }
577 
578  // Copy mesh1 internal faces (always uncoupled)
579  for (label facei = 0; facei < mesh1.nInternalFaces(); facei++)
580  {
581  allFaces[allFacei] = renumber(from1ToAllPoints, mesh1.faces()[facei]);
582  allOwner[allFacei] = mesh1.faceOwner()[facei] + mesh0.nCells();
583  allNeighbour[allFacei] = mesh1.faceNeighbour()[facei] + mesh0.nCells();
584  from1ToAllFaces[facei] = allFacei++;
585  }
586 
587  nInternalFaces = allFacei;
588 
589  // Copy (unmarked/uncoupled) external faces in new order.
590  for (label allPatchi = 0; allPatchi < nAllPatches; allPatchi++)
591  {
592  if (allPatchi < patches0.size())
593  {
594  // Patch is present in mesh0
595  const polyPatch& pp = patches0[allPatchi];
596 
597  nFacesPerPatch[allPatchi] += pp.size();
598 
599  label facei = pp.start();
600 
601  forAll(pp, i)
602  {
603  if (from0ToAllFaces[facei] == -1)
604  {
605  // Is uncoupled face since has not yet been dealt with
606  allFaces[allFacei] = renumber
607  (
608  from0ToAllPoints,
609  mesh0.faces()[facei]
610  );
611  allOwner[allFacei] = mesh0.faceOwner()[facei];
612  allNeighbour[allFacei] = -1;
613 
614  from0ToAllFaces[facei] = allFacei++;
615  }
616  facei++;
617  }
618  }
619  if (fromAllTo1Patches[allPatchi] != -1)
620  {
621  // Patch is present in mesh1
622  const polyPatch& pp = patches1[fromAllTo1Patches[allPatchi]];
623 
624  nFacesPerPatch[allPatchi] += pp.size();
625 
626  label facei = pp.start();
627 
628  forAll(pp, i)
629  {
630  if (from1ToAllFaces[facei] == -1)
631  {
632  // Is uncoupled face
633  allFaces[allFacei] = renumber
634  (
635  from1ToAllPoints,
636  mesh1.faces()[facei]
637  );
638  allOwner[allFacei] =
639  mesh1.faceOwner()[facei]
640  + mesh0.nCells();
641  allNeighbour[allFacei] = -1;
642 
643  from1ToAllFaces[facei] = allFacei++;
644  }
645  facei++;
646  }
647  }
648  }
649  allFaces.setSize(allFacei);
650  allOwner.setSize(allFacei);
651  allNeighbour.setSize(allFacei);
652 
653  // Cells
654  // ~~~~~
655 
656  from1ToAllCells.setSize(mesh1.nCells());
657  from1ToAllCells = -1;
658 
659  forAll(mesh1.cells(), i)
660  {
661  from1ToAllCells[i] = i + mesh0.nCells();
662  }
663 
664  // Make cells (= cell-face addressing)
665  nCells = mesh0.nCells() + mesh1.nCells();
666  cellList allCells(nCells);
667  primitiveMesh::calcCells(allCells, allOwner, allNeighbour, nCells);
668 
669  // Reorder faces for upper-triangular order.
670  labelList oldToNew
671  (
672  getFaceOrder
673  (
674  allCells,
675  nInternalFaces,
676  allOwner,
677  allNeighbour
678  )
679  );
680 
681  inplaceReorder(oldToNew, allFaces);
682  inplaceReorder(oldToNew, allOwner);
683  inplaceReorder(oldToNew, allNeighbour);
684  inplaceRenumber(oldToNew, from0ToAllFaces);
685  inplaceRenumber(oldToNew, from1ToAllFaces);
686 }
687 
688 
689 void Foam::polyMeshAdder::mergePointZones
690 (
691  const label nAllPoints,
692  const meshPointZones& pz0,
693  const meshPointZones& pz1,
694  const labelList& from0ToAllPoints,
695  const labelList& from1ToAllPoints,
696 
697  DynamicList<word>& zoneNames,
698  labelList& from1ToAll,
699  List<DynamicList<label>>& pzPoints
700 )
701 {
702  zoneNames.setCapacity(pz0.size() + pz1.size());
703  zoneNames.append(pz0.names());
704 
705  from1ToAll.setSize(pz1.size());
706 
707  forAll(pz1, zoneI)
708  {
709  from1ToAll[zoneI] = zoneIndex(pz1[zoneI].name(), zoneNames);
710  }
711  zoneNames.shrink();
712 
713 
714 
715  // Zone(s) per point. Two levels: if only one zone
716  // stored in pointToZone. Any extra stored in additionalPointToZones.
717  // This is so we only allocate labelLists per point if absolutely
718  // necessary.
719  labelList pointToZone(nAllPoints, -1);
720  labelListList addPointToZones(nAllPoints);
721 
722  // mesh0 zones kept
723  forAll(pz0, zoneI)
724  {
725  const pointZone& pz = pz0[zoneI];
726 
727  forAll(pz, i)
728  {
729  label point0 = pz[i];
730  label allPointi = from0ToAllPoints[point0];
731 
732  if (pointToZone[allPointi] == -1)
733  {
734  pointToZone[allPointi] = zoneI;
735  }
736  else if (pointToZone[allPointi] != zoneI)
737  {
738  labelList& pZones = addPointToZones[allPointi];
739  if (findIndex(pZones, zoneI) == -1)
740  {
741  pZones.append(zoneI);
742  }
743  }
744  }
745  }
746 
747  // mesh1 zones renumbered
748  forAll(pz1, zoneI)
749  {
750  const pointZone& pz = pz1[zoneI];
751  const label allZoneI = from1ToAll[zoneI];
752 
753  forAll(pz, i)
754  {
755  label point1 = pz[i];
756  label allPointi = from1ToAllPoints[point1];
757 
758  if (pointToZone[allPointi] == -1)
759  {
760  pointToZone[allPointi] = allZoneI;
761  }
762  else if (pointToZone[allPointi] != allZoneI)
763  {
764  labelList& pZones = addPointToZones[allPointi];
765  if (findIndex(pZones, allZoneI) == -1)
766  {
767  pZones.append(allZoneI);
768  }
769  }
770  }
771  }
772 
773  // Extract back into zones
774 
775  // 1. Count
776  labelList nPoints(zoneNames.size(), 0);
777  forAll(pointToZone, allPointi)
778  {
779  label zoneI = pointToZone[allPointi];
780  if (zoneI != -1)
781  {
782  nPoints[zoneI]++;
783  }
784  }
785  forAll(addPointToZones, allPointi)
786  {
787  const labelList& pZones = addPointToZones[allPointi];
788  forAll(pZones, i)
789  {
790  nPoints[pZones[i]]++;
791  }
792  }
793 
794  // 2. Fill
795  pzPoints.setSize(zoneNames.size());
796  forAll(pzPoints, zoneI)
797  {
798  pzPoints[zoneI].setCapacity(nPoints[zoneI]);
799  }
800  forAll(pointToZone, allPointi)
801  {
802  label zoneI = pointToZone[allPointi];
803  if (zoneI != -1)
804  {
805  pzPoints[zoneI].append(allPointi);
806  }
807  }
808  forAll(addPointToZones, allPointi)
809  {
810  const labelList& pZones = addPointToZones[allPointi];
811  forAll(pZones, i)
812  {
813  pzPoints[pZones[i]].append(allPointi);
814  }
815  }
816  forAll(pzPoints, i)
817  {
818  pzPoints[i].shrink();
819  stableSort(pzPoints[i]);
820  }
821 }
822 
823 
824 void Foam::polyMeshAdder::mergeFaceZones
825 (
826  const labelList& allOwner,
827 
828  const polyMesh& mesh0,
829  const polyMesh& mesh1,
830 
831  const labelList& from0ToAllFaces,
832  const labelList& from1ToAllFaces,
833 
834  const labelList& from1ToAllCells,
835 
836  DynamicList<word>& zoneNames,
837  labelList& from1ToAll,
838  List<DynamicList<label>>& fzFaces,
839  List<DynamicList<bool>>& fzFlips
840 )
841 {
842  const meshFaceZones& fz0 = mesh0.faceZones();
843  const labelList& owner0 = mesh0.faceOwner();
844  const meshFaceZones& fz1 = mesh1.faceZones();
845  const labelList& owner1 = mesh1.faceOwner();
846 
847 
848  zoneNames.setCapacity(fz0.size() + fz1.size());
849  zoneNames.append(fz0.names());
850 
851  from1ToAll.setSize(fz1.size());
852 
853  forAll(fz1, zoneI)
854  {
855  from1ToAll[zoneI] = zoneIndex(fz1[zoneI].name(), zoneNames);
856  }
857  zoneNames.shrink();
858 
859 
860  // Zone(s) per face
861  labelList faceToZone(allOwner.size(), -1);
862  labelListList addFaceToZones(allOwner.size());
863  boolList faceToFlip(allOwner.size(), false);
864  boolListList addFaceToFlips(allOwner.size());
865 
866  // mesh0 zones kept
867  forAll(fz0, zoneI)
868  {
869  const labelList& addressing = fz0[zoneI];
870  const boolList& flipMap = fz0[zoneI].flipMap();
871 
872  forAll(addressing, i)
873  {
874  label face0 = addressing[i];
875  bool flip0 = flipMap[i];
876 
877  label allFacei = from0ToAllFaces[face0];
878  if (allFacei != -1)
879  {
880  // Check if orientation same
881  label allCell0 = owner0[face0];
882  if (allOwner[allFacei] != allCell0)
883  {
884  flip0 = !flip0;
885  }
886 
887  if (faceToZone[allFacei] == -1)
888  {
889  faceToZone[allFacei] = zoneI;
890  faceToFlip[allFacei] = flip0;
891  }
892  else if (faceToZone[allFacei] != zoneI)
893  {
894  labelList& fZones = addFaceToZones[allFacei];
895  boolList& flipZones = addFaceToFlips[allFacei];
896 
897  if (findIndex(fZones, zoneI) == -1)
898  {
899  fZones.append(zoneI);
900  flipZones.append(flip0);
901  }
902  }
903  }
904  }
905  }
906 
907  // mesh1 zones renumbered
908  forAll(fz1, zoneI)
909  {
910  const labelList& addressing = fz1[zoneI];
911  const boolList& flipMap = fz1[zoneI].flipMap();
912 
913  const label allZoneI = from1ToAll[zoneI];
914 
915  forAll(addressing, i)
916  {
917  label face1 = addressing[i];
918  bool flip1 = flipMap[i];
919 
920  label allFacei = from1ToAllFaces[face1];
921  if (allFacei != -1)
922  {
923  // Check if orientation same
924  label allCell1 = from1ToAllCells[owner1[face1]];
925  if (allOwner[allFacei] != allCell1)
926  {
927  flip1 = !flip1;
928  }
929 
930  if (faceToZone[allFacei] == -1)
931  {
932  faceToZone[allFacei] = allZoneI;
933  faceToFlip[allFacei] = flip1;
934  }
935  else if (faceToZone[allFacei] != allZoneI)
936  {
937  labelList& fZones = addFaceToZones[allFacei];
938  boolList& flipZones = addFaceToFlips[allFacei];
939 
940  if (findIndex(fZones, allZoneI) == -1)
941  {
942  fZones.append(allZoneI);
943  flipZones.append(flip1);
944  }
945  }
946  }
947  }
948  }
949 
950 
951  // Extract back into zones
952 
953  // 1. Count
954  labelList nFaces(zoneNames.size(), 0);
955  forAll(faceToZone, allFacei)
956  {
957  label zoneI = faceToZone[allFacei];
958  if (zoneI != -1)
959  {
960  nFaces[zoneI]++;
961  }
962  }
963  forAll(addFaceToZones, allFacei)
964  {
965  const labelList& fZones = addFaceToZones[allFacei];
966  forAll(fZones, i)
967  {
968  nFaces[fZones[i]]++;
969  }
970  }
971 
972  // 2. Fill
973  fzFaces.setSize(zoneNames.size());
974  fzFlips.setSize(zoneNames.size());
975  forAll(fzFaces, zoneI)
976  {
977  fzFaces[zoneI].setCapacity(nFaces[zoneI]);
978  fzFlips[zoneI].setCapacity(nFaces[zoneI]);
979  }
980  forAll(faceToZone, allFacei)
981  {
982  label zoneI = faceToZone[allFacei];
983  bool flip = faceToFlip[allFacei];
984  if (zoneI != -1)
985  {
986  fzFaces[zoneI].append(allFacei);
987  fzFlips[zoneI].append(flip);
988  }
989  }
990  forAll(addFaceToZones, allFacei)
991  {
992  const labelList& fZones = addFaceToZones[allFacei];
993  const boolList& flipZones = addFaceToFlips[allFacei];
994 
995  forAll(fZones, i)
996  {
997  label zoneI = fZones[i];
998  fzFaces[zoneI].append(allFacei);
999  fzFlips[zoneI].append(flipZones[i]);
1000  }
1001  }
1002 
1003  forAll(fzFaces, i)
1004  {
1005  fzFaces[i].shrink();
1006  fzFlips[i].shrink();
1007 
1008  labelList order;
1009  sortedOrder(fzFaces[i], order);
1010  fzFaces[i] = UIndirectList<label>(fzFaces[i], order)();
1011  fzFlips[i] = UIndirectList<bool>(fzFlips[i], order)();
1012  }
1013 }
1014 
1015 
1016 void Foam::polyMeshAdder::mergeCellZones
1017 (
1018  const label nAllCells,
1019 
1020  const meshCellZones& cz0,
1021  const meshCellZones& cz1,
1022  const labelList& from1ToAllCells,
1023 
1024  DynamicList<word>& zoneNames,
1025  labelList& from1ToAll,
1026  List<DynamicList<label>>& czCells
1027 )
1028 {
1029  zoneNames.setCapacity(cz0.size() + cz1.size());
1030  zoneNames.append(cz0.names());
1031 
1032  from1ToAll.setSize(cz1.size());
1033  forAll(cz1, zoneI)
1034  {
1035  from1ToAll[zoneI] = zoneIndex(cz1[zoneI].name(), zoneNames);
1036  }
1037  zoneNames.shrink();
1038 
1039 
1040  // Zone(s) per cell. Two levels: if only one zone
1041  // stored in cellToZone. Any extra stored in additionalCellToZones.
1042  // This is so we only allocate labelLists per cell if absolutely
1043  // necessary.
1044  labelList cellToZone(nAllCells, -1);
1045  labelListList addCellToZones(nAllCells);
1046 
1047  // mesh0 zones kept
1048  forAll(cz0, zoneI)
1049  {
1050  const cellZone& cz = cz0[zoneI];
1051  forAll(cz, i)
1052  {
1053  label cell0 = cz[i];
1054 
1055  if (cellToZone[cell0] == -1)
1056  {
1057  cellToZone[cell0] = zoneI;
1058  }
1059  else if (cellToZone[cell0] != zoneI)
1060  {
1061  labelList& cZones = addCellToZones[cell0];
1062  if (findIndex(cZones, zoneI) == -1)
1063  {
1064  cZones.append(zoneI);
1065  }
1066  }
1067  }
1068  }
1069 
1070  // mesh1 zones renumbered
1071  forAll(cz1, zoneI)
1072  {
1073  const cellZone& cz = cz1[zoneI];
1074  const label allZoneI = from1ToAll[zoneI];
1075  forAll(cz, i)
1076  {
1077  label cell1 = cz[i];
1078  label allCelli = from1ToAllCells[cell1];
1079 
1080  if (cellToZone[allCelli] == -1)
1081  {
1082  cellToZone[allCelli] = allZoneI;
1083  }
1084  else if (cellToZone[allCelli] != allZoneI)
1085  {
1086  labelList& cZones = addCellToZones[allCelli];
1087  if (findIndex(cZones, allZoneI) == -1)
1088  {
1089  cZones.append(allZoneI);
1090  }
1091  }
1092  }
1093  }
1094 
1095  // Extract back into zones
1096 
1097  // 1. Count
1098  labelList nCells(zoneNames.size(), 0);
1099  forAll(cellToZone, allCelli)
1100  {
1101  label zoneI = cellToZone[allCelli];
1102  if (zoneI != -1)
1103  {
1104  nCells[zoneI]++;
1105  }
1106  }
1107  forAll(addCellToZones, allCelli)
1108  {
1109  const labelList& cZones = addCellToZones[allCelli];
1110  forAll(cZones, i)
1111  {
1112  nCells[cZones[i]]++;
1113  }
1114  }
1115 
1116  // 2. Fill
1117  czCells.setSize(zoneNames.size());
1118  forAll(czCells, zoneI)
1119  {
1120  czCells[zoneI].setCapacity(nCells[zoneI]);
1121  }
1122  forAll(cellToZone, allCelli)
1123  {
1124  label zoneI = cellToZone[allCelli];
1125  if (zoneI != -1)
1126  {
1127  czCells[zoneI].append(allCelli);
1128  }
1129  }
1130  forAll(addCellToZones, allCelli)
1131  {
1132  const labelList& cZones = addCellToZones[allCelli];
1133  forAll(cZones, i)
1134  {
1135  czCells[cZones[i]].append(allCelli);
1136  }
1137  }
1138  forAll(czCells, i)
1139  {
1140  czCells[i].shrink();
1141  stableSort(czCells[i]);
1142  }
1143 }
1144 
1145 
1146 void Foam::polyMeshAdder::mergeZones
1147 (
1148  const label nAllPoints,
1149  const labelList& allOwner,
1150  const label nAllCells,
1151 
1152  const polyMesh& mesh0,
1153  const polyMesh& mesh1,
1154  const labelList& from0ToAllPoints,
1155  const labelList& from0ToAllFaces,
1156 
1157  const labelList& from1ToAllPoints,
1158  const labelList& from1ToAllFaces,
1159  const labelList& from1ToAllCells,
1160 
1161  DynamicList<word>& pointZoneNames,
1162  List<DynamicList<label>>& pzPoints,
1163 
1164  DynamicList<word>& faceZoneNames,
1165  List<DynamicList<label>>& fzFaces,
1166  List<DynamicList<bool>>& fzFlips,
1167 
1168  DynamicList<word>& cellZoneNames,
1169  List<DynamicList<label>>& czCells
1170 )
1171 {
1172  labelList from1ToAllPZones;
1173  mergePointZones
1174  (
1175  nAllPoints,
1176  mesh0.pointZones(),
1177  mesh1.pointZones(),
1178  from0ToAllPoints,
1179  from1ToAllPoints,
1180 
1181  pointZoneNames,
1182  from1ToAllPZones,
1183  pzPoints
1184  );
1185 
1186  labelList from1ToAllFZones;
1187  mergeFaceZones
1188  (
1189  allOwner,
1190  mesh0,
1191  mesh1,
1192  from0ToAllFaces,
1193  from1ToAllFaces,
1194  from1ToAllCells,
1195 
1196  faceZoneNames,
1197  from1ToAllFZones,
1198  fzFaces,
1199  fzFlips
1200  );
1201 
1202  labelList from1ToAllCZones;
1203  mergeCellZones
1204  (
1205  nAllCells,
1206  mesh0.cellZones(),
1207  mesh1.cellZones(),
1208  from1ToAllCells,
1209 
1210  cellZoneNames,
1211  from1ToAllCZones,
1212  czCells
1213  );
1214 }
1215 
1216 
1217 void Foam::polyMeshAdder::addZones
1218 (
1219  const DynamicList<word>& pointZoneNames,
1220  const List<DynamicList<label>>& pzPoints,
1221 
1222  const DynamicList<word>& faceZoneNames,
1223  const List<DynamicList<label>>& fzFaces,
1224  const List<DynamicList<bool>>& fzFlips,
1225 
1226  const DynamicList<word>& cellZoneNames,
1227  const List<DynamicList<label>>& czCells,
1228 
1229  polyMesh& mesh
1230 )
1231 {
1232  List<pointZone*> pZones(pzPoints.size());
1233  forAll(pZones, i)
1234  {
1235  pZones[i] = new pointZone
1236  (
1237  pointZoneNames[i],
1238  pzPoints[i],
1239  i,
1240  mesh.pointZones()
1241  );
1242  }
1243 
1244  List<faceZone*> fZones(fzFaces.size());
1245  forAll(fZones, i)
1246  {
1247  fZones[i] = new faceZone
1248  (
1249  faceZoneNames[i],
1250  fzFaces[i],
1251  fzFlips[i],
1252  i,
1253  mesh.faceZones()
1254  );
1255  }
1256 
1257  List<cellZone*> cZones(czCells.size());
1258  forAll(cZones, i)
1259  {
1260  cZones[i] = new cellZone
1261  (
1262  cellZoneNames[i],
1263  czCells[i],
1264  i,
1265  mesh.cellZones()
1266  );
1267  }
1268 
1269  mesh.addZones
1270  (
1271  pZones,
1272  fZones,
1273  cZones
1274  );
1275 }
1276 
1277 
1278 
1279 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1280 
1281 // Returns new mesh and sets
1282 // - map from new cell/face/point/patch to either mesh0 or mesh1
1283 //
1284 // mesh0Faces/mesh1Faces: corresponding faces on both meshes.
1287  const IOobject& io,
1288  const polyMesh& mesh0,
1289  const polyMesh& mesh1,
1290  const faceCoupleInfo& coupleInfo,
1292 )
1293 {
1294  const polyBoundaryMesh& patches0 = mesh0.boundaryMesh();
1295  const polyBoundaryMesh& patches1 = mesh1.boundaryMesh();
1296 
1297 
1298  DynamicList<word> allPatchNames(patches0.size() + patches1.size());
1299  DynamicList<word> allPatchTypes(allPatchNames.size());
1300 
1301  // Patch maps
1302  labelList from1ToAllPatches(patches1.size());
1303  labelList fromAllTo1Patches(allPatchNames.size(), -1);
1304 
1305  mergePatchNames
1306  (
1307  patches0,
1308  patches1,
1309  allPatchNames,
1310  allPatchTypes,
1311  from1ToAllPatches,
1312  fromAllTo1Patches
1313  );
1314 
1315 
1316  // New points
1318 
1319  // Map from mesh0/1 points to allPoints.
1320  labelList from0ToAllPoints(mesh0.nPoints(), -1);
1321  labelList from1ToAllPoints(mesh1.nPoints(), -1);
1322 
1323  // New faces
1324  faceList allFaces;
1325  label nInternalFaces;
1326 
1327  // New cells
1328  labelList allOwner;
1329  labelList allNeighbour;
1330  label nCells;
1331 
1332  // Sizes per patch
1333  labelList nFaces(allPatchNames.size(), 0);
1334 
1335  // Maps
1336  labelList from0ToAllFaces(mesh0.nFaces(), -1);
1337  labelList from1ToAllFaces(mesh1.nFaces(), -1);
1338 
1339  // Map
1340  labelList from1ToAllCells(mesh1.nCells(), -1);
1341 
1342  mergePrimitives
1343  (
1344  mesh0,
1345  mesh1,
1346  coupleInfo,
1347 
1348  allPatchNames.size(),
1349  fromAllTo1Patches,
1350  from1ToAllPatches,
1351 
1352  allPoints,
1353  from0ToAllPoints,
1354  from1ToAllPoints,
1355 
1356  allFaces,
1357  allOwner,
1358  allNeighbour,
1359  nInternalFaces,
1360  nFaces,
1361  nCells,
1362 
1363  from0ToAllFaces,
1364  from1ToAllFaces,
1365  from1ToAllCells
1366  );
1367 
1368 
1369  // Zones
1370  // ~~~~~
1371 
1372  DynamicList<word> pointZoneNames;
1373  List<DynamicList<label>> pzPoints;
1374 
1375  DynamicList<word> faceZoneNames;
1376  List<DynamicList<label>> fzFaces;
1377  List<DynamicList<bool>> fzFlips;
1378 
1379  DynamicList<word> cellZoneNames;
1380  List<DynamicList<label>> czCells;
1381 
1382  mergeZones
1383  (
1384  allPoints.size(),
1385  allOwner,
1386  nCells,
1387 
1388  mesh0,
1389  mesh1,
1390 
1391  from0ToAllPoints,
1392  from0ToAllFaces,
1393 
1394  from1ToAllPoints,
1395  from1ToAllFaces,
1396  from1ToAllCells,
1397 
1398  pointZoneNames,
1399  pzPoints,
1400 
1401  faceZoneNames,
1402  fzFaces,
1403  fzFlips,
1404 
1405  cellZoneNames,
1406  czCells
1407  );
1408 
1409 
1410  // Patches
1411  // ~~~~~~~
1412 
1413  // Map from 0 to all patches (since gets compacted)
1414  labelList from0ToAllPatches(patches0.size(), -1);
1415 
1416  List<polyPatch*> allPatches
1417  (
1418  combinePatches
1419  (
1420  mesh0,
1421  mesh1,
1422  patches0, // Should be boundaryMesh() on new mesh.
1423  allPatchNames.size(),
1424  fromAllTo1Patches,
1425  mesh0.nInternalFaces()
1426  + mesh1.nInternalFaces()
1427  + coupleInfo.masterPatch().size(),
1428  nFaces,
1429 
1430  from0ToAllPatches,
1431  from1ToAllPatches
1432  )
1433  );
1434 
1435 
1436  // Map information
1437  // ~~~~~~~~~~~~~~~
1438 
1439  mapPtr.reset
1440  (
1441  new mapAddedPolyMesh
1442  (
1443  mesh0.nPoints(),
1444  mesh0.nFaces(),
1445  mesh0.nCells(),
1446 
1447  mesh1.nPoints(),
1448  mesh1.nFaces(),
1449  mesh1.nCells(),
1450 
1451  from0ToAllPoints,
1452  from0ToAllFaces,
1453  identity(mesh0.nCells()),
1454 
1455  from1ToAllPoints,
1456  from1ToAllFaces,
1457  from1ToAllCells,
1458 
1459  from0ToAllPatches,
1460  from1ToAllPatches,
1461  getPatchSizes(patches0),
1462  getPatchStarts(patches0)
1463  )
1464  );
1465 
1466 
1467 
1468  // Now we have extracted all information from all meshes.
1469  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1470 
1471  // Construct mesh
1472  autoPtr<polyMesh> tmesh
1473  (
1474  new polyMesh
1475  (
1476  io,
1477  move(allPoints),
1478  move(allFaces),
1479  move(allOwner),
1480  move(allNeighbour)
1481  )
1482  );
1483  polyMesh& mesh = tmesh();
1484 
1485  // Add zones to new mesh.
1486  addZones
1487  (
1488  pointZoneNames,
1489  pzPoints,
1490 
1491  faceZoneNames,
1492  fzFaces,
1493  fzFlips,
1494 
1495  cellZoneNames,
1496  czCells,
1497  mesh
1498  );
1499 
1500  // Add patches to new mesh
1501  mesh.addPatches(allPatches);
1502 
1503  return tmesh;
1504 }
1505 
1506 
1507 // Inplace add mesh1 to mesh0
1510  polyMesh& mesh0,
1511  const polyMesh& mesh1,
1512  const faceCoupleInfo& coupleInfo,
1513  const bool validBoundary
1514 )
1515 {
1516  const polyBoundaryMesh& patches0 = mesh0.boundaryMesh();
1517  const polyBoundaryMesh& patches1 = mesh1.boundaryMesh();
1518 
1519  DynamicList<word> allPatchNames(patches0.size() + patches1.size());
1520  DynamicList<word> allPatchTypes(allPatchNames.size());
1521 
1522  // Patch maps
1523  labelList from1ToAllPatches(patches1.size());
1524  labelList fromAllTo1Patches(allPatchNames.size(), -1);
1525 
1526  mergePatchNames
1527  (
1528  patches0,
1529  patches1,
1530  allPatchNames,
1531  allPatchTypes,
1532  from1ToAllPatches,
1533  fromAllTo1Patches
1534  );
1535 
1536 
1537  // New points
1539 
1540  // Map from mesh0/1 points to allPoints.
1541  labelList from0ToAllPoints(mesh0.nPoints(), -1);
1542  labelList from1ToAllPoints(mesh1.nPoints(), -1);
1543 
1544  // New faces
1545  faceList allFaces;
1546  labelList allOwner;
1547  labelList allNeighbour;
1548  label nInternalFaces;
1549  // Sizes per patch
1550  labelList nFaces(allPatchNames.size(), 0);
1551  label nCells;
1552 
1553  // Maps
1554  labelList from0ToAllFaces(mesh0.nFaces(), -1);
1555  labelList from1ToAllFaces(mesh1.nFaces(), -1);
1556  // Map
1557  labelList from1ToAllCells(mesh1.nCells(), -1);
1558 
1559  mergePrimitives
1560  (
1561  mesh0,
1562  mesh1,
1563  coupleInfo,
1564 
1565  allPatchNames.size(),
1566  fromAllTo1Patches,
1567  from1ToAllPatches,
1568 
1569  allPoints,
1570  from0ToAllPoints,
1571  from1ToAllPoints,
1572 
1573  allFaces,
1574  allOwner,
1575  allNeighbour,
1576  nInternalFaces,
1577  nFaces,
1578  nCells,
1579 
1580  from0ToAllFaces,
1581  from1ToAllFaces,
1582  from1ToAllCells
1583  );
1584 
1585 
1586  // Zones
1587  // ~~~~~
1588 
1589  DynamicList<word> pointZoneNames;
1590  List<DynamicList<label>> pzPoints;
1591 
1592  DynamicList<word> faceZoneNames;
1593  List<DynamicList<label>> fzFaces;
1594  List<DynamicList<bool>> fzFlips;
1595 
1596  DynamicList<word> cellZoneNames;
1597  List<DynamicList<label>> czCells;
1598 
1599  mergeZones
1600  (
1601  allPoints.size(),
1602  allOwner,
1603  nCells,
1604 
1605  mesh0,
1606  mesh1,
1607 
1608  from0ToAllPoints,
1609  from0ToAllFaces,
1610 
1611  from1ToAllPoints,
1612  from1ToAllFaces,
1613  from1ToAllCells,
1614 
1615  pointZoneNames,
1616  pzPoints,
1617 
1618  faceZoneNames,
1619  fzFaces,
1620  fzFlips,
1621 
1622  cellZoneNames,
1623  czCells
1624  );
1625 
1626 
1627  // Patches
1628  // ~~~~~~~
1629 
1630 
1631  // Store mesh0 patch info before modifying patches0.
1632  labelList mesh0PatchSizes(getPatchSizes(patches0));
1633  labelList mesh0PatchStarts(getPatchStarts(patches0));
1634 
1635  // Map from 0 to all patches (since gets compacted)
1636  labelList from0ToAllPatches(patches0.size(), -1);
1637 
1638  // Inplace extend mesh0 patches (note that patches0.size() now also
1639  // has changed)
1640  labelList patchSizes(allPatchNames.size());
1641  labelList patchStarts(allPatchNames.size());
1642 
1643  label startFacei = nInternalFaces;
1644 
1645  // Copy patches0 with new sizes. First patches always come from
1646  // mesh0 and will always be present.
1647  label allPatchi = 0;
1648 
1649  forAll(from0ToAllPatches, patch0)
1650  {
1651  // Originates from mesh0. Clone with new size & filter out empty
1652  // patch.
1653 
1654  if (nFaces[patch0] == 0 && isA<processorPolyPatch>(patches0[patch0]))
1655  {
1656  // Pout<< "Removing zero sized mesh0 patch "
1657  // << allPatchNames[patch0]
1658  // << endl;
1659  from0ToAllPatches[patch0] = -1;
1660 
1661  // Check if patch was also in mesh1 and update its addressing if so.
1662  if (fromAllTo1Patches[patch0] != -1)
1663  {
1664  from1ToAllPatches[fromAllTo1Patches[patch0]] = -1;
1665  }
1666  }
1667  else
1668  {
1669  patchSizes[allPatchi] = nFaces[patch0];
1670  patchStarts[allPatchi] = startFacei;
1671 
1672  // Record new index in allPatches
1673  from0ToAllPatches[patch0] = allPatchi;
1674 
1675  // Check if patch was also in mesh1 and update its addressing if so.
1676  if (fromAllTo1Patches[patch0] != -1)
1677  {
1678  from1ToAllPatches[fromAllTo1Patches[patch0]] = allPatchi;
1679  }
1680 
1681  startFacei += nFaces[patch0];
1682 
1683  allPatchi++;
1684  }
1685  }
1686 
1687  // Trim the existing patches
1688  {
1689  const label sz0 = from0ToAllPatches.size();
1690  labelList newToOld(sz0, sz0-1);
1691  label nNew = 0;
1692  forAll(from0ToAllPatches, patchi)
1693  {
1694  if (from0ToAllPatches[patchi] != -1)
1695  {
1696  newToOld[nNew++] = patchi;
1697  }
1698  }
1699  newToOld.setSize(nNew);
1700  mesh0.reorderPatches(newToOld, false);
1701  }
1702 
1703  // Copy unique patches of mesh1.
1704  forAll(from1ToAllPatches, patch1)
1705  {
1706  label uncompactAllPatchi = from1ToAllPatches[patch1];
1707 
1708  if (uncompactAllPatchi >= from0ToAllPatches.size())
1709  {
1710  // Patch has not been merged with any mesh0 patch.
1711 
1712  if
1713  (
1714  nFaces[uncompactAllPatchi] == 0
1715  && isA<processorPolyPatch>(patches1[patch1])
1716  )
1717  {
1718  // Pout<< "Removing zero sized mesh1 patch "
1719  // << allPatchNames[uncompactAllPatchi] << endl;
1720  from1ToAllPatches[patch1] = -1;
1721  }
1722  else
1723  {
1724  patchSizes[allPatchi] = nFaces[uncompactAllPatchi];
1725  patchStarts[allPatchi] = startFacei;
1726 
1727  // Clone. Note dummy size and start. Gets overwritten later in
1728  // resetPrimitives. This avoids getting temporarily illegal
1729  // SubList construction in polyPatch.
1730  mesh0.addPatch
1731  (
1732  allPatchi,
1733  patches1[patch1],
1734  dictionary(),
1735  "calculated",
1736  false
1737  );
1738 
1739  // Record new index in allPatches
1740  from1ToAllPatches[patch1] = allPatchi;
1741 
1742  startFacei += nFaces[uncompactAllPatchi];
1743 
1744  allPatchi++;
1745  }
1746  }
1747  }
1748  patchSizes.setSize(allPatchi);
1749  patchStarts.setSize(allPatchi);
1750 
1751 
1752  // Construct map information before changing mesh0 primitives
1753  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1754 
1756  (
1757  new mapAddedPolyMesh
1758  (
1759  mesh0.nPoints(),
1760  mesh0.nFaces(),
1761  mesh0.nCells(),
1762 
1763  mesh1.nPoints(),
1764  mesh1.nFaces(),
1765  mesh1.nCells(),
1766 
1767  from0ToAllPoints,
1768  from0ToAllFaces,
1769  identity(mesh0.nCells()),
1770 
1771  from1ToAllPoints,
1772  from1ToAllFaces,
1773  from1ToAllCells,
1774 
1775  from0ToAllPatches,
1776  from1ToAllPatches,
1777 
1778  mesh0PatchSizes,
1779  mesh0PatchStarts
1780  )
1781  );
1782 
1783 
1784 
1785  // Now we have extracted all information from all meshes.
1786  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1787 
1788  mesh0.resetMotion(); // delete any oldPoints.
1789  mesh0.resetPrimitives
1790  (
1791  move(allPoints),
1792  move(allFaces),
1793  move(allOwner),
1794  move(allNeighbour),
1795  patchSizes, // size
1796  patchStarts, // patchstarts
1797  validBoundary // boundary valid?
1798  );
1799 
1800  // Add zones to new mesh.
1801  mesh0.pointZones().clear();
1802  mesh0.faceZones().clear();
1803  mesh0.cellZones().clear();
1804  addZones
1805  (
1806  pointZoneNames,
1807  pzPoints,
1808 
1809  faceZoneNames,
1810  fzFaces,
1811  fzFlips,
1812 
1813  cellZoneNames,
1814  czCells,
1815  mesh0
1816  );
1817 
1818  return mapPtr;
1819 }
1820 
1821 
1822 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:453
ListType renumber(const labelUList &oldToNew, const ListType &)
Renumber the values (not the indices) of a list.
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
void inplaceReorder(const labelUList &oldToNew, ListType &)
Inplace reorder the elements of a list.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
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.
void reset(T *=nullptr)
If object pointer already set, delete object and set to given.
Definition: autoPtrI.H:114
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
label nFaces() const
const indirectPrimitivePatch & masterPatch() const
Addressing engine for coupled faces on mesh0.
List< face > faceList
Definition: faceListFwd.H:43
label nCells() const
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
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
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
Foam::indirectPrimitivePatch.
Container for information needed to couple to meshes. When constructed from two meshes and a list of ...
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
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
label nPoints
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
List< label > labelList
A List of labels.
Definition: labelList.H:56
MeshZones< pointZone, polyMesh > meshPointZones
A MeshZones with the type pointZone.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
MeshZones< cellZone, polyMesh > meshCellZones
A MeshZones with the type cellZone.
void addPatches(const List< polyPatch *> &, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:989
Foam::polyBoundaryMesh.
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
MeshZones< faceZone, polyMesh > meshFaceZones
A MeshZones with the type faceZone.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
label nPoints() const
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:76
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.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
List< cell > cellList
list of cells
Definition: cellList.H:42
IOporosityModelList pZones(mesh)