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