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-2024 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 "faceCoupleInfo.H"
29 #include "processorPolyPatch.H"
30 #include "Time.H"
31 
32 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33 
34 // Get index of patch in new set of patchnames/types
35 Foam::label Foam::polyMeshAdder::patchIndex
36 (
37  const polyPatch& p,
38  DynamicList<word>& allPatchNames,
39  DynamicList<word>& allPatchTypes
40 )
41 {
42  // Find the patch name on the list. If the patch is already there
43  // and patch types match, return index
44  const word& pType = p.type();
45  const word& pName = p.name();
46 
47  label patchi = findIndex(allPatchNames, pName);
48 
49  if (patchi == -1)
50  {
51  // Patch not found. Append to the list
52  allPatchNames.append(pName);
53  allPatchTypes.append(pType);
54 
55  return allPatchNames.size() - 1;
56  }
57  else if (allPatchTypes[patchi] == pType)
58  {
59  // Found name and types match
60  return patchi;
61  }
62  else
63  {
64  // Found the name, but type is different
65 
66  // Duplicate name is not allowed. Create a composite name from the
67  // patch name and case name
68  const word& caseName = p.boundaryMesh().mesh().time().caseName();
69 
70  allPatchNames.append(pName + "_" + caseName);
71  allPatchTypes.append(pType);
72 
73  Pout<< "label patchIndex(const polyPatch& p) : "
74  << "Patch " << p.index() << " named "
75  << pName << " in mesh " << caseName
76  << " already exists, but patch types"
77  << " do not match.\nCreating a composite name as "
78  << allPatchNames.last() << endl;
79 
80  return allPatchNames.size() - 1;
81  }
82 }
83 
84 
85 // Get index of zone in new set of zone names
86 Foam::label Foam::polyMeshAdder::zoneIndex
87 (
88  const word& curName,
89  DynamicList<word>& names
90 )
91 {
92  label zoneI = findIndex(names, curName);
93 
94  if (zoneI != -1)
95  {
96  return zoneI;
97  }
98  else
99  {
100  // Not found. Add new name to the list
101  names.append(curName);
102 
103  return names.size() - 1;
104  }
105 }
106 
107 
108 void Foam::polyMeshAdder::mergePatchNames
109 (
110  const polyBoundaryMesh& patches0,
111  const polyBoundaryMesh& patches1,
112 
113  DynamicList<word>& allPatchNames,
114  DynamicList<word>& allPatchTypes,
115 
116  labelList& from1ToAllPatches,
117  labelList& fromAllTo1Patches
118 )
119 {
120  // Insert the mesh0 patches and zones
121  allPatchNames.append(patches0.names());
122  allPatchTypes.append(patches0.types());
123 
124 
125  // Patches
126  // ~~~~~~~
127  // Patches from 0 are taken over as is; those from 1 get either merged
128  // (if they share name and type) or appended.
129  // Empty patches are filtered out much much later on.
130 
131  // Add mesh1 patches and build map both ways.
132  from1ToAllPatches.setSize(patches1.size());
133 
134  forAll(patches1, patchi)
135  {
136  from1ToAllPatches[patchi] = patchIndex
137  (
138  patches1[patchi],
139  allPatchNames,
140  allPatchTypes
141  );
142  }
143  allPatchTypes.shrink();
144  allPatchNames.shrink();
145 
146  // Invert 1 to all patch map
147  fromAllTo1Patches.setSize(allPatchNames.size());
148  fromAllTo1Patches = -1;
149 
150  forAll(from1ToAllPatches, i)
151  {
152  fromAllTo1Patches[from1ToAllPatches[i]] = i;
153  }
154 }
155 
156 
157 Foam::labelList Foam::polyMeshAdder::getPatchStarts
158 (
159  const polyBoundaryMesh& patches
160 )
161 {
162  labelList patchStarts(patches.size());
164  {
165  patchStarts[patchi] = patches[patchi].start();
166  }
167  return patchStarts;
168 }
169 
170 
171 Foam::labelList Foam::polyMeshAdder::getPatchSizes
172 (
173  const polyBoundaryMesh& patches
174 )
175 {
176  labelList patchSizes(patches.size());
178  {
179  patchSizes[patchi] = patches[patchi].size();
180  }
181  return patchSizes;
182 }
183 
184 
185 Foam::labelList Foam::polyMeshAdder::getFaceOrder
186 (
187  const cellList& cells,
188  const label nInternalFaces,
189  const labelList& owner,
190  const labelList& neighbour
191 )
192 {
193  labelList oldToNew(owner.size(), -1);
194 
195  // Leave boundary faces in order
196  for (label facei = nInternalFaces; facei < owner.size(); ++facei)
197  {
198  oldToNew[facei] = facei;
199  }
200 
201  // First unassigned face
202  label newFacei = 0;
203 
204  forAll(cells, celli)
205  {
206  const labelList& cFaces = cells[celli];
207 
208  SortableList<label> nbr(cFaces.size());
209 
210  forAll(cFaces, i)
211  {
212  label facei = cFaces[i];
213 
214  label nbrCelli = neighbour[facei];
215 
216  if (nbrCelli != -1)
217  {
218  // Internal face. Get cell on other side.
219  if (nbrCelli == celli)
220  {
221  nbrCelli = owner[facei];
222  }
223 
224  if (celli < nbrCelli)
225  {
226  // Celli is master
227  nbr[i] = nbrCelli;
228  }
229  else
230  {
231  // nbrCell is master. Let it handle this face.
232  nbr[i] = -1;
233  }
234  }
235  else
236  {
237  // External face. Do later.
238  nbr[i] = -1;
239  }
240  }
241 
242  nbr.sort();
243 
244  forAll(nbr, i)
245  {
246  if (nbr[i] != -1)
247  {
248  oldToNew[cFaces[nbr.indices()[i]]] = newFacei++;
249  }
250  }
251  }
252 
253 
254  // Check done all faces.
255  forAll(oldToNew, facei)
256  {
257  if (oldToNew[facei] == -1)
258  {
260  << "Did not determine new position"
261  << " for face " << facei
262  << abort(FatalError);
263  }
264  }
265 
266  return oldToNew;
267 }
268 
269 
270 // Adds primitives (cells, faces, points)
271 // Cells:
272 // - all of mesh0
273 // - all of mesh1
274 // Faces:
275 // - all uncoupled of mesh0
276 // - all coupled faces
277 // - all uncoupled of mesh1
278 // Points:
279 // - all coupled
280 // - all uncoupled of mesh0
281 // - all uncoupled of mesh1
282 void Foam::polyMeshAdder::mergePrimitives
283 (
284  const polyMesh& mesh0,
285  const polyMesh& mesh1,
286  const faceCoupleInfo& coupleInfo,
287 
288  const label nAllPatches, // number of patches in the new mesh
289  const labelList& fromAllTo1Patches,
290  const labelList& from1ToAllPatches,
291 
292  pointField& allPoints,
293  labelList& from0ToAllPoints,
294  labelList& from1ToAllPoints,
295 
296  faceList& allFaces,
297  labelList& allOwner,
298  labelList& allNeighbour,
299  label& nInternalFaces,
300  labelList& nFacesPerPatch,
301  label& nCells,
302 
303  labelList& from0ToAllFaces,
304  labelList& from1ToAllFaces,
305  labelList& from1ToAllCells
306 )
307 {
308  const polyBoundaryMesh& patches0 = mesh0.boundaryMesh();
309  const polyBoundaryMesh& patches1 = mesh1.boundaryMesh();
310 
311  const indirectPrimitivePatch& masterPatch = coupleInfo.masterPatch();
312  const indirectPrimitivePatch& slavePatch = coupleInfo.slavePatch();
313 
314 
315  // Points
316  // ~~~~~~
317 
318  // Storage for new points
319  allPoints.setSize(mesh0.nPoints() + mesh1.nPoints());
320  label allPointi = 0;
321 
322  from0ToAllPoints.setSize(mesh0.nPoints());
323  from0ToAllPoints = -1;
324  from1ToAllPoints.setSize(mesh1.nPoints());
325  from1ToAllPoints = -1;
326 
327  // Copy coupled points
328  {
329  const labelListList coupleToMasterPoints
330  (
331  coupleInfo.coupleToMasterPoints()
332  );
333  const labelListList coupleToSlavePoints
334  (
335  coupleInfo.coupleToSlavePoints()
336  );
337 
338  forAll(coupleToMasterPoints, couplePointi)
339  {
340  const labelList& masterPoints = coupleToMasterPoints[couplePointi];
341 
342  forAll(masterPoints, j)
343  {
344  label mesh0Pointi = masterPatch.meshPoints()[masterPoints[j]];
345  from0ToAllPoints[mesh0Pointi] = allPointi;
346  allPoints[allPointi] = mesh0.points()[mesh0Pointi];
347  }
348 
349  const labelList& slavePoints = coupleToSlavePoints[couplePointi];
350 
351  forAll(slavePoints, j)
352  {
353  label mesh1Pointi = slavePatch.meshPoints()[slavePoints[j]];
354  from1ToAllPoints[mesh1Pointi] = allPointi;
355  allPoints[allPointi] = mesh1.points()[mesh1Pointi];
356  }
357 
358  allPointi++;
359  }
360  }
361 
362  // Add uncoupled mesh0 points
363  forAll(mesh0.points(), pointi)
364  {
365  if (from0ToAllPoints[pointi] == -1)
366  {
367  allPoints[allPointi] = mesh0.points()[pointi];
368  from0ToAllPoints[pointi] = allPointi;
369  allPointi++;
370  }
371  }
372 
373  // Add uncoupled mesh1 points
374  forAll(mesh1.points(), pointi)
375  {
376  if (from1ToAllPoints[pointi] == -1)
377  {
378  allPoints[allPointi] = mesh1.points()[pointi];
379  from1ToAllPoints[pointi] = allPointi;
380  allPointi++;
381  }
382  }
383 
384  allPoints.setSize(allPointi);
385 
386 
387  // Faces
388  // ~~~~~
389 
390  // Sizes per patch
391  nFacesPerPatch.setSize(nAllPatches);
392  nFacesPerPatch = 0;
393 
394  // Storage for faces and owner/neighbour
395  allFaces.setSize(mesh0.nFaces() + mesh1.nFaces());
396  allOwner.setSize(allFaces.size());
397  allOwner = -1;
398  allNeighbour.setSize(allFaces.size());
399  allNeighbour = -1;
400  label allFacei = 0;
401 
402  from0ToAllFaces.setSize(mesh0.nFaces());
403  from0ToAllFaces = -1;
404  from1ToAllFaces.setSize(mesh1.nFaces());
405  from1ToAllFaces = -1;
406 
407  // Copy mesh0 internal faces (always uncoupled)
408  for (label facei = 0; facei < mesh0.nInternalFaces(); facei++)
409  {
410  allFaces[allFacei] = renumber(from0ToAllPoints, mesh0.faces()[facei]);
411  allOwner[allFacei] = mesh0.faceOwner()[facei];
412  allNeighbour[allFacei] = mesh0.faceNeighbour()[facei];
413  from0ToAllFaces[facei] = allFacei++;
414  }
415 
416  // Copy coupled faces. Every coupled face has an equivalent master and
417  // slave. Also uncount as boundary faces all the newly coupled faces.
418  forAll(masterPatch, coupleFacei)
419  {
420  const label mesh0Facei = masterPatch.addressing()[coupleFacei];
421 
422  if (from0ToAllFaces[mesh0Facei] == -1)
423  {
424  // First occurrence of face
425  from0ToAllFaces[mesh0Facei] = allFacei;
426 
427  // External face becomes internal so uncount
428  label patch0 = patches0.whichPatch(mesh0Facei);
429  nFacesPerPatch[patch0]--;
430  }
431 
432  const label mesh1Facei = slavePatch.addressing()[coupleFacei];
433 
434  if (from1ToAllFaces[mesh1Facei] == -1)
435  {
436  from1ToAllFaces[mesh1Facei] = allFacei;
437 
438  label patch1 = patches1.whichPatch(mesh1Facei);
439  nFacesPerPatch[from1ToAllPatches[patch1]]--;
440  }
441 
442  // Copy cut face (since cutPoints are copied first no renumbering
443  // necessary)
444  allFaces[allFacei] = coupleInfo.coupleFace(coupleFacei);
445  allOwner[allFacei] = mesh0.faceOwner()[mesh0Facei];
446  allNeighbour[allFacei] = mesh1.faceOwner()[mesh1Facei] + mesh0.nCells();
447 
448  allFacei++;
449  }
450 
451  // Copy mesh1 internal faces (always uncoupled)
452  for (label facei = 0; facei < mesh1.nInternalFaces(); facei++)
453  {
454  allFaces[allFacei] = renumber(from1ToAllPoints, mesh1.faces()[facei]);
455  allOwner[allFacei] = mesh1.faceOwner()[facei] + mesh0.nCells();
456  allNeighbour[allFacei] = mesh1.faceNeighbour()[facei] + mesh0.nCells();
457  from1ToAllFaces[facei] = allFacei++;
458  }
459 
460  nInternalFaces = allFacei;
461 
462  // Copy (unmarked/uncoupled) external faces in new order.
463  for (label allPatchi = 0; allPatchi < nAllPatches; allPatchi++)
464  {
465  if (allPatchi < patches0.size())
466  {
467  // Patch is present in mesh0
468  const polyPatch& pp = patches0[allPatchi];
469 
470  nFacesPerPatch[allPatchi] += pp.size();
471 
472  label facei = pp.start();
473 
474  forAll(pp, i)
475  {
476  if (from0ToAllFaces[facei] == -1)
477  {
478  // Is uncoupled face since has not yet been dealt with
479  allFaces[allFacei] = renumber
480  (
481  from0ToAllPoints,
482  mesh0.faces()[facei]
483  );
484  allOwner[allFacei] = mesh0.faceOwner()[facei];
485  allNeighbour[allFacei] = -1;
486 
487  from0ToAllFaces[facei] = allFacei++;
488  }
489  facei++;
490  }
491  }
492  if (fromAllTo1Patches[allPatchi] != -1)
493  {
494  // Patch is present in mesh1
495  const polyPatch& pp = patches1[fromAllTo1Patches[allPatchi]];
496 
497  nFacesPerPatch[allPatchi] += pp.size();
498 
499  label facei = pp.start();
500 
501  forAll(pp, i)
502  {
503  if (from1ToAllFaces[facei] == -1)
504  {
505  // Is uncoupled face
506  allFaces[allFacei] = renumber
507  (
508  from1ToAllPoints,
509  mesh1.faces()[facei]
510  );
511  allOwner[allFacei] =
512  mesh1.faceOwner()[facei]
513  + mesh0.nCells();
514  allNeighbour[allFacei] = -1;
515 
516  from1ToAllFaces[facei] = allFacei++;
517  }
518  facei++;
519  }
520  }
521  }
522  allFaces.setSize(allFacei);
523  allOwner.setSize(allFacei);
524  allNeighbour.setSize(allFacei);
525 
526  // Cells
527  // ~~~~~
528 
529  from1ToAllCells.setSize(mesh1.nCells());
530  from1ToAllCells = -1;
531 
532  forAll(mesh1.cells(), i)
533  {
534  from1ToAllCells[i] = i + mesh0.nCells();
535  }
536 
537  // Make cells (= cell-face addressing)
538  nCells = mesh0.nCells() + mesh1.nCells();
539  cellList allCells(nCells);
540  primitiveMesh::calcCells(allCells, allOwner, allNeighbour, nCells);
541 
542  // Reorder faces for upper-triangular order.
543  labelList oldToNew
544  (
545  getFaceOrder
546  (
547  allCells,
548  nInternalFaces,
549  allOwner,
550  allNeighbour
551  )
552  );
553 
554  inplaceReorder(oldToNew, allFaces);
555  inplaceReorder(oldToNew, allOwner);
556  inplaceReorder(oldToNew, allNeighbour);
557  inplaceRenumber(oldToNew, from0ToAllFaces);
558  inplaceRenumber(oldToNew, from1ToAllFaces);
559 }
560 
561 
562 void Foam::polyMeshAdder::mergePointZones
563 (
564  const label nAllPoints,
565  const pointZoneList& pz0,
566  const pointZoneList& pz1,
567  const labelList& from0ToAllPoints,
568  const labelList& from1ToAllPoints,
569 
570  DynamicList<word>& zoneNames,
571  labelList& from1ToAll,
572  List<DynamicList<label>>& pzPoints
573 )
574 {
575  zoneNames.setCapacity(pz0.size() + pz1.size());
576  zoneNames.append(pz0.toc());
577 
578  from1ToAll.setSize(pz1.size());
579 
580  forAll(pz1, zoneI)
581  {
582  from1ToAll[zoneI] = zoneIndex(pz1[zoneI].name(), zoneNames);
583  }
584  zoneNames.shrink();
585 
586 
587 
588  // Zone(s) per point. Two levels: if only one zone
589  // stored in pointToZone. Any extra stored in additionalPointToZones.
590  // This is so we only allocate labelLists per point if absolutely
591  // necessary.
592  labelList pointToZone(nAllPoints, -1);
593  labelListList addPointToZones(nAllPoints);
594 
595  // mesh0 zones kept
596  forAll(pz0, zoneI)
597  {
598  const pointZone& pz = pz0[zoneI];
599 
600  forAll(pz, i)
601  {
602  label point0 = pz[i];
603  label allPointi = from0ToAllPoints[point0];
604 
605  if (pointToZone[allPointi] == -1)
606  {
607  pointToZone[allPointi] = zoneI;
608  }
609  else if (pointToZone[allPointi] != zoneI)
610  {
611  labelList& pZones = addPointToZones[allPointi];
612  if (findIndex(pZones, zoneI) == -1)
613  {
614  pZones.append(zoneI);
615  }
616  }
617  }
618  }
619 
620  // mesh1 zones renumbered
621  forAll(pz1, zoneI)
622  {
623  const pointZone& pz = pz1[zoneI];
624  const label allZoneI = from1ToAll[zoneI];
625 
626  forAll(pz, i)
627  {
628  label point1 = pz[i];
629  label allPointi = from1ToAllPoints[point1];
630 
631  if (pointToZone[allPointi] == -1)
632  {
633  pointToZone[allPointi] = allZoneI;
634  }
635  else if (pointToZone[allPointi] != allZoneI)
636  {
637  labelList& pZones = addPointToZones[allPointi];
638  if (findIndex(pZones, allZoneI) == -1)
639  {
640  pZones.append(allZoneI);
641  }
642  }
643  }
644  }
645 
646  // Extract back into zones
647 
648  // 1. Count
649  labelList nPoints(zoneNames.size(), 0);
650  forAll(pointToZone, allPointi)
651  {
652  label zoneI = pointToZone[allPointi];
653  if (zoneI != -1)
654  {
655  nPoints[zoneI]++;
656  }
657  }
658  forAll(addPointToZones, allPointi)
659  {
660  const labelList& pZones = addPointToZones[allPointi];
661  forAll(pZones, i)
662  {
663  nPoints[pZones[i]]++;
664  }
665  }
666 
667  // 2. Fill
668  pzPoints.setSize(zoneNames.size());
669  forAll(pzPoints, zoneI)
670  {
671  pzPoints[zoneI].setCapacity(nPoints[zoneI]);
672  }
673  forAll(pointToZone, allPointi)
674  {
675  label zoneI = pointToZone[allPointi];
676  if (zoneI != -1)
677  {
678  pzPoints[zoneI].append(allPointi);
679  }
680  }
681  forAll(addPointToZones, allPointi)
682  {
683  const labelList& pZones = addPointToZones[allPointi];
684  forAll(pZones, i)
685  {
686  pzPoints[pZones[i]].append(allPointi);
687  }
688  }
689  forAll(pzPoints, i)
690  {
691  pzPoints[i].shrink();
692  stableSort(pzPoints[i]);
693  }
694 }
695 
696 
697 void Foam::polyMeshAdder::mergeFaceZones
698 (
699  const labelList& allOwner,
700 
701  const polyMesh& mesh0,
702  const polyMesh& mesh1,
703 
704  const labelList& from0ToAllFaces,
705  const labelList& from1ToAllFaces,
706 
707  const labelList& from1ToAllCells,
708 
709  DynamicList<word>& zoneNames,
710  labelList& from1ToAll,
711  List<DynamicList<label>>& fzFaces,
712  List<DynamicList<bool>>& fzFlips
713 )
714 {
715  const faceZoneList& fz0 = mesh0.faceZones();
716  const labelList& owner0 = mesh0.faceOwner();
717  const faceZoneList& fz1 = mesh1.faceZones();
718  const labelList& owner1 = mesh1.faceOwner();
719 
720 
721  zoneNames.setCapacity(fz0.size() + fz1.size());
722  zoneNames.append(fz0.toc());
723 
724  from1ToAll.setSize(fz1.size());
725 
726  forAll(fz1, zoneI)
727  {
728  from1ToAll[zoneI] = zoneIndex(fz1[zoneI].name(), zoneNames);
729  }
730  zoneNames.shrink();
731 
732 
733  // Zone(s) per face
734  labelList faceToZone(allOwner.size(), -1);
735  labelListList addFaceToZones(allOwner.size());
736  boolList faceToFlip(allOwner.size(), false);
737  boolListList addFaceToFlips(allOwner.size());
738 
739  // mesh0 zones kept
740  forAll(fz0, zoneI)
741  {
742  const labelList& addressing = fz0[zoneI];
743  const boolList& flipMap = fz0[zoneI].flipMap();
744 
745  forAll(addressing, i)
746  {
747  label face0 = addressing[i];
748  bool flip0 = flipMap[i];
749 
750  label allFacei = from0ToAllFaces[face0];
751  if (allFacei != -1)
752  {
753  // Check if orientation same
754  label allCell0 = owner0[face0];
755  if (allOwner[allFacei] != allCell0)
756  {
757  flip0 = !flip0;
758  }
759 
760  if (faceToZone[allFacei] == -1)
761  {
762  faceToZone[allFacei] = zoneI;
763  faceToFlip[allFacei] = flip0;
764  }
765  else if (faceToZone[allFacei] != zoneI)
766  {
767  labelList& fZones = addFaceToZones[allFacei];
768  boolList& flipZones = addFaceToFlips[allFacei];
769 
770  if (findIndex(fZones, zoneI) == -1)
771  {
772  fZones.append(zoneI);
773  flipZones.append(flip0);
774  }
775  }
776  }
777  }
778  }
779 
780  // mesh1 zones renumbered
781  forAll(fz1, zoneI)
782  {
783  const labelList& addressing = fz1[zoneI];
784  const boolList& flipMap = fz1[zoneI].flipMap();
785 
786  const label allZoneI = from1ToAll[zoneI];
787 
788  forAll(addressing, i)
789  {
790  label face1 = addressing[i];
791  bool flip1 = flipMap[i];
792 
793  label allFacei = from1ToAllFaces[face1];
794  if (allFacei != -1)
795  {
796  // Check if orientation same
797  label allCell1 = from1ToAllCells[owner1[face1]];
798  if (allOwner[allFacei] != allCell1)
799  {
800  flip1 = !flip1;
801  }
802 
803  if (faceToZone[allFacei] == -1)
804  {
805  faceToZone[allFacei] = allZoneI;
806  faceToFlip[allFacei] = flip1;
807  }
808  else if (faceToZone[allFacei] != allZoneI)
809  {
810  labelList& fZones = addFaceToZones[allFacei];
811  boolList& flipZones = addFaceToFlips[allFacei];
812 
813  if (findIndex(fZones, allZoneI) == -1)
814  {
815  fZones.append(allZoneI);
816  flipZones.append(flip1);
817  }
818  }
819  }
820  }
821  }
822 
823 
824  // Extract back into zones
825 
826  // 1. Count
827  labelList nFaces(zoneNames.size(), 0);
828  forAll(faceToZone, allFacei)
829  {
830  label zoneI = faceToZone[allFacei];
831  if (zoneI != -1)
832  {
833  nFaces[zoneI]++;
834  }
835  }
836  forAll(addFaceToZones, allFacei)
837  {
838  const labelList& fZones = addFaceToZones[allFacei];
839  forAll(fZones, i)
840  {
841  nFaces[fZones[i]]++;
842  }
843  }
844 
845  // 2. Fill
846  fzFaces.setSize(zoneNames.size());
847  fzFlips.setSize(zoneNames.size());
848  forAll(fzFaces, zoneI)
849  {
850  fzFaces[zoneI].setCapacity(nFaces[zoneI]);
851  fzFlips[zoneI].setCapacity(nFaces[zoneI]);
852  }
853  forAll(faceToZone, allFacei)
854  {
855  label zoneI = faceToZone[allFacei];
856  bool flip = faceToFlip[allFacei];
857  if (zoneI != -1)
858  {
859  fzFaces[zoneI].append(allFacei);
860  fzFlips[zoneI].append(flip);
861  }
862  }
863  forAll(addFaceToZones, allFacei)
864  {
865  const labelList& fZones = addFaceToZones[allFacei];
866  const boolList& flipZones = addFaceToFlips[allFacei];
867 
868  forAll(fZones, i)
869  {
870  label zoneI = fZones[i];
871  fzFaces[zoneI].append(allFacei);
872  fzFlips[zoneI].append(flipZones[i]);
873  }
874  }
875 
876  forAll(fzFaces, i)
877  {
878  fzFaces[i].shrink();
879  fzFlips[i].shrink();
880 
882  sortedOrder(fzFaces[i], order);
883  fzFaces[i] = UIndirectList<label>(fzFaces[i], order)();
884  fzFlips[i] = UIndirectList<bool>(fzFlips[i], order)();
885  }
886 }
887 
888 
889 void Foam::polyMeshAdder::mergeCellZones
890 (
891  const label nAllCells,
892 
893  const cellZoneList& cz0,
894  const cellZoneList& cz1,
895  const labelList& from1ToAllCells,
896 
897  DynamicList<word>& zoneNames,
898  labelList& from1ToAll,
899  List<DynamicList<label>>& czCells
900 )
901 {
902  zoneNames.setCapacity(cz0.size() + cz1.size());
903  zoneNames.append(cz0.toc());
904 
905  from1ToAll.setSize(cz1.size());
906  forAll(cz1, zoneI)
907  {
908  from1ToAll[zoneI] = zoneIndex(cz1[zoneI].name(), zoneNames);
909  }
910  zoneNames.shrink();
911 
912 
913  // Zone(s) per cell. Two levels: if only one zone
914  // stored in cellToZone. Any extra stored in additionalCellToZones.
915  // This is so we only allocate labelLists per cell if absolutely
916  // necessary.
917  labelList cellToZone(nAllCells, -1);
918  labelListList addCellToZones(nAllCells);
919 
920  // mesh0 zones kept
921  forAll(cz0, zoneI)
922  {
923  const cellZone& cz = cz0[zoneI];
924  forAll(cz, i)
925  {
926  label cell0 = cz[i];
927 
928  if (cellToZone[cell0] == -1)
929  {
930  cellToZone[cell0] = zoneI;
931  }
932  else if (cellToZone[cell0] != zoneI)
933  {
934  labelList& cZones = addCellToZones[cell0];
935  if (findIndex(cZones, zoneI) == -1)
936  {
937  cZones.append(zoneI);
938  }
939  }
940  }
941  }
942 
943  // mesh1 zones renumbered
944  forAll(cz1, zoneI)
945  {
946  const cellZone& cz = cz1[zoneI];
947  const label allZoneI = from1ToAll[zoneI];
948  forAll(cz, i)
949  {
950  label cell1 = cz[i];
951  label allCelli = from1ToAllCells[cell1];
952 
953  if (cellToZone[allCelli] == -1)
954  {
955  cellToZone[allCelli] = allZoneI;
956  }
957  else if (cellToZone[allCelli] != allZoneI)
958  {
959  labelList& cZones = addCellToZones[allCelli];
960  if (findIndex(cZones, allZoneI) == -1)
961  {
962  cZones.append(allZoneI);
963  }
964  }
965  }
966  }
967 
968  // Extract back into zones
969 
970  // 1. Count
971  labelList nCells(zoneNames.size(), 0);
972  forAll(cellToZone, allCelli)
973  {
974  label zoneI = cellToZone[allCelli];
975  if (zoneI != -1)
976  {
977  nCells[zoneI]++;
978  }
979  }
980  forAll(addCellToZones, allCelli)
981  {
982  const labelList& cZones = addCellToZones[allCelli];
983  forAll(cZones, i)
984  {
985  nCells[cZones[i]]++;
986  }
987  }
988 
989  // 2. Fill
990  czCells.setSize(zoneNames.size());
991  forAll(czCells, zoneI)
992  {
993  czCells[zoneI].setCapacity(nCells[zoneI]);
994  }
995  forAll(cellToZone, allCelli)
996  {
997  label zoneI = cellToZone[allCelli];
998  if (zoneI != -1)
999  {
1000  czCells[zoneI].append(allCelli);
1001  }
1002  }
1003  forAll(addCellToZones, allCelli)
1004  {
1005  const labelList& cZones = addCellToZones[allCelli];
1006  forAll(cZones, i)
1007  {
1008  czCells[cZones[i]].append(allCelli);
1009  }
1010  }
1011  forAll(czCells, i)
1012  {
1013  czCells[i].shrink();
1014  stableSort(czCells[i]);
1015  }
1016 }
1017 
1018 
1019 void Foam::polyMeshAdder::mergeZones
1020 (
1021  const label nAllPoints,
1022  const labelList& allOwner,
1023  const label nAllCells,
1024 
1025  const polyMesh& mesh0,
1026  const polyMesh& mesh1,
1027  const labelList& from0ToAllPoints,
1028  const labelList& from0ToAllFaces,
1029 
1030  const labelList& from1ToAllPoints,
1031  const labelList& from1ToAllFaces,
1032  const labelList& from1ToAllCells,
1033 
1034  DynamicList<word>& pointZoneNames,
1035  List<DynamicList<label>>& pzPoints,
1036 
1037  DynamicList<word>& faceZoneNames,
1038  List<DynamicList<label>>& fzFaces,
1039  List<DynamicList<bool>>& fzFlips,
1040 
1041  DynamicList<word>& cellZoneNames,
1042  List<DynamicList<label>>& czCells
1043 )
1044 {
1045  labelList from1ToAllPZones;
1046  mergePointZones
1047  (
1048  nAllPoints,
1049  mesh0.pointZones(),
1050  mesh1.pointZones(),
1051  from0ToAllPoints,
1052  from1ToAllPoints,
1053 
1054  pointZoneNames,
1055  from1ToAllPZones,
1056  pzPoints
1057  );
1058 
1059  labelList from1ToAllFZones;
1060  mergeFaceZones
1061  (
1062  allOwner,
1063  mesh0,
1064  mesh1,
1065  from0ToAllFaces,
1066  from1ToAllFaces,
1067  from1ToAllCells,
1068 
1069  faceZoneNames,
1070  from1ToAllFZones,
1071  fzFaces,
1072  fzFlips
1073  );
1074 
1075  labelList from1ToAllCZones;
1076  mergeCellZones
1077  (
1078  nAllCells,
1079  mesh0.cellZones(),
1080  mesh1.cellZones(),
1081  from1ToAllCells,
1082 
1083  cellZoneNames,
1084  from1ToAllCZones,
1085  czCells
1086  );
1087 }
1088 
1089 
1090 void Foam::polyMeshAdder::addZones
1091 (
1092  const DynamicList<word>& pointZoneNames,
1093  const List<DynamicList<label>>& pzPoints,
1094 
1095  const DynamicList<word>& faceZoneNames,
1096  const List<DynamicList<label>>& fzFaces,
1097  const List<DynamicList<bool>>& fzFlips,
1098 
1099  const DynamicList<word>& cellZoneNames,
1100  const List<DynamicList<label>>& czCells,
1101 
1102  polyMesh& mesh
1103 )
1104 {
1105  List<pointZone*> pZones(pzPoints.size());
1106  forAll(pZones, i)
1107  {
1108  pZones[i] = new pointZone
1109  (
1110  pointZoneNames[i],
1111  pzPoints[i],
1112  mesh.pointZones()
1113  );
1114  }
1115 
1116  List<faceZone*> fZones(fzFaces.size());
1117  forAll(fZones, i)
1118  {
1119  fZones[i] = new faceZone
1120  (
1121  faceZoneNames[i],
1122  fzFaces[i],
1123  fzFlips[i],
1124  mesh.faceZones()
1125  );
1126  }
1127 
1128  List<cellZone*> cZones(czCells.size());
1129  forAll(cZones, i)
1130  {
1131  cZones[i] = new cellZone
1132  (
1133  cellZoneNames[i],
1134  czCells[i],
1135  mesh.cellZones()
1136  );
1137  }
1138 
1139  mesh.addZones
1140  (
1141  pZones,
1142  fZones,
1143  cZones
1144  );
1145 }
1146 
1147 
1148 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1149 
1151 (
1152  polyMesh& mesh0,
1153  const polyMesh& mesh1,
1154  const faceCoupleInfo& coupleInfo,
1155  const bool validBoundary
1156 )
1157 {
1158  const polyBoundaryMesh& patches0 = mesh0.boundaryMesh();
1159  const polyBoundaryMesh& patches1 = mesh1.boundaryMesh();
1160 
1161  DynamicList<word> allPatchNames(patches0.size() + patches1.size());
1162  DynamicList<word> allPatchTypes(allPatchNames.size());
1163 
1164  // Patch maps
1165  labelList from1ToAllPatches(patches1.size());
1166  labelList fromAllTo1Patches(allPatchNames.size(), -1);
1167 
1168  mergePatchNames
1169  (
1170  patches0,
1171  patches1,
1172  allPatchNames,
1173  allPatchTypes,
1174  from1ToAllPatches,
1175  fromAllTo1Patches
1176  );
1177 
1178 
1179  // New points
1180  pointField allPoints;
1181 
1182  // Map from mesh0/1 points to allPoints.
1183  labelList from0ToAllPoints(mesh0.nPoints(), -1);
1184  labelList from1ToAllPoints(mesh1.nPoints(), -1);
1185 
1186  // New faces
1187  faceList allFaces;
1188  labelList allOwner;
1189  labelList allNeighbour;
1190  label nInternalFaces;
1191  // Sizes per patch
1192  labelList nFaces(allPatchNames.size(), 0);
1193  label nCells;
1194 
1195  // Maps
1196  labelList from0ToAllFaces(mesh0.nFaces(), -1);
1197  labelList from1ToAllFaces(mesh1.nFaces(), -1);
1198  // Map
1199  labelList from1ToAllCells(mesh1.nCells(), -1);
1200 
1201  mergePrimitives
1202  (
1203  mesh0,
1204  mesh1,
1205  coupleInfo,
1206 
1207  allPatchNames.size(),
1208  fromAllTo1Patches,
1209  from1ToAllPatches,
1210 
1211  allPoints,
1212  from0ToAllPoints,
1213  from1ToAllPoints,
1214 
1215  allFaces,
1216  allOwner,
1217  allNeighbour,
1218  nInternalFaces,
1219  nFaces,
1220  nCells,
1221 
1222  from0ToAllFaces,
1223  from1ToAllFaces,
1224  from1ToAllCells
1225  );
1226 
1227 
1228  // Zones
1229  // ~~~~~
1230 
1231  DynamicList<word> pointZoneNames;
1232  List<DynamicList<label>> pzPoints;
1233 
1234  DynamicList<word> faceZoneNames;
1235  List<DynamicList<label>> fzFaces;
1236  List<DynamicList<bool>> fzFlips;
1237 
1238  DynamicList<word> cellZoneNames;
1239  List<DynamicList<label>> czCells;
1240 
1241  mergeZones
1242  (
1243  allPoints.size(),
1244  allOwner,
1245  nCells,
1246 
1247  mesh0,
1248  mesh1,
1249 
1250  from0ToAllPoints,
1251  from0ToAllFaces,
1252 
1253  from1ToAllPoints,
1254  from1ToAllFaces,
1255  from1ToAllCells,
1256 
1257  pointZoneNames,
1258  pzPoints,
1259 
1260  faceZoneNames,
1261  fzFaces,
1262  fzFlips,
1263 
1264  cellZoneNames,
1265  czCells
1266  );
1267 
1268 
1269  // Patches
1270  // ~~~~~~~
1271 
1272 
1273  // Store mesh0 patch info before modifying patches0.
1274  labelList mesh0PatchSizes(getPatchSizes(patches0));
1275  labelList mesh0PatchStarts(getPatchStarts(patches0));
1276 
1277  // Map from 0 to all patches (since gets compacted)
1278  labelList from0ToAllPatches(patches0.size(), -1);
1279 
1280  // Inplace extend mesh0 patches (note that patches0.size() now also
1281  // has changed)
1282  labelList patchSizes(allPatchNames.size());
1283  labelList patchStarts(allPatchNames.size());
1284 
1285  label startFacei = nInternalFaces;
1286 
1287  // Copy patches0 with new sizes. First patches always come from
1288  // mesh0 and will always be present.
1289  label allPatchi = 0;
1290 
1291  forAll(from0ToAllPatches, patch0)
1292  {
1293  // Originates from mesh0. Clone with new size & filter out empty
1294  // patch.
1295 
1296  if (nFaces[patch0] == 0 && isA<processorPolyPatch>(patches0[patch0]))
1297  {
1298  // Pout<< "Removing zero sized mesh0 patch "
1299  // << allPatchNames[patch0]
1300  // << endl;
1301  from0ToAllPatches[patch0] = -1;
1302 
1303  // Check if patch was also in mesh1 and update its addressing if so.
1304  if (fromAllTo1Patches[patch0] != -1)
1305  {
1306  from1ToAllPatches[fromAllTo1Patches[patch0]] = -1;
1307  }
1308  }
1309  else
1310  {
1311  patchSizes[allPatchi] = nFaces[patch0];
1312  patchStarts[allPatchi] = startFacei;
1313 
1314  // Record new index in allPatches
1315  from0ToAllPatches[patch0] = allPatchi;
1316 
1317  // Check if patch was also in mesh1 and update its addressing if so.
1318  if (fromAllTo1Patches[patch0] != -1)
1319  {
1320  from1ToAllPatches[fromAllTo1Patches[patch0]] = allPatchi;
1321  }
1322 
1323  startFacei += nFaces[patch0];
1324 
1325  allPatchi++;
1326  }
1327  }
1328 
1329  // Trim the existing patches
1330  {
1331  const label sz0 = from0ToAllPatches.size();
1332  labelList newToOld(sz0, sz0-1);
1333  label nNew = 0;
1334  forAll(from0ToAllPatches, patchi)
1335  {
1336  if (from0ToAllPatches[patchi] != -1)
1337  {
1338  newToOld[nNew++] = patchi;
1339  }
1340  }
1341  newToOld.setSize(nNew);
1342  mesh0.reorderPatches(newToOld, false);
1343  }
1344 
1345  // Copy unique patches of mesh1.
1346  forAll(from1ToAllPatches, patch1)
1347  {
1348  label uncompactAllPatchi = from1ToAllPatches[patch1];
1349 
1350  if (uncompactAllPatchi >= from0ToAllPatches.size())
1351  {
1352  // Patch has not been merged with any mesh0 patch.
1353 
1354  if
1355  (
1356  nFaces[uncompactAllPatchi] == 0
1357  && isA<processorPolyPatch>(patches1[patch1])
1358  )
1359  {
1360  // Pout<< "Removing zero sized mesh1 patch "
1361  // << allPatchNames[uncompactAllPatchi] << endl;
1362  from1ToAllPatches[patch1] = -1;
1363  }
1364  else
1365  {
1366  patchSizes[allPatchi] = nFaces[uncompactAllPatchi];
1367  patchStarts[allPatchi] = startFacei;
1368 
1369  // Clone. Note dummy size and start. Gets overwritten later in
1370  // resetPrimitives. This avoids getting temporarily illegal
1371  // SubList construction in polyPatch.
1372  mesh0.addPatch
1373  (
1374  allPatchi,
1375  patches1[patch1],
1376  dictionary(),
1377  "calculated",
1378  false
1379  );
1380 
1381  // Record new index in allPatches
1382  from1ToAllPatches[patch1] = allPatchi;
1383 
1384  startFacei += nFaces[uncompactAllPatchi];
1385 
1386  allPatchi++;
1387  }
1388  }
1389  }
1390  patchSizes.setSize(allPatchi);
1391  patchStarts.setSize(allPatchi);
1392 
1393 
1394  // Construct map information before changing mesh0 primitives
1395  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1396 
1398  (
1399  new mapAddedPolyMesh
1400  (
1401  mesh0.nPoints(),
1402  mesh0.nFaces(),
1403  mesh0.nCells(),
1404 
1405  mesh1.nPoints(),
1406  mesh1.nFaces(),
1407  mesh1.nCells(),
1408 
1409  from0ToAllPoints,
1410  from0ToAllFaces,
1411  identityMap(mesh0.nCells()),
1412 
1413  from1ToAllPoints,
1414  from1ToAllFaces,
1415  from1ToAllCells,
1416 
1417  from0ToAllPatches,
1418  from1ToAllPatches,
1419 
1420  mesh0PatchSizes,
1421  mesh0PatchStarts
1422  )
1423  );
1424 
1425 
1426 
1427  // Now we have extracted all information from all meshes.
1428  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1429 
1430  mesh0.resetMotion(); // delete any oldPoints.
1431  mesh0.resetPrimitives
1432  (
1433  move(allPoints),
1434  move(allFaces),
1435  move(allOwner),
1436  move(allNeighbour),
1437  patchSizes, // size
1438  patchStarts, // patchstarts
1439  validBoundary // boundary valid?
1440  );
1441 
1442  // Add zones to new mesh.
1443  mesh0.pointZones().clear();
1444  mesh0.faceZones().clear();
1445  mesh0.cellZones().clear();
1446  addZones
1447  (
1448  pointZoneNames,
1449  pzPoints,
1450 
1451  faceZoneNames,
1452  fzFaces,
1453  fzFlips,
1454 
1455  cellZoneNames,
1456  czCells,
1457  mesh0
1458  );
1459 
1460  return mapPtr;
1461 }
1462 
1463 
1464 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
void append(T *)
Append an element at the end of the list.
Definition: PtrListI.H:39
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
void clear()
Clear the zones.
Definition: ZoneList.C:236
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Container for information needed to couple to meshes. When constructed from two meshes and a list of ...
Class containing mesh-to-mesh mapping information after a mesh addition where we add a mesh ('added m...
Foam::polyBoundaryMesh.
static autoPtr< mapAddedPolyMesh > add(polyMesh &mesh0, const polyMesh &mesh1, const faceCoupleInfo &coupleInfo, const bool validBoundary=true)
Inplace add mesh to polyMesh. Returns map construct.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const pointZoneList & pointZones() const
Return point zones.
Definition: polyMesh.H:440
const cellZoneList & cellZones() const
Return cell zones.
Definition: polyMesh.H:452
virtual void addPatch(const label insertPatchi, const polyPatch &patch, const dictionary &patchFieldDict, const word &defaultPatchFieldType, const bool validBoundary)
Add/insert single patch. If validBoundary the new situation.
Definition: polyMesh.C:1231
void resetPrimitives(pointField &&points, faceList &&faces, labelList &&owner, labelList &&neighbour, const labelList &patchSizes, const labelList &patchStarts, const bool validBoundary=true)
Reset mesh primitive data. Assumes all patch info correct.
Definition: polyMesh.C:703
void resetMotion() const
Reset motion.
Definition: polyMesh.C:1507
const faceZoneList & faceZones() const
Return face zones.
Definition: polyMesh.H:446
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
virtual void reorderPatches(const labelUList &newToOld, const bool validBoundary)
Reorder and trim existing patches. If validBoundary the new.
Definition: polyMesh.C:1186
label nCells() const
label nPoints() const
label nFaces() const
IOporosityModelList pZones(mesh)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
label nPoints
const cellShapeList & cells
const fvPatchList & patches
List< label > labelList
A List of labels.
Definition: labelList.H:56
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< cell > cellList
list of cells
Definition: cellList.H:42
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
int order(const scalar s)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
ListType renumber(const labelUList &oldToNew, const ListType &)
Renumber the values (not the indices) of a list.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
List< List< bool > > boolListList
Definition: boolList.H:51
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
error FatalError
void stableSort(UList< T > &)
Definition: UList.C:129
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
void sortedOrder(const UList< T > &, labelList &order)
Generate the (stable) sort order for the list.
List< face > faceList
Definition: faceListFwd.H:41
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
Foam::indirectPrimitivePatch.
void inplaceReorder(const labelUList &oldToNew, ListType &)
Inplace reorder the elements of a list.
volScalarField & p