meshRefinementBaffles.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-2022 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 "meshRefinement.H"
27 #include "refinementSurfaces.H"
28 #include "faceSet.H"
29 #include "polyTopoChange.H"
30 #include "meshTools.H"
31 #include "polyModifyFace.H"
32 #include "polyModifyCell.H"
33 #include "polyAddFace.H"
34 #include "polyRemoveFace.H"
35 #include "localPointRegion.H"
36 #include "duplicatePoints.H"
37 #include "regionSplit.H"
38 #include "removeCells.H"
39 #include "unitConversion.H"
40 #include "OBJstream.H"
41 #include "patchFaceOrientation.H"
42 #include "PatchEdgeFaceWave.H"
43 #include "patchEdgeFaceRegion.H"
44 #include "OSspecific.H"
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
48 Foam::label Foam::meshRefinement::createBaffle
49 (
50  const label facei,
51  const label ownPatch,
52  const label nbrPatch,
53  polyTopoChange& meshMod
54 ) const
55 {
56  const face& f = mesh_.faces()[facei];
57  const label zoneID = mesh_.faceZones().whichZone(facei);
58  bool zoneFlip = false;
59 
60  if (zoneID >= 0)
61  {
62  const faceZone& fZone = mesh_.faceZones()[zoneID];
63  zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
64  }
65 
66  meshMod.setAction
67  (
68  polyModifyFace
69  (
70  f, // modified face
71  facei, // label of face
72  mesh_.faceOwner()[facei], // owner
73  -1, // neighbour
74  false, // face flip
75  ownPatch, // patch for face
76  false, // remove from zone
77  zoneID, // zone for face
78  zoneFlip // face flip in zone
79  )
80  );
81 
82 
83  label dupFacei = -1;
84 
85  if (mesh_.isInternalFace(facei))
86  {
87  if (nbrPatch == -1)
88  {
90  << "No neighbour patch for internal face " << facei
91  << " fc:" << mesh_.faceCentres()[facei]
92  << " ownPatch:" << ownPatch << abort(FatalError);
93  }
94 
95  bool reverseFlip = false;
96  if (zoneID >= 0)
97  {
98  reverseFlip = !zoneFlip;
99  }
100 
101  dupFacei = meshMod.setAction
102  (
103  polyAddFace
104  (
105  f.reverseFace(), // modified face
106  mesh_.faceNeighbour()[facei],// owner
107  -1, // neighbour
108  -1, // masterPointID
109  -1, // masterEdgeID
110  facei, // masterFaceID,
111  true, // face flip
112  nbrPatch, // patch for face
113  zoneID, // zone for face
114  reverseFlip // face flip in zone
115  )
116  );
117  }
118  return dupFacei;
119 }
120 
121 
122 void Foam::meshRefinement::getBafflePatches
123 (
124  const labelList& globalToMasterPatch,
125  const labelList& neiLevel,
126  const pointField& neiCc,
127  labelList& ownPatch,
128  labelList& nbrPatch
129 ) const
130 {
131  const pointField& cellCentres = mesh_.cellCentres();
132 
133  // Surfaces that need to be baffled
134  const labelList surfacesToBaffle
135  (
137  );
138 
139  ownPatch.setSize(mesh_.nFaces());
140  ownPatch = -1;
141  nbrPatch.setSize(mesh_.nFaces());
142  nbrPatch = -1;
143 
144 
145  // Collect candidate faces
146  // ~~~~~~~~~~~~~~~~~~~~~~~
147 
148  labelList testFaces(intersectedFaces());
149 
150  // Collect segments
151  // ~~~~~~~~~~~~~~~~
152 
153  pointField start(testFaces.size());
154  pointField end(testFaces.size());
155 
156  forAll(testFaces, i)
157  {
158  const label facei = testFaces[i];
159  const label own = mesh_.faceOwner()[facei];
160 
161  if (mesh_.isInternalFace(facei))
162  {
163  start[i] = cellCentres[own];
164  end[i] = cellCentres[mesh_.faceNeighbour()[facei]];
165  }
166  else
167  {
168  start[i] = cellCentres[own];
169  end[i] = neiCc[facei-mesh_.nInternalFaces()];
170  }
171  }
172 
173  // Extend segments a bit
174  {
175  const vectorField smallVec(rootSmall*(end-start));
176  start -= smallVec;
177  end += smallVec;
178  }
179 
180 
181  // Do test for intersections
182  // ~~~~~~~~~~~~~~~~~~~~~~~~~
183  labelList surface1;
184  List<pointIndexHit> hit1;
185  labelList region1;
186  labelList surface2;
187  List<pointIndexHit> hit2;
188  labelList region2;
189  surfaces_.findNearestIntersection
190  (
191  surfacesToBaffle,
192  start,
193  end,
194 
195  surface1,
196  hit1,
197  region1,
198 
199  surface2,
200  hit2,
201  region2
202  );
203 
204  forAll(testFaces, i)
205  {
206  const label facei = testFaces[i];
207 
208  if (hit1[i].hit() && hit2[i].hit())
209  {
210  // Pick up the patches
211  ownPatch[facei] = globalToMasterPatch
212  [
213  surfaces_.globalRegion(surface1[i], region1[i])
214  ];
215  nbrPatch[facei] = globalToMasterPatch
216  [
217  surfaces_.globalRegion(surface2[i], region2[i])
218  ];
219 
220  if (ownPatch[facei] == -1 || nbrPatch[facei] == -1)
221  {
223  << "problem." << abort(FatalError);
224  }
225  }
226  }
227 
228  // No need to parallel sync since intersection data (surfaceIndex_ etc.)
229  // already guaranteed to be synced...
230  // However:
231  // - owncc and neicc are reversed on different procs so might pick
232  // up different regions reversed? No problem. Neighbour on one processor
233  // might not be owner on the other processor but the neighbour is
234  // not used when creating baffles from proc faces.
235  // - tolerances issues occasionally crop up.
236  syncTools::syncFaceList(mesh_, ownPatch, maxEqOp<label>());
237  syncTools::syncFaceList(mesh_, nbrPatch, maxEqOp<label>());
238 }
239 
240 
241 Foam::Map<Foam::labelPair> Foam::meshRefinement::getZoneBafflePatches
242 (
243  const bool allowBoundary,
244  const labelList& globalToMasterPatch,
245  const labelList& globalToSlavePatch
246 ) const
247 {
248  Map<labelPair> bafflePatch(mesh_.nFaces()/1000);
249 
250  const PtrList<surfaceZonesInfo>& surfZones = surfaces_.surfZones();
251  const meshFaceZones& fZones = mesh_.faceZones();
252 
253  forAll(surfZones, surfi)
254  {
255  const word& faceZoneName = surfZones[surfi].faceZoneName();
256 
257  if (faceZoneName.size())
258  {
259  // Get zone
260  const label zonei = fZones.findZoneID(faceZoneName);
261  const faceZone& fZone = fZones[zonei];
262 
263  // Get patch allocated for zone
264  const label globalRegioni = surfaces_.globalRegion(surfi, 0);
265  const labelPair zPatches
266  (
267  globalToMasterPatch[globalRegioni],
268  globalToSlavePatch[globalRegioni]
269  );
270 
271  Info<< "For zone " << fZone.name() << " found patches "
272  << mesh_.boundaryMesh()[zPatches[0]].name() << " and "
273  << mesh_.boundaryMesh()[zPatches[1]].name()
274  << endl;
275 
276  forAll(fZone, i)
277  {
278  const label facei = fZone[i];
279 
280  if (allowBoundary || mesh_.isInternalFace(facei))
281  {
282  labelPair patches = zPatches;
283  if (fZone.flipMap()[i])
284  {
285  patches = reverse(patches);
286  }
287 
288  if (!bafflePatch.insert(facei, patches))
289  {
291  << "Face " << facei
292  << " fc:" << mesh_.faceCentres()[facei]
293  << " in zone " << fZone.name()
294  << " is in multiple zones!"
295  << abort(FatalError);
296  }
297  }
298  }
299  }
300  }
301 
302  return bafflePatch;
303 }
304 
305 
307 (
308  const labelList& ownPatch,
309  const labelList& nbrPatch
310 )
311 {
312  if
313  (
314  ownPatch.size() != mesh_.nFaces()
315  || nbrPatch.size() != mesh_.nFaces()
316  )
317  {
319  << "Illegal size :"
320  << " ownPatch:" << ownPatch.size()
321  << " nbrPatch:" << nbrPatch.size()
322  << ". Should be number of faces:" << mesh_.nFaces()
323  << abort(FatalError);
324  }
325 
326  if (debug)
327  {
328  labelList syncedOwnPatch(ownPatch);
329  syncTools::syncFaceList(mesh_, syncedOwnPatch, maxEqOp<label>());
330  labelList syncedNeiPatch(nbrPatch);
331  syncTools::syncFaceList(mesh_, syncedNeiPatch, maxEqOp<label>());
332 
333  forAll(syncedOwnPatch, facei)
334  {
335  if
336  (
337  (ownPatch[facei] == -1 && syncedOwnPatch[facei] != -1)
338  || (nbrPatch[facei] == -1 && syncedNeiPatch[facei] != -1)
339  )
340  {
342  << "Non synchronised at face:" << facei
343  << " on patch:" << mesh_.boundaryMesh().whichPatch(facei)
344  << " fc:" << mesh_.faceCentres()[facei] << endl
345  << "ownPatch:" << ownPatch[facei]
346  << " syncedOwnPatch:" << syncedOwnPatch[facei]
347  << " nbrPatch:" << nbrPatch[facei]
348  << " syncedNeiPatch:" << syncedNeiPatch[facei]
349  << abort(FatalError);
350  }
351  }
352  }
353 
354  // Topochange container
355  polyTopoChange meshMod(mesh_);
356 
357  label nBaffles = 0;
358 
359  forAll(ownPatch, facei)
360  {
361  if (ownPatch[facei] != -1)
362  {
363  // Create baffle or repatch face. Return label of inserted baffle
364  // face.
365  createBaffle
366  (
367  facei,
368  ownPatch[facei], // owner side patch
369  nbrPatch[facei], // neighbour side patch
370  meshMod
371  );
372  nBaffles++;
373  }
374  }
375 
376  mesh_.clearOut();
377 
378  // Change the mesh (no inflation, parallel sync)
379  autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(mesh_, false, true);
380 
381  // Update fields
382  mesh_.topoChange(map);
383 
384  // Move mesh if in inflation mode
385  if (map().hasMotionPoints())
386  {
387  mesh_.movePoints(map().preMotionPoints());
388  }
389  else
390  {
391  // Delete mesh volumes.
392  mesh_.clearOut();
393  }
394 
395 
396  // Reset the instance for if in overwrite mode
397  mesh_.setInstance(timeName());
398 
399  //- Redo the intersections on the newly create baffle faces. Note that
400  // this changes also the cell centre positions.
401  faceSet baffledFacesSet(mesh_, "baffledFacesSet", 2*nBaffles);
402 
403  const labelList& reverseFaceMap = map().reverseFaceMap();
404  const labelList& faceMap = map().faceMap();
405 
406  // Pick up owner side of baffle
407  forAll(ownPatch, oldFacei)
408  {
409  const label facei = reverseFaceMap[oldFacei];
410 
411  if (ownPatch[oldFacei] != -1 && facei >= 0)
412  {
413  const cell& ownFaces = mesh_.cells()[mesh_.faceOwner()[facei]];
414 
415  forAll(ownFaces, i)
416  {
417  baffledFacesSet.insert(ownFaces[i]);
418  }
419  }
420  }
421  // Pick up neighbour side of baffle (added faces)
422  forAll(faceMap, facei)
423  {
424  const label oldFacei = faceMap[facei];
425 
426  if (oldFacei >= 0 && reverseFaceMap[oldFacei] != facei)
427  {
428  const cell& ownFaces = mesh_.cells()[mesh_.faceOwner()[facei]];
429 
430  forAll(ownFaces, i)
431  {
432  baffledFacesSet.insert(ownFaces[i]);
433  }
434  }
435  }
436  baffledFacesSet.sync(mesh_);
437 
438  topoChange(map, baffledFacesSet.toc());
439 
440  return map;
441 }
442 
443 
445 {
446  const meshFaceZones& fZones = mesh_.faceZones();
447 
448  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
449 
450  forAll(pbm, patchi)
451  {
452  const polyPatch& pp = pbm[patchi];
453 
454  if (isA<processorPolyPatch>(pp))
455  {
456  forAll(pp, i)
457  {
458  const label facei = pp.start() + i;
459  const label zonei = fZones.whichZone(facei);
460 
461  if (zonei != -1)
462  {
464  << "face:" << facei << " on patch " << pp.name()
465  << " is in zone " << fZones[zonei].name()
466  << exit(FatalError);
467  }
468  }
469  }
470  }
471 }
472 
473 
475 (
476  const labelList& globalToMasterPatch,
477  const labelList& globalToSlavePatch,
478  List<labelPair>& baffles
479 )
480 {
481  const labelList zonedSurfaces
482  (
484  );
485 
487 
488  // No need to sync; all processors will have all same zonedSurfaces.
489  if (zonedSurfaces.size())
490  {
491  // Split internal faces on interface surfaces
492  Info<< "Converting zoned faces into baffles ..." << endl;
493 
494  // Get faces (internal only) to be baffled. Map from face to patch
495  // label.
496  Map<labelPair> faceToPatch
497  (
498  getZoneBafflePatches
499  (
500  false,
501  globalToMasterPatch,
502  globalToSlavePatch
503  )
504  );
505 
506  const label nZoneFaces =
507  returnReduce(faceToPatch.size(), sumOp<label>());
508 
509  if (nZoneFaces > 0)
510  {
511  // Convert into labelLists
512  labelList ownPatch(mesh_.nFaces(), -1);
513  labelList nbrPatch(mesh_.nFaces(), -1);
514  forAllConstIter(Map<labelPair>, faceToPatch, iter)
515  {
516  ownPatch[iter.key()] = iter().first();
517  nbrPatch[iter.key()] = iter().second();
518  }
519 
520  // Create baffles. both sides same patch.
521  map = createBaffles(ownPatch, nbrPatch);
522 
523  // Get pairs of faces created.
524  // Just loop over faceMap and store baffle if we encounter a slave
525  // face.
526 
527  baffles.setSize(faceToPatch.size());
528  label baffleI = 0;
529 
530  const labelList& faceMap = map().faceMap();
531  const labelList& reverseFaceMap = map().reverseFaceMap();
532 
533  forAll(faceMap, facei)
534  {
535  const label oldFacei = faceMap[facei];
536 
537  // Does face originate from face-to-patch
538  Map<labelPair>::const_iterator iter = faceToPatch.find
539  (
540  oldFacei
541  );
542 
543  if (iter != faceToPatch.end())
544  {
545  const label masterFacei = reverseFaceMap[oldFacei];
546  if (facei != masterFacei)
547  {
548  baffles[baffleI++] = labelPair(masterFacei, facei);
549  }
550  }
551  }
552 
553  if (baffleI != faceToPatch.size())
554  {
556  << "Had " << faceToPatch.size() << " patches to create "
557  << " but encountered " << baffleI
558  << " slave faces originating from patcheable faces."
559  << abort(FatalError);
560  }
561 
562  if (debug&MESH)
563  {
564  const_cast<Time&>(mesh_.time())++;
565  Pout<< "Writing zone-baffled mesh to time " << timeName()
566  << endl;
567  write
568  (
569  debugType(debug),
571  mesh_.time().path()/"baffles"
572  );
573  }
574  }
575  Info<< "Created " << nZoneFaces << " baffles in = "
576  << mesh_.time().cpuTimeIncrement() << " s\n" << nl << endl;
577  }
578  return map;
579 }
580 
581 
582 Foam::List<Foam::labelPair> Foam::meshRefinement::freeStandingBaffles
583 (
584  const List<labelPair>& couples,
585  const scalar planarAngle
586 ) const
587 {
588  // Done by counting the number of baffles faces per mesh edge. If edge
589  // has 2 boundary faces and both are baffle faces it is the edge of a baffle
590  // region.
591 
592  // All duplicate faces on edge of the patch are to be merged.
593  // So we count for all edges of duplicate faces how many duplicate
594  // faces use them.
595  labelList nBafflesPerEdge(mesh_.nEdges(), 0);
596 
597 
598  // This algorithm is quite tricky. We don't want to use edgeFaces and
599  // also want it to run in parallel so it is now an algorithm over
600  // all (boundary) faces instead.
601  // We want to pick up any edges that are only used by the baffle
602  // or internal faces but not by any other boundary faces. So
603  // - increment count on an edge by 1 if it is used by any (uncoupled)
604  // boundary face.
605  // - increment count on an edge by 1000000 if it is used by a baffle face
606  // - sum in parallel
607  //
608  // So now any edge that is used by baffle faces only will have the
609  // value 2*1000000+2*1.
610 
611 
612  const label baffleValue = 1000000;
613 
614 
615  // Count number of boundary faces per edge
616  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
617 
618  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
619 
620  forAll(patches, patchi)
621  {
622  const polyPatch& pp = patches[patchi];
623 
624  // Count number of boundary faces. Discard coupled boundary faces.
625  if (!pp.coupled())
626  {
627  label facei = pp.start();
628 
629  forAll(pp, i)
630  {
631  const labelList& fEdges = mesh_.faceEdges(facei);
632 
633  forAll(fEdges, fEdgei)
634  {
635  nBafflesPerEdge[fEdges[fEdgei]]++;
636  }
637  facei++;
638  }
639  }
640  }
641 
642 
643  DynamicList<label> fe0;
644  DynamicList<label> fe1;
645 
646 
647  // Count number of duplicate boundary faces per edge
648  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
649 
650  forAll(couples, i)
651  {
652  {
653  const label f0 = couples[i].first();
654  const labelList& fEdges0 = mesh_.faceEdges(f0, fe0);
655  forAll(fEdges0, fEdgei)
656  {
657  nBafflesPerEdge[fEdges0[fEdgei]] += baffleValue;
658  }
659  }
660 
661  {
662  const label f1 = couples[i].second();
663  const labelList& fEdges1 = mesh_.faceEdges(f1, fe1);
664  forAll(fEdges1, fEdgei)
665  {
666  nBafflesPerEdge[fEdges1[fEdgei]] += baffleValue;
667  }
668  }
669  }
670 
671  // Add nBaffles on shared edges
673  (
674  mesh_,
675  nBafflesPerEdge,
676  plusEqOp<label>(), // in-place add
677  label(0) // initial value
678  );
679 
680 
681  // Baffles which are not next to other boundaries and baffles will have
682  // nBafflesPerEdge value 2*baffleValue+2*1 (from 2 boundary faces which
683  // are both baffle faces)
684 
685  List<labelPair> filteredCouples(couples.size());
686  label filterI = 0;
687 
688  forAll(couples, i)
689  {
690  const labelPair& couple = couples[i];
691 
692  if
693  (
694  patches.whichPatch(couple.first())
695  == patches.whichPatch(couple.second())
696  )
697  {
698  const labelList& fEdges = mesh_.faceEdges(couple.first());
699 
700  forAll(fEdges, fEdgei)
701  {
702  const label edgei = fEdges[fEdgei];
703 
704  if (nBafflesPerEdge[edgei] == 2*baffleValue+2*1)
705  {
706  filteredCouples[filterI++] = couple;
707  break;
708  }
709  }
710  }
711  }
712  filteredCouples.setSize(filterI);
713 
714 
715  const label nFiltered =
716  returnReduce(filteredCouples.size(), sumOp<label>());
717 
718  Info<< "freeStandingBaffles : detected "
719  << nFiltered
720  << " free-standing baffles out of "
721  << returnReduce(couples.size(), sumOp<label>())
722  << nl << endl;
723 
724 
725  if (nFiltered > 0)
726  {
727  // Collect segments
728  // ~~~~~~~~~~~~~~~~
729 
730  pointField start(filteredCouples.size());
731  pointField end(filteredCouples.size());
732 
733  const pointField& cellCentres = mesh_.cellCentres();
734 
735  forAll(filteredCouples, i)
736  {
737  const labelPair& couple = filteredCouples[i];
738  start[i] = cellCentres[mesh_.faceOwner()[couple.first()]];
739  end[i] = cellCentres[mesh_.faceOwner()[couple.second()]];
740  }
741 
742  // Extend segments a bit
743  {
744  const vectorField smallVec(rootSmall*(end-start));
745  start -= smallVec;
746  end += smallVec;
747  }
748 
749 
750  // Do test for intersections
751  // ~~~~~~~~~~~~~~~~~~~~~~~~~
752  labelList surface1;
753  List<pointIndexHit> hit1;
754  labelList region1;
755  vectorField normal1;
756 
757  labelList surface2;
758  List<pointIndexHit> hit2;
759  labelList region2;
760  vectorField normal2;
761 
762  surfaces_.findNearestIntersection
763  (
764  identity(surfaces_.surfaces().size()),
765  start,
766  end,
767 
768  surface1,
769  hit1,
770  region1,
771  normal1,
772 
773  surface2,
774  hit2,
775  region2,
776  normal2
777  );
778 
779  const scalar planarAngleCos = Foam::cos(degToRad(planarAngle));
780 
781  label filterI = 0;
782  forAll(filteredCouples, i)
783  {
784  const labelPair& couple = filteredCouples[i];
785 
786  if
787  (
788  hit1[i].hit()
789  && hit2[i].hit()
790  && (
791  surface1[i] != surface2[i]
792  || hit1[i].index() != hit2[i].index()
793  )
794  )
795  {
796  if ((normal1[i ]& normal2[i]) > planarAngleCos)
797  {
798  // Both normals aligned
799  const vector n = end[i] - start[i];
800  const scalar magN = mag(n);
801  if (magN > vSmall)
802  {
803  filteredCouples[filterI++] = couple;
804  }
805  }
806  }
807  else if (hit1[i].hit() || hit2[i].hit())
808  {
809  // Single hit. Do not include in freestanding baffles.
810  }
811  }
812 
813  filteredCouples.setSize(filterI);
814 
815  Info<< "freeStandingBaffles : detected "
816  << returnReduce(filterI, sumOp<label>())
817  << " planar (within " << planarAngle
818  << " degrees) free-standing baffles out of "
819  << nFiltered
820  << nl << endl;
821  }
822 
823  return filteredCouples;
824 }
825 
826 
828 (
829  const List<labelPair>& couples
830 )
831 {
832  // Mesh change engine
833  polyTopoChange meshMod(mesh_);
834 
835  const faceList& faces = mesh_.faces();
836  const labelList& faceOwner = mesh_.faceOwner();
837  const meshFaceZones& faceZones = mesh_.faceZones();
838 
839  forAll(couples, i)
840  {
841  const label face0 = couples[i].first();
842  const label face1 = couples[i].second();
843 
844  // face1 < 0 signals a coupled face that has been converted to baffle.
845 
846  const label own0 = faceOwner[face0];
847  const label own1 = faceOwner[face1];
848 
849  if (face1 < 0 || own0 < own1)
850  {
851  // Use face0 as the new internal face.
852  const label zoneID = faceZones.whichZone(face0);
853  bool zoneFlip = false;
854 
855  if (zoneID >= 0)
856  {
857  const faceZone& fZone = faceZones[zoneID];
858  zoneFlip = fZone.flipMap()[fZone.whichFace(face0)];
859  }
860 
861  const label nei = (face1 < 0 ? -1 : own1);
862 
863  meshMod.setAction(polyRemoveFace(face1));
864  meshMod.setAction
865  (
867  (
868  faces[face0], // modified face
869  face0, // label of face being modified
870  own0, // owner
871  nei, // neighbour
872  false, // face flip
873  -1, // patch for face
874  false, // remove from zone
875  zoneID, // zone for face
876  zoneFlip // face flip in zone
877  )
878  );
879  }
880  else
881  {
882  // Use face1 as the new internal face.
883  const label zoneID = faceZones.whichZone(face1);
884  bool zoneFlip = false;
885 
886  if (zoneID >= 0)
887  {
888  const faceZone& fZone = faceZones[zoneID];
889  zoneFlip = fZone.flipMap()[fZone.whichFace(face1)];
890  }
891 
892  meshMod.setAction(polyRemoveFace(face0));
893  meshMod.setAction
894  (
896  (
897  faces[face1], // modified face
898  face1, // label of face being modified
899  own1, // owner
900  own0, // neighbour
901  false, // face flip
902  -1, // patch for face
903  false, // remove from zone
904  zoneID, // zone for face
905  zoneFlip // face flip in zone
906  )
907  );
908  }
909  }
910 
911  mesh_.clearOut();
912 
913  // Change the mesh (no inflation)
914  autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(mesh_, false, true);
915 
916  // Update fields
917  mesh_.topoChange(map);
918 
919  // Move mesh (since morphing does not do this)
920  if (map().hasMotionPoints())
921  {
922  mesh_.movePoints(map().preMotionPoints());
923  }
924  else
925  {
926  // Delete mesh volumes.
927  mesh_.clearOut();
928  }
929 
930  // Reset the instance for if in overwrite mode
931  mesh_.setInstance(timeName());
932 
933  // Update intersections. Recalculate intersections on merged faces since
934  // this seems to give problems? Note: should not be necessary since
935  // baffles preserve intersections from when they were created.
936  labelList newExposedFaces(2*couples.size());
937  label newI = 0;
938 
939  forAll(couples, i)
940  {
941  const label newFace0 = map().reverseFaceMap()[couples[i].first()];
942  if (newFace0 != -1)
943  {
944  newExposedFaces[newI++] = newFace0;
945  }
946 
947  const label newFace1 = map().reverseFaceMap()[couples[i].second()];
948  if (newFace1 != -1)
949  {
950  newExposedFaces[newI++] = newFace1;
951  }
952  }
953  newExposedFaces.setSize(newI);
954  topoChange(map, newExposedFaces);
955 
956  return map;
957 }
958 
959 
960 // Finds region per cell for cells inside closed named surfaces
961 void Foam::meshRefinement::findCellZoneGeometric
962 (
963  const pointField& neiCc,
964  const labelList& closedNamedSurfaces, // indices of closed surfaces
965  labelList& namedSurfaceIndex, // per face index of named surface
966  const labelList& surfaceToCellZone, // cell zone index per surface
967 
968  labelList& cellToZone
969 ) const
970 {
971  const pointField& cellCentres = mesh_.cellCentres();
972  const labelList& faceOwner = mesh_.faceOwner();
973  const labelList& faceNeighbour = mesh_.faceNeighbour();
974 
975  // Check if cell centre is inside
976  labelList insideSurfaces;
977  surfaces_.findInside
978  (
979  closedNamedSurfaces,
980  cellCentres,
981  insideSurfaces
982  );
983 
984  forAll(insideSurfaces, celli)
985  {
986  if (cellToZone[celli] == -2)
987  {
988  label surfi = insideSurfaces[celli];
989 
990  if (surfi != -1)
991  {
992  cellToZone[celli] = surfaceToCellZone[surfi];
993  }
994  }
995  }
996 
997 
998  // Some cells with cell centres close to surface might have
999  // had been put into wrong surface. Recheck with perturbed cell centre.
1000 
1001 
1002  // 1. Collect points
1003 
1004  // Count points to test.
1005  label nCandidates = 0;
1006  forAll(namedSurfaceIndex, facei)
1007  {
1008  const label surfi = namedSurfaceIndex[facei];
1009 
1010  if (surfi != -1)
1011  {
1012  if (mesh_.isInternalFace(facei))
1013  {
1014  nCandidates += 2;
1015  }
1016  else
1017  {
1018  nCandidates += 1;
1019  }
1020  }
1021  }
1022 
1023  // Collect points.
1024  pointField candidatePoints(nCandidates);
1025  nCandidates = 0;
1026  forAll(namedSurfaceIndex, facei)
1027  {
1028  const label surfi = namedSurfaceIndex[facei];
1029 
1030  if (surfi != -1)
1031  {
1032  const label own = faceOwner[facei];
1033  const point& ownCc = cellCentres[own];
1034 
1035  if (mesh_.isInternalFace(facei))
1036  {
1037  const label nei = faceNeighbour[facei];
1038  const point& neiCc = cellCentres[nei];
1039 
1040  // Perturbed cc
1041  const vector d = 1e-4*(neiCc - ownCc);
1042  candidatePoints[nCandidates++] = ownCc-d;
1043  candidatePoints[nCandidates++] = neiCc+d;
1044  }
1045  else
1046  {
1047  const point& neiFc = neiCc[facei - mesh_.nInternalFaces()];
1048 
1049  // Perturbed cc
1050  const vector d = 1e-4*(neiFc - ownCc);
1051  candidatePoints[nCandidates++] = ownCc-d;
1052  }
1053  }
1054  }
1055 
1056 
1057  // 2. Test points for inside
1058 
1059  surfaces_.findInside
1060  (
1061  closedNamedSurfaces,
1062  candidatePoints,
1063  insideSurfaces
1064  );
1065 
1066 
1067  // 3. Update zone information
1068 
1069  nCandidates = 0;
1070  forAll(namedSurfaceIndex, facei)
1071  {
1072  const label surfi = namedSurfaceIndex[facei];
1073 
1074  if (surfi != -1)
1075  {
1076  const label own = faceOwner[facei];
1077 
1078  if (mesh_.isInternalFace(facei))
1079  {
1080  const label ownSurfI = insideSurfaces[nCandidates++];
1081  if (ownSurfI != -1)
1082  {
1083  cellToZone[own] = surfaceToCellZone[ownSurfI];
1084  }
1085 
1086  const label neiSurfI = insideSurfaces[nCandidates++];
1087  if (neiSurfI != -1)
1088  {
1089  label nei = faceNeighbour[facei];
1090 
1091  cellToZone[nei] = surfaceToCellZone[neiSurfI];
1092  }
1093  }
1094  else
1095  {
1096  const label ownSurfI = insideSurfaces[nCandidates++];
1097  if (ownSurfI != -1)
1098  {
1099  cellToZone[own] = surfaceToCellZone[ownSurfI];
1100  }
1101  }
1102  }
1103  }
1104 
1105 
1106  // Adapt the namedSurfaceIndex
1107  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
1108  // for if any cells were not completely covered.
1109 
1110  for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
1111  {
1112  const label ownZone = cellToZone[mesh_.faceOwner()[facei]];
1113  const label neiZone = cellToZone[mesh_.faceNeighbour()[facei]];
1114 
1115  if (namedSurfaceIndex[facei] == -1 && (ownZone != neiZone))
1116  {
1117  // Give face the zone of max cell zone
1118  namedSurfaceIndex[facei] = findIndex
1119  (
1120  surfaceToCellZone,
1121  max(ownZone, neiZone)
1122  );
1123  }
1124  }
1125 
1126  labelList neiCellZone(mesh_.nFaces() - mesh_.nInternalFaces());
1127  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
1128 
1129  forAll(patches, patchi)
1130  {
1131  const polyPatch& pp = patches[patchi];
1132 
1133  if (pp.coupled())
1134  {
1135  forAll(pp, i)
1136  {
1137  const label facei = pp.start() + i;
1138  const label ownZone = cellToZone[mesh_.faceOwner()[facei]];
1139  neiCellZone[facei-mesh_.nInternalFaces()] = ownZone;
1140  }
1141  }
1142  }
1143  syncTools::swapBoundaryFaceList(mesh_, neiCellZone);
1144 
1145  forAll(patches, patchi)
1146  {
1147  const polyPatch& pp = patches[patchi];
1148 
1149  if (pp.coupled())
1150  {
1151  forAll(pp, i)
1152  {
1153  const label facei = pp.start() + i;
1154  const label ownZone = cellToZone[mesh_.faceOwner()[facei]];
1155  const label neiZone = neiCellZone[facei-mesh_.nInternalFaces()];
1156 
1157  if (namedSurfaceIndex[facei] == -1 && (ownZone != neiZone))
1158  {
1159  // Give face the max cell zone
1160  namedSurfaceIndex[facei] = findIndex
1161  (
1162  surfaceToCellZone,
1163  max(ownZone, neiZone)
1164  );
1165  }
1166  }
1167  }
1168  }
1169 
1170  // Sync
1171  syncTools::syncFaceList(mesh_, namedSurfaceIndex, maxEqOp<label>());
1172 }
1173 
1174 
1175 void Foam::meshRefinement::findCellZoneInsideWalk
1176 (
1177  const labelList& locationSurfaces, // indices of surfaces with inside point
1178  const labelList& namedSurfaceIndex, // per face index of named surface
1179  const labelList& surfaceToCellZone, // cell zone index per surface
1180 
1181  labelList& cellToZone
1182 ) const
1183 {
1184  // Analyse regions. Reuse regionsplit
1185  boolList blockedFace(mesh_.nFaces());
1186  // selectSeparatedCoupledFaces(blockedFace);
1187 
1188  forAll(namedSurfaceIndex, facei)
1189  {
1190  if (namedSurfaceIndex[facei] == -1)
1191  {
1192  blockedFace[facei] = false;
1193  }
1194  else
1195  {
1196  blockedFace[facei] = true;
1197  }
1198  }
1199  // No need to sync since namedSurfaceIndex already is synced
1200 
1201  // Set region per cell based on walking
1202  regionSplit cellRegion(mesh_, blockedFace);
1203  blockedFace.clear();
1204 
1205 
1206  // Force calculation of face decomposition (used in findCell)
1207  (void)mesh_.tetBasePtIs();
1208 
1209  const PtrList<surfaceZonesInfo>& surfZones = surfaces_.surfZones();
1210 
1211  // For all locationSurface find the cell
1212  forAll(locationSurfaces, i)
1213  {
1214  const label surfi = locationSurfaces[i];
1215 
1216  const point& insidePoint = surfZones[surfi].zoneInsidePoint();
1217 
1218  Info<< "For surface " << surfaces_.names()[surfi]
1219  << " finding inside point " << insidePoint
1220  << endl;
1221 
1222  // Find the region containing the insidePoint
1223  label regioninMeshi = findRegion
1224  (
1225  mesh_,
1226  cellRegion,
1227  mergeDistance_*vector::one,
1228  insidePoint
1229  );
1230 
1231  Info<< "For surface " << surfaces_.names()[surfi]
1232  << " found point " << insidePoint
1233  << " in global region " << regioninMeshi
1234  << " out of " << cellRegion.nRegions() << " regions." << endl;
1235 
1236  if (regioninMeshi == -1)
1237  {
1239  << "Point " << insidePoint
1240  << " is not inside the mesh." << nl
1241  << "Bounding box of the mesh:" << mesh_.bounds()
1242  << exit(FatalError);
1243  }
1244 
1245  // Set all cells with this region
1246  forAll(cellRegion, celli)
1247  {
1248  if (cellRegion[celli] == regioninMeshi)
1249  {
1250  if (cellToZone[celli] == -2)
1251  {
1252  cellToZone[celli] = surfaceToCellZone[surfi];
1253  }
1254  else if (cellToZone[celli] != surfaceToCellZone[surfi])
1255  {
1257  << "Cell " << celli
1258  << " at " << mesh_.cellCentres()[celli]
1259  << " is inside surface " << surfaces_.names()[surfi]
1260  << " but already marked as being in zone "
1261  << cellToZone[celli] << endl
1262  << "This can happen if your surfaces are not"
1263  << " (sufficiently) closed."
1264  << endl;
1265  }
1266  }
1267  }
1268  }
1269 }
1270 
1271 
1272 bool Foam::meshRefinement::calcRegionToZone
1273 (
1274  const label surfZoneI,
1275  const label ownRegion,
1276  const label neiRegion,
1277 
1278  labelList& regionToCellZone
1279 ) const
1280 {
1281  bool changed = false;
1282 
1283  // Check whether in between different regions
1284  if (ownRegion != neiRegion)
1285  {
1286  // Jump. Change one of the sides to my type.
1287 
1288  // 1. Interface between my type and unset region.
1289  // Set region to regioninMesh
1290 
1291  if (regionToCellZone[ownRegion] == -2)
1292  {
1293  if (regionToCellZone[neiRegion] == surfZoneI)
1294  {
1295  // Face between unset and my region. Put unset
1296  // region into regioninMesh
1297  regionToCellZone[ownRegion] = -1;
1298  changed = true;
1299  }
1300  else if (regionToCellZone[neiRegion] != -2)
1301  {
1302  // Face between unset and other region.
1303  // Put unset region into my region
1304  regionToCellZone[ownRegion] = surfZoneI;
1305  changed = true;
1306  }
1307  }
1308  else if (regionToCellZone[neiRegion] == -2)
1309  {
1310  if (regionToCellZone[ownRegion] == surfZoneI)
1311  {
1312  // Face between unset and my region. Put unset
1313  // region into regioninMesh
1314  regionToCellZone[neiRegion] = -1;
1315  changed = true;
1316  }
1317  else if (regionToCellZone[ownRegion] != -2)
1318  {
1319  // Face between unset and other region.
1320  // Put unset region into my region
1321  regionToCellZone[neiRegion] = surfZoneI;
1322  changed = true;
1323  }
1324  }
1325  }
1326  return changed;
1327 }
1328 
1329 
1330 void Foam::meshRefinement::findCellZoneTopo
1331 (
1332  const List<point>& insidePoints,
1333  const labelList& namedSurfaceIndex,
1334  const labelList& surfaceToCellZone,
1335  labelList& cellToZone
1336 ) const
1337 {
1338  // Assumes:
1339  // - region containing insidePoints does not go into a cellZone
1340  // - all other regions can be found by crossing faces marked in
1341  // namedSurfaceIndex.
1342 
1343  // Analyse regions. Reuse regionsplit
1344  boolList blockedFace(mesh_.nFaces());
1345 
1346  forAll(namedSurfaceIndex, facei)
1347  {
1348  if (namedSurfaceIndex[facei] == -1)
1349  {
1350  blockedFace[facei] = false;
1351  }
1352  else
1353  {
1354  blockedFace[facei] = true;
1355  }
1356  }
1357  // No need to sync since namedSurfaceIndex already is synced
1358 
1359  // Set region per cell based on walking
1360  regionSplit cellRegion(mesh_, blockedFace);
1361  blockedFace.clear();
1362 
1363  // Per mesh region the zone the cell should be put in.
1364  // -2 : not analysed yet
1365  // -1 : insidePoint region. Not put into any cellzone.
1366  // >= 0 : index of cellZone
1367  labelList regionToCellZone(cellRegion.nRegions(), -2);
1368 
1369  // See which cells already are set in the cellToZone (from geometric
1370  // searching) and use these to take over their zones.
1371  // Note: could be improved to count number of cells per region.
1372  forAll(cellToZone, celli)
1373  {
1374  if (cellToZone[celli] != -2)
1375  {
1376  regionToCellZone[cellRegion[celli]] = cellToZone[celli];
1377  }
1378  }
1379 
1380 
1381  // Find the regions containing the insidePoints
1382  forAll(insidePoints, i)
1383  {
1384  const label regioninMeshi = findRegion
1385  (
1386  mesh_,
1387  cellRegion,
1388  mergeDistance_*vector::one,
1389  insidePoints[i]
1390  );
1391 
1392  Info<< "Found point " << insidePoints[i]
1393  << " in global region " << regioninMeshi
1394  << " out of " << cellRegion.nRegions() << " regions." << endl;
1395 
1396  if (regioninMeshi == -1)
1397  {
1399  << "Point " << insidePoints[i]
1400  << " is not inside the mesh." << nl
1401  << "Bounding box of the mesh:" << mesh_.bounds()
1402  << exit(FatalError);
1403  }
1404 
1405  // Mark default region with zone -1.
1406  if (regionToCellZone[regioninMeshi] == -2)
1407  {
1408  regionToCellZone[regioninMeshi] = -1;
1409  }
1410  }
1411 
1412 
1413  // Find correspondence between cell zone and surface
1414  // by changing cell zone every time we cross a surface.
1415  while (true)
1416  {
1417  // Synchronise regionToCellZone.
1418  // Note:
1419  // - region numbers are identical on all processors
1420  // - regioninMesh is identical ,,
1421  // - cellZones are identical ,,
1422  // This done at top of loop to account for geometric matching
1423  // not being synchronised.
1424  Pstream::listCombineGather(regionToCellZone, maxEqOp<label>());
1425  Pstream::listCombineScatter(regionToCellZone);
1426 
1427 
1428  bool changed = false;
1429 
1430  // Internal faces
1431 
1432  for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
1433  {
1434  const label surfi = namedSurfaceIndex[facei];
1435 
1436  // Connected even if no cellZone defined for surface
1437  if (surfi != -1)
1438  {
1439  // Calculate region to zone from cellRegions on either side
1440  // of internal face.
1441  const bool changedCell = calcRegionToZone
1442  (
1443  surfaceToCellZone[surfi],
1444  cellRegion[mesh_.faceOwner()[facei]],
1445  cellRegion[mesh_.faceNeighbour()[facei]],
1446  regionToCellZone
1447  );
1448 
1449  changed = changed | changedCell;
1450  }
1451  }
1452 
1453  // Boundary faces
1454 
1455  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
1456 
1457  // Get coupled neighbour cellRegion
1458  labelList neiCellRegion(mesh_.nFaces() - mesh_.nInternalFaces());
1459  forAll(patches, patchi)
1460  {
1461  const polyPatch& pp = patches[patchi];
1462 
1463  if (pp.coupled())
1464  {
1465  forAll(pp, i)
1466  {
1467  const label facei = pp.start() + i;
1468  neiCellRegion[facei-mesh_.nInternalFaces()] =
1469  cellRegion[mesh_.faceOwner()[facei]];
1470  }
1471  }
1472  }
1473  syncTools::swapBoundaryFaceList(mesh_, neiCellRegion);
1474 
1475  // Calculate region to zone from cellRegions on either side of coupled
1476  // face.
1477  forAll(patches, patchi)
1478  {
1479  const polyPatch& pp = patches[patchi];
1480 
1481  if (pp.coupled())
1482  {
1483  forAll(pp, i)
1484  {
1485  const label facei = pp.start() + i;
1486  const label surfi = namedSurfaceIndex[facei];
1487 
1488  // Connected even if no cellZone defined for surface
1489  if (surfi != -1)
1490  {
1491  const bool changedCell = calcRegionToZone
1492  (
1493  surfaceToCellZone[surfi],
1494  cellRegion[mesh_.faceOwner()[facei]],
1495  neiCellRegion[facei-mesh_.nInternalFaces()],
1496  regionToCellZone
1497  );
1498 
1499  changed = changed | changedCell;
1500  }
1501  }
1502  }
1503  }
1504 
1505  if (!returnReduce(changed, orOp<bool>()))
1506  {
1507  break;
1508  }
1509  }
1510 
1511 
1512  forAll(regionToCellZone, regioni)
1513  {
1514  const label zonei = regionToCellZone[regioni];
1515 
1516  if (zonei == -2)
1517  {
1519  << "For region " << regioni << " haven't set cell zone."
1520  << exit(FatalError);
1521  }
1522  }
1523 
1524  if (debug)
1525  {
1526  forAll(regionToCellZone, regioni)
1527  {
1528  Pout<< "Region " << regioni
1529  << " becomes cellZone:" << regionToCellZone[regioni]
1530  << endl;
1531  }
1532  }
1533 
1534  // Rework into cellToZone
1535  forAll(cellToZone, celli)
1536  {
1537  cellToZone[celli] = regionToCellZone[cellRegion[celli]];
1538  }
1539 }
1540 
1541 
1542 void Foam::meshRefinement::makeConsistentFaceIndex
1543 (
1544  const labelList& cellToZone,
1545  labelList& namedSurfaceIndex
1546 ) const
1547 {
1548  const labelList& faceOwner = mesh_.faceOwner();
1549  const labelList& faceNeighbour = mesh_.faceNeighbour();
1550 
1551  for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
1552  {
1553  const label ownZone = cellToZone[faceOwner[facei]];
1554  const label neiZone = cellToZone[faceNeighbour[facei]];
1555 
1556  if (ownZone == neiZone && namedSurfaceIndex[facei] != -1)
1557  {
1558  namedSurfaceIndex[facei] = -1;
1559  }
1560  else if (ownZone != neiZone && namedSurfaceIndex[facei] == -1)
1561  {
1563  << "Different cell zones on either side of face " << facei
1564  << " at " << mesh_.faceCentres()[facei]
1565  << " but face not marked with a surface."
1566  << abort(FatalError);
1567  }
1568  }
1569 
1570  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
1571 
1572  // Get coupled neighbour cellZone
1573  labelList neiCellZone(mesh_.nFaces() - mesh_.nInternalFaces());
1574  forAll(patches, patchi)
1575  {
1576  const polyPatch& pp = patches[patchi];
1577 
1578  if (pp.coupled())
1579  {
1580  forAll(pp, i)
1581  {
1582  const label facei = pp.start() + i;
1583  neiCellZone[facei-mesh_.nInternalFaces()] =
1584  cellToZone[mesh_.faceOwner()[facei]];
1585  }
1586  }
1587  }
1588  syncTools::swapBoundaryFaceList(mesh_, neiCellZone);
1589 
1590  // Use coupled cellZone to do check
1591  forAll(patches, patchi)
1592  {
1593  const polyPatch& pp = patches[patchi];
1594 
1595  if (pp.coupled())
1596  {
1597  forAll(pp, i)
1598  {
1599  const label facei = pp.start() + i;
1600 
1601  const label ownZone = cellToZone[faceOwner[facei]];
1602  const label neiZone = neiCellZone[facei-mesh_.nInternalFaces()];
1603 
1604  if (ownZone == neiZone && namedSurfaceIndex[facei] != -1)
1605  {
1606  namedSurfaceIndex[facei] = -1;
1607  }
1608  else if (ownZone != neiZone && namedSurfaceIndex[facei] == -1)
1609  {
1611  << "Different cell zones on either side of face "
1612  << facei << " at " << mesh_.faceCentres()[facei]
1613  << " but face not marked with a surface."
1614  << abort(FatalError);
1615  }
1616  }
1617  }
1618  else
1619  {
1620  // Unzonify boundary faces
1621  forAll(pp, i)
1622  {
1623  const label facei = pp.start() + i;
1624  namedSurfaceIndex[facei] = -1;
1625  }
1626  }
1627  }
1628 }
1629 
1630 
1631 void Foam::meshRefinement::handleSnapProblems
1632 (
1633  const snapParameters& snapParams,
1634  const bool useTopologicalSnapDetection,
1635  const bool removeEdgeConnectedCells,
1636  const scalarField& perpendicularAngle,
1637  const dictionary& motionDict,
1638  Time& runTime,
1639  const labelList& globalToMasterPatch,
1640  const labelList& globalToSlavePatch
1641 )
1642 {
1643  Info<< nl
1644  << "Introducing baffles to block off problem cells" << nl
1645  << "----------------------------------------------" << nl
1646  << endl;
1647 
1648  labelList facePatch;
1649  if (useTopologicalSnapDetection)
1650  {
1651  facePatch = markFacesOnProblemCells
1652  (
1653  motionDict,
1654  removeEdgeConnectedCells,
1655  perpendicularAngle,
1656  globalToMasterPatch
1657  );
1658  }
1659  else
1660  {
1661  facePatch = markFacesOnProblemCellsGeometric(snapParams, motionDict);
1662  }
1663  Info<< "Analyzed problem cells in = "
1664  << runTime.cpuTimeIncrement() << " s\n" << nl << endl;
1665 
1666  if (debug&MESH)
1667  {
1668  faceSet problemFaces(mesh_, "problemFaces", mesh_.nFaces()/100);
1669 
1670  forAll(facePatch, facei)
1671  {
1672  if (facePatch[facei] != -1)
1673  {
1674  problemFaces.insert(facei);
1675  }
1676  }
1677  problemFaces.instance() = timeName();
1678  Pout<< "Dumping " << problemFaces.size()
1679  << " problem faces to " << problemFaces.objectPath() << endl;
1680  problemFaces.write();
1681  }
1682 
1683  Info<< "Introducing baffles to delete problem cells." << nl << endl;
1684 
1685  if (debug)
1686  {
1687  runTime++;
1688  }
1689 
1690  // Create baffles with same owner and neighbour for now.
1691  createBaffles(facePatch, facePatch);
1692 
1693  if (debug)
1694  {
1695  // Debug:test all is still synced across proc patches
1696  checkData();
1697  }
1698  Info<< "Created baffles in = "
1699  << runTime.cpuTimeIncrement() << " s\n" << nl << endl;
1700 
1701  printMeshInfo(debug, "After introducing baffles");
1702 
1703  if (debug&MESH)
1704  {
1705  Pout<< "Writing extra baffled mesh to time "
1706  << timeName() << endl;
1707  write
1708  (
1709  debugType(debug),
1711  runTime.path()/"extraBaffles"
1712  );
1713  Pout<< "Dumped debug data in = "
1714  << runTime.cpuTimeIncrement() << " s\n" << nl << endl;
1715  }
1716 }
1717 
1718 
1719 Foam::labelList Foam::meshRefinement::freeStandingBaffleFaces
1720 (
1721  const labelList& faceToZone,
1722  const labelList& cellToZone,
1723  const labelList& neiCellZone
1724 ) const
1725 {
1726  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
1727  const labelList& faceOwner = mesh_.faceOwner();
1728  const labelList& faceNeighbour = mesh_.faceNeighbour();
1729 
1730 
1731  // We want to pick up the faces to orient. These faces come in
1732  // two variants:
1733  // - faces originating from stand-alone faceZones
1734  // (these will most likely have no cellZone on either side so
1735  // ownZone and neiZone both -1)
1736  // - sticky-up faces originating from a 'bulge' in a outside of
1737  // a cellZone. These will have the same cellZone on either side.
1738  // How to orient these is not really clearly defined so do them
1739  // same as stand-alone faceZone faces for now. (Normally these will
1740  // already have been removed by the 'allowFreeStandingZoneFaces=false'
1741  // default setting)
1742 
1743  // Note that argument neiCellZone will have -1 on uncoupled boundaries.
1744 
1745  DynamicList<label> faceLabels(mesh_.nFaces()/100);
1746 
1747  for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
1748  {
1749  if (faceToZone[facei] != -1)
1750  {
1751  // Free standing baffle?
1752  const label ownZone = cellToZone[faceOwner[facei]];
1753  const label neiZone = cellToZone[faceNeighbour[facei]];
1754  if (ownZone == neiZone)
1755  {
1756  faceLabels.append(facei);
1757  }
1758  }
1759  }
1760  forAll(patches, patchi)
1761  {
1762  const polyPatch& pp = patches[patchi];
1763 
1764  forAll(pp, i)
1765  {
1766  const label facei = pp.start() + i;
1767  if (faceToZone[facei] != -1)
1768  {
1769  // Free standing baffle?
1770  const label ownZone = cellToZone[faceOwner[facei]];
1771  const label neiZone =
1772  neiCellZone[facei - mesh_.nInternalFaces()];
1773 
1774  if (ownZone == neiZone)
1775  {
1776  faceLabels.append(facei);
1777  }
1778  }
1779  }
1780  }
1781  return faceLabels.shrink();
1782 }
1783 
1784 
1785 void Foam::meshRefinement::calcPatchNumMasterFaces
1786 (
1787  const PackedBoolList& isMasterFace,
1788  const indirectPrimitivePatch& patch,
1789  labelList& nMasterFacesPerEdge
1790 ) const
1791 {
1792  // Number of (master)faces per edge
1793  nMasterFacesPerEdge.setSize(patch.nEdges());
1794  nMasterFacesPerEdge = 0;
1795 
1796  forAll(patch.addressing(), facei)
1797  {
1798  const label meshFacei = patch.addressing()[facei];
1799 
1800  if (isMasterFace[meshFacei])
1801  {
1802  const labelList& fEdges = patch.faceEdges()[facei];
1803  forAll(fEdges, fEdgei)
1804  {
1805  nMasterFacesPerEdge[fEdges[fEdgei]]++;
1806  }
1807  }
1808  }
1809 
1811  (
1812  mesh_,
1813  patch.meshEdges(mesh_.edges(), mesh_.pointEdges()),
1814  nMasterFacesPerEdge,
1815  plusEqOp<label>(),
1816  label(0)
1817  );
1818 }
1819 
1820 
1821 Foam::label Foam::meshRefinement::markPatchZones
1822 (
1823  const indirectPrimitivePatch& patch,
1824  const labelList& nMasterFacesPerEdge,
1825  labelList& faceToZone
1826 ) const
1827 {
1828  List<patchEdgeFaceRegion> allEdgeInfo(patch.nEdges());
1829  List<patchEdgeFaceRegion> allFaceInfo(patch.size());
1830 
1831 
1832  // Protect all non-manifold edges
1833  {
1834  label nProtected = 0;
1835 
1836  forAll(nMasterFacesPerEdge, edgei)
1837  {
1838  if (nMasterFacesPerEdge[edgei] > 2)
1839  {
1840  allEdgeInfo[edgei] = -2;
1841  nProtected++;
1842  }
1843  }
1844  }
1845 
1846 
1847  // Hand out zones
1848 
1849  DynamicList<label> changedEdges;
1851 
1852  const scalar tol = PatchEdgeFaceWave
1853  <
1856  >::propagationTol();
1857 
1858  int dummyTrackData;
1859 
1860  const globalIndex globalFaces(patch.size());
1861 
1862  label facei = 0;
1863 
1864  label currentZoneI = 0;
1865 
1866  while (true)
1867  {
1868  // Pick an unset face
1869  label globalSeed = labelMax;
1870  for (; facei < allFaceInfo.size(); facei++)
1871  {
1872  if (!allFaceInfo[facei].valid(dummyTrackData))
1873  {
1874  globalSeed = globalFaces.toGlobal(facei);
1875  break;
1876  }
1877  }
1878 
1879  reduce(globalSeed, minOp<label>());
1880 
1881  if (globalSeed == labelMax)
1882  {
1883  break;
1884  }
1885 
1886  const label proci = globalFaces.whichProcID(globalSeed);
1887  const label seedFacei = globalFaces.toLocal(proci, globalSeed);
1888 
1889  // Info<< "Seeding zone " << currentZoneI
1890  // << " from processor " << proci << " face " << seedFacei
1891  // << endl;
1892 
1893  if (proci == Pstream::myProcNo())
1894  {
1895  patchEdgeFaceRegion& faceInfo = allFaceInfo[seedFacei];
1896 
1897  // Set face
1898  faceInfo = currentZoneI;
1899 
1900  // .. and seed its edges
1901  const labelList& fEdges = patch.faceEdges()[seedFacei];
1902  forAll(fEdges, fEdgei)
1903  {
1904  const label edgei = fEdges[fEdgei];
1905 
1906  patchEdgeFaceRegion& edgeinfo = allEdgeInfo[edgei];
1907 
1908  if
1909  (
1910  edgeinfo.updateEdge<int>
1911  (
1912  mesh_,
1913  patch,
1914  edgei,
1915  seedFacei,
1916  faceInfo,
1917  tol,
1918  dummyTrackData
1919  )
1920  )
1921  {
1922  changedEdges.append(edgei);
1923  changedInfo.append(edgeinfo);
1924  }
1925  }
1926  }
1927 
1928 
1929  if (returnReduce(changedEdges.size(), sumOp<label>()) == 0)
1930  {
1931  break;
1932  }
1933 
1934 
1935  // Walk
1937  <
1940  > calc
1941  (
1942  mesh_,
1943  patch,
1944  changedEdges,
1945  changedInfo,
1946  allEdgeInfo,
1947  allFaceInfo,
1948  returnReduce(patch.nEdges(), sumOp<label>())
1949  );
1950 
1951  currentZoneI++;
1952  }
1953 
1954 
1955  faceToZone.setSize(patch.size());
1956  forAll(allFaceInfo, facei)
1957  {
1958  if (!allFaceInfo[facei].valid(dummyTrackData))
1959  {
1961  << "Problem: unvisited face " << facei
1962  << " at " << patch.faceCentres()[facei]
1963  << exit(FatalError);
1964  }
1965  faceToZone[facei] = allFaceInfo[facei].region();
1966  }
1967 
1968  return currentZoneI;
1969 }
1970 
1971 
1972 void Foam::meshRefinement::consistentOrientation
1973 (
1974  const PackedBoolList& isMasterFace,
1975  const indirectPrimitivePatch& patch,
1976  const labelList& nMasterFacesPerEdge,
1977  const labelList& faceToZone,
1978  const Map<label>& zoneToOrientation,
1979  boolList& meshFlipMap
1980 ) const
1981 {
1982  const polyBoundaryMesh& bm = mesh_.boundaryMesh();
1983 
1984  // Data on all edges and faces
1985  List<patchFaceOrientation> allEdgeInfo(patch.nEdges());
1986  List<patchFaceOrientation> allFaceInfo(patch.size());
1987 
1988  // Make sure we don't walk through
1989  // - slaves of coupled faces
1990  // - non-manifold edges
1991  {
1992  label nProtected = 0;
1993 
1994  forAll(patch.addressing(), facei)
1995  {
1996  const label meshFacei = patch.addressing()[facei];
1997  const label patchi = bm.whichPatch(meshFacei);
1998 
1999  if
2000  (
2001  patchi != -1
2002  && bm[patchi].coupled()
2003  && !isMasterFace[meshFacei]
2004  )
2005  {
2006  // Slave side. Mark so doesn't get visited.
2007  allFaceInfo[facei] = orientedSurface::NOFLIP;
2008  nProtected++;
2009  }
2010  }
2011  }
2012  {
2013  label nProtected = 0;
2014 
2015  forAll(nMasterFacesPerEdge, edgei)
2016  {
2017  if (nMasterFacesPerEdge[edgei] > 2)
2018  {
2019  allEdgeInfo[edgei] = orientedSurface::NOFLIP;
2020  nProtected++;
2021  }
2022  }
2023 
2024  Info<< "Protected from visiting "
2025  << returnReduce(nProtected, sumOp<label>())
2026  << " non-manifold edges" << nl << endl;
2027  }
2028 
2029 
2030 
2031  DynamicList<label> changedEdges;
2033 
2034  const scalar tol = PatchEdgeFaceWave
2035  <
2038  >::propagationTol();
2039 
2040  int dummyTrackData;
2041 
2042  globalIndex globalFaces(patch.size());
2043 
2044  while (true)
2045  {
2046  // Pick an unset face
2047  label globalSeed = labelMax;
2048  forAll(allFaceInfo, facei)
2049  {
2050  if (allFaceInfo[facei] == orientedSurface::UNVISITED)
2051  {
2052  globalSeed = globalFaces.toGlobal(facei);
2053  break;
2054  }
2055  }
2056 
2057  reduce(globalSeed, minOp<label>());
2058 
2059  if (globalSeed == labelMax)
2060  {
2061  break;
2062  }
2063 
2064  const label proci = globalFaces.whichProcID(globalSeed);
2065  const label seedFacei = globalFaces.toLocal(proci, globalSeed);
2066 
2067  // Info<< "Seeding from processor " << proci << " face " << seedFacei
2068  // << endl;
2069 
2070  if (proci == Pstream::myProcNo())
2071  {
2072  // Determine orientation of seedFace
2073 
2074  patchFaceOrientation& faceInfo = allFaceInfo[seedFacei];
2075 
2076  // Start off with correct orientation
2077  faceInfo = orientedSurface::NOFLIP;
2078 
2079  if (zoneToOrientation[faceToZone[seedFacei]] < 0)
2080  {
2081  faceInfo.flip();
2082  }
2083 
2084 
2085  const labelList& fEdges = patch.faceEdges()[seedFacei];
2086  forAll(fEdges, fEdgei)
2087  {
2088  const label edgei = fEdges[fEdgei];
2089 
2090  patchFaceOrientation& edgeinfo = allEdgeInfo[edgei];
2091 
2092  if
2093  (
2094  edgeinfo.updateEdge<int>
2095  (
2096  mesh_,
2097  patch,
2098  edgei,
2099  seedFacei,
2100  faceInfo,
2101  tol,
2102  dummyTrackData
2103  )
2104  )
2105  {
2106  changedEdges.append(edgei);
2107  changedInfo.append(edgeinfo);
2108  }
2109  }
2110  }
2111 
2112 
2113  if (returnReduce(changedEdges.size(), sumOp<label>()) == 0)
2114  {
2115  break;
2116  }
2117 
2118 
2119 
2120  // Walk
2122  <
2125  > calc
2126  (
2127  mesh_,
2128  patch,
2129  changedEdges,
2130  changedInfo,
2131  allEdgeInfo,
2132  allFaceInfo,
2133  returnReduce(patch.nEdges(), sumOp<label>())
2134  );
2135  }
2136 
2137 
2138  // Push master zone info over to slave (since slave faces never visited)
2139  {
2140  labelList neiStatus
2141  (
2142  mesh_.nFaces() - mesh_.nInternalFaces(),
2144  );
2145 
2146  forAll(patch.addressing(), i)
2147  {
2148  const label meshFacei = patch.addressing()[i];
2149  if (!mesh_.isInternalFace(meshFacei))
2150  {
2151  neiStatus[meshFacei-mesh_.nInternalFaces()] =
2152  allFaceInfo[i].flipStatus();
2153  }
2154  }
2155  syncTools::swapBoundaryFaceList(mesh_, neiStatus);
2156 
2157  forAll(patch.addressing(), i)
2158  {
2159  const label meshFacei = patch.addressing()[i];
2160  const label patchi = bm.whichPatch(meshFacei);
2161 
2162  if
2163  (
2164  patchi != -1
2165  && bm[patchi].coupled()
2166  && !isMasterFace[meshFacei]
2167  )
2168  {
2169  // Slave side. Take flipped from neighbour
2170  label bFacei = meshFacei-mesh_.nInternalFaces();
2171 
2172  if (neiStatus[bFacei] == orientedSurface::NOFLIP)
2173  {
2174  allFaceInfo[i] = orientedSurface::FLIP;
2175  }
2176  else if (neiStatus[bFacei] == orientedSurface::FLIP)
2177  {
2178  allFaceInfo[i] = orientedSurface::NOFLIP;
2179  }
2180  else
2181  {
2183  << "Incorrect status for face " << meshFacei
2184  << abort(FatalError);
2185  }
2186  }
2187  }
2188  }
2189 
2190 
2191  // Convert to meshFlipMap and adapt faceZones
2192 
2193  meshFlipMap.setSize(mesh_.nFaces());
2194  meshFlipMap = false;
2195 
2196  forAll(allFaceInfo, facei)
2197  {
2198  label meshFacei = patch.addressing()[facei];
2199 
2200  if (allFaceInfo[facei] == orientedSurface::NOFLIP)
2201  {
2202  meshFlipMap[meshFacei] = false;
2203  }
2204  else if (allFaceInfo[facei] == orientedSurface::FLIP)
2205  {
2206  meshFlipMap[meshFacei] = true;
2207  }
2208  else
2209  {
2211  << "Problem : unvisited face " << facei
2212  << " centre:" << mesh_.faceCentres()[meshFacei]
2213  << abort(FatalError);
2214  }
2215  }
2216 }
2217 
2218 
2219 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
2220 
2223  const bool doHandleSnapProblems,
2224  const snapParameters& snapParams,
2225  const bool useTopologicalSnapDetection,
2226  const bool removeEdgeConnectedCells,
2227  const scalarField& perpendicularAngle,
2228  const bool mergeFreeStanding,
2229  const scalar planarAngle,
2230  const dictionary& motionDict,
2231  Time& runTime,
2232  const labelList& globalToMasterPatch,
2233  const labelList& globalToSlavePatch,
2234  const refinementParameters::cellSelectionPoints& selectionPoints
2235 )
2236 {
2237  // Introduce baffles
2238  // ~~~~~~~~~~~~~~~~~
2239 
2240  // Split the mesh along internal faces wherever there is a pierce between
2241  // two cell centres.
2242 
2243  Info<< "Introducing baffles for "
2245  << " faces that are intersected by the surface." << nl << endl;
2246 
2247  // Swap neighbouring cell centres and cell level
2248  labelList neiLevel(mesh_.nFaces() - mesh_.nInternalFaces());
2249  pointField neiCc(mesh_.nFaces() - mesh_.nInternalFaces());
2250  calcNeighbourData(neiLevel, neiCc);
2251 
2252  labelList ownPatch, nbrPatch;
2253  getBafflePatches
2254  (
2255  globalToMasterPatch,
2256  neiLevel,
2257  neiCc,
2258 
2259  ownPatch,
2260  nbrPatch
2261  );
2262 
2263  createBaffles(ownPatch, nbrPatch);
2264 
2265  if (debug)
2266  {
2267  // Debug:test all is still synced across proc patches
2268  checkData();
2269  }
2270 
2271  Info<< "Created baffles in = "
2272  << runTime.cpuTimeIncrement() << " s\n" << nl << endl;
2273 
2274  printMeshInfo(debug, "After introducing baffles");
2275 
2276  if (debug&MESH)
2277  {
2278  Pout<< "Writing baffled mesh to time " << timeName()
2279  << endl;
2280  write
2281  (
2282  debugType(debug),
2284  runTime.path()/"baffles"
2285  );
2286  Pout<< "Dumped debug data in = "
2287  << runTime.cpuTimeIncrement() << " s\n" << nl << endl;
2288  }
2289 
2290 
2291  // Introduce baffles to delete problem cells
2292  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2293 
2294  // Create some additional baffles where we want surface cells removed.
2295 
2296  if (doHandleSnapProblems)
2297  {
2298  handleSnapProblems
2299  (
2300  snapParams,
2301  useTopologicalSnapDetection,
2302  removeEdgeConnectedCells,
2303  perpendicularAngle,
2304  motionDict,
2305  runTime,
2306  globalToMasterPatch,
2307  globalToSlavePatch
2308  );
2309  }
2310 
2311 
2312  // Select part of mesh
2313  // ~~~~~~~~~~~~~~~~~~~
2314 
2315  Info<< nl
2316  << "Remove unreachable sections of mesh" << nl
2317  << "-----------------------------------" << nl
2318  << endl;
2319 
2320  if (debug)
2321  {
2322  runTime++;
2323  }
2324 
2325  splitMeshRegions(globalToMasterPatch, globalToSlavePatch, selectionPoints);
2326 
2327  if (debug)
2328  {
2329  // Debug:test all is still synced across proc patches
2330  checkData();
2331  }
2332  Info<< "Split mesh in = "
2333  << runTime.cpuTimeIncrement() << " s\n" << nl << endl;
2334 
2335  printMeshInfo(debug, "After subsetting");
2336 
2337  if (debug&MESH)
2338  {
2339  Pout<< "Writing subsetted mesh to time " << timeName()
2340  << endl;
2341  write
2342  (
2343  debugType(debug),
2345  runTime.path()/timeName()
2346  );
2347  Pout<< "Dumped debug data in = "
2348  << runTime.cpuTimeIncrement() << " s\n" << nl << endl;
2349  }
2350 
2351 
2352  // Merge baffles
2353  // ~~~~~~~~~~~~~
2354 
2355  if (mergeFreeStanding)
2356  {
2357  Info<< nl
2358  << "Merge free-standing baffles" << nl
2359  << "---------------------------" << nl
2360  << endl;
2361 
2362 
2363  // List of pairs of freestanding baffle faces.
2364  List<labelPair> couples
2365  (
2366  freeStandingBaffles // filter out freestanding baffles
2367  (
2369  planarAngle
2370  )
2371  );
2372 
2373  label nCouples = couples.size();
2374  reduce(nCouples, sumOp<label>());
2375 
2376  Info<< "Detected free-standing baffles : " << nCouples << endl;
2377 
2378  if (nCouples > 0)
2379  {
2380  // Actually merge baffles. Note: not exactly parallelised. Should
2381  // convert baffle faces into processor faces if they resulted
2382  // from them.
2383  mergeBaffles(couples);
2384 
2385  // Detect any problem cells resulting from merging of baffles
2386  // and delete them
2387  handleSnapProblems
2388  (
2389  snapParams,
2390  useTopologicalSnapDetection,
2391  removeEdgeConnectedCells,
2392  perpendicularAngle,
2393  motionDict,
2394  runTime,
2395  globalToMasterPatch,
2396  globalToSlavePatch
2397  );
2398 
2399  if (debug)
2400  {
2401  // Debug:test all is still synced across proc patches
2402  checkData();
2403  }
2404  }
2405  Info<< "Merged free-standing baffles in = "
2406  << runTime.cpuTimeIncrement() << " s\n" << nl << endl;
2407  }
2408 }
2409 
2410 
2413  const label nBufferLayers,
2414  const labelList& globalToMasterPatch,
2415  const labelList& globalToSlavePatch,
2416  const refinementParameters::cellSelectionPoints& selectionPoints
2417 )
2418 {
2419  // Determine patches to put intersections into
2420  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2421 
2422  // Swap neighbouring cell centres and cell level
2423  labelList neiLevel(mesh_.nFaces() - mesh_.nInternalFaces());
2424  pointField neiCc(mesh_.nFaces() - mesh_.nInternalFaces());
2425  calcNeighbourData(neiLevel, neiCc);
2426 
2427  labelList ownPatch;
2428  labelList nbrPatch;
2429  getBafflePatches
2430  (
2431  globalToMasterPatch,
2432  neiLevel,
2433  neiCc,
2434 
2435  ownPatch,
2436  nbrPatch
2437  );
2438 
2439  // Analyse regions. Reuse regionsplit
2440  boolList blockedFace(mesh_.nFaces(), false);
2441 
2442  forAll(ownPatch, facei)
2443  {
2444  if (ownPatch[facei] != -1 || nbrPatch[facei] != -1)
2445  {
2446  blockedFace[facei] = true;
2447  }
2448  }
2449  syncTools::syncFaceList(mesh_, blockedFace, orEqOp<bool>());
2450 
2451  // Set region per cell based on walking
2452  regionSplit cellRegion(mesh_, blockedFace);
2453  blockedFace.clear();
2454 
2455  // Find the regions containing any of the points in selectionPoints
2456  findRegions
2457  (
2458  mesh_,
2459  cellRegion,
2460  mergeDistance_*vector::one,
2461  selectionPoints
2462  );
2463 
2464 
2465  // Walk out nBufferlayers from region split
2466  // (modifies cellRegion, ownPatch)
2467  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2468  // Takes over face patch onto points and then back to faces and cells
2469  // (so cell-face-point walk)
2470 
2471  const labelList& faceOwner = mesh_.faceOwner();
2472  const labelList& faceNeighbour = mesh_.faceNeighbour();
2473 
2474  // Patch for exposed faces for lack of anything sensible.
2475  label defaultPatch = 0;
2476  if (globalToMasterPatch.size())
2477  {
2478  defaultPatch = globalToMasterPatch[0];
2479  }
2480 
2481  for (label i = 0; i < nBufferLayers; i++)
2482  {
2483  // 1. From cells (via faces) to points
2484 
2485  labelList pointBaffle(mesh_.nPoints(), -1);
2486 
2487  forAll(faceNeighbour, facei)
2488  {
2489  const face& f = mesh_.faces()[facei];
2490  const label ownRegion = cellRegion[faceOwner[facei]];
2491  const label neiRegion = cellRegion[faceNeighbour[facei]];
2492 
2493  if (ownRegion == -1 && neiRegion != -1)
2494  {
2495  // Note max(..) since possibly regionSplit might have split
2496  // off extra unreachable parts of mesh. Note: or can this only
2497  // happen for boundary faces?
2498  forAll(f, fp)
2499  {
2500  pointBaffle[f[fp]] = max(defaultPatch, ownPatch[facei]);
2501  }
2502  }
2503  else if (ownRegion != -1 && neiRegion == -1)
2504  {
2505  label newPatchi = nbrPatch[facei];
2506  if (newPatchi == -1)
2507  {
2508  newPatchi = max(defaultPatch, ownPatch[facei]);
2509  }
2510  forAll(f, fp)
2511  {
2512  pointBaffle[f[fp]] = newPatchi;
2513  }
2514  }
2515  }
2516 
2517  for
2518  (
2519  label facei = mesh_.nInternalFaces();
2520  facei < mesh_.nFaces();
2521  facei++
2522  )
2523  {
2524  const face& f = mesh_.faces()[facei];
2525  const label ownRegion = cellRegion[faceOwner[facei]];
2526 
2527  if (ownRegion == -1)
2528  {
2529  forAll(f, fp)
2530  {
2531  pointBaffle[f[fp]] = max(defaultPatch, ownPatch[facei]);
2532  }
2533  }
2534  }
2535 
2536  // Sync
2538  (
2539  mesh_,
2540  pointBaffle,
2541  maxEqOp<label>(),
2542  label(-1) // null value
2543  );
2544 
2545 
2546  // 2. From points back to faces
2547 
2548  const labelListList& pointFaces = mesh_.pointFaces();
2549 
2550  forAll(pointFaces, pointi)
2551  {
2552  if (pointBaffle[pointi] != -1)
2553  {
2554  const labelList& pFaces = pointFaces[pointi];
2555 
2556  forAll(pFaces, pFacei)
2557  {
2558  const label facei = pFaces[pFacei];
2559 
2560  if (ownPatch[facei] == -1)
2561  {
2562  ownPatch[facei] = pointBaffle[pointi];
2563  }
2564  }
2565  }
2566  }
2567  syncTools::syncFaceList(mesh_, ownPatch, maxEqOp<label>());
2568 
2569 
2570  // 3. From faces to cells (cellRegion) and back to faces (ownPatch)
2571 
2572  labelList newOwnPatch(ownPatch);
2573 
2574  forAll(ownPatch, facei)
2575  {
2576  if (ownPatch[facei] != -1)
2577  {
2578  const label own = faceOwner[facei];
2579 
2580  if (cellRegion[own] == -1)
2581  {
2582  cellRegion[own] = labelMax;
2583 
2584  const cell& ownFaces = mesh_.cells()[own];
2585  forAll(ownFaces, j)
2586  {
2587  if (ownPatch[ownFaces[j]] == -1)
2588  {
2589  newOwnPatch[ownFaces[j]] = ownPatch[facei];
2590  }
2591  }
2592  }
2593  if (mesh_.isInternalFace(facei))
2594  {
2595  const label nei = faceNeighbour[facei];
2596 
2597  if (cellRegion[nei] == -1)
2598  {
2599  cellRegion[nei] = labelMax;
2600 
2601  const cell& neiFaces = mesh_.cells()[nei];
2602  forAll(neiFaces, j)
2603  {
2604  if (ownPatch[neiFaces[j]] == -1)
2605  {
2606  newOwnPatch[neiFaces[j]] = ownPatch[facei];
2607  }
2608  }
2609  }
2610  }
2611  }
2612  }
2613 
2614  ownPatch.transfer(newOwnPatch);
2615 
2616  syncTools::syncFaceList(mesh_, ownPatch, maxEqOp<label>());
2617  }
2618 
2619 
2620 
2621  // Subset
2622  // ~~~~~~
2623 
2624  // Get cells to remove
2625  DynamicList<label> cellsToRemove(mesh_.nCells());
2626  forAll(cellRegion, celli)
2627  {
2628  if (cellRegion[celli] == -1)
2629  {
2630  cellsToRemove.append(celli);
2631  }
2632  }
2633  cellsToRemove.shrink();
2634 
2635  label nCellsInMesh = mesh_.nCells() - cellsToRemove.size();
2636  reduce(nCellsInMesh, sumOp<label>());
2637 
2638  Info<< "Selecting all cells in regions containing any of the points in "
2639  << selectionPoints.inside() << endl
2640  << "Selected: " << nCellsInMesh << " cells." << endl;
2641 
2642 
2643  // Remove cells
2644  removeCells cellRemover(mesh_);
2645 
2646  // Pick up patches for exposed faces
2647  labelList exposedFaces(cellRemover.getExposedFaces(cellsToRemove));
2648  labelList exposedPatches(exposedFaces.size());
2649 
2650  forAll(exposedFaces, i)
2651  {
2652  const label facei = exposedFaces[i];
2653 
2654  if (ownPatch[facei] != -1)
2655  {
2656  exposedPatches[i] = ownPatch[facei];
2657  }
2658  else
2659  {
2661  << "For exposed face " << facei
2662  << " fc:" << mesh_.faceCentres()[facei]
2663  << " found no patch." << endl
2664  << " Taking patch " << defaultPatch
2665  << " instead." << endl;
2666  exposedPatches[i] = defaultPatch;
2667  }
2668  }
2669 
2670  return doRemoveCells
2671  (
2672  cellsToRemove,
2673  exposedFaces,
2674  exposedPatches,
2675  cellRemover
2676  );
2677 }
2678 
2679 
2683  const localPointRegion& regionSide
2684 )
2685 {
2686  // Topochange container
2687  polyTopoChange meshMod(mesh_);
2688 
2689  const label nNonManifPoints = returnReduce
2690  (
2691  regionSide.meshPointMap().size(),
2692  sumOp<label>()
2693  );
2694 
2695  Info<< "dupNonManifoldPoints : Found : " << nNonManifPoints
2696  << " non-manifold points (out of "
2697  << mesh_.globalData().nTotalPoints()
2698  << ')' << endl;
2699 
2700  // Topo change engine
2701  duplicatePoints pointDuplicator(mesh_);
2702 
2703  // Insert changes into meshMod
2704  pointDuplicator.setRefinement(regionSide, meshMod);
2705 
2706  mesh_.clearOut();
2707 
2708  // Change the mesh (no inflation, parallel sync)
2709  autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(mesh_, false, true);
2710 
2711  // Update fields
2712  mesh_.topoChange(map);
2713 
2714  // Move mesh if in inflation mode
2715  if (map().hasMotionPoints())
2716  {
2717  mesh_.movePoints(map().preMotionPoints());
2718  }
2719  else
2720  {
2721  // Delete mesh volumes.
2722  mesh_.clearOut();
2723  }
2724 
2725  // Reset the instance for if in overwrite mode
2726  mesh_.setInstance(timeName());
2727 
2728  // Update intersections. Is mapping only (no faces created, positions stay
2729  // same) so no need to recalculate intersections.
2730  topoChange(map, labelList(0));
2731 
2732  return map;
2733 }
2734 
2735 
2738 {
2739  // Analyse which points need to be duplicated
2740  localPointRegion regionSide(mesh_);
2741 
2742  return dupNonManifoldPoints(regionSide);
2743 }
2744 
2745 
2748  const List<point>& insidePoints,
2749  const bool allowFreeStandingZoneFaces
2750 )
2751 {
2752  const PtrList<surfaceZonesInfo>& surfZones = surfaces_.surfZones();
2753 
2754  const labelList namedSurfaces
2755  (
2757  );
2758 
2759  forAll(namedSurfaces, i)
2760  {
2761  label surfi = namedSurfaces[i];
2762 
2763  Info<< "Surface : " << surfaces_.names()[surfi] << nl
2764  << " faceZone : " << surfZones[surfi].faceZoneName() << nl
2765  << " cellZone : " << surfZones[surfi].cellZoneName() << endl;
2766  }
2767 
2768 
2769  // Add zones to mesh
2770  const labelList surfaceToFaceZone =
2772  (
2773  surfZones,
2774  namedSurfaces,
2775  mesh_
2776  );
2777 
2778  const labelList surfaceToCellZone =
2780  (
2781  surfZones,
2782  namedSurfaces,
2783  mesh_
2784  );
2785 
2786 
2787  const pointField& cellCentres = mesh_.cellCentres();
2788  const labelList& faceOwner = mesh_.faceOwner();
2789  const labelList& faceNeighbour = mesh_.faceNeighbour();
2790  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
2791 
2792 
2793  // Swap neighbouring cell centres and cell level
2794  labelList neiLevel(mesh_.nFaces() - mesh_.nInternalFaces());
2795  pointField neiCc(mesh_.nFaces() - mesh_.nInternalFaces());
2796  calcNeighbourData(neiLevel, neiCc);
2797 
2798 
2799  // Mark faces intersecting zoned surfaces
2800  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2801 
2802 
2803  // Like surfaceIndex_ but only for named surfaces.
2804  labelList namedSurfaceIndex(mesh_.nFaces(), -1);
2805  PackedBoolList posOrientation(mesh_.nFaces());
2806 
2807  {
2808  // Statistics: number of faces per faceZone
2809  labelList nSurfFaces(surfZones.size(), 0);
2810 
2811  // Note: for all internal faces? internal + coupled?
2812  // Since zonify is run after baffling the surfaceIndex_ on baffles is
2813  // not synchronised across both baffle faces. Fortunately we don't
2814  // do zonify baffle faces anyway (they are normal boundary faces).
2815 
2816  // Collect candidate faces
2817  // ~~~~~~~~~~~~~~~~~~~~~~~
2818 
2819  labelList testFaces(intersectedFaces());
2820 
2821  // Collect segments
2822  // ~~~~~~~~~~~~~~~~
2823 
2824  pointField start(testFaces.size());
2825  pointField end(testFaces.size());
2826 
2827  forAll(testFaces, i)
2828  {
2829  label facei = testFaces[i];
2830 
2831  if (mesh_.isInternalFace(facei))
2832  {
2833  start[i] = cellCentres[faceOwner[facei]];
2834  end[i] = cellCentres[faceNeighbour[facei]];
2835  }
2836  else
2837  {
2838  start[i] = cellCentres[faceOwner[facei]];
2839  end[i] = neiCc[facei-mesh_.nInternalFaces()];
2840  }
2841  }
2842 
2843  // Extend segments a bit
2844  {
2845  const vectorField smallVec(rootSmall*(end-start));
2846  start -= smallVec;
2847  end += smallVec;
2848  }
2849 
2850 
2851  // Do test for intersections
2852  // ~~~~~~~~~~~~~~~~~~~~~~~~~
2853  // Note that we intersect all intersected faces again. Could reuse
2854  // the information already in surfaceIndex_.
2855 
2856  labelList surface1;
2857  List<pointIndexHit> hit1;
2858  vectorField normal1;
2859  labelList surface2;
2860  List<pointIndexHit> hit2;
2861  vectorField normal2;
2862  {
2863  labelList region1;
2864  labelList region2;
2865  surfaces_.findNearestIntersection
2866  (
2867  namedSurfaces,
2868  start,
2869  end,
2870 
2871  surface1,
2872  hit1,
2873  region1,
2874  normal1,
2875 
2876  surface2,
2877  hit2,
2878  region2,
2879  normal2
2880  );
2881  }
2882 
2883  forAll(testFaces, i)
2884  {
2885  label facei = testFaces[i];
2886  const vector& area = mesh_.faceAreas()[facei];
2887 
2888  if (surface1[i] != -1)
2889  {
2890  // If both hit should probably choose 'nearest'
2891  if
2892  (
2893  surface2[i] != -1
2894  && (
2895  magSqr(hit2[i].hitPoint())
2896  < magSqr(hit1[i].hitPoint())
2897  )
2898  )
2899  {
2900  namedSurfaceIndex[facei] = surface2[i];
2901  posOrientation[facei] = ((area&normal2[i]) > 0);
2902  nSurfFaces[surface2[i]]++;
2903  }
2904  else
2905  {
2906  namedSurfaceIndex[facei] = surface1[i];
2907  posOrientation[facei] = ((area&normal1[i]) > 0);
2908  nSurfFaces[surface1[i]]++;
2909  }
2910  }
2911  else if (surface2[i] != -1)
2912  {
2913  namedSurfaceIndex[facei] = surface2[i];
2914  posOrientation[facei] = ((area&normal2[i]) > 0);
2915  nSurfFaces[surface2[i]]++;
2916  }
2917  }
2918 
2919 
2920  // surfaceIndex might have different surfaces on both sides if
2921  // there happen to be a (obviously thin) surface with different
2922  // regions between the cell centres. If one is on a named surface
2923  // and the other is not this might give problems so sync.
2925  (
2926  mesh_,
2927  namedSurfaceIndex,
2928  maxEqOp<label>()
2929  );
2930 
2931  // Print a bit
2932  if (debug)
2933  {
2934  forAll(nSurfFaces, surfi)
2935  {
2936  Pout<< "Surface:"
2937  << surfaces_.names()[surfi]
2938  << " nZoneFaces:" << nSurfFaces[surfi] << nl;
2939  }
2940  Pout<< endl;
2941  }
2942  }
2943 
2944 
2945  // Put the cells into the correct zone
2946  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2947 
2948  // Zone per cell:
2949  // -2 : unset
2950  // -1 : not in any zone
2951  // >=0: zoneID
2952  labelList cellToZone(mesh_.nCells(), -2);
2953 
2954 
2955  // Set using geometric test
2956  // ~~~~~~~~~~~~~~~~~~~~~~~~
2957 
2958  // Closed surfaces with cellZone specified.
2959  labelList closedNamedSurfaces
2960  (
2962  (
2963  surfZones,
2964  surfaces_.geometry(),
2965  surfaces_.surfaces()
2966  )
2967  );
2968 
2969  if (closedNamedSurfaces.size())
2970  {
2971  Info<< "Found " << closedNamedSurfaces.size()
2972  << " closed, named surfaces. Assigning cells in/outside"
2973  << " these surfaces to the corresponding cellZone."
2974  << nl << endl;
2975 
2976  findCellZoneGeometric
2977  (
2978  neiCc,
2979  closedNamedSurfaces, // indices of closed surfaces
2980  namedSurfaceIndex, // per face index of named surface
2981  surfaceToCellZone, // cell zone index per surface
2982 
2983  cellToZone
2984  );
2985  }
2986 
2987 
2988  // Set using provided locations
2989  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2990 
2991  labelList locationSurfaces
2992  (
2994  );
2995 
2996  if (locationSurfaces.size())
2997  {
2998  Info<< "Found " << locationSurfaces.size()
2999  << " named surfaces with a provided inside point."
3000  << " Assigning cells inside these surfaces"
3001  << " to the corresponding cellZone."
3002  << nl << endl;
3003 
3004  findCellZoneInsideWalk
3005  (
3006  locationSurfaces, // indices of closed surfaces
3007  namedSurfaceIndex, // per face index of named surface
3008  surfaceToCellZone, // cell zone index per surface
3009 
3010  cellToZone
3011  );
3012  }
3013 
3014 
3015  // Set using walking
3016  // ~~~~~~~~~~~~~~~~~
3017 
3018  {
3019  Info<< "Walking from locations-in-mesh " << insidePoints
3020  << " to assign cellZones "
3021  << "- crossing a faceZone face changes cellZone" << nl << endl;
3022 
3023  // Topological walk
3024  findCellZoneTopo
3025  (
3026  insidePoints,
3027  namedSurfaceIndex,
3028  surfaceToCellZone,
3029 
3030  cellToZone
3031  );
3032  }
3033 
3034 
3035  // Make sure namedSurfaceIndex is unset in between same cell cell zones.
3036  if (!allowFreeStandingZoneFaces)
3037  {
3038  Info<< "Only selecting zone faces in between different cellZones."
3039  << nl << endl;
3040 
3041  makeConsistentFaceIndex(cellToZone, namedSurfaceIndex);
3042  }
3043 
3044 
3045  //- Per face index of faceZone or -1
3046  labelList faceToZone(mesh_.nFaces(), -1);
3047 
3048  // Convert namedSurfaceIndex (index of named surfaces) to
3049  // actual faceZone index
3050 
3051  forAll(namedSurfaceIndex, facei)
3052  {
3053  const label surfi = namedSurfaceIndex[facei];
3054  if (surfi != -1)
3055  {
3056  faceToZone[facei] = surfaceToFaceZone[surfi];
3057  }
3058  }
3059 
3060 
3061  // Topochange container
3062  polyTopoChange meshMod(mesh_);
3063 
3064 
3065 
3066  // Get coupled neighbour cellZone. Set to -1 on non-coupled patches.
3067  labelList neiCellZone;
3068  syncTools::swapBoundaryCellList(mesh_, cellToZone, neiCellZone);
3069  forAll(patches, patchi)
3070  {
3071  const polyPatch& pp = patches[patchi];
3072 
3073  if (!pp.coupled())
3074  {
3075  label bFacei = pp.start() - mesh_.nInternalFaces();
3076  forAll(pp, i)
3077  {
3078  neiCellZone[bFacei++] = -1;
3079  }
3080  }
3081  }
3082 
3083 
3084 
3085  // Get per face whether is it master (of a coupled set of faces)
3086  const PackedBoolList isMasterFace(syncTools::getMasterFaces(mesh_));
3087 
3088 
3089 
3090  // faceZones
3091  // ~~~~~~~~~
3092  // Faces on faceZones come in two variants:
3093  // - faces on the outside of a cellZone. They will be oriented to
3094  // point out of the maximum cellZone.
3095  // - free-standing faces. These will be oriented according to the
3096  // local surface normal. We do this in a two step algorithm:
3097  // - do a consistent orientation
3098  // - check number of faces with consistent orientation
3099  // - if <0 flip the whole patch
3100  boolList meshFlipMap(mesh_.nFaces(), false);
3101  {
3102  // Collect all data on zone faces without cellZones on either side.
3103  const indirectPrimitivePatch patch
3104  (
3106  (
3107  mesh_.faces(),
3108  freeStandingBaffleFaces
3109  (
3110  faceToZone,
3111  cellToZone,
3112  neiCellZone
3113  )
3114  ),
3115  mesh_.points()
3116  );
3117 
3118  label nFreeStanding = returnReduce(patch.size(), sumOp<label>());
3119  if (nFreeStanding > 0)
3120  {
3121  Info<< "Detected " << nFreeStanding << " free-standing zone faces"
3122  << endl;
3123 
3124  if (debug)
3125  {
3126  OBJstream str(mesh_.time().path()/"freeStanding.obj");
3127  str.write(patch.localFaces(), patch.localPoints(), false);
3128  }
3129 
3130 
3131  // Detect non-manifold edges
3132  labelList nMasterFacesPerEdge;
3133  calcPatchNumMasterFaces(isMasterFace, patch, nMasterFacesPerEdge);
3134 
3135  // Mark zones. Even a single original surface might create multiple
3136  // disconnected/non-manifold-connected zones
3137  labelList faceToConnectedZone;
3138  const label nZones = markPatchZones
3139  (
3140  patch,
3141  nMasterFacesPerEdge,
3142  faceToConnectedZone
3143  );
3144 
3145  Map<label> nPosOrientation(2*nZones);
3146  for (label zonei = 0; zonei < nZones; zonei++)
3147  {
3148  nPosOrientation.insert(zonei, 0);
3149  }
3150 
3151  // Make orientations consistent in a topological way. This just
3152  // checks the first face per zone for whether nPosOrientation
3153  // is negative (which it never is at this point)
3154  consistentOrientation
3155  (
3156  isMasterFace,
3157  patch,
3158  nMasterFacesPerEdge,
3159  faceToConnectedZone,
3160  nPosOrientation,
3161 
3162  meshFlipMap
3163  );
3164 
3165  // Count per region the number of orientations (taking the new
3166  // flipMap into account)
3167  forAll(patch.addressing(), facei)
3168  {
3169  label meshFacei = patch.addressing()[facei];
3170 
3171  if (isMasterFace[meshFacei])
3172  {
3173  label n = 1;
3174  if
3175  (
3176  bool(posOrientation[meshFacei])
3177  == meshFlipMap[meshFacei]
3178  )
3179  {
3180  n = -1;
3181  }
3182 
3183  nPosOrientation.find(faceToConnectedZone[facei])() += n;
3184  }
3185  }
3186  Pstream::mapCombineGather(nPosOrientation, plusEqOp<label>());
3187  Pstream::mapCombineScatter(nPosOrientation);
3188 
3189 
3190  Info<< "Split " << nFreeStanding << " free-standing zone faces"
3191  << " into " << nZones << " disconnected regions with size"
3192  << " (negative denotes wrong orientation) :"
3193  << endl;
3194 
3195  for (label zonei = 0; zonei < nZones; zonei++)
3196  {
3197  Info<< " " << zonei << "\t" << nPosOrientation[zonei]
3198  << endl;
3199  }
3200  Info<< endl;
3201 
3202 
3203  // Reapply with new counts (in nPosOrientation). This will cause
3204  // zones with a negative count to be flipped.
3205  consistentOrientation
3206  (
3207  isMasterFace,
3208  patch,
3209  nMasterFacesPerEdge,
3210  faceToConnectedZone,
3211  nPosOrientation,
3212 
3213  meshFlipMap
3214  );
3215  }
3216  }
3217 
3218 
3219  // Put the faces into the correct zone
3220  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
3221 
3222  for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
3223  {
3224  label faceZoneI = faceToZone[facei];
3225 
3226  if (faceZoneI != -1)
3227  {
3228  // Orient face zone to have slave cells in max cell zone.
3229  // Note: logic to use flipMap should be consistent with logic
3230  // to pick up the freeStandingBaffleFaces!
3231 
3232  const label ownZone = cellToZone[faceOwner[facei]];
3233  const label neiZone = cellToZone[faceNeighbour[facei]];
3234 
3235  bool flip;
3236 
3237  if (ownZone == neiZone)
3238  {
3239  // free-standing face. Use geometrically derived orientation
3240  flip = meshFlipMap[facei];
3241  }
3242  else
3243  {
3244  flip =
3245  (
3246  ownZone == -1
3247  || (neiZone != -1 && ownZone > neiZone)
3248  );
3249  }
3250 
3251  meshMod.setAction
3252  (
3254  (
3255  mesh_.faces()[facei], // modified face
3256  facei, // label of face
3257  faceOwner[facei], // owner
3258  faceNeighbour[facei], // neighbour
3259  false, // face flip
3260  -1, // patch for face
3261  false, // remove from zone
3262  faceZoneI, // zone for face
3263  flip // face flip in zone
3264  )
3265  );
3266  }
3267  }
3268 
3269 
3270  // Set owner as no-flip
3271  forAll(patches, patchi)
3272  {
3273  const polyPatch& pp = patches[patchi];
3274 
3275  label facei = pp.start();
3276 
3277  forAll(pp, i)
3278  {
3279  const label faceZoneI = faceToZone[facei];
3280 
3281  if (faceZoneI != -1)
3282  {
3283  const label ownZone = cellToZone[faceOwner[facei]];
3284  const label neiZone = neiCellZone[facei-mesh_.nInternalFaces()];
3285 
3286  bool flip;
3287 
3288  if (ownZone == neiZone)
3289  {
3290  // free-standing face. Use geometrically derived orientation
3291  flip = meshFlipMap[facei];
3292  }
3293  else
3294  {
3295  flip =
3296  (
3297  ownZone == -1
3298  || (neiZone != -1 && ownZone > neiZone)
3299  );
3300  }
3301 
3302  meshMod.setAction
3303  (
3305  (
3306  mesh_.faces()[facei], // modified face
3307  facei, // label of face
3308  faceOwner[facei], // owner
3309  -1, // neighbour
3310  false, // face flip
3311  patchi, // patch for face
3312  false, // remove from zone
3313  faceZoneI, // zone for face
3314  flip // face flip in zone
3315  )
3316  );
3317  }
3318  facei++;
3319  }
3320  }
3321 
3322 
3323  // Put the cells into the correct zone
3324  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
3325 
3326  forAll(cellToZone, celli)
3327  {
3328  const label zonei = cellToZone[celli];
3329 
3330  if (zonei >= 0)
3331  {
3332  meshMod.setAction
3333  (
3335  (
3336  celli,
3337  false, // removeFromZone
3338  zonei
3339  )
3340  );
3341  }
3342  }
3343 
3344  mesh_.clearOut();
3345 
3346  // Change the mesh (no inflation, parallel sync)
3347  autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(mesh_, false, true);
3348 
3349  // Update fields
3350  mesh_.topoChange(map);
3351 
3352  // Move mesh if in inflation mode
3353  if (map().hasMotionPoints())
3354  {
3355  mesh_.movePoints(map().preMotionPoints());
3356  }
3357  else
3358  {
3359  // Delete mesh volumes.
3360  mesh_.clearOut();
3361  }
3362 
3363  // Reset the instance for if in overwrite mode
3364  mesh_.setInstance(timeName());
3365 
3366  // Print some stats (note: zones are synchronised)
3367  if (mesh_.cellZones().size() > 0)
3368  {
3369  Info<< "CellZones:" << endl;
3370  forAll(mesh_.cellZones(), zonei)
3371  {
3372  const cellZone& cz = mesh_.cellZones()[zonei];
3373  Info<< " " << cz.name()
3374  << "\tsize:" << returnReduce(cz.size(), sumOp<label>())
3375  << endl;
3376  }
3377  Info<< endl;
3378  }
3379  if (mesh_.faceZones().size() > 0)
3380  {
3381  Info<< "FaceZones:" << endl;
3382  forAll(mesh_.faceZones(), zonei)
3383  {
3384  const faceZone& fz = mesh_.faceZones()[zonei];
3385  Info<< " " << fz.name()
3386  << "\tsize:" << returnReduce(fz.size(), sumOp<label>())
3387  << endl;
3388  }
3389  Info<< endl;
3390  }
3391 
3392  // None of the faces has changed, only the zones. Still...
3393  topoChange(map, labelList());
3394 
3395  return map;
3396 }
3397 
3398 
3399 // ************************************************************************* //
const fvPatchList & patches
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:453
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Definition: regionSplit.H:116
Class containing data for face removal.
void clearOut()
Clear all geometry and addressing.
Definition: fvMesh.C:258
void printMeshInfo(const bool, const string &) const
Print some mesh stats.
bool updateEdge(const polyMesh &mesh, const indirectPrimitivePatch &patch, const label edgeI, const label facei, const patchEdgeFaceRegion &faceInfo, const scalar tol, TrackingData &td)
Influence of face on edge.
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
autoPtr< polyTopoChangeMap > changeMesh(polyMesh &mesh, const bool inflate, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Inplace changes mesh without change of patches.
void checkZoneFaces() const
Debug helper: check faceZones are not on processor patches.
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
virtual Ostream & write(const char)=0
Write character.
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const word & name() const
Return name.
const word & name() const
Return name.
Definition: IOobject.H:315
const labelList & surfaces() const
A list of face labels.
Definition: faceSet.H:48
An STL-conforming const_iterator.
Definition: HashTable.H:481
autoPtr< polyTopoChangeMap > createBaffles(const labelList &ownPatch, const labelList &nbrPatch)
Create baffle for every internal face where ownPatch != -1.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Class describing modification of a face.
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
const meshCellZones & cellZones() const
Return cell zones.
Definition: polyMesh.H:501
label countHits() const
Count number of intersections (local)
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
const labelListList & faceEdges() const
autoPtr< polyTopoChangeMap > zonify(const List< point > &insidePoints, const bool allowFreeStandingZoneFaces)
Put faces/cells into zones according to surface specification.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
const Field< PointType > & faceCentres() const
Return face centres for patch.
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
label nInternalFaces() const
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1243
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:928
const labelListList & pointEdges() const
label nFaces() const
Wave propagation of information along patch. Every iteration information goes through one layer of fa...
Unit conversion functions.
label nTotalPoints() const
Return total number of points in decomposed mesh. Not.
void topoChange(const polyTopoChangeMap &, const labelList &changedFaces)
Update for external change to mesh. changedFaces are in new mesh.
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:249
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
Given list of cells to remove insert all the topology changes.
Definition: removeCells.H:59
label nCells() const
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
autoPtr< polyTopoChangeMap > dupNonManifoldPoints()
Find boundary points that connect to more than one cell.
labelList intersectedFaces() const
Get faces with intersection.
static writeType writeLevel()
Get/set write level.
Holds information (coordinate and normal) regarding nearest wall point.
autoPtr< polyTopoChangeMap > createZoneBaffles(const labelList &globalToMasterPatch, const labelList &globalToSlavePatch, List< labelPair > &)
Create baffles for faces straddling zoned surfaces. Return.
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
Transport of orientation for use in PatchEdgeFaceWave.
static labelList addCellZonesToMesh(const PtrList< surfaceZonesInfo > &surfList, const labelList &namedSurfaces, polyMesh &mesh)
const cellList & cells() const
label whichFace(const label globalCellID) const
Helper function to re-direct to zone::localID(...)
Definition: faceZone.C:307
const Field< PointType > & localPoints() const
Return pointField of points in patch.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:372
virtual Ostream & write(const char)
Write character.
Definition: OBJstream.C:81
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
Foam::indirectPrimitivePatch.
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:111
T & first()
Return the first element of the list.
Definition: UListI.H:114
Class describing modification of a cell.
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
const Map< label > & meshPointMap() const
Per point that is to be duplicated the local index.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
static PackedBoolList getMasterFaces(const polyMesh &)
Get per face whether it is uncoupled or a master of a.
Definition: syncTools.C:153
scalar f1
Definition: createFields.H:15
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
virtual void topoChange(const polyTopoChangeMap &map)
Update mesh corresponding to the given map.
Definition: fvMesh.C:1151
autoPtr< polyTopoChangeMap > splitMesh(const label nBufferLayers, const labelList &globalToMasterPatch, const labelList &globalToSlavePatch, const refinementParameters::cellSelectionPoints &selectionPoints)
Split off (with optional buffer layers) unreachable areas.
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Takes mesh with &#39;baffles&#39; (= boundary faces sharing points). Determines for selected points on bounda...
static labelPairList findDuplicateFacePairs(const polyMesh &)
Helper routine to find all baffles (two boundary faces.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1211
static void mapCombineScatter(const List< commsStruct > &comms, Container &Values, const int tag, const label comm)
Scatter data. Reverse of combineGather.
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
iterator find(const label &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:142
const searchableSurfaces & geometry() const
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
A list of faces which address into the list of points.
const List< point > & inside() const
Return the points inside the surface regions to selected cells.
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Duplicate points.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:316
dimensionedScalar cos(const dimensionedScalar &ds)
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
static void syncEdgeList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronise values on all mesh edges.
Simple container to keep together snap specific information.
double cpuTimeIncrement() const
Return CPU time (in seconds) since last call to cpuTimeIncrement()
Definition: cpuTime.C:60
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
Transport of region for use in PatchEdgeFaceWave.
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1237
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
void flip()
Reverse orientation.
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1435
static const label labelMax
Definition: label.H:62
List< label > labelList
A List of labels.
Definition: labelList.H:56
const word & name() const
Return name.
Definition: zone.H:147
static void mapCombineGather(const List< commsStruct > &comms, Container &Values, const CombineOp &cop, const int tag, const label comm)
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1224
static labelList getInsidePointNamedSurfaces(const PtrList< surfaceZonesInfo > &surfList)
Get indices of surfaces with a cellZone that have &#39;insidePoint&#39;.
const PtrList< surfaceZonesInfo > & surfZones() const
const vectorField & cellCentres() const
static labelList getUnnamedSurfaces(const PtrList< surfaceZonesInfo > &surfList)
Get indices of unnamed surfaces (surfaces without faceZoneName)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
word timeName() const
Replacement for Time::timeName() : return oldInstance (if.
Foam::polyBoundaryMesh.
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=*(cellModeller::lookup("hex"));labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells]=cellShape(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< small) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};const label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &mergedCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< SLList< face > > pFaces[nBCs]
Definition: readKivaGrid.H:235
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
label nEdges() const
Return number of edges in patch.
void reverse(UList< T > &, const label n)
Definition: UListI.H:334
OFstream which keeps track of vertices.
Definition: OBJstream.H:53
static const char nl
Definition: Ostream.H:260
void checkData()
Debugging: check that all faces still obey start()>end()
label globalRegion(const label surfi, const label regioni) const
From surface and region on surface to global region.
void baffleAndSplitMesh(const bool handleSnapProblems, const snapParameters &snapParams, const bool useTopologicalSnapDetection, const bool removeEdgeConnectedCells, const scalarField &perpendicularAngle, const bool mergeFreeStanding, const scalar freeStandingAngle, const dictionary &motionDict, Time &runTime, const labelList &globalToMasterPatch, const labelList &globalToSlavePatch, const refinementParameters::cellSelectionPoints &selectionPoints)
Split off unreachable areas of mesh.
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.
A subset of mesh cells.
Definition: cellZone.H:61
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:77
const Type & second() const
Return second.
Definition: Pair.H:110
label nEdges() const
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: fvMesh.C:1065
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData)
Swap to obtain neighbour cell values for all boundary faces.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
const vectorField & faceCentres() const
static void syncFaceList(const polyMesh &mesh, UList< T > &l, const CombineOp &cop)
Synchronise values on all mesh faces.
Definition: syncTools.H:387
void setSize(const label)
Reset size of List.
Definition: List.C:281
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: MeshZones.C:221
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:459
A bit-packed bool list.
const wordList & names() const
Names of surfaces.
label patchi
autoPtr< polyTopoChangeMap > splitMeshRegions(const labelList &globalToMasterPatch, const labelList &globalToSlavePatch, const refinementParameters::cellSelectionPoints &selectionPoints)
Split mesh according to selectionPoints.
#define WarningInFunction
Report a warning using Foam::Warning.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
labelList meshEdges(const edgeList &allEdges, const labelListList &cellEdges, const labelList &faceCells) const
Return labels of patch edges in the global edge list using.
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
const meshFaceZones & faceZones() const
Return face zones.
Definition: polyMesh.H:495
const List< FaceType > & localFaces() const
Return patch faces addressing into local point list.
const vectorField & faceAreas() const
Direct mesh changes based on v1.3 polyTopoChange syntax.
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:306
static labelList addFaceZonesToMesh(const PtrList< surfaceZonesInfo > &surfList, const labelList &namedSurfaces, polyMesh &mesh)
void setRefinement(const localPointRegion &regionSide, polyTopoChange &)
Play commands into polyTopoChange to duplicate points. Gets.
static void syncPointList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronise values on all mesh points.
label nRegions() const
Return total number of regions.
Definition: regionSplit.H:198
static labelList getClosedNamedSurfaces(const PtrList< surfaceZonesInfo > &surfList, const searchableSurfaces &allGeometry, const labelList &surfaces)
Get indices of surfaces with a cellZone that are closed and.
static labelList getNamedSurfaces(const PtrList< surfaceZonesInfo > &surfList)
Get indices of named surfaces (surfaces with faceZoneName)
messageStream Info
const labelListList & pointFaces() const
const labelListList & faceEdges() const
Return face-edge addressing.
dimensioned< scalar > mag(const dimensioned< Type > &)
static void findRegions(const polyMesh &, labelList &cellRegion, const vector &perturbVec, const refinementParameters::cellSelectionPoints &selectionPoints)
Find regions points are in and update cellRegion.
label n
label nPoints() const
Field< vector > vectorField
Specialisation of Field<T> for vector.
void setSize(const label)
Dummy setSize function.
Definition: FixedListI.H:183
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
labelList getExposedFaces(const labelList &cellsToRemove) const
Get labels of exposed faces.
Definition: removeCells.C:73
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &l)
Swap coupled boundary face values.
Definition: syncTools.H:436
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:65
List< Key > toc() const
Return the table of contents.
Definition: HashTable.C:202
virtual void sync(const polyMesh &mesh)
Sync faceSet across coupled patches.
Definition: faceSet.C:113
bool write() const
Write mesh and all data.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
void findNearestIntersection(const labelList &surfacesToTest, const pointField &start, const pointField &end, labelList &surface1, List< pointIndexHit > &hit1, labelList &region1, labelList &surface2, List< pointIndexHit > &hit2, labelList &region2) const
Find intersection nearest to the endpoints. surface1,2 are.
autoPtr< polyTopoChangeMap > mergeBaffles(const List< labelPair > &)
Merge baffles. Gets pairs of faces.
fileName path() const
Return path.
Definition: TimePaths.H:139
static label findRegion(const polyMesh &, const labelList &cellRegion, const vector &perturbVec, const point &location)
Find region point is in. Uses optional perturbation to re-test.
const Type & first() const
Return first.
Definition: Pair.H:98
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
bool updateEdge(const polyMesh &mesh, const indirectPrimitivePatch &patch, const label edgei, const label facei, const patchFaceOrientation &faceInfo, const scalar tol, TrackingData &td)
Influence of face on edge.
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
A HashTable to objects of type <T> with a label key.
Definition: Map.H:49
Class to hold the points to select cells inside and outside.
void findInside(const labelList &surfacesToTest, const pointField &pt, labelList &insideSurfaces) const
Detect if a point is &#39;inside&#39; (closed) surfaces.
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.