meshRefinementRefine.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2015-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "meshRefinement.H"
30 #include "trackedParticle.H"
31 #include "syncTools.H"
32 #include "Time.H"
33 #include "refinementSurfaces.H"
34 #include "refinementFeatures.H"
35 #include "shellSurfaces.H"
36 #include "faceSet.H"
37 #include "decompositionMethod.H"
38 #include "fvMeshDistribute.H"
39 #include "polyTopoChange.H"
40 #include "mapDistributePolyMesh.H"
41 #include "Cloud.H"
42 //#include "globalIndex.H"
43 #include "cellSet.H"
44 #include "treeDataCell.H"
45 
46 #include "cellCuts.H"
47 #include "refineCell.H"
48 #include "hexCellLooper.H"
49 #include "meshCutter.H"
50 
51 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
52 
53 namespace Foam
54 {
55  //- To compare normals
56  static bool less(const vector& x, const vector& y)
57  {
58  for (direction i = 0; i < vector::nComponents; i++)
59  {
60  if (x[i] < y[i])
61  {
62  return true;
63  }
64  else if (x[i] > y[i])
65  {
66  return false;
67  }
68  }
69  // All components the same
70  return false;
71  }
72 
73 
74  //- To compare normals
75  class normalLess
76  {
77  const vectorList& values_;
78 
79  public:
80 
82  :
83  values_(values)
84  {}
85 
86  bool operator()(const label a, const label b) const
87  {
88  return less(values_[a], values_[b]);
89  }
90  };
91 
92 
93  //- Template specialization for pTraits<labelList> so we can have fields
94  template<>
95  class pTraits<labelList>
96  {
97 
98  public:
99 
100  //- Component type
101  typedef labelList cmptType;
102  };
103 
104  //- Template specialization for pTraits<labelList> so we can have fields
105  template<>
106  class pTraits<vectorList>
107  {
108 
109  public:
110 
111  //- Component type
112  typedef vectorList cmptType;
113  };
114 }
115 
116 
117 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
118 
119 // Get faces (on the new mesh) that have in some way been affected by the
120 // mesh change. Picks up all faces but those that are between two
121 // unrefined faces. (Note that of an unchanged face the edge still might be
122 // split but that does not change any face centre or cell centre.
123 Foam::labelList Foam::meshRefinement::getChangedFaces
124 (
125  const mapPolyMesh& map,
126  const labelList& oldCellsToRefine
127 )
128 {
129  const polyMesh& mesh = map.mesh();
130 
131  labelList changedFaces;
132 
133  // For reporting: number of masterFaces changed
134  label nMasterChanged = 0;
135  bitSet isMasterFace(syncTools::getMasterFaces(mesh));
136 
137  {
138  // Mark any face on a cell which has been added or changed
139  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
140  // Note that refining a face changes the face centre (for a warped face)
141  // which changes the cell centre. This again changes the cellcentre-
142  // cellcentre edge across all faces of the cell.
143  // Note: this does not happen for unwarped faces but unfortunately
144  // we don't have this information.
145 
146  const labelList& faceOwner = mesh.faceOwner();
147  const labelList& faceNeighbour = mesh.faceNeighbour();
148  const cellList& cells = mesh.cells();
149  const label nInternalFaces = mesh.nInternalFaces();
150 
151  // Mark refined cells on old mesh
152  bitSet oldRefineCell(map.nOldCells(), oldCellsToRefine);
153 
154  // Mark refined faces
155  bitSet refinedInternalFace(nInternalFaces);
156 
157  // 1. Internal faces
158 
159  for (label faceI = 0; faceI < nInternalFaces; faceI++)
160  {
161  label oldOwn = map.cellMap()[faceOwner[faceI]];
162  label oldNei = map.cellMap()[faceNeighbour[faceI]];
163 
164  if
165  (
166  oldOwn >= 0 && !oldRefineCell.test(oldOwn)
167  && oldNei >= 0 && !oldRefineCell.test(oldNei)
168  )
169  {
170  // Unaffected face since both neighbours were not refined.
171  }
172  else
173  {
174  refinedInternalFace.set(faceI);
175  }
176  }
177 
178 
179  // 2. Boundary faces
180 
181  boolList refinedBoundaryFace(mesh.nBoundaryFaces(), false);
182 
183  forAll(mesh.boundaryMesh(), patchI)
184  {
185  const polyPatch& pp = mesh.boundaryMesh()[patchI];
186 
187  label faceI = pp.start();
188 
189  forAll(pp, i)
190  {
191  label oldOwn = map.cellMap()[faceOwner[faceI]];
192 
193  if (oldOwn >= 0 && !oldRefineCell.test(oldOwn))
194  {
195  // owner did exist and wasn't refined.
196  }
197  else
198  {
199  refinedBoundaryFace[faceI-nInternalFaces] = true;
200  }
201  faceI++;
202  }
203  }
204 
205  // Synchronise refined face status
207  (
208  mesh,
209  refinedBoundaryFace,
210  orEqOp<bool>()
211  );
212 
213 
214  // 3. Mark all faces affected by refinement. Refined faces are in
215  // - refinedInternalFace
216  // - refinedBoundaryFace
217  boolList changedFace(mesh.nFaces(), false);
218 
219  forAll(refinedInternalFace, faceI)
220  {
221  if (refinedInternalFace.test(faceI))
222  {
223  const cell& ownFaces = cells[faceOwner[faceI]];
224  forAll(ownFaces, ownI)
225  {
226  changedFace[ownFaces[ownI]] = true;
227  }
228  const cell& neiFaces = cells[faceNeighbour[faceI]];
229  forAll(neiFaces, neiI)
230  {
231  changedFace[neiFaces[neiI]] = true;
232  }
233  }
234  }
235 
236  forAll(refinedBoundaryFace, i)
237  {
238  if (refinedBoundaryFace[i])
239  {
240  const cell& ownFaces = cells[faceOwner[i+nInternalFaces]];
241  forAll(ownFaces, ownI)
242  {
243  changedFace[ownFaces[ownI]] = true;
244  }
245  }
246  }
247 
249  (
250  mesh,
251  changedFace,
252  orEqOp<bool>()
253  );
254 
255 
256  // Now we have in changedFace marked all affected faces. Pack.
257  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
258 
259  changedFaces = findIndices(changedFace, true);
260 
261  // Count changed master faces.
262  nMasterChanged = 0;
263 
264  forAll(changedFace, faceI)
265  {
266  if (changedFace[faceI] && isMasterFace.test(faceI))
267  {
268  nMasterChanged++;
269  }
270  }
271 
272  }
273 
275  {
276  Pout<< "getChangedFaces : Detected "
277  << " local:" << changedFaces.size()
278  << " global:" << returnReduce(nMasterChanged, sumOp<label>())
279  << " changed faces out of " << mesh.globalData().nTotalFaces()
280  << endl;
281 
282  faceSet changedFacesSet(mesh, "changedFaces", changedFaces);
283  Pout<< "getChangedFaces : Writing " << changedFaces.size()
284  << " changed faces to faceSet " << changedFacesSet.name()
285  << endl;
286  changedFacesSet.write();
287  }
288 
289  return changedFaces;
290 }
291 
292 
293 // Mark cell for refinement (if not already marked). Return false if
294 // refinelimit hit. Keeps running count (in nRefine) of cells marked for
295 // refinement
296 bool Foam::meshRefinement::markForRefine
297 (
298  const label markValue,
299  const label nAllowRefine,
300 
301  label& cellValue,
302  label& nRefine
303 )
304 {
305  if (cellValue == -1)
306  {
307  cellValue = markValue;
308  nRefine++;
309  }
310 
311  return nRefine <= nAllowRefine;
312 }
313 
314 
315 void Foam::meshRefinement::markFeatureCellLevel
316 (
317  const pointField& keepPoints,
318  labelList& maxFeatureLevel
319 ) const
320 {
321  // We want to refine all cells containing a feature edge.
322  // - don't want to search over whole mesh
323  // - don't want to build octree for whole mesh
324  // - so use tracking to follow the feature edges
325  //
326  // 1. find non-manifold points on feature edges (i.e. start of feature edge
327  // or 'knot')
328  // 2. seed particle starting at keepPoint going to this non-manifold point
329  // 3. track particles to their non-manifold point
330  //
331  // 4. track particles across their connected feature edges, marking all
332  // visited cells with their level (through trackingData)
333  // 5. do 4 until all edges have been visited.
334 
335 
336  // Find all start cells of features. Is done by tracking from keepPoint.
337  Cloud<trackedParticle> startPointCloud
338  (
339  mesh_,
340  "startPointCloud",
341  IDLList<trackedParticle>()
342  );
343 
344 
345  // Features are identical on all processors. Number them so we know
346  // what to seed. Do this on only the processor that
347  // holds the keepPoint.
348 
349  for (const point& keepPoint : keepPoints)
350  {
351  const label celli = mesh_.cellTree().findInside(keepPoint);
352 
353  if (celli != -1)
354  {
355  // I am the processor that holds the keepPoint
356 
357  forAll(features_, feati)
358  {
359  const edgeMesh& featureMesh = features_[feati];
360  const label featureLevel = features_.levels()[feati][0];
361  const labelListList& pointEdges = featureMesh.pointEdges();
362 
363  // Find regions on edgeMesh
364  labelList edgeRegion;
365  label nRegions = featureMesh.regions(edgeRegion);
366 
367 
368  bitSet regionVisited(nRegions);
369 
370 
371  // 1. Seed all 'knots' in edgeMesh
372 
373 
374  forAll(pointEdges, pointi)
375  {
376  if (pointEdges[pointi].size() != 2)
377  {
379  {
380  Pout<< "Adding particle from point:" << pointi
381  << " coord:" << featureMesh.points()[pointi]
382  << " since number of emanating edges:"
383  << pointEdges[pointi].size()
384  << endl;
385  }
386 
387  // Non-manifold point. Create particle.
388  startPointCloud.addParticle
389  (
390  new trackedParticle
391  (
392  mesh_,
393  keepPoint,
394  celli,
395  featureMesh.points()[pointi], // endpos
396  featureLevel, // level
397  feati, // featureMesh
398  pointi, // end point
399  -1 // feature edge
400  )
401  );
402 
403  // Mark
404  if (pointEdges[pointi].size() > 0)
405  {
406  label e0 = pointEdges[pointi][0];
407  label regioni = edgeRegion[e0];
408  regionVisited.set(regioni);
409  }
410  }
411  }
412 
413 
414  // 2. Any regions that have not been visited at all? These can
415  // only be circular regions!
416  forAll(featureMesh.edges(), edgei)
417  {
418  if (regionVisited.set(edgeRegion[edgei]))
419  {
420  const edge& e = featureMesh.edges()[edgei];
421  label pointi = e.start();
423  {
424  Pout<< "Adding particle from point:" << pointi
425  << " coord:" << featureMesh.points()[pointi]
426  << " on circular region:" << edgeRegion[edgei]
427  << endl;
428  }
429 
430  // Non-manifold point. Create particle.
431  startPointCloud.addParticle
432  (
433  new trackedParticle
434  (
435  mesh_,
436  keepPoint,
437  celli,
438  featureMesh.points()[pointi], // endpos
439  featureLevel, // level
440  feati, // featureMesh
441  pointi, // end point
442  -1 // feature edge
443  )
444  );
445  }
446  }
447  }
448  }
449  }
450 
451 
452  // Largest refinement level of any feature passed through
453  maxFeatureLevel = labelList(mesh_.nCells(), -1);
454 
455  // Whether edge has been visited.
456  List<bitSet> featureEdgeVisited(features_.size());
457 
458  forAll(features_, featI)
459  {
460  featureEdgeVisited[featI].setSize(features_[featI].edges().size());
461  featureEdgeVisited[featI] = false;
462  }
463 
464  // Database to pass into trackedParticle::move
465  trackedParticle::trackingData td
466  (
467  startPointCloud,
468  maxFeatureLevel,
469  featureEdgeVisited
470  );
471 
472 
473  // Track all particles to their end position (= starting feature point)
474  // Note that the particle might have started on a different processor
475  // so this will transfer across nicely until we can start tracking proper.
476  scalar maxTrackLen = 2.0*mesh_.bounds().mag();
477 
479  {
480  Pout<< "Tracking " << startPointCloud.size()
481  << " particles over distance " << maxTrackLen
482  << " to find the starting cell" << endl;
483  }
484  startPointCloud.move(startPointCloud, td, maxTrackLen);
485 
486 
487  // Reset levels
488  maxFeatureLevel = -1;
489  forAll(features_, featI)
490  {
491  featureEdgeVisited[featI] = false;
492  }
493 
494 
495  Cloud<trackedParticle> cloud
496  (
497  mesh_,
498  "featureCloud",
499  IDLList<trackedParticle>()
500  );
501 
503  {
504  Pout<< "Constructing cloud for cell marking" << endl;
505  }
506 
507  for (const trackedParticle& startTp : startPointCloud)
508  {
509  const label featI = startTp.i();
510  const label pointI = startTp.j();
511 
512  const edgeMesh& featureMesh = features_[featI];
513  const labelList& pEdges = featureMesh.pointEdges()[pointI];
514 
515  // Now shoot particles down all pEdges.
516  forAll(pEdges, pEdgeI)
517  {
518  label edgeI = pEdges[pEdgeI];
519 
520  if (featureEdgeVisited[featI].set(edgeI))
521  {
522  // Unvisited edge. Make the particle go to the other point
523  // on the edge.
524 
525  const edge& e = featureMesh.edges()[edgeI];
526  label otherPointi = e.otherVertex(pointI);
527 
528  trackedParticle* tp(new trackedParticle(startTp));
529  tp->start() = tp->position();
530  tp->end() = featureMesh.points()[otherPointi];
531  tp->j() = otherPointi;
532  tp->k() = edgeI;
533 
535  {
536  Pout<< "Adding particle for point:" << pointI
537  << " coord:" << tp->position()
538  << " feature:" << featI
539  << " to track to:" << tp->end()
540  << endl;
541  }
542 
543  cloud.addParticle(tp);
544  }
545  }
546  }
547 
548  startPointCloud.clear();
549 
550 
551  while (true)
552  {
553  // Track all particles to their end position.
555  {
556  Pout<< "Tracking " << cloud.size()
557  << " particles over distance " << maxTrackLen
558  << " to mark cells" << endl;
559  }
560  cloud.move(cloud, td, maxTrackLen);
561 
562  // Make particle follow edge.
563  for (trackedParticle& tp : cloud)
564  {
565  const label featI = tp.i();
566  const label pointI = tp.j();
567 
568  const edgeMesh& featureMesh = features_[featI];
569  const labelList& pEdges = featureMesh.pointEdges()[pointI];
570 
571  // Particle now at pointI. Check connected edges to see which one
572  // we have to visit now.
573 
574  bool keepParticle = false;
575 
576  forAll(pEdges, i)
577  {
578  label edgeI = pEdges[i];
579 
580  if (featureEdgeVisited[featI].set(edgeI))
581  {
582  // Unvisited edge. Make the particle go to the other point
583  // on the edge.
584 
585  const edge& e = featureMesh.edges()[edgeI];
586  label otherPointi = e.otherVertex(pointI);
587 
588  tp.start() = tp.position();
589  tp.end() = featureMesh.points()[otherPointi];
590  tp.j() = otherPointi;
591  tp.k() = edgeI;
592  keepParticle = true;
593  break;
594  }
595  }
596 
597  if (!keepParticle)
598  {
599  // Particle at 'knot' where another particle already has been
600  // seeded. Delete particle.
601  cloud.deleteParticle(tp);
602  }
603  }
604 
605 
607  {
608  Pout<< "Remaining particles " << cloud.size() << endl;
609  }
610 
611  if (!returnReduceOr(cloud.size()))
612  {
613  break;
614  }
615  }
616 
617 
618 
619  //if (debug&meshRefinement::FEATURESEEDS)
620  //{
621  // forAll(maxFeatureLevel, cellI)
622  // {
623  // if (maxFeatureLevel[cellI] != -1)
624  // {
625  // Pout<< "Feature went through cell:" << cellI
626  // << " coord:" << mesh_.cellCentres()[cellI]
627  // << " level:" << maxFeatureLevel[cellI]
628  // << endl;
629  // }
630  // }
631  //}
632 }
633 
634 
635 // Calculates list of cells to refine based on intersection with feature edge.
636 Foam::label Foam::meshRefinement::markFeatureRefinement
637 (
638  const pointField& keepPoints,
639  const label nAllowRefine,
640 
641  labelList& refineCell,
642  label& nRefine
643 ) const
644 {
645  // Largest refinement level of any feature passed through
646  labelList maxFeatureLevel;
647  markFeatureCellLevel(keepPoints, maxFeatureLevel);
648 
649  // See which cells to refine. maxFeatureLevel will hold highest level
650  // of any feature edge that passed through.
651 
652  const labelList& cellLevel = meshCutter_.cellLevel();
653 
654  label oldNRefine = nRefine;
655 
656  forAll(maxFeatureLevel, cellI)
657  {
658  if (maxFeatureLevel[cellI] > cellLevel[cellI])
659  {
660  // Mark
661  if
662  (
663  !markForRefine
664  (
665  0, // surface (n/a)
666  nAllowRefine,
667  refineCell[cellI],
668  nRefine
669  )
670  )
671  {
672  // Reached limit
673  break;
674  }
675  }
676  }
677 
678  if
679  (
680  returnReduce(nRefine, sumOp<label>())
681  > returnReduce(nAllowRefine, sumOp<label>())
682  )
683  {
684  Info<< "Reached refinement limit." << endl;
685  }
686 
687  return returnReduce(nRefine-oldNRefine, sumOp<label>());
688 }
689 
690 
691 // Mark cells for distance-to-feature based refinement.
692 Foam::label Foam::meshRefinement::markInternalDistanceToFeatureRefinement
693 (
694  const label nAllowRefine,
695 
696  labelList& refineCell,
697  label& nRefine
698 ) const
699 {
700  const labelList& cellLevel = meshCutter_.cellLevel();
701  const pointField& cellCentres = mesh_.cellCentres();
702 
703  // Detect if there are any distance shells
704  if (features_.maxDistance() <= 0.0)
705  {
706  return 0;
707  }
708 
709  label oldNRefine = nRefine;
710 
711  // Collect cells to test
712  pointField testCc(cellLevel.size()-nRefine);
713  labelList testLevels(cellLevel.size()-nRefine);
714  label testI = 0;
715 
716  forAll(cellLevel, cellI)
717  {
718  if (refineCell[cellI] == -1)
719  {
720  testCc[testI] = cellCentres[cellI];
721  testLevels[testI] = cellLevel[cellI];
722  testI++;
723  }
724  }
725 
726  // Do test to see whether cells are inside/outside shell with higher level
727  labelList maxLevel;
728  features_.findHigherLevel(testCc, testLevels, maxLevel);
729 
730  // Mark for refinement. Note that we didn't store the original cellID so
731  // now just reloop in same order.
732  testI = 0;
733  forAll(cellLevel, cellI)
734  {
735  if (refineCell[cellI] == -1)
736  {
737  if (maxLevel[testI] > testLevels[testI])
738  {
739  bool reachedLimit = !markForRefine
740  (
741  maxLevel[testI], // mark with any positive value
742  nAllowRefine,
743  refineCell[cellI],
744  nRefine
745  );
746 
747  if (reachedLimit)
748  {
749  if (debug)
750  {
751  Pout<< "Stopped refining internal cells"
752  << " since reaching my cell limit of "
753  << mesh_.nCells()+7*nRefine << endl;
754  }
755  break;
756  }
757  }
758  testI++;
759  }
760  }
761 
762  if
763  (
764  returnReduce(nRefine, sumOp<label>())
765  > returnReduce(nAllowRefine, sumOp<label>())
766  )
767  {
768  Info<< "Reached refinement limit." << endl;
769  }
770 
771  return returnReduce(nRefine-oldNRefine, sumOp<label>());
772 }
773 
774 
775 // Mark cells for non-surface intersection based refinement.
776 Foam::label Foam::meshRefinement::markInternalRefinement
777 (
778  const label nAllowRefine,
779 
780  labelList& refineCell,
781  label& nRefine
782 ) const
783 {
784  const labelList& cellLevel = meshCutter_.cellLevel();
785  const pointField& cellCentres = mesh_.cellCentres();
786 
787  label oldNRefine = nRefine;
788 
789  // Collect cells to test
790  pointField testCc(cellLevel.size()-nRefine);
791  labelList testLevels(cellLevel.size()-nRefine);
792  label testI = 0;
793 
794  forAll(cellLevel, cellI)
795  {
796  if (refineCell[cellI] == -1)
797  {
798  testCc[testI] = cellCentres[cellI];
799  testLevels[testI] = cellLevel[cellI];
800  testI++;
801  }
802  }
803 
804  // Do test to see whether cells are inside/outside shell with higher level
805  labelList maxLevel;
806  shells_.findHigherLevel(testCc, testLevels, maxLevel);
807 
808  // Mark for refinement. Note that we didn't store the original cellID so
809  // now just reloop in same order.
810  testI = 0;
811  forAll(cellLevel, cellI)
812  {
813  if (refineCell[cellI] == -1)
814  {
815  if (maxLevel[testI] > testLevels[testI])
816  {
817  bool reachedLimit = !markForRefine
818  (
819  maxLevel[testI], // mark with any positive value
820  nAllowRefine,
821  refineCell[cellI],
822  nRefine
823  );
824 
825  if (reachedLimit)
826  {
827  if (debug)
828  {
829  Pout<< "Stopped refining internal cells"
830  << " since reaching my cell limit of "
831  << mesh_.nCells()+7*nRefine << endl;
832  }
833  break;
834  }
835  }
836  testI++;
837  }
838  }
839 
840  if
841  (
842  returnReduce(nRefine, sumOp<label>())
843  > returnReduce(nAllowRefine, sumOp<label>())
844  )
845  {
846  Info<< "Reached refinement limit." << endl;
847  }
848 
849  return returnReduce(nRefine-oldNRefine, sumOp<label>());
850 }
851 
852 
853 Foam::label Foam::meshRefinement::unmarkInternalRefinement
854 (
855  labelList& refineCell,
856  label& nRefine
857 ) const
858 {
859  const labelList& cellLevel = meshCutter_.cellLevel();
860  const pointField& cellCentres = mesh_.cellCentres();
861 
862  label oldNRefine = nRefine;
863 
864  // Collect cells to test
865  pointField testCc(nRefine);
866  labelList testLevels(nRefine);
867  label testI = 0;
868 
869  forAll(cellLevel, cellI)
870  {
871  if (refineCell[cellI] >= 0)
872  {
873  testCc[testI] = cellCentres[cellI];
874  testLevels[testI] = cellLevel[cellI];
875  testI++;
876  }
877  }
878 
879  // Do test to see whether cells are inside/outside shell with higher level
880  labelList levelShell;
881  limitShells_.findLevel(testCc, testLevels, levelShell);
882 
883  // Mark for refinement. Note that we didn't store the original cellID so
884  // now just reloop in same order.
885  testI = 0;
886  forAll(cellLevel, cellI)
887  {
888  if (refineCell[cellI] >= 0)
889  {
890  if (levelShell[testI] != -1)
891  {
892  refineCell[cellI] = -1;
893  nRefine--;
894  }
895  testI++;
896  }
897  }
898 
899  return returnReduce(oldNRefine-nRefine, sumOp<label>());
900 }
901 
902 
903 // Collect faces that are intersected and whose neighbours aren't yet marked
904 // for refinement.
905 Foam::labelList Foam::meshRefinement::getRefineCandidateFaces
906 (
907  const labelList& refineCell
908 ) const
909 {
910  labelList testFaces(mesh_.nFaces());
911 
912  label nTest = 0;
913 
914  const labelList& surfIndex = surfaceIndex();
915 
916  labelList boundaryRefineCell;
917  syncTools::swapBoundaryCellList(mesh_, refineCell, boundaryRefineCell);
918 
919  forAll(surfIndex, faceI)
920  {
921  if (surfIndex[faceI] != -1)
922  {
923  label own = mesh_.faceOwner()[faceI];
924 
925  if (mesh_.isInternalFace(faceI))
926  {
927  label nei = mesh_.faceNeighbour()[faceI];
928 
929  if (refineCell[own] == -1 || refineCell[nei] == -1)
930  {
931  testFaces[nTest++] = faceI;
932  }
933  }
934  else
935  {
936  const label bFacei = faceI - mesh_.nInternalFaces();
937  if (refineCell[own] == -1 || boundaryRefineCell[bFacei] == -1)
938  {
939  testFaces[nTest++] = faceI;
940  }
941  }
942  }
943  }
944  testFaces.setSize(nTest);
945 
946  return testFaces;
947 }
948 
949 
950 // Mark cells for surface intersection based refinement.
951 Foam::label Foam::meshRefinement::markSurfaceRefinement
952 (
953  const label nAllowRefine,
954  const labelList& neiLevel,
955  const pointField& neiCc,
956 
957  labelList& refineCell,
958  label& nRefine
959 ) const
960 {
961  const labelList& cellLevel = meshCutter_.cellLevel();
962 
963  label oldNRefine = nRefine;
964 
965  // Use cached surfaceIndex_ to detect if any intersection. If so
966  // re-intersect to determine level wanted.
967 
968  // Collect candidate faces
969  // ~~~~~~~~~~~~~~~~~~~~~~~
970 
971  labelList testFaces(getRefineCandidateFaces(refineCell));
972 
973  // Collect segments
974  // ~~~~~~~~~~~~~~~~
975 
976  pointField start(testFaces.size());
977  pointField end(testFaces.size());
978  labelList minLevel(testFaces.size());
979 
980  calcCellCellRays
981  (
982  neiCc,
983  neiLevel,
984  testFaces,
985  start,
986  end,
987  minLevel
988  );
989 
990 
991  // Do test for higher intersections
992  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
993 
994  labelList surfaceHit;
995  labelList surfaceMinLevel;
996  surfaces_.findHigherIntersection
997  (
998  shells_,
999  start,
1000  end,
1001  minLevel,
1002 
1003  surfaceHit,
1004  surfaceMinLevel
1005  );
1006 
1007 
1008  // Mark cells for refinement
1009  // ~~~~~~~~~~~~~~~~~~~~~~~~~
1010 
1011  forAll(testFaces, i)
1012  {
1013  label faceI = testFaces[i];
1014 
1015  label surfI = surfaceHit[i];
1016 
1017  if (surfI != -1)
1018  {
1019  // Found intersection with surface with higher wanted
1020  // refinement. Check if the level field on that surface
1021  // specifies an even higher level. Note:this is weird. Should
1022  // do the check with the surfaceMinLevel whilst intersecting the
1023  // surfaces?
1024 
1025  label own = mesh_.faceOwner()[faceI];
1026 
1027  if (surfaceMinLevel[i] > cellLevel[own])
1028  {
1029  // Owner needs refining
1030  if
1031  (
1032  !markForRefine
1033  (
1034  surfI,
1035  nAllowRefine,
1036  refineCell[own],
1037  nRefine
1038  )
1039  )
1040  {
1041  break;
1042  }
1043  }
1044 
1045  if (mesh_.isInternalFace(faceI))
1046  {
1047  label nei = mesh_.faceNeighbour()[faceI];
1048  if (surfaceMinLevel[i] > cellLevel[nei])
1049  {
1050  // Neighbour needs refining
1051  if
1052  (
1053  !markForRefine
1054  (
1055  surfI,
1056  nAllowRefine,
1057  refineCell[nei],
1058  nRefine
1059  )
1060  )
1061  {
1062  break;
1063  }
1064  }
1065  }
1066  }
1067  }
1068 
1069  if
1070  (
1071  returnReduce(nRefine, sumOp<label>())
1072  > returnReduce(nAllowRefine, sumOp<label>())
1073  )
1074  {
1075  Info<< "Reached refinement limit." << endl;
1076  }
1077 
1078  return returnReduce(nRefine-oldNRefine, sumOp<label>());
1079 }
1080 
1081 
1082 // Count number of matches of first argument in second argument
1083 Foam::label Foam::meshRefinement::countMatches
1084 (
1085  const List<point>& normals1,
1086  const List<point>& normals2,
1087  const scalar tol
1088 ) const
1089 {
1090  label nMatches = 0;
1091 
1092  forAll(normals1, i)
1093  {
1094  const vector& n1 = normals1[i];
1095 
1096  forAll(normals2, j)
1097  {
1098  const vector& n2 = normals2[j];
1099 
1100  if (magSqr(n1-n2) < tol)
1101  {
1102  nMatches++;
1103  break;
1104  }
1105  }
1106  }
1107 
1108  return nMatches;
1109 }
1110 
1111 
1112 //bool Foam::meshRefinement::highCurvature
1113 //(
1114 // const scalar minCosAngle,
1115 // const point& p0,
1116 // const vector& n0,
1117 // const point& p1,
1118 // const vector& n1
1119 //) const
1120 //{
1121 // return ((n0&n1) < minCosAngle);
1122 //}
1123 bool Foam::meshRefinement::highCurvature
1124 (
1125  const scalar minCosAngle,
1126  const scalar lengthScale,
1127  const point& p0,
1128  const vector& n0,
1129  const point& p1,
1130  const vector& n1
1131 ) const
1132 {
1133  const scalar cosAngle = (n0&n1);
1134 
1135  if (cosAngle < minCosAngle)
1136  {
1137  // Sharp feature, independent of intersection points
1138  return true;
1139  }
1140  else if (cosAngle > 1-1e-6)
1141  {
1142  // Co-planar
1143  return false;
1144  }
1145  else if (lengthScale > SMALL)
1146  {
1147  // Calculate radius of curvature
1148 
1149  vector axis(n0 ^ n1);
1150  const plane pl0(p0, (n0 ^ axis));
1151  const scalar r1 = pl0.normalIntersect(p1, n1);
1152 
1153  //const point radiusPoint(p1+r1*n1);
1154  //DebugVar(radiusPoint);
1155  const plane pl1(p1, (n1 ^ axis));
1156  const scalar r0 = pl1.normalIntersect(p0, n0);
1157 
1158  //const point radiusPoint(p0+r0*n0);
1159  //DebugVar(radiusPoint);
1160  //- Take average radius. Bit ad hoc
1161  const scalar radius = 0.5*(mag(r1)+mag(r0));
1162 
1163  if (radius < lengthScale)
1164  {
1165  return true;
1166  }
1167  else
1168  {
1169  return false;
1170  }
1171  }
1172  else
1173  {
1174  return false;
1175  }
1176 }
1177 //XXXXX
1178 
1179 
1180 // Mark cells for surface curvature based refinement.
1181 Foam::label Foam::meshRefinement::markSurfaceCurvatureRefinement
1182 (
1183  const scalar curvature,
1184  const label nAllowRefine,
1185  const labelList& neiLevel,
1186  const pointField& neiCc,
1187 
1188  labelList& refineCell,
1189  label& nRefine
1190 ) const
1191 {
1192  const labelList& cellLevel = meshCutter_.cellLevel();
1193  const pointField& cellCentres = mesh_.cellCentres();
1194 
1195  label oldNRefine = nRefine;
1196 
1197  // 1. local test: any cell on more than one surface gets refined
1198  // (if its current level is < max of the surface max level)
1199 
1200  // 2. 'global' test: any cell on only one surface with a neighbour
1201  // on a different surface gets refined (if its current level etc.)
1202 
1203 
1204  const bitSet isMasterFace(syncTools::getMasterFaces(mesh_));
1205 
1206 
1207  // Collect candidate faces (i.e. intersecting any surface and
1208  // owner/neighbour not yet refined.
1209  labelList testFaces(getRefineCandidateFaces(refineCell));
1210 
1211  // Collect segments
1212  pointField start(testFaces.size());
1213  pointField end(testFaces.size());
1214  labelList minLevel(testFaces.size());
1215 
1216  // Note: uses isMasterFace otherwise could be call to calcCellCellRays
1217  forAll(testFaces, i)
1218  {
1219  label faceI = testFaces[i];
1220 
1221  label own = mesh_.faceOwner()[faceI];
1222 
1223  if (mesh_.isInternalFace(faceI))
1224  {
1225  label nei = mesh_.faceNeighbour()[faceI];
1226 
1227  start[i] = cellCentres[own];
1228  end[i] = cellCentres[nei];
1229  minLevel[i] = min(cellLevel[own], cellLevel[nei]);
1230  }
1231  else
1232  {
1233  label bFaceI = faceI - mesh_.nInternalFaces();
1234 
1235  start[i] = cellCentres[own];
1236  end[i] = neiCc[bFaceI];
1237 
1238  if (!isMasterFace[faceI])
1239  {
1240  std::swap(start[i], end[i]);
1241  }
1242 
1243  minLevel[i] = min(cellLevel[own], neiLevel[bFaceI]);
1244  }
1245  }
1246 
1247  // Extend segments a bit
1248  {
1249  const vectorField smallVec(ROOTSMALL*(end-start));
1250  start -= smallVec;
1251  end += smallVec;
1252  }
1253 
1254 
1255  // If no curvature supplied behave as before
1256  const bool hasCurvatureLevels = (max(surfaces_.maxCurvatureLevel()) > 0);
1257 
1258 
1259  // Test for all intersections (with surfaces of higher max level than
1260  // minLevel) and cache per cell the interesting inter
1261  labelListList cellSurfLevels(mesh_.nCells());
1262  List<pointList> cellSurfLocations(mesh_.nCells());
1263  List<vectorList> cellSurfNormals(mesh_.nCells());
1264 
1265  {
1266  // Per segment the intersection point of the surfaces
1267  List<pointList> surfaceLocation;
1268  // Per segment the normals of the surfaces
1269  List<vectorList> surfaceNormal;
1270  // Per segment the list of levels of the surfaces
1271  labelListList surfaceLevel;
1272 
1273  surfaces_.findAllIntersections
1274  (
1275  start,
1276  end,
1277  minLevel, // max level of surface has to be bigger
1278  // than min level of neighbouring cells
1279 
1280  labelList(surfaces_.maxLevel().size(), 0), // min level
1281  surfaces_.maxLevel(),
1282 
1283  surfaceLocation,
1284  surfaceNormal,
1285  surfaceLevel
1286  );
1287 
1288 
1289  // Sort the data according to surface location. This will guarantee
1290  // that on coupled faces both sides visit the intersections in
1291  // the same order so will decide the same
1292  labelList visitOrder;
1293  forAll(surfaceNormal, pointI)
1294  {
1295  pointList& pLocations = surfaceLocation[pointI];
1296  vectorList& pNormals = surfaceNormal[pointI];
1297  labelList& pLevel = surfaceLevel[pointI];
1298 
1299  sortedOrder(pNormals, visitOrder, normalLess(pLocations));
1300 
1301  pLocations = List<point>(pLocations, visitOrder);
1302  pNormals = List<vector>(pNormals, visitOrder);
1303  pLevel = labelUIndList(pLevel, visitOrder);
1304  }
1305 
1306 
1307  //- At some point could just return the intersected surface+region
1308  // and derive all the surface information (maxLevel, curvatureLevel)
1309  // from that - we're now doing it inside refinementSurfaces itself.
1311  //List<labelList> hitSurface;
1312  //List<pointList> hitLocation;
1313  //List<labelList> hitRegion;
1314  //List<vectorField> hitNormal;
1315  //surfaces_.findAllIntersections
1316  //(
1317  // identity(surfaces_.surfaces()), // all refinement geometries
1318  // start,
1319  // end,
1320  //
1321  // hitSurface,
1322  // hitLocation,
1323  // hitRegion,
1324  // hitNormal
1325  //);
1326  //
1329  //const auto& maxLevel = surfaces_.maxLevel();
1330  //labelList visitOrder;
1331  //DynamicList<label> valid;
1332  //forAll(hitSurface, segmenti)
1333  //{
1334  // const label meshLevel = minLevel[segmenti];
1335  //
1336  // auto& fSurface = hitSurface[segmenti];
1337  // auto& fLocation = hitLocation[segmenti];
1338  // auto& fRegion = hitRegion[segmenti];
1339  // auto& fNormal = hitNormal[segmenti];
1340  //
1341  // // Sort the data according to intersection location. This will
1342  // // guarantee
1343  // // that on coupled faces both sides visit the intersections in
1344  // // the same order so will decide the same
1345  // sortedOrder(fLocation, visitOrder, normalLess(hfLocation));
1346  // fLocation = List<point>(fLocation, visitOrder);
1347  // fSurface = labelUIndList(fSurface, visitOrder);
1348  // fRegion = labelUIndList(fRegion, visitOrder);
1349  // fNormal = List<vector>(fNormal, visitOrder);
1350  //
1351  // // Filter out any intersections with surfaces outside cell level.
1352  // // Note that min refinement level of surfaces is ignored.
1353  // valid.clear();
1354  // forAll(fSurface, hiti)
1355  // {
1356  // const label regioni =
1357  // surfaces_.globalRegion(fSurface[hiti], fRegion[hiti]);
1358  // if (meshLevel < maxLevel[regioni]) //&& >= minLevel(regioni)
1359  // {
1360  // valid.append(hiti);
1361  // }
1362  // }
1363  // fLocation = List<point>(fLocation, valid);
1364  // fSurface = labelUIndList(fSurface, valid);
1365  // fRegion = labelUIndList(fRegion, valid);
1366  // fNormal = List<vector>(fNormal, valid);
1367  //}
1368 
1369  // Clear out unnecessary data
1370  start.clear();
1371  end.clear();
1372  minLevel.clear();
1373 
1374  // Convert face-wise data to cell.
1375  forAll(testFaces, i)
1376  {
1377  label faceI = testFaces[i];
1378  label own = mesh_.faceOwner()[faceI];
1379 
1380  const pointList& fPoints = surfaceLocation[i];
1381  const vectorList& fNormals = surfaceNormal[i];
1382  const labelList& fLevels = surfaceLevel[i];
1383 
1384  forAll(fNormals, hitI)
1385  {
1386  if (fLevels[hitI] > cellLevel[own])
1387  {
1388  cellSurfLevels[own].append(fLevels[hitI]);
1389  cellSurfLocations[own].append(fPoints[hitI]);
1390  cellSurfNormals[own].append(fNormals[hitI]);
1391  }
1392 
1393  if (mesh_.isInternalFace(faceI))
1394  {
1395  label nei = mesh_.faceNeighbour()[faceI];
1396  if (fLevels[hitI] > cellLevel[nei])
1397  {
1398  cellSurfLevels[nei].append(fLevels[hitI]);
1399  cellSurfLocations[nei].append(fPoints[hitI]);
1400  cellSurfNormals[nei].append(fNormals[hitI]);
1401  }
1402  }
1403  }
1404  }
1405  }
1406 
1407 
1408 
1409  // Bit of statistics
1410  if (debug)
1411  {
1412  label nSet = 0;
1413  label nNormals = 0;
1414  forAll(cellSurfNormals, cellI)
1415  {
1416  const vectorList& normals = cellSurfNormals[cellI];
1417  if (normals.size())
1418  {
1419  nSet++;
1420  nNormals += normals.size();
1421  }
1422  }
1423  reduce(nSet, sumOp<label>());
1424  reduce(nNormals, sumOp<label>());
1425  Info<< "markSurfaceCurvatureRefinement :"
1426  << " cells:" << mesh_.globalData().nTotalCells()
1427  << " of which with normals:" << nSet
1428  << " ; total normals stored:" << nNormals
1429  << endl;
1430  }
1431 
1432 
1433 
1434  bool reachedLimit = false;
1435 
1436 
1437  // 1. Check normals per cell
1438  // ~~~~~~~~~~~~~~~~~~~~~~~~~
1439 
1440  for
1441  (
1442  label cellI = 0;
1443  !reachedLimit && cellI < cellSurfLocations.size();
1444  cellI++
1445  )
1446  {
1447  const pointList& points = cellSurfLocations[cellI];
1448  const vectorList& normals = cellSurfNormals[cellI];
1449  const labelList& levels = cellSurfLevels[cellI];
1450 
1451  // Current cell size
1452  const scalar cellSize =
1453  meshCutter_.level0EdgeLength()/pow(2.0, cellLevel[cellI]);
1454 
1455  // n^2 comparison of all normals in a cell
1456  for (label i = 0; !reachedLimit && i < normals.size(); i++)
1457  {
1458  for (label j = i+1; !reachedLimit && j < normals.size(); j++)
1459  {
1460  // TBD: calculate curvature size (if curvatureLevel specified)
1461  // and pass in instead of cellSize
1462 
1463  //if ((normals[i] & normals[j]) < curvature)
1464  if
1465  (
1466  highCurvature
1467  (
1468  curvature,
1469  (hasCurvatureLevels ? cellSize : scalar(0)),
1470  points[i],
1471  normals[i],
1472  points[j],
1473  normals[j]
1474  )
1475  )
1476  {
1477  const label maxLevel = max(levels[i], levels[j]);
1478 
1479  if (cellLevel[cellI] < maxLevel)
1480  {
1481  if
1482  (
1483  !markForRefine
1484  (
1485  maxLevel,
1486  nAllowRefine,
1487  refineCell[cellI],
1488  nRefine
1489  )
1490  )
1491  {
1492  if (debug)
1493  {
1494  Pout<< "Stopped refining since reaching my cell"
1495  << " limit of " << mesh_.nCells()+7*nRefine
1496  << endl;
1497  }
1498  reachedLimit = true;
1499  break;
1500  }
1501  }
1502  }
1503  }
1504  }
1505  }
1506 
1507 
1508 
1509  // 2. Find out a measure of surface curvature
1510  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1511  // Look at normals between neighbouring surfaces
1512  // Loop over all faces. Could only be checkFaces, except if they're coupled
1513 
1514  // Internal faces
1515  for
1516  (
1517  label faceI = 0;
1518  !reachedLimit && faceI < mesh_.nInternalFaces();
1519  faceI++
1520  )
1521  {
1522  label own = mesh_.faceOwner()[faceI];
1523  label nei = mesh_.faceNeighbour()[faceI];
1524 
1525  const pointList& ownPoints = cellSurfLocations[own];
1526  const vectorList& ownNormals = cellSurfNormals[own];
1527  const labelList& ownLevels = cellSurfLevels[own];
1528  const pointList& neiPoints = cellSurfLocations[nei];
1529  const vectorList& neiNormals = cellSurfNormals[nei];
1530  const labelList& neiLevels = cellSurfLevels[nei];
1531 
1532 
1533  // Special case: owner normals set is a subset of the neighbour
1534  // normals. Do not do curvature refinement since other cell's normals
1535  // do not add any information. Avoids weird corner extensions of
1536  // refinement regions.
1537  bool ownIsSubset =
1538  countMatches(ownNormals, neiNormals)
1539  == ownNormals.size();
1540 
1541  bool neiIsSubset =
1542  countMatches(neiNormals, ownNormals)
1543  == neiNormals.size();
1544 
1545 
1546  if (!ownIsSubset && !neiIsSubset)
1547  {
1548  // Current average cell size. Is min? max? average?
1549  const scalar cellSize =
1550  meshCutter_.level0EdgeLength()
1551  / pow(2.0, min(cellLevel[own], cellLevel[nei]));
1552 
1553  // n^2 comparison of between ownNormals and neiNormals
1554  for (label i = 0; !reachedLimit && i < ownNormals.size(); i++)
1555  {
1556  for (label j = 0; !reachedLimit && j < neiNormals.size(); j++)
1557  {
1558  // Have valid data on both sides. Check curvature.
1559  //if ((ownNormals[i] & neiNormals[j]) < curvature)
1560  if
1561  (
1562  highCurvature
1563  (
1564  curvature,
1565  (hasCurvatureLevels ? cellSize : scalar(0)),
1566  ownPoints[i],
1567  ownNormals[i],
1568  neiPoints[j],
1569  neiNormals[j]
1570  )
1571  )
1572  {
1573  // See which side to refine.
1574  if (cellLevel[own] < ownLevels[i])
1575  {
1576  if
1577  (
1578  !markForRefine
1579  (
1580  ownLevels[i],
1581  nAllowRefine,
1582  refineCell[own],
1583  nRefine
1584  )
1585  )
1586  {
1587  if (debug)
1588  {
1589  Pout<< "Stopped refining since reaching"
1590  << " my cell limit of "
1591  << mesh_.nCells()+7*nRefine << endl;
1592  }
1593  reachedLimit = true;
1594  break;
1595  }
1596  }
1597  if (cellLevel[nei] < neiLevels[j])
1598  {
1599  if
1600  (
1601  !markForRefine
1602  (
1603  neiLevels[j],
1604  nAllowRefine,
1605  refineCell[nei],
1606  nRefine
1607  )
1608  )
1609  {
1610  if (debug)
1611  {
1612  Pout<< "Stopped refining since reaching"
1613  << " my cell limit of "
1614  << mesh_.nCells()+7*nRefine << endl;
1615  }
1616  reachedLimit = true;
1617  break;
1618  }
1619  }
1620  }
1621  }
1622  }
1623  }
1624  }
1625 
1626 
1627  // Send over surface point/normal to neighbour cell.
1628 // labelListList neiSurfaceLevel;
1629 // syncTools::swapBoundaryCellList(mesh_, cellSurfLevels, neiSurfaceLevel);
1630  List<vectorList> neiSurfacePoints;
1631  syncTools::swapBoundaryCellList(mesh_, cellSurfLocations, neiSurfacePoints);
1632  List<vectorList> neiSurfaceNormals;
1633  syncTools::swapBoundaryCellList(mesh_, cellSurfNormals, neiSurfaceNormals);
1634 
1635  // Boundary faces
1636  for
1637  (
1638  label faceI = mesh_.nInternalFaces();
1639  !reachedLimit && faceI < mesh_.nFaces();
1640  faceI++
1641  )
1642  {
1643  label own = mesh_.faceOwner()[faceI];
1644  label bFaceI = faceI - mesh_.nInternalFaces();
1645 
1646  const pointList& ownPoints = cellSurfLocations[own];
1647  const vectorList& ownNormals = cellSurfNormals[own];
1648  const labelList& ownLevels = cellSurfLevels[own];
1649 
1650  const pointList& neiPoints = neiSurfacePoints[bFaceI];
1651  const vectorList& neiNormals = neiSurfaceNormals[bFaceI];
1652  //const labelList& neiLevels = neiSurfaceLevel[bFacei];
1653 
1654  // Special case: owner normals set is a subset of the neighbour
1655  // normals. Do not do curvature refinement since other cell's normals
1656  // do not add any information. Avoids weird corner extensions of
1657  // refinement regions.
1658  bool ownIsSubset =
1659  countMatches(ownNormals, neiNormals)
1660  == ownNormals.size();
1661 
1662  bool neiIsSubset =
1663  countMatches(neiNormals, ownNormals)
1664  == neiNormals.size();
1665 
1666 
1667  if (!ownIsSubset && !neiIsSubset)
1668  {
1669  // Current average cell size. Is min? max? average?
1670  const scalar cellSize =
1671  meshCutter_.level0EdgeLength()
1672  / pow(2.0, min(cellLevel[own], neiLevel[bFaceI]));
1673 
1674  // n^2 comparison of between ownNormals and neiNormals
1675  for (label i = 0; !reachedLimit && i < ownNormals.size(); i++)
1676  {
1677  for (label j = 0; !reachedLimit && j < neiNormals.size(); j++)
1678  {
1679  // Have valid data on both sides. Check curvature.
1680  //if ((ownNormals[i] & neiNormals[j]) < curvature)
1681  if
1682  (
1683  highCurvature
1684  (
1685  curvature,
1686  (hasCurvatureLevels ? cellSize : scalar(0)),
1687  ownPoints[i],
1688  ownNormals[i],
1689  neiPoints[j],
1690  neiNormals[j]
1691  )
1692  )
1693  {
1694  if (cellLevel[own] < ownLevels[i])
1695  {
1696  if
1697  (
1698  !markForRefine
1699  (
1700  ownLevels[i],
1701  nAllowRefine,
1702  refineCell[own],
1703  nRefine
1704  )
1705  )
1706  {
1707  if (debug)
1708  {
1709  Pout<< "Stopped refining since reaching"
1710  << " my cell limit of "
1711  << mesh_.nCells()+7*nRefine
1712  << endl;
1713  }
1714  reachedLimit = true;
1715  break;
1716  }
1717  }
1718  }
1719  }
1720  }
1721  }
1722  }
1723 
1724 
1725  if
1726  (
1727  returnReduce(nRefine, sumOp<label>())
1728  > returnReduce(nAllowRefine, sumOp<label>())
1729  )
1730  {
1731  Info<< "Reached refinement limit." << endl;
1732  }
1733 
1734  return returnReduce(nRefine-oldNRefine, sumOp<label>());
1735 }
1736 
1737 
1739 (
1740  const scalar planarCos,
1741  const vector& point0,
1742  const vector& normal0,
1743 
1744  const vector& point1,
1745  const vector& normal1
1746 ) const
1747 {
1748  //- Hits differ and angles are oppositeish and
1749  // hits have a normal distance
1750  vector d = point1-point0;
1751  scalar magD = mag(d);
1752 
1753  if (magD > mergeDistance())
1754  {
1755  scalar cosAngle = (normal0 & normal1);
1756 
1757  vector avg = Zero;
1758  if (cosAngle < (-1+planarCos))
1759  {
1760  // Opposite normals
1761  avg = 0.5*(normal0-normal1);
1762  }
1763  else if (cosAngle > (1-planarCos))
1764  {
1765  avg = 0.5*(normal0+normal1);
1766  }
1767 
1768  if (avg != vector::zero)
1769  {
1770  avg /= mag(avg);
1771 
1772  // Check normal distance of intersection locations
1773  if (mag(avg&d) > mergeDistance())
1774  {
1775  return true;
1776  }
1777  }
1778  }
1779 
1780  return false;
1781 }
1782 
1783 
1784 // Mark small gaps
1786 (
1787  const scalar planarCos,
1788  const label level0,
1789  const vector& point0,
1790  const vector& normal0,
1791 
1792  const label level1,
1793  const vector& point1,
1794  const vector& normal1
1795 ) const
1796 {
1797  //const scalar edge0Len = meshCutter_.level0EdgeLength();
1798 
1799  //- Hits differ and angles are oppositeish and
1800  // hits have a normal distance
1801  vector d = point1-point0;
1802  scalar magD = mag(d);
1803 
1804  if (magD > mergeDistance())
1805  {
1806  scalar cosAngle = (normal0 & normal1);
1807 
1808  vector avgNormal = Zero;
1809  if (cosAngle < (-1+planarCos))
1810  {
1811  // Opposite normals
1812  avgNormal = 0.5*(normal0-normal1);
1813  }
1814  else if (cosAngle > (1-planarCos))
1815  {
1816  avgNormal = 0.5*(normal0+normal1);
1817  }
1818 
1819  if (avgNormal != vector::zero)
1820  {
1821  avgNormal /= mag(avgNormal);
1822 
1823  //scalar normalDist = mag(d&avgNormal);
1824  //const scalar avgCellSize = edge0Len/pow(2.0, 0.5*(level0+level1));
1825  //if (normalDist > 1*avgCellSize)
1826  //{
1827  // Pout<< "** skipping large distance " << endl;
1828  // return false;
1829  //}
1830 
1831  // Check average normal with respect to intersection locations
1832  if (mag(avgNormal&d/magD) > Foam::cos(degToRad(45.0)))
1833  {
1834  return true;
1835  }
1836  }
1837  }
1838 
1839  return false;
1840 }
1841 
1842 
1843 bool Foam::meshRefinement::checkProximity
1844 (
1845  const scalar planarCos,
1846  const label nAllowRefine,
1847 
1848  const label surfaceLevel, // current intersection max level
1849  const vector& surfaceLocation, // current intersection location
1850  const vector& surfaceNormal, // current intersection normal
1851 
1852  const label cellI,
1853 
1854  label& cellMaxLevel, // cached max surface level for this cell
1855  vector& cellMaxLocation, // cached surface location for this cell
1856  vector& cellMaxNormal, // cached surface normal for this cell
1857 
1858  labelList& refineCell,
1859  label& nRefine
1860 ) const
1861 {
1862  const labelList& cellLevel = meshCutter_.cellLevel();
1863 
1864  // Test if surface applicable
1865  if (surfaceLevel > cellLevel[cellI])
1866  {
1867  if (cellMaxLevel == -1)
1868  {
1869  // First visit of cell. Store
1870  cellMaxLevel = surfaceLevel;
1871  cellMaxLocation = surfaceLocation;
1872  cellMaxNormal = surfaceNormal;
1873  }
1874  else
1875  {
1876  // Second or more visit.
1877  // Check if
1878  // - different location
1879  // - opposite surface
1880 
1881  bool closeSurfaces = isNormalGap
1882  (
1883  planarCos,
1884  cellLevel[cellI], //cellMaxLevel,
1885  cellMaxLocation,
1886  cellMaxNormal,
1887  cellLevel[cellI], //surfaceLevel,
1888  surfaceLocation,
1889  surfaceNormal
1890  );
1891 
1892  // Set normal to that of highest surface. Not really necessary
1893  // over here but we reuse cellMax info when doing coupled faces.
1894  if (surfaceLevel > cellMaxLevel)
1895  {
1896  cellMaxLevel = surfaceLevel;
1897  cellMaxLocation = surfaceLocation;
1898  cellMaxNormal = surfaceNormal;
1899  }
1900 
1901 
1902  if (closeSurfaces)
1903  {
1904  //Pout<< "Found gap:" << nl
1905  // << " location:" << surfaceLocation
1906  // << "\tnormal :" << surfaceNormal << nl
1908  // << "\tnormal :" << cellMaxNormal << nl
1909  // << "\tcos :" << (surfaceNormal&cellMaxNormal) << nl
1910  // << endl;
1911 
1912  return markForRefine
1913  (
1914  surfaceLevel, // mark with any non-neg number.
1915  nAllowRefine,
1916  refineCell[cellI],
1917  nRefine
1918  );
1919  }
1920  }
1921  }
1922 
1923  // Did not reach refinement limit.
1924  return true;
1925 }
1926 
1927 
1928 Foam::label Foam::meshRefinement::markProximityRefinement
1929 (
1930  const scalar planarCos,
1931 
1932  const labelList& surfaceMinLevel,
1933  const labelList& surfaceMaxLevel,
1934 
1935  const label nAllowRefine,
1936  const labelList& neiLevel,
1937  const pointField& neiCc,
1938 
1939  labelList& refineCell,
1940  label& nRefine
1941 ) const
1942 {
1943  const labelList& cellLevel = meshCutter_.cellLevel();
1944 
1945  label oldNRefine = nRefine;
1946 
1947  // 1. local test: any cell on more than one surface gets refined
1948  // (if its current level is < max of the surface max level)
1949 
1950  // 2. 'global' test: any cell on only one surface with a neighbour
1951  // on a different surface gets refined (if its current level etc.)
1952 
1953 
1954  // Collect candidate faces (i.e. intersecting any surface and
1955  // owner/neighbour not yet refined.
1956  const labelList testFaces(getRefineCandidateFaces(refineCell));
1957 
1958  // Collect segments
1959  pointField start(testFaces.size());
1960  pointField end(testFaces.size());
1961  labelList minLevel(testFaces.size());
1962 
1963  calcCellCellRays
1964  (
1965  neiCc,
1966  neiLevel,
1967  testFaces,
1968  start,
1969  end,
1970  minLevel
1971  );
1972 
1973 
1974  // Test for all intersections (with surfaces of higher gap level than
1975  // minLevel) and cache per cell the max surface level and the local normal
1976  // on that surface.
1977  labelList cellMaxLevel(mesh_.nCells(), -1);
1978  vectorField cellMaxNormal(mesh_.nCells(), Zero);
1979  pointField cellMaxLocation(mesh_.nCells(), Zero);
1980 
1981  {
1982  // Per segment the normals of the surfaces
1983  List<vectorList> surfaceLocation;
1984  List<vectorList> surfaceNormal;
1985  // Per segment the list of levels of the surfaces
1986  labelListList surfaceLevel;
1987 
1988  surfaces_.findAllIntersections
1989  (
1990  start,
1991  end,
1992  minLevel, // gap level of surface has to be bigger
1993  // than min level of neighbouring cells
1994 
1995  surfaceMinLevel,
1996  surfaceMaxLevel,
1997 
1998  surfaceLocation,
1999  surfaceNormal,
2000  surfaceLevel
2001  );
2002  // Clear out unnecessary data
2003  start.clear();
2004  end.clear();
2005  minLevel.clear();
2006 
2007 
2009 
2010  forAll(testFaces, i)
2011  {
2012  label faceI = testFaces[i];
2013  label own = mesh_.faceOwner()[faceI];
2014 
2015  const labelList& fLevels = surfaceLevel[i];
2016  const vectorList& fPoints = surfaceLocation[i];
2017  const vectorList& fNormals = surfaceNormal[i];
2018 
2019  forAll(fLevels, hitI)
2020  {
2021  checkProximity
2022  (
2023  planarCos,
2024  nAllowRefine,
2025 
2026  fLevels[hitI],
2027  fPoints[hitI],
2028  fNormals[hitI],
2029 
2030  own,
2031  cellMaxLevel[own],
2032  cellMaxLocation[own],
2033  cellMaxNormal[own],
2034 
2035  refineCell,
2036  nRefine
2037  );
2038  }
2039 
2040  if (mesh_.isInternalFace(faceI))
2041  {
2042  label nei = mesh_.faceNeighbour()[faceI];
2043 
2044  forAll(fLevels, hitI)
2045  {
2046  checkProximity
2047  (
2048  planarCos,
2049  nAllowRefine,
2050 
2051  fLevels[hitI],
2052  fPoints[hitI],
2053  fNormals[hitI],
2054 
2055  nei,
2056  cellMaxLevel[nei],
2057  cellMaxLocation[nei],
2058  cellMaxNormal[nei],
2059 
2060  refineCell,
2061  nRefine
2062  );
2063  }
2064  }
2065  }
2066  }
2067 
2068  // 2. Find out a measure of surface curvature and region edges.
2069  // Send over surface region and surface normal to neighbour cell.
2070 
2071  labelList neiBndMaxLevel(mesh_.nBoundaryFaces());
2072  pointField neiBndMaxLocation(mesh_.nBoundaryFaces());
2073  vectorField neiBndMaxNormal(mesh_.nBoundaryFaces());
2074 
2075  for (label faceI = mesh_.nInternalFaces(); faceI < mesh_.nFaces(); faceI++)
2076  {
2077  label bFaceI = faceI-mesh_.nInternalFaces();
2078  label own = mesh_.faceOwner()[faceI];
2079 
2080  neiBndMaxLevel[bFaceI] = cellMaxLevel[own];
2081  neiBndMaxLocation[bFaceI] = cellMaxLocation[own];
2082  neiBndMaxNormal[bFaceI] = cellMaxNormal[own];
2083  }
2084  syncTools::swapBoundaryFaceList(mesh_, neiBndMaxLevel);
2085  syncTools::swapBoundaryFaceList(mesh_, neiBndMaxLocation);
2086  syncTools::swapBoundaryFaceList(mesh_, neiBndMaxNormal);
2087 
2088  // Loop over all faces. Could only be checkFaces.. except if they're coupled
2089 
2090  // Internal faces
2091  for (label faceI = 0; faceI < mesh_.nInternalFaces(); faceI++)
2092  {
2093  label own = mesh_.faceOwner()[faceI];
2094  label nei = mesh_.faceNeighbour()[faceI];
2095 
2096  if (cellMaxLevel[own] != -1 && cellMaxLevel[nei] != -1)
2097  {
2098  // Have valid data on both sides. Check planarCos.
2099  if
2100  (
2101  isNormalGap
2102  (
2103  planarCos,
2104  cellLevel[own], //cellMaxLevel[own],
2105  cellMaxLocation[own],
2106  cellMaxNormal[own],
2107  cellLevel[nei], //cellMaxLevel[nei],
2108  cellMaxLocation[nei],
2109  cellMaxNormal[nei]
2110  )
2111  )
2112  {
2113  // See which side to refine
2114  if (cellLevel[own] < cellMaxLevel[own])
2115  {
2116  if
2117  (
2118  !markForRefine
2119  (
2120  cellMaxLevel[own],
2121  nAllowRefine,
2122  refineCell[own],
2123  nRefine
2124  )
2125  )
2126  {
2127  if (debug)
2128  {
2129  Pout<< "Stopped refining since reaching my cell"
2130  << " limit of " << mesh_.nCells()+7*nRefine
2131  << endl;
2132  }
2133  break;
2134  }
2135  }
2136 
2137  if (cellLevel[nei] < cellMaxLevel[nei])
2138  {
2139  if
2140  (
2141  !markForRefine
2142  (
2143  cellMaxLevel[nei],
2144  nAllowRefine,
2145  refineCell[nei],
2146  nRefine
2147  )
2148  )
2149  {
2150  if (debug)
2151  {
2152  Pout<< "Stopped refining since reaching my cell"
2153  << " limit of " << mesh_.nCells()+7*nRefine
2154  << endl;
2155  }
2156  break;
2157  }
2158  }
2159  }
2160  }
2161  }
2162  // Boundary faces
2163  for (label faceI = mesh_.nInternalFaces(); faceI < mesh_.nFaces(); faceI++)
2164  {
2165  label own = mesh_.faceOwner()[faceI];
2166  label bFaceI = faceI - mesh_.nInternalFaces();
2167 
2168  if (cellLevel[own] < cellMaxLevel[own] && neiBndMaxLevel[bFaceI] != -1)
2169  {
2170  // Have valid data on both sides. Check planarCos.
2171  if
2172  (
2173  isNormalGap
2174  (
2175  planarCos,
2176  cellLevel[own], //cellMaxLevel[own],
2177  cellMaxLocation[own],
2178  cellMaxNormal[own],
2179  neiLevel[bFaceI], //neiBndMaxLevel[bFaceI],
2180  neiBndMaxLocation[bFaceI],
2181  neiBndMaxNormal[bFaceI]
2182  )
2183  )
2184  {
2185  if
2186  (
2187  !markForRefine
2188  (
2189  cellMaxLevel[own],
2190  nAllowRefine,
2191  refineCell[own],
2192  nRefine
2193  )
2194  )
2195  {
2196  if (debug)
2197  {
2198  Pout<< "Stopped refining since reaching my cell"
2199  << " limit of " << mesh_.nCells()+7*nRefine
2200  << endl;
2201  }
2202  break;
2203  }
2204  }
2205  }
2206  }
2207 
2208  if
2209  (
2210  returnReduce(nRefine, sumOp<label>())
2211  > returnReduce(nAllowRefine, sumOp<label>())
2212  )
2213  {
2214  Info<< "Reached refinement limit." << endl;
2215  }
2216 
2217  return returnReduce(nRefine-oldNRefine, sumOp<label>());
2218 }
2219 
2220 
2221 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
2222 
2223 // Calculate list of cells to refine. Gets for any edge (start - end)
2224 // whether it intersects the surface. Does more accurate test and checks
2225 // the wanted level on the surface intersected.
2226 // Does approximate precalculation of how many cells can be refined before
2227 // hitting overall limit maxGlobalCells.
2229 (
2230  const pointField& keepPoints,
2231  const scalar curvature,
2232  const scalar planarAngle,
2233 
2234  const bool featureRefinement,
2235  const bool featureDistanceRefinement,
2236  const bool internalRefinement,
2237  const bool surfaceRefinement,
2238  const bool curvatureRefinement,
2239  const bool smallFeatureRefinement,
2240  const bool gapRefinement,
2241  const bool bigGapRefinement,
2242  const bool spreadGapSize,
2243  const label maxGlobalCells,
2244  const label maxLocalCells
2245 ) const
2246 {
2247  label totNCells = mesh_.globalData().nTotalCells();
2248 
2249  labelList cellsToRefine;
2250 
2251  if (totNCells >= maxGlobalCells)
2252  {
2253  Info<< "No cells marked for refinement since reached limit "
2254  << maxGlobalCells << '.' << endl;
2255  }
2256  else
2257  {
2258  // Every cell I refine adds 7 cells. Estimate the number of cells
2259  // I am allowed to refine.
2260  // Assume perfect distribution so can only refine as much the fraction
2261  // of the mesh I hold. This prediction step prevents us having to do
2262  // lots of reduces to keep count of the total number of cells selected
2263  // for refinement.
2264 
2265  //scalar fraction = scalar(mesh_.nCells())/totNCells;
2266  //label myMaxCells = label(maxGlobalCells*fraction);
2267  //label nAllowRefine = (myMaxCells - mesh_.nCells())/7;
2269  //
2270  //Pout<< "refineCandidates:" << nl
2271  // << " total cells:" << totNCells << nl
2272  // << " local cells:" << mesh_.nCells() << nl
2273  // << " local fraction:" << fraction << nl
2274  // << " local allowable cells:" << myMaxCells << nl
2275  // << " local allowable refinement:" << nAllowRefine << nl
2276  // << endl;
2277 
2278  //- Disable refinement shortcut. nAllowRefine is per processor limit.
2279  label nAllowRefine = labelMax / Pstream::nProcs();
2280 
2281  // Marked for refinement (>= 0) or not (-1). Actual value is the
2282  // index of the surface it intersects / shell it is inside.
2283  labelList refineCell(mesh_.nCells(), -1);
2284  label nRefine = 0;
2285 
2286 
2287  // Swap neighbouring cell centres and cell level
2288  labelList neiLevel(mesh_.nBoundaryFaces());
2289  pointField neiCc(mesh_.nBoundaryFaces());
2290  calcNeighbourData(neiLevel, neiCc);
2291 
2292 
2293  const scalar planarCos = Foam::cos(degToRad(planarAngle));
2294 
2295 
2296  // Cells pierced by feature lines
2297  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2298 
2299  if (featureRefinement)
2300  {
2301  label nFeatures = markFeatureRefinement
2302  (
2303  keepPoints,
2304  nAllowRefine,
2305 
2306  refineCell,
2307  nRefine
2308  );
2309 
2310  Info<< "Marked for refinement due to explicit features "
2311  << ": " << nFeatures << " cells." << endl;
2312  }
2313 
2314  // Inside distance-to-feature shells
2315  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2316 
2317  if (featureDistanceRefinement)
2318  {
2319  label nShell = markInternalDistanceToFeatureRefinement
2320  (
2321  nAllowRefine,
2322 
2323  refineCell,
2324  nRefine
2325  );
2326  Info<< "Marked for refinement due to distance to explicit features "
2327  ": " << nShell << " cells." << endl;
2328  }
2329 
2330  // Inside refinement shells
2331  // ~~~~~~~~~~~~~~~~~~~~~~~~
2332 
2333  if (internalRefinement)
2334  {
2335  label nShell = markInternalRefinement
2336  (
2337  nAllowRefine,
2338 
2339  refineCell,
2340  nRefine
2341  );
2342  Info<< "Marked for refinement due to refinement shells "
2343  << ": " << nShell << " cells." << endl;
2344  }
2345 
2346  // Refinement based on intersection of surface
2347  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2348 
2349  if (surfaceRefinement)
2350  {
2351  label nSurf = markSurfaceRefinement
2352  (
2353  nAllowRefine,
2354  neiLevel,
2355  neiCc,
2356 
2357  refineCell,
2358  nRefine
2359  );
2360  Info<< "Marked for refinement due to surface intersection "
2361  << ": " << nSurf << " cells." << endl;
2362 
2363  // Refine intersected-cells only inside gaps. See
2364  // markInternalGapRefinement to refine all cells inside gaps.
2365  if
2366  (
2367  planarCos >= -1
2368  && planarCos <= 1
2369  && max(shells_.maxGapLevel()) > 0
2370  )
2371  {
2372  label nGapSurf = markSurfaceGapRefinement
2373  (
2374  planarCos,
2375  nAllowRefine,
2376  neiLevel,
2377  neiCc,
2378 
2379  refineCell,
2380  nRefine
2381  );
2382  Info<< "Marked for refinement due to surface intersection"
2383  << " (at gaps)"
2384  << ": " << nGapSurf << " cells." << endl;
2385  }
2386  }
2387 
2388  // Refinement based on curvature of surface
2389  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2390 
2391  if
2392  (
2393  curvatureRefinement
2394  && (curvature >= -1 && curvature <= 1)
2395  && (surfaces_.minLevel() != surfaces_.maxLevel())
2396  )
2397  {
2398  label nCurv = markSurfaceCurvatureRefinement
2399  (
2400  curvature,
2401  nAllowRefine,
2402  neiLevel,
2403  neiCc,
2404 
2405  refineCell,
2406  nRefine
2407  );
2408  Info<< "Marked for refinement due to curvature/regions "
2409  << ": " << nCurv << " cells." << endl;
2410  }
2411 
2412 
2413  // Refinement based on features smaller than cell size
2414  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2415 
2416  if (smallFeatureRefinement)
2417  {
2418  label nGap = 0;
2419  if
2420  (
2421  planarCos >= -1
2422  && planarCos <= 1
2423  && max(shells_.maxGapLevel()) > 0
2424  )
2425  {
2426  nGap = markSmallFeatureRefinement
2427  (
2428  planarCos,
2429  nAllowRefine,
2430  neiLevel,
2431  neiCc,
2432 
2433  refineCell,
2434  nRefine
2435  );
2436  }
2437  Info<< "Marked for refinement due to close opposite surfaces "
2438  << ": " << nGap << " cells." << endl;
2439 
2440  label nCurv = 0;
2441  if (max(surfaces_.maxCurvatureLevel()) > 0)
2442  {
2443  nCurv = markSurfaceFieldRefinement
2444  (
2445  nAllowRefine,
2446  neiLevel,
2447  neiCc,
2448 
2449  refineCell,
2450  nRefine
2451  );
2452  Info<< "Marked for refinement"
2453  << " due to curvature "
2454  << ": " << nCurv << " cells." << endl;
2455  }
2456  }
2457 
2458 
2459  // Refinement based on gap (only neighbouring cells)
2460  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2461 
2462  const labelList& surfaceGapLevel = surfaces_.gapLevel();
2463 
2464  if
2465  (
2466  gapRefinement
2467  && (planarCos >= -1 && planarCos <= 1)
2468  && (max(surfaceGapLevel) > -1)
2469  )
2470  {
2471  Info<< "Specified gap level : " << max(surfaceGapLevel)
2472  << ", planar angle " << planarAngle << endl;
2473 
2474  label nGap = markProximityRefinement
2475  (
2476  planarCos,
2477 
2478  labelList(surfaceGapLevel.size(), -1), // surface min level
2479  surfaceGapLevel, // surface max level
2480 
2481  nAllowRefine,
2482  neiLevel,
2483  neiCc,
2484 
2485  refineCell,
2486  nRefine
2487  );
2488  Info<< "Marked for refinement due to close opposite surfaces "
2489  << ": " << nGap << " cells." << endl;
2490  }
2491 
2492 
2493  // Refinement based on gaps larger than cell size
2494  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2495 
2496  if
2497  (
2498  bigGapRefinement
2499  && (planarCos >= -1 && planarCos <= 1)
2500  && max(shells_.maxGapLevel()) > 0
2501  )
2502  {
2503  // Refine based on gap information provided by shell and nearest
2504  // surface
2505  labelList numGapCells;
2506  scalarField gapSize;
2507  label nGap = markInternalGapRefinement
2508  (
2509  planarCos,
2510  spreadGapSize,
2511  nAllowRefine,
2512 
2513  refineCell,
2514  nRefine,
2515  numGapCells,
2516  gapSize
2517  );
2518  Info<< "Marked for refinement due to opposite surfaces"
2519  << " "
2520  << ": " << nGap << " cells." << endl;
2521  }
2522 
2523 
2524  // Limit refinement
2525  // ~~~~~~~~~~~~~~~~
2526 
2527  if (limitShells_.shells().size())
2528  {
2529  label nUnmarked = unmarkInternalRefinement(refineCell, nRefine);
2530  if (nUnmarked > 0)
2531  {
2532  Info<< "Unmarked for refinement due to limit shells"
2533  << " : " << nUnmarked << " cells." << endl;
2534  }
2535  }
2536 
2537 
2538 
2539  // Pack cells-to-refine
2540  // ~~~~~~~~~~~~~~~~~~~~
2541 
2542  cellsToRefine.setSize(nRefine);
2543  nRefine = 0;
2544 
2545  forAll(refineCell, cellI)
2546  {
2547  if (refineCell[cellI] != -1)
2548  {
2549  cellsToRefine[nRefine++] = cellI;
2550  }
2551  }
2552  }
2553 
2554  return cellsToRefine;
2555 }
2556 
2557 
2559 (
2560  const labelList& cellsToRefine
2561 )
2562 {
2563  // Mesh changing engine.
2564  polyTopoChange meshMod(mesh_);
2565 
2566  // Play refinement commands into mesh changer.
2567  meshCutter_.setRefinement(cellsToRefine, meshMod);
2568 
2569  // Remove any unnecessary fields
2570  mesh_.clearOut();
2571  mesh_.moving(false);
2572 
2573  // Create mesh (no inflation), return map from old to new mesh.
2574  autoPtr<mapPolyMesh> mapPtr = meshMod.changeMesh(mesh_, false);
2575  mapPolyMesh& map = *mapPtr;
2576 
2577  // Update fields
2578  mesh_.updateMesh(map);
2579 
2580  // Optionally inflate mesh
2581  if (map.hasMotionPoints())
2582  {
2583  mesh_.movePoints(map.preMotionPoints());
2584  }
2585  else
2586  {
2587  // Delete mesh volumes.
2588  mesh_.clearOut();
2589  }
2590 
2591  // Reset the instance for if in overwrite mode
2592  mesh_.setInstance(timeName());
2593 
2594  // Update intersection info
2595  updateMesh(map, getChangedFaces(map, cellsToRefine));
2596 
2597  return mapPtr;
2598 }
2599 
2600 
2601 // Do refinement of consistent set of cells followed by truncation and
2602 // load balancing.
2605 (
2606  const string& msg,
2607  decompositionMethod& decomposer,
2608  fvMeshDistribute& distributor,
2609  const labelList& cellsToRefine,
2610  const scalar maxLoadUnbalance
2611 )
2612 {
2613  // Do all refinement
2614  refine(cellsToRefine);
2615 
2617  {
2618  Pout<< "Writing refined but unbalanced " << msg
2619  << " mesh to time " << timeName() << endl;
2620  write
2621  (
2622  debugType(debug),
2623  writeType(writeLevel() | WRITEMESH),
2624  mesh_.time().path()/timeName()
2625  );
2626  Pout<< "Dumped debug data in = "
2627  << mesh_.time().cpuTimeIncrement() << " s" << endl;
2628 
2629  // test all is still synced across proc patches
2630  checkData();
2631  }
2632 
2633  Info<< "Refined mesh in = "
2634  << mesh_.time().cpuTimeIncrement() << " s" << endl;
2635  printMeshInfo(debug, "After refinement " + msg);
2636 
2637 
2638  // Load balancing
2639  // ~~~~~~~~~~~~~~
2640 
2642 
2643  if (Pstream::nProcs() > 1)
2644  {
2645  scalar nIdealCells =
2646  mesh_.globalData().nTotalCells()
2647  / Pstream::nProcs();
2648 
2649  scalar unbalance = returnReduce
2650  (
2651  mag(1.0-mesh_.nCells()/nIdealCells),
2652  maxOp<scalar>()
2653  );
2654 
2655  if (unbalance <= maxLoadUnbalance)
2656  {
2657  Info<< "Skipping balancing since max unbalance " << unbalance
2658  << " is less than allowable " << maxLoadUnbalance
2659  << endl;
2660  }
2661  else
2662  {
2663  scalarField cellWeights(mesh_.nCells(), 1);
2664 
2665  distMap = balance
2666  (
2667  false, //keepZoneFaces
2668  false, //keepBaffles
2669  cellWeights,
2670  decomposer,
2671  distributor
2672  );
2673 
2674  Info<< "Balanced mesh in = "
2675  << mesh_.time().cpuTimeIncrement() << " s" << endl;
2676 
2677  printMeshInfo(debug, "After balancing " + msg);
2678 
2679 
2681  {
2682  Pout<< "Writing balanced " << msg
2683  << " mesh to time " << timeName() << endl;
2684  write
2685  (
2686  debugType(debug),
2687  writeType(writeLevel() | WRITEMESH),
2688  mesh_.time().path()/timeName()
2689  );
2690  Pout<< "Dumped debug data in = "
2691  << mesh_.time().cpuTimeIncrement() << " s" << endl;
2692 
2693  // test all is still synced across proc patches
2694  checkData();
2695  }
2696  }
2697  }
2698 
2699  return distMap;
2700 }
2701 
2702 
2703 // Do load balancing followed by refinement of consistent set of cells.
2706 (
2707  const string& msg,
2708  decompositionMethod& decomposer,
2709  fvMeshDistribute& distributor,
2710  const labelList& initCellsToRefine,
2711  const scalar maxLoadUnbalance
2712 )
2713 {
2714  labelList cellsToRefine(initCellsToRefine);
2715 
2716  //{
2717  // globalIndex globalCells(mesh_.nCells());
2718  //
2719  // Info<< "** Distribution before balancing/refining:" << endl;
2720  // for (const int procI : Pstream::allProcs())
2721  // {
2722  // Info<< " " << procI << '\t'
2723  // << globalCells.localSize(procI) << endl;
2724  // }
2725  // Info<< endl;
2726  //}
2727  //{
2728  // globalIndex globalCells(cellsToRefine.size());
2729  //
2730  // Info<< "** Cells to be refined:" << endl;
2731  // for (const int procI : Pstream::allProcs())
2732  // {
2733  // Info<< " " << procI << '\t'
2734  // << globalCells.localSize(procI) << endl;
2735  // }
2736  // Info<< endl;
2737  //}
2738 
2739 
2740  // Load balancing
2741  // ~~~~~~~~~~~~~~
2742 
2744 
2745  if (Pstream::nProcs() > 1)
2746  {
2747  // First check if we need to balance at all. Precalculate number of
2748  // cells after refinement and see what maximum difference is.
2749  scalar nNewCells = scalar(mesh_.nCells() + 7*cellsToRefine.size());
2750  scalar nIdealNewCells =
2751  returnReduce(nNewCells, sumOp<scalar>())
2752  / Pstream::nProcs();
2753  scalar unbalance = returnReduce
2754  (
2755  mag(1.0-nNewCells/nIdealNewCells),
2756  maxOp<scalar>()
2757  );
2758 
2759  if (unbalance <= maxLoadUnbalance)
2760  {
2761  Info<< "Skipping balancing since max unbalance " << unbalance
2762  << " is less than allowable " << maxLoadUnbalance
2763  << endl;
2764  }
2765  else
2766  {
2767  scalarField cellWeights(mesh_.nCells(), 1);
2768  forAll(cellsToRefine, i)
2769  {
2770  cellWeights[cellsToRefine[i]] += 7;
2771  }
2772 
2773  distMap = balance
2774  (
2775  false, //keepZoneFaces
2776  false, //keepBaffles
2777  cellWeights,
2778  decomposer,
2779  distributor
2780  );
2781 
2782  // Update cells to refine
2783  distMap().distributeCellIndices(cellsToRefine);
2784 
2785  Info<< "Balanced mesh in = "
2786  << mesh_.time().cpuTimeIncrement() << " s" << endl;
2787  }
2788 
2789  //{
2790  // globalIndex globalCells(mesh_.nCells());
2791  //
2792  // Info<< "** Distribution after balancing:" << endl;
2793  // for (const int procI : Pstream::allProcs())
2794  // {
2795  // Info<< " " << procI << '\t'
2796  // << globalCells.localSize(procI) << endl;
2797  // }
2798  // Info<< endl;
2799  //}
2800 
2801  printMeshInfo(debug, "After balancing " + msg);
2802 
2804  {
2805  Pout<< "Writing balanced " << msg
2806  << " mesh to time " << timeName() << endl;
2807  write
2808  (
2809  debugType(debug),
2810  writeType(writeLevel() | WRITEMESH),
2811  mesh_.time().path()/timeName()
2812  );
2813  Pout<< "Dumped debug data in = "
2814  << mesh_.time().cpuTimeIncrement() << " s" << endl;
2815 
2816  // test all is still synced across proc patches
2817  checkData();
2818  }
2819  }
2820 
2821 
2822  // Refinement
2823  // ~~~~~~~~~~
2824 
2825  refine(cellsToRefine);
2826 
2828  {
2829  Pout<< "Writing refined " << msg
2830  << " mesh to time " << timeName() << endl;
2831  write
2832  (
2833  debugType(debug),
2834  writeType(writeLevel() | WRITEMESH),
2835  mesh_.time().path()/timeName()
2836  );
2837  Pout<< "Dumped debug data in = "
2838  << mesh_.time().cpuTimeIncrement() << " s" << endl;
2839 
2840  // test all is still synced across proc patches
2841  checkData();
2842  }
2843 
2844  Info<< "Refined mesh in = "
2845  << mesh_.time().cpuTimeIncrement() << " s" << endl;
2846 
2847  //{
2848  // globalIndex globalCells(mesh_.nCells());
2849  //
2850  // Info<< "** After refinement distribution:" << endl;
2851  // for (const int procI : Pstream::allProcs())
2852  // {
2853  // Info<< " " << procI << '\t'
2854  // << globalCells.localSize(procI) << endl;
2855  // }
2856  // Info<< endl;
2857  //}
2858 
2859  printMeshInfo(debug, "After refinement " + msg);
2860 
2861  return distMap;
2862 }
2863 
2864 
2866 (
2867  const label maxGlobalCells,
2868  const label maxLocalCells,
2869  const labelList& currentLevel,
2870  const direction dir
2871 ) const
2872 {
2873  const labelList& cellLevel = meshCutter_.cellLevel();
2874  const pointField& cellCentres = mesh_.cellCentres();
2875 
2876  label totNCells = mesh_.globalData().nTotalCells();
2877 
2878  labelList cellsToRefine;
2879 
2880  if (totNCells >= maxGlobalCells)
2881  {
2882  Info<< "No cells marked for refinement since reached limit "
2883  << maxGlobalCells << '.' << endl;
2884  }
2885  else
2886  {
2887  // Disable refinement shortcut. nAllowRefine is per processor limit.
2888  label nAllowRefine = labelMax / Pstream::nProcs();
2889 
2890  // Marked for refinement (>= 0) or not (-1). Actual value is the
2891  // index of the surface it intersects / shell it is inside
2892  labelList refineCell(mesh_.nCells(), -1);
2893  label nRefine = 0;
2894 
2895  // Find cells inside the shells with directional levels
2896  labelList insideShell;
2897  shells_.findDirectionalLevel
2898  (
2899  cellCentres,
2900  cellLevel,
2901  currentLevel, // current directional level
2902  dir,
2903  insideShell
2904  );
2905 
2906  // Mark for refinement
2907  forAll(insideShell, celli)
2908  {
2909  if (insideShell[celli] >= 0)
2910  {
2911  bool reachedLimit = !markForRefine
2912  (
2913  insideShell[celli], // mark with any positive value
2914  nAllowRefine,
2915  refineCell[celli],
2916  nRefine
2917  );
2918 
2919  if (reachedLimit)
2920  {
2921  if (debug)
2922  {
2923  Pout<< "Stopped refining cells"
2924  << " since reaching my cell limit of "
2925  << mesh_.nCells()+nAllowRefine << endl;
2926  }
2927  break;
2928  }
2929  }
2930  }
2931 
2932  // Limit refinement
2933  // ~~~~~~~~~~~~~~~~
2934 
2935  {
2936  label nUnmarked = unmarkInternalRefinement(refineCell, nRefine);
2937  if (nUnmarked > 0)
2938  {
2939  Info<< "Unmarked for refinement due to limit shells"
2940  << " : " << nUnmarked << " cells." << endl;
2941  }
2942  }
2943 
2944 
2945 
2946  // Pack cells-to-refine
2947  // ~~~~~~~~~~~~~~~~~~~~
2948 
2949  cellsToRefine.setSize(nRefine);
2950  nRefine = 0;
2951 
2952  forAll(refineCell, cellI)
2953  {
2954  if (refineCell[cellI] != -1)
2955  {
2956  cellsToRefine[nRefine++] = cellI;
2957  }
2958  }
2959  }
2960 
2961  return cellsToRefine;
2962 }
2963 
2964 
2966 (
2967  const string& msg,
2968  const direction cmpt,
2969  const labelList& cellsToRefine
2970 )
2971 {
2972  // Set splitting direction
2973  vector refDir(Zero);
2974  refDir[cmpt] = 1;
2975  List<refineCell> refCells(cellsToRefine.size());
2976  forAll(cellsToRefine, i)
2977  {
2978  refCells[i] = refineCell(cellsToRefine[i], refDir);
2979  }
2980 
2981  // How to walk circumference of cells
2982  hexCellLooper cellWalker(mesh_);
2983 
2984  // Analyse cuts
2985  cellCuts cuts(mesh_, cellWalker, refCells);
2986 
2987  // Cell cutter
2988  Foam::meshCutter meshRefiner(mesh_);
2989 
2990  polyTopoChange meshMod(mesh_);
2991 
2992  // Insert mesh refinement into polyTopoChange.
2993  meshRefiner.setRefinement(cuts, meshMod);
2994 
2995  // Remove any unnecessary fields
2996  mesh_.clearOut();
2997  mesh_.moving(false);
2998 
2999  autoPtr<mapPolyMesh> morphMap = meshMod.changeMesh(mesh_, false);
3000 
3001  // Update fields
3002  mesh_.updateMesh(*morphMap);
3003 
3004  // Move mesh (since morphing does not do this)
3005  if (morphMap().hasMotionPoints())
3006  {
3007  mesh_.movePoints(morphMap().preMotionPoints());
3008  }
3009  else
3010  {
3011  // Delete mesh volumes.
3012  mesh_.clearOut();
3013  }
3014 
3015  // Reset the instance for if in overwrite mode
3016  mesh_.setInstance(timeName());
3017 
3018  // Update stored refinement pattern
3019  meshRefiner.updateMesh(*morphMap);
3020 
3021  // Update intersection info
3022  updateMesh(*morphMap, getChangedFaces(*morphMap, cellsToRefine));
3023 
3024  return morphMap;
3025 }
3026 
3027 
3028 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:584
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:51
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:118
uint8_t direction
Definition: direction.H:46
A list of face labels.
Definition: faceSet.H:47
label nOldCells() const
Number of old cells.
Definition: mapPolyMesh.H:469
const labelList & cellMap() const
Old cell map.
Definition: mapPolyMesh.H:534
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static bitSet getMasterFaces(const polyMesh &mesh)
Get per face whether it is uncoupled or a master of a coupled set of faces.
Definition: syncTools.C:119
autoPtr< mapDistributePolyMesh > balanceAndRefine(const string &msg, decompositionMethod &decomposer, fvMeshDistribute &distributor, const labelList &cellsToRefine, const scalar maxLoadUnbalance)
Balance before refining some cells.
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1110
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:150
bool isNormalGap(const scalar planarCos, const label level0, const vector &point0, const vector &normal0, const label level1, const vector &point1, const vector &normal1) const
Is local topology a small gap normal to the test vector.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:463
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
A traits class, which is primarily used for primitives.
Definition: pTraits.H:50
static void syncFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop)
Synchronize values on all mesh faces.
Definition: syncTools.H:432
labelList findIndices(const ListType &input, typename ListType::const_reference val, label start=0)
Linear search to find all occurrences of given element.
const cellList & cells() const
labelList refineCandidates(const pointField &keepPoints, const scalar curvature, const scalar planarAngle, const bool featureRefinement, const bool featureDistanceRefinement, const bool internalRefinement, const bool surfaceRefinement, const bool curvatureRefinement, const bool smallFeatureRefinement, const bool gapRefinement, const bool bigGapRefinement, const bool spreadGapSize, const label maxGlobalCells, const label maxLocalCells) const
Calculate list of cells to refine.
List< point > pointList
A List of points.
Definition: pointList.H:61
label nFaces() const noexcept
Number of mesh faces.
autoPtr< mapPolyMesh > directionalRefine(const string &msg, const direction cmpt, const labelList &cellsToRefine)
Directionally refine in direction cmpt.
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Description of cuts across cells.
Definition: cellCuts.H:106
static bool less(const vector &x, const vector &y)
To compare normals.
bool hasMotionPoints() const
Has valid preMotionPoints?
Definition: mapPolyMesh.H:764
List< vector > vectorList
A List of vectors.
Definition: vectorList.H:61
static void syncBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop, const TransformOp &top, const bool parRun=Pstream::parRun())
Synchronize values on boundary faces only.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:157
const fvMesh & mesh() const
Reference to mesh.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:164
word timeName
Definition: getTimeIndex.H:3
scalar y
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator) is 1 for serial run.
Definition: UPstream.H:656
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
void setSize(const label n)
Alias for resize()
Definition: List.H:289
dimensionedScalar cos(const dimensionedScalar &ds)
const cellShapeList & cells
const pointField & points
autoPtr< mapDistributePolyMesh > refineAndBalance(const string &msg, decompositionMethod &decomposer, fvMeshDistribute &distributor, const labelList &cellsToRefine, const scalar maxLoadUnbalance)
Refine some cells and rebalance.
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:109
label checkData(const fvMesh &mesh, const instantList &timeDirs, wordList &objectNames)
Check if fields are good to use (available at all times)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Cuts (splits) cells.
Definition: meshCutter.H:134
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1104
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1293
label nInternalFaces() const noexcept
Number of internal faces.
Vector< scalar > vector
Definition: vector.H:57
Container with cells to refine. Refinement given as single direction.
Definition: refineCell.H:52
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
labelList directionalRefineCandidates(const label maxGlobalCells, const label maxLocalCells, const labelList &currentLevel, const direction dir) const
Calculate list of cells to directionally refine.
bool test(const label pos) const
Test value at specified position, never auto-vivify entries.
Definition: bitSetI.H:514
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
int debug
Static debugging option.
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:193
autoPtr< mapPolyMesh > refine(const labelList &cellsToRefine)
Refine some cells.
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData)
Swap to obtain neighbour cell values for all boundary faces.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
To compare normals.
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: IndirectList.H:65
vector point
Point is a vector.
Definition: point.H:37
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:53
Direct mesh changes based on v1.3 polyTopoChange syntax.
const polyMesh & mesh() const
Return polyMesh.
Definition: mapPolyMesh.H:437
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:441
normalLess(const vectorList &values)
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
messageStream Info
Information stream (stdout output on master, null elsewhere)
label nTotalFaces() const noexcept
Return total number of faces in decomposed mesh. Not.
constexpr label labelMax
Definition: label.H:55
writeType
Enumeration for what to write. Used as a bit-pattern.
Field< vector > vectorField
Specialisation of Field<T> for vector.
debugType
Enumeration for what to debug. Used as a bit-pattern.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:73
bool isGap(const scalar, const vector &, const vector &, const vector &, const vector &) const
Is local topology a small gap?
List< label > labelList
A List of labels.
Definition: List.H:62
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
virtual bool write(const bool valid=true) const
Write using setting from DB.
label nBoundaryFaces() const noexcept
Number of boundary faces (== nFaces - nInternalFaces)
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
const volScalarField & p0
Definition: EEqn.H:36
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
const pointField & preMotionPoints() const
Pre-motion point positions.
Definition: mapPolyMesh.H:756
Namespace for OpenFOAM.
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues)
Swap coupled boundary face values. Uses eqOp.
Definition: syncTools.H:485
bool operator()(const label a, const label b) const
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157