meshRefinement.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "meshRefinement.H"
27 #include "volMesh.H"
28 #include "volFields.H"
29 #include "surfaceMesh.H"
30 #include "syncTools.H"
31 #include "Time.H"
32 #include "refinementHistory.H"
33 #include "refinementSurfaces.H"
34 #include "refinementFeatures.H"
35 #include "decompositionMethod.H"
36 #include "regionSplit.H"
37 #include "fvMeshDistribute.H"
38 #include "indirectPrimitivePatch.H"
39 #include "polyTopoChange.H"
40 #include "removeCells.H"
41 #include "polyDistributionMap.H"
42 #include "localPointRegion.H"
43 #include "pointMesh.H"
44 #include "pointFields.H"
45 #include "slipPointPatchFields.H"
49 #include "processorPointPatch.H"
50 #include "globalIndex.H"
51 #include "meshTools.H"
52 #include "OFstream.H"
53 #include "geometric.H"
54 #include "randomGenerator.H"
55 #include "searchableSurfaces.H"
56 #include "treeBoundBox.H"
57 #include "fvMeshTools.H"
58 
59 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
60 
61 namespace Foam
62 {
64 
65  template<>
66  const char* Foam::NamedEnum
67  <
69  5
70  >::names[] =
71  {
72  "mesh",
73  "intersections",
74  "featureSeeds",
75  "attraction",
76  "layerInfo"
77  };
78 
79  template<>
80  const char* Foam::NamedEnum
81  <
83  1
84  >::names[] =
85  {
86  "layerInfo"
87  };
88 
89  template<>
90  const char* Foam::NamedEnum
91  <
93  5
94  >::names[] =
95  {
96  "mesh",
97  "noRefinement",
98  "scalarLevels",
99  "layerSets",
100  "layerFields"
101  };
102 
103 }
104 
107 
110 
113 
114 
115 Foam::meshRefinement::writeType Foam::meshRefinement::writeLevel_;
116 
117 Foam::meshRefinement::outputType Foam::meshRefinement::outputLevel_;
118 
119 
120 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
121 
122 void Foam::meshRefinement::calcNeighbourData
123 (
124  labelList& neiLevel,
125  pointField& neiCc
126 ) const
127 {
128  const labelList& cellLevel = meshCutter_.cellLevel();
129  const pointField& cellCentres = mesh_.cellCentres();
130 
131  const label nBoundaryFaces = mesh_.nFaces() - mesh_.nInternalFaces();
132 
133  if (neiLevel.size() != nBoundaryFaces || neiCc.size() != nBoundaryFaces)
134  {
136  << nBoundaryFaces << " neiLevel:" << neiLevel.size()
137  << abort(FatalError);
138  }
139 
140  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
141 
142  labelHashSet addedPatchIDSet(meshedPatches());
143 
145  {
146  const polyPatch& pp = patches[patchi];
147 
148  const labelUList& faceCells = pp.faceCells();
149  const vectorField::subField faceCentres = pp.faceCentres();
150  const vectorField::subField faceAreas = pp.faceAreas();
151 
152  label bFacei = pp.start() - mesh_.nInternalFaces();
153 
154  if (pp.coupled())
155  {
156  forAll(faceCells, i)
157  {
158  neiLevel[bFacei] = cellLevel[faceCells[i]];
159  neiCc[bFacei] = cellCentres[faceCells[i]];
160  bFacei++;
161  }
162  }
163  else if (addedPatchIDSet.found(patchi))
164  {
165  // Face was introduced from cell-cell intersection. Try to
166  // reconstruct other side cell(centre). Three possibilities:
167  // - cells same size.
168  // - preserved cell smaller. Not handled.
169  // - preserved cell larger.
170  forAll(faceCells, i)
171  {
172  // Extrapolate the face centre.
173  const vector fn = faceAreas[i]/(mag(faceAreas[i]) + vSmall);
174 
175  const label own = faceCells[i];
176  const label ownLevel = cellLevel[own];
177  const label faceLevel = meshCutter_.faceLevel(pp.start() + i);
178 
179  // Normal distance from face centre to cell centre
180  scalar d = ((faceCentres[i] - cellCentres[own]) & fn);
181  if (faceLevel > ownLevel)
182  {
183  // Other cell more refined. Adjust normal distance
184  d *= 0.5;
185  }
186  neiLevel[bFacei] = faceLevel;
187  // Calculate other cell centre by extrapolation
188  neiCc[bFacei] = faceCentres[i] + d*fn;
189  bFacei++;
190  }
191  }
192  else
193  {
194  forAll(faceCells, i)
195  {
196  neiLevel[bFacei] = cellLevel[faceCells[i]];
197  neiCc[bFacei] = faceCentres[i];
198  bFacei++;
199  }
200  }
201  }
202 
203  // Swap coupled boundaries. Apply separation to cc since is coordinate.
205  syncTools::swapBoundaryFaceList(mesh_, neiLevel);
206 }
207 
208 
209 void Foam::meshRefinement::updateIntersections(const labelList& changedFaces)
210 {
211  const pointField& cellCentres = mesh_.cellCentres();
212 
213  // Stats on edges to test. Count proc faces only once.
214  PackedBoolList isMasterFace(syncTools::getMasterFaces(mesh_));
215 
216  {
217  label nMasterFaces = 0;
218  forAll(isMasterFace, facei)
219  {
220  if (isMasterFace.get(facei) == 1)
221  {
222  nMasterFaces++;
223  }
224  }
225  reduce(nMasterFaces, sumOp<label>());
226 
227  label nChangedFaces = 0;
228  forAll(changedFaces, i)
229  {
230  if (isMasterFace.get(changedFaces[i]) == 1)
231  {
232  nChangedFaces++;
233  }
234  }
235  reduce(nChangedFaces, sumOp<label>());
236 
237  Info<< "Edge intersection testing:" << nl
238  << " Number of edges : " << nMasterFaces << nl
239  << " Number of edges to retest : " << nChangedFaces
240  << endl;
241  }
242 
243 
244  // Get boundary face centre and level. Coupled aware.
245  labelList neiLevel(mesh_.nFaces() - mesh_.nInternalFaces());
246  pointField neiCc(mesh_.nFaces() - mesh_.nInternalFaces());
247  calcNeighbourData(neiLevel, neiCc);
248 
249  // Collect segments we want to test for
250  pointField start(changedFaces.size());
251  pointField end(changedFaces.size());
252 
253  forAll(changedFaces, i)
254  {
255  const label facei = changedFaces[i];
256  const label own = mesh_.faceOwner()[facei];
257 
258  start[i] = cellCentres[own];
259  if (mesh_.isInternalFace(facei))
260  {
261  end[i] = cellCentres[mesh_.faceNeighbour()[facei]];
262  }
263  else
264  {
265  end[i] = neiCc[facei-mesh_.nInternalFaces()];
266  }
267  }
268 
269  // Extend segments a bit
270  {
271  const vectorField smallVec(rootSmall*(end-start));
272  start -= smallVec;
273  end += smallVec;
274  }
275 
276 
277  // Do tests in one go
278  labelList surfaceHit;
279  {
280  labelList surfaceLevel;
281  surfaces_.findHigherIntersection
282  (
283  start,
284  end,
285  labelList(start.size(), -1), // accept any intersection
286  surfaceHit,
287  surfaceLevel
288  );
289  }
290 
291  // Keep just surface hit
292  forAll(surfaceHit, i)
293  {
294  surfaceIndex_[changedFaces[i]] = surfaceHit[i];
295  }
296 
297  // Make sure both sides have same information. This should be
298  // case in general since same vectors but just to make sure.
299  syncTools::syncFaceList(mesh_, surfaceIndex_, maxEqOp<label>());
300 
301  const label nHits = countHits();
302  const label nTotHits = returnReduce(nHits, sumOp<label>());
303 
304  Info<< " Number of intersected edges : " << nTotHits << endl;
305 
306  // Set files to same time as mesh
307  setInstance(mesh_.facesInstance());
308 }
309 
310 
312 (
313  const string& msg,
314  const polyMesh& mesh,
315  const List<scalar>& fld
316 )
317 {
318  if (fld.size() != mesh.nPoints())
319  {
321  << msg << endl
322  << "fld size:" << fld.size() << " mesh points:" << mesh.nPoints()
323  << abort(FatalError);
324  }
325 
326  Pout<< "Checking field " << msg << endl;
327  scalarField minFld(fld);
329  (
330  mesh,
331  minFld,
332  minEqOp<scalar>(),
333  great
334  );
335  scalarField maxFld(fld);
337  (
338  mesh,
339  maxFld,
340  maxEqOp<scalar>(),
341  -great
342  );
343  forAll(minFld, pointi)
344  {
345  const scalar& minVal = minFld[pointi];
346  const scalar& maxVal = maxFld[pointi];
347  if (mag(minVal-maxVal) > small)
348  {
349  Pout<< msg << " at:" << mesh.points()[pointi] << nl
350  << " minFld:" << minVal << nl
351  << " maxFld:" << maxVal << nl
352  << endl;
353  }
354  }
355 }
356 
357 
359 (
360  const string& msg,
361  const polyMesh& mesh,
362  const List<point>& fld
363 )
364 {
365  if (fld.size() != mesh.nPoints())
366  {
368  << msg << endl
369  << "fld size:" << fld.size() << " mesh points:" << mesh.nPoints()
370  << abort(FatalError);
371  }
372 
373  Pout<< "Checking field " << msg << endl;
374  pointField minFld(fld);
376  (
377  mesh,
378  minFld,
380  point(great, great, great)
381  );
382  pointField maxFld(fld);
384  (
385  mesh,
386  maxFld,
389  );
390  forAll(minFld, pointi)
391  {
392  const point& minVal = minFld[pointi];
393  const point& maxVal = maxFld[pointi];
394  if (mag(minVal-maxVal) > small)
395  {
396  Pout<< msg << " at:" << mesh.points()[pointi] << nl
397  << " minFld:" << minVal << nl
398  << " maxFld:" << maxVal << nl
399  << endl;
400  }
401  }
402 }
403 
404 
406 {
407  Pout<< "meshRefinement::checkData() : Checking refinement structure."
408  << endl;
409  meshCutter_.checkMesh();
410 
411  Pout<< "meshRefinement::checkData() : Checking refinement levels."
412  << endl;
413  meshCutter_.checkRefinementLevels(1, labelList(0));
414 
415 
416  label nBnd = mesh_.nFaces() - mesh_.nInternalFaces();
417 
418  Pout<< "meshRefinement::checkData() : Checking synchronisation."
419  << endl;
420 
421  // Check face centres
422  {
423  // Boundary face centres
424  pointField::subList boundaryFc
425  (
426  mesh_.faceCentres(),
427  nBnd,
428  mesh_.nInternalFaces()
429  );
430 
431  // Get neighbouring face centres
432  pointField neiBoundaryFc(boundaryFc);
434  (
435  mesh_,
436  neiBoundaryFc,
437  eqOp<point>()
438  );
439 
440  // Compare
441  testSyncBoundaryFaceList
442  (
443  mergeDistance_,
444  "testing faceCentres : ",
445  boundaryFc,
446  neiBoundaryFc
447  );
448  }
449  // Check meshRefinement
450  {
451  // Get boundary face centre and level. Coupled aware.
452  labelList neiLevel(nBnd);
453  pointField neiCc(nBnd);
454  calcNeighbourData(neiLevel, neiCc);
455 
456  // Collect segments we want to test for
457  pointField start(mesh_.nFaces());
458  pointField end(mesh_.nFaces());
459 
460  forAll(start, facei)
461  {
462  start[facei] = mesh_.cellCentres()[mesh_.faceOwner()[facei]];
463 
464  if (mesh_.isInternalFace(facei))
465  {
466  end[facei] = mesh_.cellCentres()[mesh_.faceNeighbour()[facei]];
467  }
468  else
469  {
470  end[facei] = neiCc[facei-mesh_.nInternalFaces()];
471  }
472  }
473 
474  // Extend segments a bit
475  {
476  const vectorField smallVec(rootSmall*(end-start));
477  start -= smallVec;
478  end += smallVec;
479  }
480 
481 
482  // Do tests in one go
483  labelList surfaceHit;
484  {
485  labelList surfaceLevel;
486  surfaces_.findHigherIntersection
487  (
488  start,
489  end,
490  labelList(start.size(), -1), // accept any intersection
491  surfaceHit,
492  surfaceLevel
493  );
494  }
495  // Get the coupled hit
496  labelList neiHit
497  (
499  (
500  surfaceHit,
501  nBnd,
502  mesh_.nInternalFaces()
503  )
504  );
505  syncTools::swapBoundaryFaceList(mesh_, neiHit);
506 
507  // Check
508  forAll(surfaceHit, facei)
509  {
510  if (surfaceIndex_[facei] != surfaceHit[facei])
511  {
512  if (mesh_.isInternalFace(facei))
513  {
515  << "Internal face:" << facei
516  << " fc:" << mesh_.faceCentres()[facei]
517  << " cached surfaceIndex_:" << surfaceIndex_[facei]
518  << " current:" << surfaceHit[facei]
519  << " ownCc:"
520  << mesh_.cellCentres()[mesh_.faceOwner()[facei]]
521  << " neiCc:"
522  << mesh_.cellCentres()[mesh_.faceNeighbour()[facei]]
523  << endl;
524  }
525  else if
526  (
527  surfaceIndex_[facei]
528  != neiHit[facei-mesh_.nInternalFaces()]
529  )
530  {
532  << "Boundary face:" << facei
533  << " fc:" << mesh_.faceCentres()[facei]
534  << " cached surfaceIndex_:" << surfaceIndex_[facei]
535  << " current:" << surfaceHit[facei]
536  << " ownCc:"
537  << mesh_.cellCentres()[mesh_.faceOwner()[facei]]
538  << " end:" << end[facei]
539  << endl;
540  }
541  }
542  }
543  }
544  {
545  labelList::subList boundarySurface
546  (
547  surfaceIndex_,
548  mesh_.nFaces() - mesh_.nInternalFaces(),
549  mesh_.nInternalFaces()
550  );
551 
552  labelList neiBoundarySurface(boundarySurface);
554  (
555  mesh_,
556  neiBoundarySurface
557  );
558 
559  // Compare
560  testSyncBoundaryFaceList
561  (
562  0, // tolerance
563  "testing surfaceIndex() : ",
564  boundarySurface,
565  neiBoundarySurface
566  );
567  }
568 
569 
570  // Find duplicate faces
571  Pout<< "meshRefinement::checkData() : Counting duplicate faces."
572  << endl;
573 
574  labelList duplicateFace
575  (
577  (
578  mesh_,
579  identityMap(mesh_.nFaces() - mesh_.nInternalFaces())
580  + mesh_.nInternalFaces()
581  )
582  );
583 
584  // Count
585  {
586  label nDup = 0;
587 
588  forAll(duplicateFace, i)
589  {
590  if (duplicateFace[i] != -1)
591  {
592  nDup++;
593  }
594  }
595  nDup /= 2; // will have counted both faces of duplicate
596  Pout<< "meshRefinement::checkData() : Found " << nDup
597  << " duplicate pairs of faces." << endl;
598  }
599 }
600 
601 
603 {
604  meshCutter_.setInstance(inst);
605  surfaceIndex_.instance() = inst;
606 }
607 
608 
609 Foam::autoPtr<Foam::polyTopoChangeMap> Foam::meshRefinement::doRemoveCells
610 (
611  const labelList& cellsToRemove,
612  const labelList& exposedFaces,
613  const labelList& exposedPatchIDs,
614  removeCells& cellRemover
615 )
616 {
617  polyTopoChange meshMod(mesh_);
618 
619  // Arbitrary: put exposed faces into last patch.
620  cellRemover.setRefinement
621  (
622  cellsToRemove,
623  exposedFaces,
624  exposedPatchIDs,
625  meshMod
626  );
627 
628  // Change the mesh
629  autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(mesh_, true);
630 
631  // Update fields
632  mesh_.topoChange(map);
633 
634  // Reset the instance for if in overwrite mode
635  mesh_.setInstance(name());
636  setInstance(mesh_.facesInstance());
637 
638  // Update local mesh data
639  cellRemover.topoChange(map);
640 
641  // Update intersections. Recalculate intersections for exposed faces.
642  labelList newExposedFaces = renumber
643  (
644  map().reverseFaceMap(),
645  exposedFaces
646  );
647 
648  // Pout<< "removeCells : updating intersections for "
649  // << newExposedFaces.size() << " newly exposed faces." << endl;
650 
651  topoChange(map, newExposedFaces);
652 
653  return map;
654 }
655 
656 
658 (
659  const labelList& splitFaces,
660  const labelPairList& splits
661 )
662 {
663  polyTopoChange meshMod(mesh_);
664 
665  forAll(splitFaces, i)
666  {
667  const label facei = splitFaces[i];
668  const face& f = mesh_.faces()[facei];
669 
670  // Split as start and end index in face
671  const labelPair& split = splits[i];
672 
673  label nVerts = split[1] - split[0];
674  if (nVerts < 0)
675  {
676  nVerts += f.size();
677  }
678  nVerts += 1;
679 
680 
681  // Split into f0, f1
682  face f0(nVerts);
683 
684  label fp = split[0];
685  forAll(f0, i)
686  {
687  f0[i] = f[fp];
688  fp = f.fcIndex(fp);
689  }
690 
691  face f1(f.size() - f0.size() + 2);
692  fp = split[1];
693  forAll(f1, i)
694  {
695  f1[i] = f[fp];
696  fp = f.fcIndex(fp);
697  }
698 
699 
700  // Determine face properties
701  const label own = mesh_.faceOwner()[facei];
702  label nei = -1;
703  label patchi = -1;
704  if (facei >= mesh_.nInternalFaces())
705  {
706  patchi = mesh_.boundaryMesh().whichPatch(facei);
707  }
708  else
709  {
710  nei = mesh_.faceNeighbour()[facei];
711  }
712 
713  if (debug)
714  {
715  Pout<< "face:" << facei << " verts:" << f
716  << " split into f0:" << f0
717  << " f1:" << f1 << endl;
718  }
719 
720  // Change/add faces
721  meshMod.modifyFace
722  (
723  f0, // modified face
724  facei, // label of face
725  own, // owner
726  nei, // neighbour
727  false, // face flip
728  patchi // patch for face
729  );
730 
731  meshMod.addFace
732  (
733  f1, // modified face
734  own, // owner
735  nei, // neighbour
736  facei, // master face
737  false, // face flip
738  patchi // patch for face
739  );
740  }
741 
742 
743  // Change the mesh (without keeping old points)
744  autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(mesh_, true);
745 
746  // Update fields
747  mesh_.topoChange(map);
748 
749  // Reset the instance for if in overwrite mode
750  mesh_.setInstance(name());
751  setInstance(mesh_.facesInstance());
752 
753  // Update local mesh data
754  const labelList& oldToNew = map().reverseFaceMap();
755  labelList newSplitFaces(renumber(oldToNew, splitFaces));
756  // Add added faces (every splitFaces becomes two faces)
757  label sz = newSplitFaces.size();
758  newSplitFaces.setSize(2*sz);
759  forAll(map().faceMap(), facei)
760  {
761  label oldFacei = map().faceMap()[facei];
762  if (oldToNew[oldFacei] != facei)
763  {
764  // Added face
765  newSplitFaces[sz++] = facei;
766  }
767  }
768 
769  topoChange(map, newSplitFaces);
770 
771  return map;
772 }
773 
774 
777 //void Foam::meshRefinement::getCoupledRegionMaster
778 //(
779 // const globalIndex& globalCells,
780 // const boolList& blockedFace,
781 // const regionSplit& globalRegion,
782 // Map<label>& regionToMaster
783 //) const
784 //{
785 // const polyBoundaryMesh& patches = mesh_.boundaryMesh();
786 //
787 // forAll(patches, patchi)
788 // {
789 // const polyPatch& pp = patches[patchi];
790 //
791 // if (isA<processorPolyPatch>(pp))
792 // {
793 // forAll(pp, i)
794 // {
795 // label facei = pp.start() + i;
796 //
797 // if (!blockedFace[facei])
798 // {
799 // // Only if there is a connection to the neighbour
800 // // will there be a multi-domain region; if not through
801 // // this face then through another.
802 //
803 // label celli = mesh_.faceOwner()[facei];
804 // label globalCelli = globalCells.toGlobal(celli);
805 //
806 // Map<label>::iterator iter =
807 // regionToMaster.find(globalRegion[celli]);
808 //
809 // if (iter != regionToMaster.end())
810 // {
811 // label master = iter();
812 // iter() = min(master, globalCelli);
813 // }
814 // else
815 // {
816 // regionToMaster.insert
817 // (
818 // globalRegion[celli],
819 // globalCelli
820 // );
821 // }
822 // }
823 // }
824 // }
825 // }
826 //
827 // // Do reduction
828 // Pstream::mapCombineGather(regionToMaster, minEqOp<label>());
829 // Pstream::mapCombineScatter(regionToMaster);
830 //}
831 //
832 //
833 //void Foam::meshRefinement::calcLocalRegions
834 //(
835 // const globalIndex& globalCells,
836 // const labelList& globalRegion,
837 // const Map<label>& coupledRegionToMaster,
838 // const scalarField& cellWeights,
839 //
840 // Map<label>& globalToLocalRegion,
841 // pointField& localPoints,
842 // scalarField& localWeights
843 //) const
844 //{
845 // globalToLocalRegion.resize(globalRegion.size());
846 // DynamicList<point> localCc(globalRegion.size()/2);
847 // DynamicList<scalar> localWts(globalRegion.size()/2);
848 //
849 // forAll(globalRegion, celli)
850 // {
851 // Map<label>::const_iterator fndMaster =
852 // coupledRegionToMaster.find(globalRegion[celli]);
853 //
854 // if (fndMaster != coupledRegionToMaster.end())
855 // {
856 // // Multi-processor region.
857 // if (globalCells.toGlobal(celli) == fndMaster())
858 // {
859 // // I am master. Allocate region for me.
860 // globalToLocalRegion.insert
861 // (
862 // globalRegion[celli],
863 // localCc.size()
864 // );
865 // localCc.append(mesh_.cellCentres()[celli]);
866 // localWts.append(cellWeights[celli]);
867 // }
868 // }
869 // else
870 // {
871 // // Single processor region.
872 // if
873 // (
874 // globalToLocalRegion.insert
875 // (
876 // globalRegion[celli],
877 // localCc.size()
878 // )
879 // )
880 // {
881 // localCc.append(mesh_.cellCentres()[celli]);
882 // localWts.append(cellWeights[celli]);
883 // }
884 // }
885 // }
886 //
887 // localPoints.transfer(localCc);
888 // localWeights.transfer(localWts);
889 //
890 // if (localPoints.size() != globalToLocalRegion.size())
891 // {
892 // FatalErrorInFunction
893 // << "localPoints:" << localPoints.size()
894 // << " globalToLocalRegion:" << globalToLocalRegion.size()
895 // << exit(FatalError);
896 // }
897 //}
898 //
899 //
900 //Foam::label Foam::meshRefinement::getShiftedRegion
901 //(
902 // const globalIndex& indexer,
903 // const Map<label>& globalToLocalRegion,
904 // const Map<label>& coupledRegionToShifted,
905 // const label globalRegion
906 //)
907 //{
908 // Map<label>::const_iterator iter =
909 // globalToLocalRegion.find(globalRegion);
910 //
911 // if (iter != globalToLocalRegion.end())
912 // {
913 // // Region is 'owned' locally. Convert local region index into global.
914 // return indexer.toGlobal(iter());
915 // }
916 // else
917 // {
918 // return coupledRegionToShifted[globalRegion];
919 // }
920 //}
921 //
922 //
924 //void Foam::meshRefinement::addUnique(const label elem, labelList& lst)
925 //{
926 // if (findIndex(lst, elem) == -1)
927 // {
928 // label sz = lst.size();
929 // lst.setSize(sz+1);
930 // lst[sz] = elem;
931 // }
932 //}
933 //
934 //
935 //void Foam::meshRefinement::calcRegionRegions
936 //(
937 // const labelList& globalRegion,
938 // const Map<label>& globalToLocalRegion,
939 // const Map<label>& coupledRegionToMaster,
940 // labelListList& regionRegions
941 //) const
942 //{
943 // // Global region indexing since we now know the shifted regions.
944 // globalIndex shiftIndexer(globalToLocalRegion.size());
945 //
946 // // Redo the coupledRegionToMaster to be in shifted region indexing.
947 // Map<label> coupledRegionToShifted(coupledRegionToMaster.size());
948 // forAllConstIter(Map<label>, coupledRegionToMaster, iter)
949 // {
950 // label region = iter.key();
951 //
952 // Map<label>::const_iterator fndRegion = globalToLocalRegion.find
953 // (region);
954 //
955 // if (fndRegion != globalToLocalRegion.end())
956 // {
957 // // A local cell is master of this region. Get its shifted region.
958 // coupledRegionToShifted.insert
959 // (
960 // region,
961 // shiftIndexer.toGlobal(fndRegion())
962 // );
963 // }
964 // Pstream::mapCombineGather(coupledRegionToShifted, minEqOp<label>());
965 // Pstream::mapCombineScatter(coupledRegionToShifted);
966 // }
967 //
968 //
969 // // Determine region-region connectivity.
970 // // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
971 // // This is for all my regions (so my local ones or the ones I am master
972 // // of) the neighbouring region indices.
973 //
974 //
975 // // Transfer lists.
976 // PtrList<HashSet<edge, Hash<edge>>> regionConnectivity
977 // (Pstream::nProcs());
978 // forAll(regionConnectivity, proci)
979 // {
980 // if (proci != Pstream::myProcNo())
981 // {
982 // regionConnectivity.set
983 // (
984 // proci,
985 // new HashSet<edge, Hash<edge>>
986 // (
987 // coupledRegionToShifted.size()
988 // / Pstream::nProcs()
989 // )
990 // );
991 // }
992 // }
993 //
994 //
995 // // Connectivity. For all my local regions the connected regions.
996 // regionRegions.setSize(globalToLocalRegion.size());
997 //
998 // // Add all local connectivity to regionRegions; add all non-local
999 // // to the transferlists.
1000 // for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
1001 // {
1002 // label ownRegion = globalRegion[mesh_.faceOwner()[facei]];
1003 // label neiRegion = globalRegion[mesh_.faceNeighbour()[facei]];
1004 //
1005 // if (ownRegion != neiRegion)
1006 // {
1007 // label shiftOwnRegion = getShiftedRegion
1008 // (
1009 // shiftIndexer,
1010 // globalToLocalRegion,
1011 // coupledRegionToShifted,
1012 // ownRegion
1013 // );
1014 // label shiftNeiRegion = getShiftedRegion
1015 // (
1016 // shiftIndexer,
1017 // globalToLocalRegion,
1018 // coupledRegionToShifted,
1019 // neiRegion
1020 // );
1021 //
1022 //
1023 // // Connection between two regions. Send to owner of region.
1024 // // - is ownRegion 'owned' by me
1025 // // - is neiRegion 'owned' by me
1026 //
1027 // if (shiftIndexer.isLocal(shiftOwnRegion))
1028 // {
1029 // label localI = shiftIndexer.toLocal(shiftOwnRegion);
1030 // addUnique(shiftNeiRegion, regionRegions[localI]);
1031 // }
1032 // else
1033 // {
1034 // label masterProc = shiftIndexer.whichProcID(shiftOwnRegion);
1035 // regionConnectivity[masterProc].insert
1036 // (
1037 // edge(shiftOwnRegion, shiftNeiRegion)
1038 // );
1039 // }
1040 //
1041 // if (shiftIndexer.isLocal(shiftNeiRegion))
1042 // {
1043 // label localI = shiftIndexer.toLocal(shiftNeiRegion);
1044 // addUnique(shiftOwnRegion, regionRegions[localI]);
1045 // }
1046 // else
1047 // {
1048 // label masterProc = shiftIndexer.whichProcID(shiftNeiRegion);
1049 // regionConnectivity[masterProc].insert
1050 // (
1051 // edge(shiftOwnRegion, shiftNeiRegion)
1052 // );
1053 // }
1054 // }
1055 // }
1056 //
1057 //
1058 // // Send
1059 // forAll(regionConnectivity, proci)
1060 // {
1061 // if (proci != Pstream::myProcNo())
1062 // {
1063 // OPstream str(Pstream::commsTypes::blocking, proci);
1064 // str << regionConnectivity[proci];
1065 // }
1066 // }
1067 // // Receive
1068 // forAll(regionConnectivity, proci)
1069 // {
1070 // if (proci != Pstream::myProcNo())
1071 // {
1072 // IPstream str(Pstream::commsTypes::blocking, proci);
1073 // str >> regionConnectivity[proci];
1074 // }
1075 // }
1076 //
1077 // // Add to addressing.
1078 // forAll(regionConnectivity, proci)
1079 // {
1080 // if (proci != Pstream::myProcNo())
1081 // {
1082 // for
1083 // (
1084 // HashSet<edge, Hash<edge>>::const_iterator iter =
1085 // regionConnectivity[proci].begin();
1086 // iter != regionConnectivity[proci].end();
1087 // ++iter
1088 // )
1089 // {
1090 // const edge& e = iter.key();
1091 //
1092 // bool someLocal = false;
1093 // if (shiftIndexer.isLocal(e[0]))
1094 // {
1095 // label localI = shiftIndexer.toLocal(e[0]);
1096 // addUnique(e[1], regionRegions[localI]);
1097 // someLocal = true;
1098 // }
1099 // if (shiftIndexer.isLocal(e[1]))
1100 // {
1101 // label localI = shiftIndexer.toLocal(e[1]);
1102 // addUnique(e[0], regionRegions[localI]);
1103 // someLocal = true;
1104 // }
1105 //
1106 // if (!someLocal)
1107 // {
1108 // FatalErrorInFunction
1109 // << "Received from processor " << proci
1110 // << " connection " << e
1111 // << " where none of the elements is local to me."
1112 // << abort(FatalError);
1113 // }
1114 // }
1115 // }
1116 // }
1117 //}
1118 
1119 
1120 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1121 
1123 (
1124  fvMesh& mesh,
1125  const dictionary& refineDict,
1126  const scalar mergeDistance,
1127  const bool overwrite,
1128  refinementSurfaces& surfaces,
1129  const refinementFeatures& features,
1130  const refinementRegions& shells
1131 )
1132 :
1133  mesh_(mesh),
1134  mergeDistance_(mergeDistance),
1135  overwrite_(overwrite),
1136  oldInstance_(mesh.pointsInstance()),
1137  surfaces_(surfaces),
1138  features_(features),
1139  shells_(shells),
1140  meshCutter_
1141  (
1142  mesh,
1143  false // do not try to read history.
1144  ),
1145  surfaceIndex_
1146  (
1147  IOobject
1148  (
1149  "surfaceIndex",
1150  mesh_.facesInstance(),
1151  fvMesh::meshSubDir,
1152  mesh_,
1153  IOobject::NO_READ,
1154  IOobject::NO_WRITE,
1155  false
1156  ),
1157  labelList(mesh_.nFaces(), -1)
1158  ),
1159  userFaceData_(0)
1160 {
1162  (
1163  shells_,
1164  meshCutter_.level0EdgeLength(),
1165  refineDict.lookupOrDefault<Switch>("extendedRefinementSpan", true)
1166  );
1167 
1168  // recalculate intersections for all faces
1169  updateIntersections(identityMap(mesh_.nFaces()));
1170 }
1171 
1172 
1173 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1174 
1176 {
1177  // Stats on edges to test. Count proc faces only once.
1178  PackedBoolList isMasterFace(syncTools::getMasterFaces(mesh_));
1179 
1180  label nHits = 0;
1181 
1182  forAll(surfaceIndex_, facei)
1183  {
1184  if (surfaceIndex_[facei] >= 0 && isMasterFace.get(facei) == 1)
1185  {
1186  nHits++;
1187  }
1188  }
1189  return nHits;
1190 }
1191 
1192 
1194 //Foam::labelList Foam::meshRefinement::decomposeCombineRegions
1195 //(
1196 // const scalarField& cellWeights,
1197 // const boolList& blockedFace,
1198 // const List<labelPair>& explicitConnections,
1199 // decompositionMethod& decomposer
1200 //) const
1201 //{
1202 // // Determine global regions, separated by blockedFaces
1203 // regionSplit globalRegion(mesh_, blockedFace, explicitConnections);
1204 //
1205 // // Now globalRegion has global region per cell. Problem is that
1206 // // the region might span multiple domains so we want to get
1207 // // a "region master" per domain. Note that multi-processor
1208 // // regions can only occur on cells on coupled patches.
1209 // // Note: since the number of regions does not change by this the
1210 // // process can be seen as just a shift of a region onto a single
1211 // // processor.
1212 //
1213 //
1214 // // Global cell numbering engine
1215 // globalIndex globalCells(mesh_.nCells());
1216 //
1217 //
1218 // // Determine per coupled region the master cell (lowest numbered cell
1219 // // on lowest numbered processor)
1220 // // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1221 // // (does not determine master for non-coupled=fully-local regions)
1222 //
1223 // Map<label> coupledRegionToMaster(mesh_.nFaces() - mesh_.nInternalFaces());
1224 // getCoupledRegionMaster
1225 // (
1226 // globalCells,
1227 // blockedFace,
1228 // globalRegion,
1229 // coupledRegionToMaster
1230 // );
1231 //
1232 // // Determine my regions
1233 // // ~~~~~~~~~~~~~~~~~~~~
1234 // // These are the fully local ones or the coupled ones of which I am master
1235 //
1236 // Map<label> globalToLocalRegion;
1237 // pointField localPoints;
1238 // scalarField localWeights;
1239 // calcLocalRegions
1240 // (
1241 // globalCells,
1242 // globalRegion,
1243 // coupledRegionToMaster,
1244 // cellWeights,
1245 //
1246 // globalToLocalRegion,
1247 // localPoints,
1248 // localWeights
1249 // );
1250 //
1251 //
1252 //
1253 // // Find distribution for regions
1254 // // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1255 //
1256 // labelList regionDistribution;
1257 //
1258 // if (isA<decompositionMethods::geometric>(decomposer))
1259 // {
1260 // regionDistribution = decomposer.decompose(localPoints, localWeights);
1261 // }
1262 // else
1263 // {
1264 // labelListList regionRegions;
1265 // calcRegionRegions
1266 // (
1267 // globalRegion,
1268 // globalToLocalRegion,
1269 // coupledRegionToMaster,
1270 //
1271 // regionRegions
1272 // );
1273 //
1274 // regionDistribution = decomposer.decompose
1275 // (
1276 // regionRegions,
1277 // localPoints,
1278 // localWeights
1279 // );
1280 // }
1281 //
1282 //
1283 //
1284 // // Convert region-based decomposition back to cell-based one
1285 // // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1286 //
1287 // // Transfer destination processor back to all. Use global reduce for now.
1288 // Map<label> regionToDist(coupledRegionToMaster.size());
1289 // forAllConstIter(Map<label>, coupledRegionToMaster, iter)
1290 // {
1291 // label region = iter.key();
1292 //
1293 // Map<label>::const_iterator regionFnd = globalToLocalRegion.find
1294 // (region);
1295 //
1296 // if (regionFnd != globalToLocalRegion.end())
1297 // {
1298 // // Master cell is local. Store my distribution.
1299 // regionToDist.insert(iter.key(), regionDistribution[regionFnd()]);
1300 // }
1301 // else
1302 // {
1303 // // Master cell is not on this processor. Make sure it is
1304 // // overridden.
1305 // regionToDist.insert(iter.key(), labelMax);
1306 // }
1307 // }
1308 // Pstream::mapCombineGather(regionToDist, minEqOp<label>());
1309 // Pstream::mapCombineScatter(regionToDist);
1310 //
1311 //
1312 // // Determine destination for all cells
1313 // labelList distribution(mesh_.nCells());
1314 //
1315 // forAll(globalRegion, celli)
1316 // {
1317 // Map<label>::const_iterator fndRegion =
1318 // regionToDist.find(globalRegion[celli]);
1319 //
1320 // if (fndRegion != regionToDist.end())
1321 // {
1322 // distribution[celli] = fndRegion();
1323 // }
1324 // else
1325 // {
1326 // // region is local to the processor.
1327 // label localRegioni = globalToLocalRegion[globalRegion[celli]];
1328 //
1329 // distribution[celli] = regionDistribution[localRegioni];
1330 // }
1331 // }
1332 //
1333 // return distribution;
1334 //}
1335 
1336 
1338 (
1339  const bool keepZoneFaces,
1340  const bool keepBaffles,
1341  const scalarField& cellWeights,
1342  decompositionMethod& decomposer,
1343  fvMeshDistribute& distributor
1344 )
1345 {
1347 
1348  if (Pstream::parRun())
1349  {
1350  // Wanted distribution
1352 
1353 
1354  // Faces where owner and neighbour are not 'connected' so can
1355  // go to different processors.
1356  boolList blockedFace;
1357  label nUnblocked = 0;
1358 
1359  // Faces that move as block onto single processor
1360  PtrList<labelList> specifiedProcessorFaces;
1361  labelList specifiedProcessor;
1362 
1363  // Pairs of baffles
1364  List<labelPair> couples;
1365 
1366  // Constraints from decomposeParDict
1367  decomposer.setConstraints
1368  (
1369  mesh_,
1370  blockedFace,
1371  specifiedProcessorFaces,
1372  specifiedProcessor,
1373  couples
1374  );
1375 
1376 
1377  if (keepZoneFaces || keepBaffles)
1378  {
1379  if (keepZoneFaces)
1380  {
1381  // Determine decomposition to keep/move surface zones
1382  // on one processor. The reason is that snapping will make these
1383  // into baffles, move and convert them back so if they were
1384  // proc boundaries after baffling&moving the points might be no
1385  // longer synchronised so recoupling will fail. To prevent this
1386  // keep owner&neighbour of such a surface zone on the same
1387  // processor.
1388 
1389  const PtrList<surfaceZonesInfo>& surfZones =
1390  surfaces().surfZones();
1391  const faceZoneList& fZones = mesh_.faceZones();
1392  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
1393 
1394  // Get faces whose owner and neighbour should stay together,
1395  // i.e. they are not 'blocked'.
1396 
1397  forAll(surfZones, surfi)
1398  {
1399  const word& fzName = surfZones[surfi].faceZoneName();
1400 
1401  if (fzName.size())
1402  {
1403  // Get zone
1404  const faceZone& fZone = fZones[fzName];
1405 
1406  forAll(fZone, i)
1407  {
1408  const label facei = fZone[i];
1409  if (blockedFace[facei])
1410  {
1411  if
1412  (
1413  mesh_.isInternalFace(facei)
1414  || pbm[pbm.whichPatch(facei)].coupled()
1415  )
1416  {
1417  blockedFace[facei] = false;
1418  nUnblocked++;
1419  }
1420  }
1421  }
1422  }
1423  }
1424 
1425 
1426  // If the faceZones are not synchronised the blockedFace
1427  // might not be synchronised. If you are sure the faceZones
1428  // are synchronised remove below check.
1430  (
1431  mesh_,
1432  blockedFace,
1433  andEqOp<bool>() // combine operator
1434  );
1435  }
1436  reduce(nUnblocked, sumOp<label>());
1437 
1438  if (keepZoneFaces)
1439  {
1440  Info<< "Found " << nUnblocked
1441  << " zoned faces to keep together." << endl;
1442  }
1443 
1444 
1445  // Extend un-blockedFaces with any cyclics
1446  {
1447  boolList separatedCoupledFace(mesh_.nFaces(), false);
1448  selectSeparatedCoupledFaces(separatedCoupledFace);
1449 
1450  label nSeparated = 0;
1451  forAll(separatedCoupledFace, facei)
1452  {
1453  if (separatedCoupledFace[facei])
1454  {
1455  if (blockedFace[facei])
1456  {
1457  blockedFace[facei] = false;
1458  nSeparated++;
1459  }
1460  }
1461  }
1462  reduce(nSeparated, sumOp<label>());
1463  Info<< "Found " << nSeparated
1464  << " separated coupled faces to keep together." << endl;
1465 
1466  nUnblocked += nSeparated;
1467  }
1468 
1469 
1470  if (keepBaffles)
1471  {
1472  const label nBnd = mesh_.nFaces() - mesh_.nInternalFaces();
1473 
1474  labelList coupledFace(mesh_.nFaces(), -1);
1475  {
1476  // Get boundary baffles that need to stay together
1477  List<labelPair> allCouples
1478  (
1480  );
1481 
1482  // Merge with any couples from
1483  // decompositionMethod::setConstraints
1484  forAll(couples, i)
1485  {
1486  const labelPair& baffle = couples[i];
1487  coupledFace[baffle.first()] = baffle.second();
1488  coupledFace[baffle.second()] = baffle.first();
1489  }
1490  forAll(allCouples, i)
1491  {
1492  const labelPair& baffle = allCouples[i];
1493  coupledFace[baffle.first()] = baffle.second();
1494  coupledFace[baffle.second()] = baffle.first();
1495  }
1496  }
1497 
1498  couples.setSize(nBnd);
1499  label nCpl = 0;
1500  forAll(coupledFace, facei)
1501  {
1502  if (coupledFace[facei] != -1 && facei < coupledFace[facei])
1503  {
1504  couples[nCpl++] = labelPair(facei, coupledFace[facei]);
1505  }
1506  }
1507  couples.setSize(nCpl);
1508  }
1509  label nCouples = returnReduce(couples.size(), sumOp<label>());
1510 
1511  if (keepBaffles)
1512  {
1513  Info<< "Found " << nCouples << " baffles to keep together."
1514  << endl;
1515  }
1516 
1517  // if (nUnblocked > 0 || nCouples > 0)
1518  //{
1519  // Info<< "Applying special decomposition to keep baffles"
1520  // << " and zoned faces together." << endl;
1521  //
1522  // distribution = decomposeCombineRegions
1523  // (
1524  // cellWeights,
1525  // blockedFace,
1526  // couples,
1527  // decomposer
1528  // );
1529  //
1530  // labelList nProcCells(distributor.countCells(distribution));
1531  // Pstream::listCombineGather(nProcCells, plusEqOp<label>());
1532  // Pstream::listCombineScatter(nProcCells);
1533  //
1534  // Info<< "Calculated decomposition:" << endl;
1535  // forAll(nProcCells, proci)
1536  // {
1537  // Info<< " " << proci << '\t' << nProcCells[proci]
1538  // << endl;
1539  // }
1540  // Info<< endl;
1541  //}
1542  // else
1543  //{
1544  // // Normal decomposition
1545  // distribution = decomposer.decompose
1546  // (
1547  // mesh_,
1548  // mesh_.cellCentres(),
1549  // cellWeights
1550  // );
1551  //}
1552  }
1553  // else
1554  //{
1555  // // Normal decomposition
1556  // distribution = decomposer.decompose
1557  // (
1558  // mesh_,
1559  // cellWeights
1560  // );
1561  //}
1562 
1563 
1564  // Make sure blockedFace not set on couples
1565  forAll(couples, i)
1566  {
1567  const labelPair& baffle = couples[i];
1568  blockedFace[baffle.first()] = false;
1569  blockedFace[baffle.second()] = false;
1570  }
1571 
1572  distribution = decomposer.decompose
1573  (
1574  mesh_,
1575  cellWeights,
1576  blockedFace,
1577  specifiedProcessorFaces,
1578  specifiedProcessor,
1579  couples // explicit connections
1580  );
1581 
1582  if (debug)
1583  {
1584  labelList nProcCells(distributor.countCells(distribution));
1585  Pout<< "Wanted distribution:" << nProcCells << endl;
1586 
1588  Pstream::listCombineScatter(nProcCells);
1589 
1590  Pout<< "Wanted resulting decomposition:" << endl;
1591  forAll(nProcCells, proci)
1592  {
1593  Pout<< " " << proci << '\t' << nProcCells[proci] << endl;
1594  }
1595  Pout<< endl;
1596  }
1597 
1598  mesh_.clearOut();
1599 
1600  // Do actual sending/receiving of mesh
1601  map = distributor.distribute(distribution);
1602 
1603  // Update numbering of meshRefiner
1604  distribute(map);
1605 
1606  // Set correct instance (for if overwrite)
1607  mesh_.setInstance(name());
1608  setInstance(mesh_.facesInstance());
1609  }
1610 
1611  return map;
1612 }
1613 
1614 
1616 {
1617  label nBoundaryFaces = 0;
1618 
1619  forAll(surfaceIndex_, facei)
1620  {
1621  if (surfaceIndex_[facei] != -1)
1622  {
1623  nBoundaryFaces++;
1624  }
1625  }
1626 
1627  labelList surfaceFaces(nBoundaryFaces);
1628  nBoundaryFaces = 0;
1629 
1630  forAll(surfaceIndex_, facei)
1631  {
1632  if (surfaceIndex_[facei] != -1)
1633  {
1634  surfaceFaces[nBoundaryFaces++] = facei;
1635  }
1636  }
1637  return surfaceFaces;
1638 }
1639 
1640 
1642 {
1643  const faceList& faces = mesh_.faces();
1644 
1645  // Mark all points on faces that will become baffles
1646  PackedBoolList isBoundaryPoint(mesh_.nPoints());
1647  label nBoundaryPoints = 0;
1648 
1649  forAll(surfaceIndex_, facei)
1650  {
1651  if (surfaceIndex_[facei] != -1)
1652  {
1653  const face& f = faces[facei];
1654 
1655  forAll(f, fp)
1656  {
1657  if (isBoundaryPoint.set(f[fp], 1u))
1658  {
1659  nBoundaryPoints++;
1660  }
1661  }
1662  }
1663  }
1664 
1666  // labelList adaptPatchIDs(meshedPatches());
1667  // forAll(adaptPatchIDs, i)
1668  //{
1669  // label patchi = adaptPatchIDs[i];
1670  //
1671  // if (patchi != -1)
1672  // {
1673  // const polyPatch& pp = mesh_.boundaryMesh()[patchi];
1674  //
1675  // label facei = pp.start();
1676  //
1677  // forAll(pp, i)
1678  // {
1679  // const face& f = faces[facei];
1680  //
1681  // forAll(f, fp)
1682  // {
1683  // if (isBoundaryPoint.set(f[fp], 1u))
1684  // nBoundaryPoints++;
1685  // }
1686  // }
1687  // facei++;
1688  // }
1689  // }
1690  //}
1691 
1692 
1693  // Pack
1694  labelList boundaryPoints(nBoundaryPoints);
1695  nBoundaryPoints = 0;
1696  forAll(isBoundaryPoint, pointi)
1697  {
1698  if (isBoundaryPoint.get(pointi) == 1u)
1699  {
1700  boundaryPoints[nBoundaryPoints++] = pointi;
1701  }
1702  }
1703 
1704  return boundaryPoints;
1705 }
1706 
1707 
1709 (
1710  const polyMesh& mesh,
1711  const labelList& patchIDs
1712 )
1713 {
1714  const polyBoundaryMesh& patches = mesh.boundaryMesh();
1715 
1716  // Count faces.
1717  label nFaces = 0;
1718 
1719  forAll(patchIDs, i)
1720  {
1721  const polyPatch& pp = patches[patchIDs[i]];
1722 
1723  nFaces += pp.size();
1724  }
1725 
1726  // Collect faces.
1727  labelList addressing(nFaces);
1728  nFaces = 0;
1729 
1730  forAll(patchIDs, i)
1731  {
1732  const polyPatch& pp = patches[patchIDs[i]];
1733 
1734  label meshFacei = pp.start();
1735 
1736  forAll(pp, i)
1737  {
1738  addressing[nFaces++] = meshFacei++;
1739  }
1740  }
1741 
1743  (
1745  (
1746  IndirectList<face>(mesh.faces(), addressing),
1747  mesh.points()
1748  )
1749  );
1750 }
1751 
1752 
1754 (
1755  const pointMesh& pMesh,
1756  const labelList& adaptPatchIDs
1757 )
1758 {
1759  // Construct displacement field.
1760  const pointBoundaryMesh& pointPatches = pMesh.boundary();
1761 
1762  wordList patchFieldTypes
1763  (
1764  pointPatches.size(),
1765  slipPointPatchVectorField::typeName
1766  );
1767 
1768  forAll(adaptPatchIDs, i)
1769  {
1770  patchFieldTypes[adaptPatchIDs[i]] =
1771  fixedValuePointPatchVectorField::typeName;
1772  }
1773 
1774  forAll(pointPatches, patchi)
1775  {
1776  if (isA<processorPointPatch>(pointPatches[patchi]))
1777  {
1778  patchFieldTypes[patchi] = calculatedPointPatchVectorField::typeName;
1779  }
1780  else if (isA<cyclicPointPatch>(pointPatches[patchi]))
1781  {
1782  patchFieldTypes[patchi] = cyclicSlipPointPatchVectorField::typeName;
1783  }
1784  }
1785 
1786  // Note: time().name() instead of meshRefinement::name() since
1787  // postprocessable field.
1789  (
1791  (
1792  "pointDisplacement",
1793  pMesh,
1795  patchFieldTypes
1796  )
1797  );
1798  return tfld;
1799 }
1800 
1801 
1803 {
1804  const faceZoneList& fZones = mesh.faceZones();
1805 
1806  // Check any zones are present anywhere and in same order
1807 
1808  {
1809  List<wordList> zoneNames(Pstream::nProcs());
1810  zoneNames[Pstream::myProcNo()] = fZones.toc();
1811  Pstream::gatherList(zoneNames);
1812  Pstream::scatterList(zoneNames);
1813  // All have same data now. Check.
1814  forAll(zoneNames, proci)
1815  {
1816  if (proci != Pstream::myProcNo())
1817  {
1818  if (zoneNames[proci] != zoneNames[Pstream::myProcNo()])
1819  {
1821  << "faceZones are not synchronised on processors." << nl
1822  << "Processor " << proci << " has faceZones "
1823  << zoneNames[proci] << nl
1824  << "Processor " << Pstream::myProcNo()
1825  << " has faceZones "
1826  << zoneNames[Pstream::myProcNo()] << nl
1827  << exit(FatalError);
1828  }
1829  }
1830  }
1831  }
1832 
1833  // Check that coupled faces are present on both sides.
1834 
1835  labelList faceToZone(mesh.nFaces() - mesh.nInternalFaces(), -1);
1836 
1837  forAll(fZones, zonei)
1838  {
1839  const faceZone& fZone = fZones[zonei];
1840 
1841  forAll(fZone, i)
1842  {
1843  const label bFacei = fZone[i] - mesh.nInternalFaces();
1844 
1845  if (bFacei >= 0)
1846  {
1847  if (faceToZone[bFacei] == -1)
1848  {
1849  faceToZone[bFacei] = zonei;
1850  }
1851  else if (faceToZone[bFacei] == zonei)
1852  {
1854  << "Face " << fZone[i] << " in zone "
1855  << fZone.name()
1856  << " is twice in zone!"
1857  << abort(FatalError);
1858  }
1859  }
1860  }
1861  }
1862 
1863  labelList neiFaceToZone(faceToZone);
1864  syncTools::swapBoundaryFaceList(mesh, neiFaceToZone);
1865 
1866  forAll(faceToZone, i)
1867  {
1868  if (faceToZone[i] != neiFaceToZone[i])
1869  {
1871  << "Face " << mesh.nInternalFaces() + i
1872  << " is in zone " << faceToZone[i]
1873  << ", its coupled face is in zone " << neiFaceToZone[i]
1874  << abort(FatalError);
1875  }
1876  }
1877 }
1878 
1879 
1881 (
1882  const polyMesh& mesh,
1883  const PackedBoolList& isMasterEdge,
1884  const labelList& meshPoints,
1885  const edgeList& edges,
1886  scalarField& edgeWeights,
1887  scalarField& invSumWeight
1888 )
1889 {
1890  const pointField& pts = mesh.points();
1891 
1892  // Calculate edgeWeights and inverse sum of edge weights
1893  edgeWeights.setSize(isMasterEdge.size());
1894  invSumWeight.setSize(meshPoints.size());
1895 
1896  forAll(edges, edgei)
1897  {
1898  const edge& e = edges[edgei];
1899  scalar eMag = max
1900  (
1901  small,
1902  mag
1903  (
1904  pts[meshPoints[e[1]]]
1905  - pts[meshPoints[e[0]]]
1906  )
1907  );
1908  edgeWeights[edgei] = 1.0/eMag;
1909  }
1910 
1911  // Sum per point all edge weights
1912  weightedSum
1913  (
1914  mesh,
1915  isMasterEdge,
1916  meshPoints,
1917  edges,
1918  edgeWeights,
1919  scalarField(meshPoints.size(), 1.0), // data
1920  invSumWeight
1921  );
1922 
1923  // Inplace invert
1924  forAll(invSumWeight, pointi)
1925  {
1926  scalar w = invSumWeight[pointi];
1927 
1928  if (w > 0.0)
1929  {
1930  invSumWeight[pointi] = 1.0/w;
1931  }
1932  }
1933 }
1934 
1935 
1937 (
1938  const word& name,
1939  const dictionary& patchInfo
1940 )
1941 {
1942  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
1943 
1944  const label meshedI = findIndex(meshedPatches_, name);
1945 
1946  if (meshedI != -1)
1947  {
1948  // Already there. Get corresponding polypatch
1949  return pbm.findIndex(name);
1950  }
1951  else
1952  {
1953  // Add patch
1954 
1955  label patchi = pbm.findIndex(name);
1956  if (patchi == -1)
1957  {
1958  patchi = pbm.size();
1959  forAll(pbm, i)
1960  {
1961  const polyPatch& pp = pbm[i];
1962 
1963  if (isA<processorPolyPatch>(pp))
1964  {
1965  patchi = i;
1966  break;
1967  }
1968  }
1969 
1970  dictionary patchDict(patchInfo);
1971  patchDict.set("nFaces", 0);
1972  patchDict.set("startFace", 0);
1973  autoPtr<polyPatch> ppPtr
1974  (
1976  (
1977  name,
1978  patchDict,
1979  0,
1980  pbm
1981  )
1982  );
1983  mesh_.addPatch
1984  (
1985  patchi,
1986  ppPtr(),
1987  dictionary(), // optional field values
1989  true // validBoundary
1990  );
1991  }
1992 
1993  // Store
1994  meshedPatches_.append(name);
1995 
1996  return patchi;
1997  }
1998 }
1999 
2000 
2002 {
2003  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
2004 
2005  DynamicList<label> patchIDs(meshedPatches_.size());
2006  forAll(meshedPatches_, i)
2007  {
2008  label patchi = patches.findIndex(meshedPatches_[i]);
2009 
2010  if (patchi == -1)
2011  {
2013  << "Problem : did not find patch " << meshedPatches_[i]
2014  << endl << "Valid patches are " << patches.names()
2015  << abort(FatalError);
2016  }
2018  {
2019  patchIDs.append(patchi);
2020  }
2021  }
2022 
2023  return patchIDs;
2024 }
2025 
2026 
2028 {
2029  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
2030 
2032  {
2033  if (isA<coupledPolyPatch>(patches[patchi]))
2034  {
2035  const coupledPolyPatch& cpp =
2036  refCast<const coupledPolyPatch>(patches[patchi]);
2037 
2038  if (cpp.transform().transformsPosition())
2039  {
2040  forAll(cpp, i)
2041  {
2042  selected[cpp.start() + i] = true;
2043  }
2044  }
2045  }
2046  }
2047 }
2048 
2049 
2051 (
2052  const polyMesh& mesh,
2053  const labelList& cellRegion,
2054  const vector& perturbVec,
2055  const point& location
2056 )
2057 {
2058  label regioni = -1;
2059  label celli = mesh.findCell(location);
2060  if (celli != -1)
2061  {
2062  regioni = cellRegion[celli];
2063  }
2064  reduce(regioni, maxOp<label>());
2065 
2066  if (regioni == -1)
2067  {
2068  // See if we can perturb a bit
2069  celli = mesh.findCell(location + perturbVec);
2070  if (celli != -1)
2071  {
2072  regioni = cellRegion[celli];
2073  }
2074  reduce(regioni, maxOp<label>());
2075  }
2076 
2077  return regioni;
2078 }
2079 
2080 
2082 (
2083  const polyMesh& mesh,
2084  labelList& cellRegion,
2085  const vector& perturbVec,
2086  const refinementParameters::cellSelectionPoints& selectionPoints
2087 )
2088 {
2089  // List of all cells selected cells
2090  PackedBoolList selectedCells(mesh.nCells());
2091 
2092  if (selectionPoints.outside().size())
2093  {
2094  // Select all cells before removing those in regions
2095  // containing locations in selectionPoints.outside()
2096  selectedCells = true;
2097 
2098  // For each of the selectionPoints.outside() find the corresponding
2099  // region and deselect the cells
2100  forAll(selectionPoints.outside(), i)
2101  {
2102  // Find the region corresponding to the selectionPoints.outside()[i]
2103  const label regioni = findRegion
2104  (
2105  mesh,
2106  cellRegion,
2107  perturbVec,
2108  selectionPoints.outside()[i]
2109  );
2110 
2111  // Deselect the cells in the region from selectedCells
2112  forAll(cellRegion, celli)
2113  {
2114  if (cellRegion[celli] == regioni)
2115  {
2116  selectedCells[celli] = false;
2117  }
2118  }
2119  }
2120  }
2121 
2122  // For each of the selectionPoints.inside() find the corresponding region
2123  // and select the cells
2124  forAll(selectionPoints.inside(), i)
2125  {
2126  // Find the region corresponding to the selectionPoints.inside()[i]
2127  const label regioni = findRegion
2128  (
2129  mesh,
2130  cellRegion,
2131  perturbVec,
2132  selectionPoints.inside()[i]
2133  );
2134 
2135  // Add all the cells in the region to selectedCells
2136  forAll(cellRegion, celli)
2137  {
2138  if (cellRegion[celli] == regioni)
2139  {
2140  selectedCells[celli] = true;
2141  }
2142  }
2143  }
2144 
2145  // For all unmarked cells set cellRegion to -1
2146  forAll(selectedCells, celli)
2147  {
2148  if (!selectedCells[celli])
2149  {
2150  cellRegion[celli] = -1;
2151  }
2152  }
2153 }
2154 
2155 
2157 (
2158  const labelList& globalToMasterPatch,
2159  const labelList& globalToSlavePatch,
2160  const refinementParameters::cellSelectionPoints& selectionPoints
2161 )
2162 {
2163  // Force calculation of face decomposition (used in findCell)
2164  (void)mesh_.tetBasePtIs();
2165 
2166  // Determine connected regions. regionSplit is the labelList with the
2167  // region per cell.
2168 
2169  boolList blockedFace(mesh_.nFaces(), false);
2170  selectSeparatedCoupledFaces(blockedFace);
2171 
2172  regionSplit cellRegion(mesh_, blockedFace);
2173 
2174  findRegions
2175  (
2176  mesh_,
2177  cellRegion,
2178  mergeDistance_*vector::one,
2179  selectionPoints
2180  );
2181 
2182  // Subset
2183  // ~~~~~~
2184 
2185  // Get cells to remove
2186  DynamicList<label> cellsToRemove(mesh_.nCells());
2187  forAll(cellRegion, celli)
2188  {
2189  if (cellRegion[celli] == -1)
2190  {
2191  cellsToRemove.append(celli);
2192  }
2193  }
2194  cellsToRemove.shrink();
2195 
2196  const label nCellsToRemove = returnReduce
2197  (
2198  cellsToRemove.size(),
2199  sumOp<label>()
2200  );
2201 
2202  if (nCellsToRemove)
2203  {
2204  label nCellsToKeep = mesh_.nCells() - cellsToRemove.size();
2205  reduce(nCellsToKeep, sumOp<label>());
2206 
2207  Info<< "Keeping all cells in regions containing any point in "
2208  << selectionPoints.inside() << endl
2209  << "Selected for keeping : " << nCellsToKeep << " cells." << endl;
2210 
2211  // Remove cells
2212  removeCells cellRemover(mesh_);
2213 
2214  const labelList exposedFaces
2215  (
2216  cellRemover.getExposedFaces(cellsToRemove)
2217  );
2218 
2219  labelList exposedPatch;
2220 
2221  const label nExposedFaces =
2222  returnReduce(exposedFaces.size(), sumOp<label>());
2223 
2224  if (nExposedFaces)
2225  {
2226  // Patch for exposed faces for lack of anything sensible.
2227  label defaultPatch = 0;
2228  if (globalToMasterPatch.size())
2229  {
2230  defaultPatch = globalToMasterPatch[0];
2231  }
2232 
2234  << "Removing non-reachable cells exposes "
2235  << nExposedFaces << " internal or coupled faces." << endl
2236  << " These get put into patch " << defaultPatch << endl;
2237 
2238  exposedPatch.setSize(exposedFaces.size(), defaultPatch);
2239  }
2240 
2241  return doRemoveCells
2242  (
2243  cellsToRemove,
2244  exposedFaces,
2245  exposedPatch,
2246  cellRemover
2247  );
2248  }
2249  else
2250  {
2251  return autoPtr<polyTopoChangeMap>();
2252  }
2253 }
2254 
2255 
2257 {
2258  // mesh_ already distributed; distribute my member data
2259 
2260  // surfaceQueries_ ok.
2261 
2262  // refinement
2263  meshCutter_.distribute(map);
2264 
2265  // surfaceIndex is face data.
2266  map.distributeFaceData(surfaceIndex_);
2267 
2268  // maintainedFaces are indices of faces.
2269  forAll(userFaceData_, i)
2270  {
2271  map.distributeFaceData(userFaceData_[i].second());
2272  }
2273 
2274  // Redistribute surface and any fields on it.
2275  {
2276  // Get local mesh bounding box. Single box for now.
2278  treeBoundBox& bb = meshBb[0];
2279  bb = treeBoundBox(mesh_.points()).extend(1e-4);
2280 
2281  // Distribute all geometry (so refinementSurfaces and refinementRegions)
2282  searchableSurfaces& geometry =
2283  const_cast<searchableSurfaces&>(surfaces_.geometry());
2284 
2285  forAll(geometry, i)
2286  {
2288  autoPtr<distributionMap> pointMap;
2289  geometry[i].distribute
2290  (
2291  meshBb,
2292  false, // do not keep outside triangles
2293  faceMap,
2294  pointMap
2295  );
2296 
2297  if (faceMap.valid())
2298  {
2299  // (ab)use the instance() to signal current modification time
2300  geometry[i].instance() = geometry[i].time().name();
2301  }
2302 
2303  faceMap.clear();
2304  pointMap.clear();
2305  }
2306  }
2307 }
2308 
2309 
2311 (
2312  const polyTopoChangeMap& map,
2313  const labelList& changedFaces
2314 )
2315 {
2316  Map<label> dummyMap(0);
2317 
2318  topoChange(map, changedFaces, dummyMap, dummyMap, dummyMap);
2319 }
2320 
2321 
2323 (
2324  const labelList& pointsToStore,
2325  const labelList& facesToStore,
2326  const labelList& cellsToStore
2327 )
2328 {
2329  // For now only meshCutter has storable/retrievable data.
2330  meshCutter_.storeData
2331  (
2332  pointsToStore,
2333  facesToStore,
2334  cellsToStore
2335  );
2336 }
2337 
2338 
2340 (
2341  const polyTopoChangeMap& map,
2342  const labelList& changedFaces,
2343  const Map<label>& pointsToRestore,
2344  const Map<label>& facesToRestore,
2345  const Map<label>& cellsToRestore
2346 )
2347 {
2348  // For now only meshCutter has storable/retrievable data.
2349 
2350  // Update numbering of cells/vertices.
2351  meshCutter_.topoChange
2352  (
2353  map,
2354  pointsToRestore,
2355  facesToRestore,
2356  cellsToRestore
2357  );
2358 
2359  // Update surfaceIndex
2360  updateList(map.faceMap(), label(-1), surfaceIndex_);
2361 
2362  // Update cached intersection information
2363  updateIntersections(changedFaces);
2364 
2365  // Update maintained faces
2366  forAll(userFaceData_, i)
2367  {
2368  labelList& data = userFaceData_[i].second();
2369 
2370  if (userFaceData_[i].first() == KEEPALL)
2371  {
2372  // extend list with face-from-face data
2373  updateList(map.faceMap(), label(-1), data);
2374  }
2375  else if (userFaceData_[i].first() == MASTERONLY)
2376  {
2377  // keep master only
2378  labelList newFaceData(map.faceMap().size(), -1);
2379 
2380  forAll(newFaceData, facei)
2381  {
2382  label oldFacei = map.faceMap()[facei];
2383 
2384  if (oldFacei >= 0 && map.reverseFaceMap()[oldFacei] == facei)
2385  {
2386  newFaceData[facei] = data[oldFacei];
2387  }
2388  }
2389  data.transfer(newFaceData);
2390  }
2391  else
2392  {
2393  // remove any face that has been refined i.e. referenced more than
2394  // once.
2395 
2396  // 1. Determine all old faces that get referenced more than once.
2397  // These get marked with -1 in reverseFaceMap
2398  labelList reverseFaceMap(map.reverseFaceMap());
2399 
2400  forAll(map.faceMap(), facei)
2401  {
2402  const label oldFacei = map.faceMap()[facei];
2403 
2404  if (oldFacei >= 0)
2405  {
2406  if (reverseFaceMap[oldFacei] != facei)
2407  {
2408  // facei is slave face. Mark old face.
2409  reverseFaceMap[oldFacei] = -1;
2410  }
2411  }
2412  }
2413 
2414  // 2. Map only faces with intact reverseFaceMap
2415  labelList newFaceData(map.faceMap().size(), -1);
2416  forAll(newFaceData, facei)
2417  {
2418  const label oldFacei = map.faceMap()[facei];
2419 
2420  if (oldFacei >= 0)
2421  {
2422  if (reverseFaceMap[oldFacei] == facei)
2423  {
2424  newFaceData[facei] = data[oldFacei];
2425  }
2426  }
2427  }
2428  data.transfer(newFaceData);
2429  }
2430  }
2431 }
2432 
2433 
2435 {
2436  // Set correct instance (for if overwrite)
2437  mesh_.setInstance(name());
2438  // setInstance(mesh_.facesInstance());
2439 
2440  bool writeOk = mesh_.write();
2441 
2442  // Make sure that any distributed surfaces (so ones which probably have
2443  // been changed) get written as well.
2444  // Note: should ideally have some 'modified' flag to say whether it
2445  // has been changed or not.
2446  searchableSurfaces& geometry =
2447  const_cast<searchableSurfaces&>(surfaces_.geometry());
2448 
2449  forAll(geometry, i)
2450  {
2451  searchableSurface& s = geometry[i];
2452 
2453  // Check if instance() of surface is not constant or system.
2454  // Is good hint that surface is distributed.
2455  if
2456  (
2457  s.instance() != s.time().system()
2458  && s.instance() != s.time().caseSystem()
2459  && s.instance() != s.time().constant()
2460  && s.instance() != s.time().caseConstant()
2461  )
2462  {
2463  // Make sure it gets written to current time, not constant.
2464  s.instance() = s.time().name();
2465  writeOk = writeOk && s.write();
2466  }
2467  }
2468 
2469  return writeOk;
2470 }
2471 
2472 
2474 (
2475  const polyMesh& mesh,
2476  const labelList& meshPoints
2477 )
2478 {
2479  const globalIndex globalPoints(meshPoints.size());
2480 
2481  labelList myPoints(meshPoints.size());
2482  forAll(meshPoints, pointi)
2483  {
2484  myPoints[pointi] = globalPoints.toGlobal(pointi);
2485  }
2486 
2488  (
2489  mesh,
2490  meshPoints,
2491  myPoints,
2492  minEqOp<label>(),
2493  labelMax
2494  );
2495 
2496 
2497  PackedBoolList isPatchMasterPoint(meshPoints.size());
2498  forAll(meshPoints, pointi)
2499  {
2500  if (myPoints[pointi] == globalPoints.toGlobal(pointi))
2501  {
2502  isPatchMasterPoint[pointi] = true;
2503  }
2504  }
2505 
2506  return isPatchMasterPoint;
2507 }
2508 
2509 
2511 (
2512  const polyMesh& mesh,
2513  const labelList& meshEdges
2514 )
2515 {
2516  const globalIndex globalEdges(meshEdges.size());
2517 
2518  labelList myEdges(meshEdges.size());
2519  forAll(meshEdges, edgei)
2520  {
2521  myEdges[edgei] = globalEdges.toGlobal(edgei);
2522  }
2523 
2525  (
2526  mesh,
2527  meshEdges,
2528  myEdges,
2529  minEqOp<label>(),
2530  labelMax
2531  );
2532 
2533 
2534  PackedBoolList isMasterEdge(meshEdges.size());
2535  forAll(meshEdges, edgei)
2536  {
2537  if (myEdges[edgei] == globalEdges.toGlobal(edgei))
2538  {
2539  isMasterEdge[edgei] = true;
2540  }
2541  }
2542 
2543  return isMasterEdge;
2544 }
2545 
2546 
2547 void Foam::meshRefinement::printMeshInfo(const bool debug, const string& msg)
2548 const
2549 {
2550  const globalMeshData& pData = mesh_.globalData();
2551 
2552  if (debug)
2553  {
2554  Pout<< msg.c_str()
2555  << " : cells(local):" << mesh_.nCells()
2556  << " faces(local):" << mesh_.nFaces()
2557  << " points(local):" << mesh_.nPoints()
2558  << endl;
2559  }
2560 
2561  {
2562  PackedBoolList isMasterFace(syncTools::getMasterFaces(mesh_));
2563  label nMasterFaces = 0;
2564  forAll(isMasterFace, i)
2565  {
2566  if (isMasterFace[i])
2567  {
2568  nMasterFaces++;
2569  }
2570  }
2571 
2572  PackedBoolList isMeshMasterPoint(syncTools::getMasterPoints(mesh_));
2573  label nMasterPoints = 0;
2574  forAll(isMeshMasterPoint, i)
2575  {
2576  if (isMeshMasterPoint[i])
2577  {
2578  nMasterPoints++;
2579  }
2580  }
2581 
2582  Info<< msg.c_str()
2583  << " : cells:" << pData.nTotalCells()
2584  << " faces:" << returnReduce(nMasterFaces, sumOp<label>())
2585  << " points:" << returnReduce(nMasterPoints, sumOp<label>())
2586  << endl;
2587  }
2588 
2589 
2590  // if (debug)
2591  {
2592  const labelList& cellLevel = meshCutter_.cellLevel();
2593 
2594  labelList nCells(gMax(cellLevel) + 1, 0);
2595 
2596  forAll(cellLevel, celli)
2597  {
2598  nCells[cellLevel[celli]]++;
2599  }
2600 
2603 
2604  Info<< "Cells per refinement level:" << endl;
2605  forAll(nCells, leveli)
2606  {
2607  Info<< " " << leveli << '\t' << nCells[leveli]
2608  << endl;
2609  }
2610  }
2611 }
2612 
2613 
2615 {
2616  if (overwrite_ && mesh_.time().timeIndex() == 0)
2617  {
2618  return oldInstance_;
2619  }
2620  else
2621  {
2622  return mesh_.time().name();
2623  }
2624 }
2625 
2626 
2628 {
2629  // Note: use time().name(), not meshRefinement::name()
2630  // so as to dump the fields to 0, not to constant.
2631  {
2632  volScalarField volRefLevel
2633  (
2634  IOobject
2635  (
2636  "cellLevel",
2637  mesh_.time().name(),
2638  mesh_,
2641  false
2642  ),
2643  mesh_,
2645  );
2646 
2647  const labelList& cellLevel = meshCutter_.cellLevel();
2648 
2649  forAll(volRefLevel, celli)
2650  {
2651  volRefLevel[celli] = cellLevel[celli];
2652  }
2653 
2654  volRefLevel.write();
2655  }
2656 
2657  // Dump pointLevel
2658  {
2659  const pointMesh& pMesh = pointMesh::New(mesh_);
2660 
2661  pointScalarField pointRefLevel
2662  (
2663  IOobject
2664  (
2665  "pointLevel",
2666  mesh_.time().name(),
2667  mesh_,
2670  false
2671  ),
2672  pMesh,
2674  );
2675 
2676  const labelList& pointLevel = meshCutter_.pointLevel();
2677 
2678  forAll(pointRefLevel, pointi)
2679  {
2680  pointRefLevel[pointi] = pointLevel[pointi];
2681  }
2682 
2683  pointRefLevel.write();
2684  }
2685 }
2686 
2687 
2689 {
2690  {
2691  const pointField& cellCentres = mesh_.cellCentres();
2692 
2693  OFstream str(prefix + "_edges.obj");
2694  label vertI = 0;
2695  Pout<< "meshRefinement::dumpIntersections :"
2696  << " Writing cellcentre-cellcentre intersections to file "
2697  << str.name() << endl;
2698 
2699 
2700  // Redo all intersections
2701  // ~~~~~~~~~~~~~~~~~~~~~~
2702 
2703  // Get boundary face centre and level. Coupled aware.
2704  labelList neiLevel(mesh_.nFaces() - mesh_.nInternalFaces());
2705  pointField neiCc(mesh_.nFaces() - mesh_.nInternalFaces());
2706  calcNeighbourData(neiLevel, neiCc);
2707 
2708  labelList intersectionFaces(intersectedFaces());
2709 
2710  // Collect segments we want to test for
2711  pointField start(intersectionFaces.size());
2712  pointField end(intersectionFaces.size());
2713 
2714  forAll(intersectionFaces, i)
2715  {
2716  const label facei = intersectionFaces[i];
2717  start[i] = cellCentres[mesh_.faceOwner()[facei]];
2718 
2719  if (mesh_.isInternalFace(facei))
2720  {
2721  end[i] = cellCentres[mesh_.faceNeighbour()[facei]];
2722  }
2723  else
2724  {
2725  end[i] = neiCc[facei-mesh_.nInternalFaces()];
2726  }
2727  }
2728 
2729  // Extend segments a bit
2730  {
2731  const vectorField smallVec(rootSmall*(end-start));
2732  start -= smallVec;
2733  end += smallVec;
2734  }
2735 
2736 
2737  // Do tests in one go
2738  labelList surfaceHit;
2739  List<pointIndexHit> surfaceHitInfo;
2740  surfaces_.findAnyIntersection
2741  (
2742  start,
2743  end,
2744  surfaceHit,
2745  surfaceHitInfo
2746  );
2747 
2748  forAll(intersectionFaces, i)
2749  {
2750  if (surfaceHit[i] != -1)
2751  {
2752  meshTools::writeOBJ(str, start[i]);
2753  vertI++;
2754  meshTools::writeOBJ(str, surfaceHitInfo[i].hitPoint());
2755  vertI++;
2756  meshTools::writeOBJ(str, end[i]);
2757  vertI++;
2758  str << "l " << vertI-2 << ' ' << vertI-1 << nl
2759  << "l " << vertI-1 << ' ' << vertI << nl;
2760  }
2761  }
2762  }
2763 
2764  Pout<< endl;
2765 }
2766 
2767 
2769 (
2770  const debugType debugFlags,
2771  const writeType writeFlags,
2772  const fileName& prefix
2773 ) const
2774 {
2775  if (writeFlags & WRITEMESH)
2776  {
2777  write();
2778  }
2779 
2780  if (writeFlags && !(writeFlags & NOWRITEREFINEMENT))
2781  {
2782  meshCutter_.write();
2783  surfaceIndex_.write();
2784  }
2785 
2786  if (writeFlags & WRITELEVELS)
2787  {
2788  dumpRefinementLevel();
2789  }
2790 
2791  if (debugFlags & OBJINTERSECTIONS && prefix.size())
2792  {
2793  dumpIntersections(prefix);
2794  }
2795 }
2796 
2797 
2799 {
2800  return writeLevel_;
2801 }
2802 
2803 
2805 {
2806  writeLevel_ = flags;
2807 }
2808 
2809 
2811 {
2812  return outputLevel_;
2813 }
2814 
2815 
2817 {
2818  outputLevel_ = flags;
2819 }
2820 
2821 
2822 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
static pointMesh & New(const word &name, const polyMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
wordList toc() const
Return the table of contents.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:252
SubField< vector > subField
Declare type of subField.
Definition: Field.H:101
Generic GeometricField class.
static tmp< GeometricField< Type, pointPatchField, pointMesh > > New(const word &name, const Internal &, const PtrList< pointPatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
A List with indirect addressing.
Definition: IndirectList.H:105
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:54
Output to file stream.
Definition: OFstream.H:86
const fileName & name() const
Return the name of the stream.
Definition: OFstream.H:120
A bit-packed bool list.
void set(const PackedList< 1 > &)
Set specified bits.
label size() const
Number of entries.
Definition: PackedListI.H:711
unsigned int get(const label) const
Get value at index I.
Definition: PackedListI.H:954
const Type & second() const
Return second.
Definition: Pair.H:110
const Type & first() const
Return first.
Definition: Pair.H:98
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
A List obtained as a section of another List.
Definition: SubList.H:56
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
label fcIndex(const label i) const
Return the forward circular index, i.e. the next index.
Definition: UListI.H:58
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
static const Form zero
Definition: VectorSpace.H:118
static const Form one
Definition: VectorSpace.H:119
const word & name() const
Return name.
Definition: Zone.H:184
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
void clear()
Delete object (if the pointer is valid) and set pointer to.
Definition: autoPtrI.H:126
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
virtual const transformer & transform() const =0
Return transformation between the coupled patches.
Abstract base class for decomposition.
virtual labelList decompose(const pointField &points, const scalarField &pointWeights)
Return for every coordinate the wanted processor number.
void setConstraints(const polyMesh &mesh, boolList &blockedFace, PtrList< labelList > &specifiedProcessorFaces, labelList &specifiedProcessor, List< labelPair > &explicitConnections)
Helper: extract constraints:
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T, if not found return the given default.
void set(entry *)
Assign a new entry, overwrite any existing entry.
Definition: dictionary.C:1152
Base class for statistical distributions.
Definition: distribution.H:76
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
Definition: ops.H:70
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:59
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
A class for handling file names.
Definition: fileName.H:82
Sends/receives parts of mesh+fvfields to neighbouring processors. Used in load balancing.
autoPtr< polyDistributionMap > distribute(const labelList &dist)
Send cells to neighbours according to distribution.
static labelList countCells(const labelList &)
Helper function: count cells per processor in wanted distribution.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:64
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
label nTotalCells() const
Return total number of cells in decomposed mesh.
Calculates points shared by more than two processor patches or cyclic patches.
Definition: globalPoints.H:101
scalar level0EdgeLength() const
Typical edge length between unrefined points.
Definition: hexRef8.H:413
const labelIOList & cellLevel() const
Definition: hexRef8.H:397
label faceLevel(const label facei) const
Gets level such that the face has four points <= level.
Definition: hexRef8.C:675
static labelList findDuplicateFaces(const primitiveMesh &, const labelList &)
Helper routine to find baffles (two boundary faces using the.
static labelPairList findDuplicateFacePairs(const polyMesh &)
Helper routine to find all baffles (two boundary faces.
Helper class which maintains intersections of (changing) mesh with (static) surfaces.
autoPtr< polyTopoChangeMap > splitMeshRegions(const labelList &globalToMasterPatch, const labelList &globalToSlavePatch, const refinementParameters::cellSelectionPoints &selectionPoints)
Split mesh according to selectionPoints.
labelList intersectedFaces() const
Get faces with intersection.
labelList intersectedPoints() const
Get points on surfaces with intersection and boundary faces.
IOwriteType
Enumeration for what to write.
meshRefinement(fvMesh &mesh, const dictionary &refineDict, const scalar mergeDistance, const bool overwrite, refinementSurfaces &, const refinementFeatures &, const refinementRegions &)
Construct from components.
static void testSyncPointList(const string &msg, const polyMesh &mesh, const List< scalar > &fld)
void checkData()
Debugging: check that all faces still obey start()>end()
const refinementSurfaces & surfaces() const
Reference to surface search engines.
void topoChange(const polyTopoChangeMap &, const labelList &changedFaces)
Update for external change to mesh. changedFaces are in new mesh.
word name() const
Replacement for Time::name() : return oldInstance (if.
static const NamedEnum< IOwriteType, 5 > IOwriteTypeNames
void printMeshInfo(const bool, const string &) const
Print some mesh stats.
void distribute(const polyDistributionMap &)
Update local numbering for mesh redistribution.
void dumpRefinementLevel() const
Write refinement level as volScalarFields for postprocessing.
static void calculateEdgeWeights(const polyMesh &mesh, const PackedBoolList &isMasterEdge, const labelList &meshPoints, const edgeList &edges, scalarField &edgeWeights, scalarField &invSumWeight)
Helper: calculate edge weights (1/length)
label countHits() const
Count number of intersections (local)
static const NamedEnum< IOdebugType, 5 > IOdebugTypeNames
labelList meshedPatches() const
Get patchIDs for patches added in addMeshedPatch.
static void checkCoupledFaceZones(const polyMesh &)
Helper function: check that face zones are synced.
label addMeshedPatch(const word &name, const dictionary &)
Add patch originating from meshing. Update meshedPatches_.
void selectSeparatedCoupledFaces(boolList &) const
Select coupled faces that are not collocated.
static tmp< pointVectorField > makeDisplacementField(const pointMesh &pMesh, const labelList &adaptPatchIDs)
Helper function to make a pointVectorField with correct.
static const NamedEnum< IOoutputType, 1 > IOoutputTypeNames
static PackedBoolList getMasterEdges(const polyMesh &mesh, const labelList &meshEdges)
Determine master edge for subset of edges. If coupled.
static PackedBoolList getMasterPoints(const polyMesh &mesh, const labelList &meshPoints)
Determine master point for subset of points. If coupled.
autoPtr< polyDistributionMap > balance(const bool keepZoneFaces, const bool keepBaffles, const scalarField &cellWeights, decompositionMethod &decomposer, fvMeshDistribute &distributor)
Redecompose according to cell count.
static autoPtr< indirectPrimitivePatch > makePatch(const polyMesh &, const labelList &)
Create patch from set of patches.
void storeData(const labelList &pointsToStore, const labelList &facesToStore, const labelList &cellsToStore)
Signal points/face/cells for which to store data.
bool write() const
Write mesh and all data.
IOoutputType
Enumeration for what to output.
void setInstance(const fileName &)
Set instance of all local IOobjects.
autoPtr< polyTopoChangeMap > splitFaces(const labelList &splitFaces, const labelPairList &splits)
Split faces into two.
static label findRegion(const polyMesh &, const labelList &cellRegion, const vector &perturbVec, const point &location)
Find region point is in. Uses optional perturbation to re-test.
static writeType writeLevel()
Get/set write level.
IOdebugType
Enumeration for what to debug.
static void findRegions(const polyMesh &, labelList &cellRegion, const vector &perturbVec, const refinementParameters::cellSelectionPoints &selectionPoints)
Find regions points are in and update cellRegion.
void dumpIntersections(const fileName &prefix) const
Debug: Write intersection information to OBJ format.
static outputType outputLevel()
Get/set output level.
Foam::pointBoundaryMesh.
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:53
const pointBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: pointMesh.H:125
Foam::polyBoundaryMesh.
label findIndex(const word &patchName) const
Find patch index given a name.
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
void distributeFaceData(List< T > &lst) const
Distribute list of face data.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1326
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1727
const faceZoneList & faceZones() const
Return face zones.
Definition: polyMesh.H:446
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1313
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
static autoPtr< polyPatch > New(const word &patchType, const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm)
Return a pointer to a new patch created on freestore from.
Definition: polyPatchNew.C:32
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:280
static bool constraintType(const word &pt)
Return true if the given type is a constraint type.
Definition: polyPatch.C:238
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const labelList & reverseFaceMap() const
Reverse face map.
const labelList & faceMap() const
Old face map.
Direct mesh changes based on v1.3 polyTopoChange syntax.
autoPtr< polyTopoChangeMap > changeMesh(polyMesh &mesh, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Inplace changes mesh without change of patches.
void modifyFace(const face &f, const label facei, const label own, const label nei, const bool flipFaceFlux, const label patchID)
Modify vertices or cell of face.
label addFace(const face &f, const label own, const label nei, const label masterFaceID, const bool flipFaceFlux, const label patchID)
Add face to cells and return new face index.
label nCells() const
label nPoints() const
const vectorField & cellCentres() const
label nInternalFaces() const
label nFaces() const
Encapsulates queries for features.
Class to hold the points to select cells inside and outside.
const List< point > & outside() const
Return the points outside the surface region to deselect cells.
const List< point > & inside() const
Return the points inside the surface regions to selected cells.
Encapsulates queries for volume refinement ('refine all cells within shell').
Container for data on surfaces used for surface-driven refinement. Contains all the data about the le...
void setMinLevelFields(const refinementRegions &shells, const scalar level0EdgeLength, const bool extendedRefinementSpan)
Calculate the refinement level for every element.
virtual bool write(const bool write=true) const
Write using setting from DB.
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Definition: regionSplit.H:118
Given list of cells to remove insert all the topology changes.
Definition: removeCells.H:60
void setRefinement(const labelList &cellsToRemove, const labelList &facesToExpose, const labelList &patchIDs, polyTopoChange &) const
Play commands into polyTopoChange to remove cells.
Definition: removeCells.C:175
void topoChange(const polyTopoChangeMap &)
Force recalculation of locally stored data on topological change.
Definition: removeCells.H:115
labelList getExposedFaces(const labelList &cellsToRemove) const
Get labels of exposed faces.
Definition: removeCells.C:69
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
Container for searchableSurfaces.
static void syncPointList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronise values on all mesh points.
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &l)
Swap coupled boundary face values.
Definition: syncTools.H:436
static PackedBoolList getMasterPoints(const polyMesh &)
Get per point whether it is uncoupled or a master of a.
Definition: syncTools.C:65
static PackedBoolList getMasterFaces(const polyMesh &)
Get per face whether it is uncoupled or a master of a.
Definition: syncTools.C:153
static void syncEdgeList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronise values on all mesh edges.
static void syncBoundaryFacePositions(const polyMesh &mesh, UList< point > &l, const CombineOp &cop)
Synchronise locations on boundary faces only.
Definition: syncTools.H:369
static void syncFaceList(const polyMesh &mesh, UList< T > &l, const CombineOp &cop)
Synchronise values on all mesh faces.
Definition: syncTools.H:387
static void swapBoundaryFacePositions(const polyMesh &mesh, UList< point > &l)
Swap coupled positions.
Definition: syncTools.H:452
A class for managing temporary objects.
Definition: tmp.H:55
bool transformsPosition() const
Return true if the transformer transforms a point.
Definition: transformerI.H:146
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:90
treeBoundBox extend(const scalar s) const
Return asymmetrically extended bounding box, with guaranteed.
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(Cloud< passiveParticle >, particles, iter) { positions[particlei++]=iter().position(mesh);} for(i=0;i< pTraits< point >::nComponents;i++) { forAll(positions, particlei) { gmvFile<< component(positions[particlei], i)<< ' ';} gmvFile<< nl;}}forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const fvPatchList & patches
#define WarningInFunction
Report a warning using Foam::Warning.
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
const doubleScalar e
Definition: doubleScalar.H:106
List< label > labelList
A List of labels.
Definition: labelList.H:56
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimless
ListType renumber(const labelUList &oldToNew, const ListType &)
Renumber the values (not the indices) of a list.
messageStream Info
const dimensionSet dimLength
labelList second(const UList< labelPair > &p)
Definition: patchToPatch.C:49
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
vector point
Point is a vector.
Definition: point.H:41
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Field< vector > vectorField
Specialisation of Field<T> for vector.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:213
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
error FatalError
static const label labelMax
Definition: label.H:62
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
UList< label > labelUList
Definition: UList.H:65
static const char nl
Definition: Ostream.H:266
Type gMax(const FieldField< Field, Type > &f)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
labelList f(nPoints)
List< treeBoundBox > meshBb(1, treeBoundBox(boundBox(coarseMesh.points(), false)).extend(1e-3))