meshRefinementRefine.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-2023 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 "trackedParticle.H"
28 #include "syncTools.H"
29 #include "Time.H"
30 #include "refinementSurfaces.H"
31 #include "refinementFeatures.H"
32 #include "refinementRegions.H"
33 #include "faceSet.H"
34 #include "decompositionMethod.H"
35 #include "fvMeshDistribute.H"
36 #include "polyTopoChange.H"
37 #include "polyDistributionMap.H"
38 #include "Cloud.H"
39 #include "OBJstream.H"
40 #include "cellSet.H"
41 #include "treeDataCell.H"
42 
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47  //- To compare normals
48  static bool less(const vector& x, const vector& y)
49  {
50  for (direction i = 0; i < vector::nComponents; i++)
51  {
52  if (x[i] < y[i])
53  {
54  return true;
55  }
56  else if (x[i] > y[i])
57  {
58  return false;
59  }
60  }
61  // All components the same
62  return false;
63  }
64 
65 
66  //- To compare normals
67  class normalLess
68  {
69  const vectorList& values_;
70 
71  public:
72 
73  normalLess(const vectorList& values)
74  :
75  values_(values)
76  {}
77 
78  bool operator()(const label a, const label b) const
79  {
80  return less(values_[a], values_[b]);
81  }
82  };
83 
84 
85  //- Template specialisation for pTraits<labelList> so we can have fields
86  template<>
88  {
89 
90  public:
91 
92  //- Component type
94  };
95 
96  //- Template specialisation for pTraits<labelList> so we can have fields
97  template<>
99  {
100 
101  public:
102 
103  //- Component type
105  };
106 }
107 
108 
109 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
110 
111 // Get faces (on the new mesh) that have in some way been affected by the
112 // mesh change. Picks up all faces but those that are between two
113 // unrefined faces. (Note that of an unchanged face the edge still might be
114 // split but that does not change any face centre or cell centre.
115 Foam::labelList Foam::meshRefinement::getChangedFaces
116 (
117  const polyTopoChangeMap& map,
118  const labelList& oldCellsToRefine
119 )
120 {
121  const polyMesh& mesh = map.mesh();
122 
123  labelList changedFaces;
124 
125  // For reporting: number of masterFaces changed
126  label nMasterChanged = 0;
127  PackedBoolList isMasterFace(syncTools::getMasterFaces(mesh));
128 
129  {
130  // Mark any face on a cell which has been added or changed
131  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
132  // Note that refining a face changes the face centre (for a warped face)
133  // which changes the cell centre. This again changes the cellcentre-
134  // cellcentre edge across all faces of the cell.
135  // Note: this does not happen for unwarped faces but unfortunately
136  // we don't have this information.
137 
138  const labelList& faceOwner = mesh.faceOwner();
139  const labelList& faceNeighbour = mesh.faceNeighbour();
140  const cellList& cells = mesh.cells();
141  const label nInternalFaces = mesh.nInternalFaces();
142 
143  // Mark refined cells on old mesh
144  PackedBoolList oldRefineCell(map.nOldCells());
145 
146  forAll(oldCellsToRefine, i)
147  {
148  oldRefineCell.set(oldCellsToRefine[i], 1u);
149  }
150 
151  // Mark refined faces
152  PackedBoolList refinedInternalFace(nInternalFaces);
153 
154  // 1. Internal faces
155 
156  for (label facei = 0; facei < nInternalFaces; facei++)
157  {
158  const label oldOwn = map.cellMap()[faceOwner[facei]];
159  const label oldNei = map.cellMap()[faceNeighbour[facei]];
160 
161  if
162  (
163  oldOwn >= 0
164  && oldRefineCell.get(oldOwn) == 0u
165  && oldNei >= 0
166  && oldRefineCell.get(oldNei) == 0u
167  )
168  {
169  // Unaffected face since both neighbours were not refined.
170  }
171  else
172  {
173  refinedInternalFace.set(facei, 1u);
174  }
175  }
176 
177 
178  // 2. Boundary faces
179 
180  boolList refinedBoundaryFace(mesh.nFaces()-nInternalFaces, false);
181 
183  {
184  const polyPatch& pp = mesh.boundaryMesh()[patchi];
185 
186  label facei = pp.start();
187 
188  forAll(pp, i)
189  {
190  label oldOwn = map.cellMap()[faceOwner[facei]];
191 
192  if (oldOwn >= 0 && oldRefineCell.get(oldOwn) == 0u)
193  {
194  // owner did exist and wasn't refined.
195  }
196  else
197  {
198  refinedBoundaryFace[facei - nInternalFaces] = true;
199  }
200  facei++;
201  }
202  }
203 
204  // Synchronise refined face status
206  (
207  mesh,
208  refinedBoundaryFace,
209  orEqOp<bool>()
210  );
211 
212 
213  // 3. Mark all faces affected by refinement. Refined faces are in
214  // - refinedInternalFace
215  // - refinedBoundaryFace
216  boolList changedFace(mesh.nFaces(), false);
217 
218  forAll(refinedInternalFace, facei)
219  {
220  if (refinedInternalFace.get(facei) == 1u)
221  {
222  const cell& ownFaces = cells[faceOwner[facei]];
223  forAll(ownFaces, owni)
224  {
225  changedFace[ownFaces[owni]] = true;
226  }
227  const cell& neiFaces = cells[faceNeighbour[facei]];
228  forAll(neiFaces, neii)
229  {
230  changedFace[neiFaces[neii]] = true;
231  }
232  }
233  }
234 
235  forAll(refinedBoundaryFace, i)
236  {
237  if (refinedBoundaryFace[i])
238  {
239  const cell& ownFaces = cells[faceOwner[i+nInternalFaces]];
240  forAll(ownFaces, owni)
241  {
242  changedFace[ownFaces[owni]] = true;
243  }
244  }
245  }
246 
248  (
249  mesh,
250  changedFace,
251  orEqOp<bool>()
252  );
253 
254 
255  // Now we have in changedFace marked all affected faces. Pack.
256  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
257 
258  changedFaces = findIndices(changedFace, true);
259 
260  // Count changed master faces.
261  nMasterChanged = 0;
262 
263  forAll(changedFace, facei)
264  {
265  if (changedFace[facei] && isMasterFace[facei])
266  {
267  nMasterChanged++;
268  }
269  }
270 
271  }
272 
273  if (debug&meshRefinement::MESH)
274  {
275  Pout<< "getChangedFaces : Detected "
276  << " local:" << changedFaces.size()
277  << " global:" << returnReduce(nMasterChanged, sumOp<label>())
278  << " changed faces out of " << mesh.globalData().nTotalFaces()
279  << endl;
280 
281  faceSet changedFacesSet(mesh, "changedFaces", changedFaces);
282  Pout<< "getChangedFaces : Writing " << changedFaces.size()
283  << " changed faces to faceSet " << changedFacesSet.name()
284  << endl;
285  changedFacesSet.write();
286  }
287 
288  return changedFaces;
289 }
290 
291 
292 // Mark cell for refinement (if not already marked). Return false if
293 // refinelimit hit. Keeps running count (in nRefine) of cells marked for
294 // refinement
295 bool Foam::meshRefinement::markForRefine
296 (
297  const label markValue,
298  const label nAllowRefine,
299 
300  label& cellValue,
301  label& nRefine
302 )
303 {
304  if (cellValue == -1)
305  {
306  cellValue = markValue;
307  nRefine++;
308  }
309 
310  return nRefine <= nAllowRefine;
311 }
312 
313 
314 void Foam::meshRefinement::markFeatureCellLevel
315 (
316  const List<point>& insidePoints,
317  labelList& maxFeatureLevel
318 ) const
319 {
320  // We want to refine all cells containing a feature edge.
321  // - don't want to search over whole mesh
322  // - don't want to build octree for whole mesh
323  // - so use tracking to follow the feature edges
324  //
325  // 1. find non-manifold points on feature edges (i.e. start of feature edge
326  // or 'knot')
327  // 2. seed particle starting at insidePoint going to this non-manifold
328  // 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.
337  // Is done by tracking from insidePoint.
338  Cloud<trackedParticle> startPointCloud
339  (
340  mesh_,
341  "startPointCloud",
342  IDLList<trackedParticle>()
343  );
344 
345 
346  // Features are identical on all processors. Number them so we know
347  // what to seed. Do this on only the processor that
348  // holds the insidePoint.
349 
350  forAll(insidePoints, i)
351  {
352  const point& insidePoint = insidePoints[i];
353 
354  const label celli = mesh_.cellTree().findInside(insidePoint);
355 
356  if (celli != -1)
357  {
358  // I am the processor that holds the insidePoint
359 
360  forAll(features_, feati)
361  {
362  const edgeMesh& featureMesh = features_[feati];
363  const label featureLevel = features_.levels()[feati][0];
364  const labelListList& pointEdges = featureMesh.pointEdges();
365 
366  // Find regions on edgeMesh
367  labelList edgeRegion;
368  const label nRegions = featureMesh.regions(edgeRegion);
369 
370 
371  PackedBoolList regionVisited(nRegions);
372 
373 
374  // 1. Seed all 'knots' in edgeMesh
375 
376 
377  forAll(pointEdges, pointi)
378  {
379  if (pointEdges[pointi].size() != 2)
380  {
382  {
383  Pout<< "Adding particle from point:" << pointi
384  << " coord:" << featureMesh.points()[pointi]
385  << " since number of emanating edges:"
386  << pointEdges[pointi].size()
387  << endl;
388  }
389 
390  // Non-manifold point. Create particle.
391  startPointCloud.addParticle
392  (
393  new trackedParticle
394  (
395  mesh_,
396  insidePoint,
397  celli,
398  featureMesh.points()[pointi], // endpos
399  featureLevel, // level
400  feati, // featureMesh
401  pointi, // end point
402  -1 // feature edge
403  )
404  );
405 
406  // Mark
407  if (pointEdges[pointi].size() > 0)
408  {
409  label e0 = pointEdges[pointi][0];
410  label regioni = edgeRegion[e0];
411  regionVisited[regioni] = 1u;
412  }
413  }
414  }
415 
416 
417  // 2. Any regions that have not been visited at all? These can
418  // only be circular regions!
419  forAll(featureMesh.edges(), edgei)
420  {
421  if (regionVisited.set(edgeRegion[edgei], 1u))
422  {
423  const edge& e = featureMesh.edges()[edgei];
424  const label pointi = e.start();
426  {
427  Pout<< "Adding particle from point:" << pointi
428  << " coord:" << featureMesh.points()[pointi]
429  << " on circular region:" << edgeRegion[edgei]
430  << endl;
431  }
432 
433  // Non-manifold point. Create particle.
434  startPointCloud.addParticle
435  (
436  new trackedParticle
437  (
438  mesh_,
439  insidePoint,
440  celli,
441  featureMesh.points()[pointi], // endpos
442  featureLevel, // level
443  feati, // featureMesh
444  pointi, // end point
445  -1 // feature edge
446  )
447  );
448  }
449  }
450  }
451  }
452  }
453 
454 
455  // Largest refinement level of any feature passed through
456  maxFeatureLevel = labelList(mesh_.nCells(), -1);
457 
458  // Whether edge has been visited.
459  List<PackedBoolList> featureEdgeVisited(features_.size());
460 
461  forAll(features_, feati)
462  {
463  featureEdgeVisited[feati].setSize(features_[feati].edges().size());
464  featureEdgeVisited[feati] = 0u;
465  }
466 
467  // Database to pass into trackedParticle::move
468  trackedParticle::trackingData td
469  (
470  startPointCloud,
471  2.0*mesh_.bounds().mag(),
472  maxFeatureLevel,
473  featureEdgeVisited
474  );
475 
476 
477  // Track all particles to their end position (= starting feature point)
478  // Note that the particle might have started on a different processor
479  // so this will transfer across nicely until we can start tracking proper.
481  {
482  Pout<< "Tracking " << startPointCloud.size()
483  << " particles over distance " << td.maxTrackLen_
484  << " to find the starting cell" << endl;
485  }
486  startPointCloud.move(startPointCloud, td);
487 
488 
489  // Reset levels
490  maxFeatureLevel = -1;
491  forAll(features_, feati)
492  {
493  featureEdgeVisited[feati] = 0u;
494  }
495 
496 
497  Cloud<trackedParticle> cloud
498  (
499  mesh_,
500  "featureCloud",
501  IDLList<trackedParticle>()
502  );
503 
505  {
506  Pout<< "Constructing cloud for cell marking" << endl;
507  }
508 
509  forAllIter(Cloud<trackedParticle>, startPointCloud, iter)
510  {
511  const trackedParticle& startTp = iter();
512 
513  const label feati = startTp.i();
514  const label pointi = startTp.j();
515 
516  const edgeMesh& featureMesh = features_[feati];
517  const labelList& pEdges = featureMesh.pointEdges()[pointi];
518 
519  // Now shoot particles down all pEdges.
520  forAll(pEdges, pEdgei)
521  {
522  const label edgei = pEdges[pEdgei];
523 
524  if (featureEdgeVisited[feati].set(edgei, 1u))
525  {
526  // Unvisited edge. Make the particle go to the other point
527  // on the edge.
528 
529  const edge& e = featureMesh.edges()[edgei];
530  label otherPointi = e.otherVertex(pointi);
531 
532  trackedParticle* tp(new trackedParticle(startTp));
533  tp->start() = tp->position(mesh_);
534  tp->end() = featureMesh.points()[otherPointi];
535  tp->j() = otherPointi;
536  tp->k() = edgei;
537 
539  {
540  Pout<< "Adding particle for point:" << pointi
541  << " coord:" << tp->position(mesh_)
542  << " feature:" << feati
543  << " to track to:" << tp->end()
544  << endl;
545  }
546 
547  cloud.addParticle(tp);
548  }
549  }
550  }
551 
552  startPointCloud.clear();
553 
554 
555  while (true)
556  {
557  // Track all particles to their end position.
559  {
560  Pout<< "Tracking " << cloud.size()
561  << " particles over distance " << td.maxTrackLen_
562  << " to mark cells" << endl;
563  }
564  cloud.move(cloud, td);
565 
566  // Make particle follow edge.
567  forAllIter(Cloud<trackedParticle>, cloud, iter)
568  {
569  trackedParticle& tp = iter();
570 
571  const label feati = tp.i();
572  const label pointi = tp.j();
573 
574  const edgeMesh& featureMesh = features_[feati];
575  const labelList& pEdges = featureMesh.pointEdges()[pointi];
576 
577  // Particle now at pointi. Check connected edges to see which one
578  // we have to visit now.
579 
580  bool keepParticle = false;
581 
582  forAll(pEdges, i)
583  {
584  const label edgei = pEdges[i];
585 
586  if (featureEdgeVisited[feati].set(edgei, 1u))
587  {
588  // Unvisited edge. Make the particle go to the other point
589  // on the edge.
590 
591  const edge& e = featureMesh.edges()[edgei];
592  const label otherPointi = e.otherVertex(pointi);
593 
594  tp.start() = tp.position(mesh_);
595  tp.end() = featureMesh.points()[otherPointi];
596  tp.j() = otherPointi;
597  tp.k() = edgei;
598  keepParticle = true;
599  break;
600  }
601  }
602 
603  if (!keepParticle)
604  {
605  // Particle at 'knot' where another particle already has been
606  // seeded. Delete particle.
607  cloud.deleteParticle(tp);
608  }
609  }
610 
611 
613  {
614  Pout<< "Remaining particles " << cloud.size() << endl;
615  }
616 
617  if (returnReduce(cloud.size(), sumOp<label>()) == 0)
618  {
619  break;
620  }
621  }
622 
623 
624 
625  // if (debug&meshRefinement::FEATURESEEDS)
626  //{
627  // forAll(maxFeatureLevel, celli)
628  // {
629  // if (maxFeatureLevel[celli] != -1)
630  // {
631  // Pout<< "Feature went through cell:" << celli
632  // << " coord:" << mesh_.cellCentres()[celli]
633  // << " level:" << maxFeatureLevel[celli]
634  // << endl;
635  // }
636  // }
637  //}
638 }
639 
640 
641 // Calculates list of cells to refine based on intersection with feature edge.
642 Foam::label Foam::meshRefinement::markFeatureRefinement
643 (
644  const List<point>& insidePoints,
645  const label nAllowRefine,
646 
647  labelList& refineCell,
648  label& nRefine
649 ) const
650 {
651  // Largest refinement level of any feature passed through
652  labelList maxFeatureLevel;
653  markFeatureCellLevel(insidePoints, maxFeatureLevel);
654 
655  // See which cells to refine. maxFeatureLevel will hold highest level
656  // of any feature edge that passed through.
657 
658  const labelList& cellLevel = meshCutter_.cellLevel();
659 
660  const label oldNRefine = nRefine;
661 
662  forAll(maxFeatureLevel, celli)
663  {
664  if (maxFeatureLevel[celli] > cellLevel[celli])
665  {
666  // Mark
667  if
668  (
669  !markForRefine
670  (
671  0, // surface (n/a)
672  nAllowRefine,
673  refineCell[celli],
674  nRefine
675  )
676  )
677  {
678  // Reached limit
679  break;
680  }
681  }
682  }
683 
684  if
685  (
686  returnReduce(nRefine, sumOp<label>())
687  > returnReduce(nAllowRefine, sumOp<label>())
688  )
689  {
690  Info<< "Reached refinement limit." << endl;
691  }
692 
693  return returnReduce(nRefine-oldNRefine, sumOp<label>());
694 }
695 
696 
697 // Mark cells for distance-to-feature based refinement.
698 Foam::label Foam::meshRefinement::markInternalDistanceToFeatureRefinement
699 (
700  const label nAllowRefine,
701 
702  labelList& refineCell,
703  label& nRefine
704 ) const
705 {
706  const labelList& cellLevel = meshCutter_.cellLevel();
707  const pointField& cellCentres = mesh_.cellCentres();
708 
709  // Detect if there are any distance features
710  if (features_.maxDistance() <= 0.0)
711  {
712  return 0;
713  }
714 
715  label oldNRefine = nRefine;
716 
717  // Collect cells to test
718  pointField testCc(cellLevel.size()-nRefine);
719  labelList testLevels(cellLevel.size()-nRefine);
720  label testi = 0;
721 
722  forAll(cellLevel, celli)
723  {
724  if (refineCell[celli] == -1)
725  {
726  testCc[testi] = cellCentres[celli];
727  testLevels[testi] = cellLevel[celli];
728  testi++;
729  }
730  }
731 
732  // Do test to see whether cells are in/near features with higher level
733  labelList maxLevel;
734  features_.findHigherLevel(testCc, testLevels, maxLevel);
735 
736  // Mark for refinement. Note that we didn't store the original cellID so
737  // now just reloop in same order.
738  testi = 0;
739  forAll(cellLevel, celli)
740  {
741  if (refineCell[celli] == -1)
742  {
743  if (maxLevel[testi] > testLevels[testi])
744  {
745  const bool reachedLimit = !markForRefine
746  (
747  maxLevel[testi], // mark with any positive value
748  nAllowRefine,
749  refineCell[celli],
750  nRefine
751  );
752 
753  if (reachedLimit)
754  {
755  if (debug)
756  {
757  Pout<< "Stopped refining internal cells"
758  << " since reaching my cell limit of "
759  << mesh_.nCells()+7*nRefine << endl;
760  }
761  break;
762  }
763  }
764  testi++;
765  }
766  }
767 
768  if
769  (
770  returnReduce(nRefine, sumOp<label>())
771  > returnReduce(nAllowRefine, sumOp<label>())
772  )
773  {
774  Info<< "Reached refinement limit." << endl;
775  }
776 
777  return returnReduce(nRefine-oldNRefine, sumOp<label>());
778 }
779 
780 
781 Foam::label Foam::meshRefinement::markInternalRefinement
782 (
783  const label nAllowRefine,
784 
785  labelList& refineCell,
786  label& nRefine
787 ) const
788 {
789  const labelList& cellLevel = meshCutter_.cellLevel();
790  const pointField& cellCentres = mesh_.cellCentres();
791 
792  const label oldNRefine = nRefine;
793 
794  // Collect cells to test
795  pointField testCc(cellLevel.size()-nRefine);
796  labelList testLevels(cellLevel.size()-nRefine);
797  label testi = 0;
798 
799  forAll(cellLevel, celli)
800  {
801  if (refineCell[celli] == -1)
802  {
803  testCc[testi] = cellCentres[celli];
804  testLevels[testi] = cellLevel[celli];
805  testi++;
806  }
807  }
808 
809  // Do test to see whether cells is inside/outside shell with higher level
810  labelList maxLevel;
811  shells_.findHigherLevel
812  (
813  testCc,
814  testLevels,
815  meshCutter_.level0EdgeLength(),
816  maxLevel
817  );
818 
819  // Mark for refinement. Note that we didn't store the original cellID so
820  // now just reloop in same order.
821  testi = 0;
822  forAll(cellLevel, celli)
823  {
824  if (refineCell[celli] == -1)
825  {
826  if (maxLevel[testi] > testLevels[testi])
827  {
828  bool reachedLimit = !markForRefine
829  (
830  maxLevel[testi], // mark with any positive value
831  nAllowRefine,
832  refineCell[celli],
833  nRefine
834  );
835 
836  if (reachedLimit)
837  {
838  if (debug)
839  {
840  Pout<< "Stopped refining internal cells"
841  << " since reaching my cell limit of "
842  << mesh_.nCells()+7*nRefine << endl;
843  }
844  break;
845  }
846  }
847  testi++;
848  }
849  }
850 
851  if
852  (
853  returnReduce(nRefine, sumOp<label>())
854  > returnReduce(nAllowRefine, sumOp<label>())
855  )
856  {
857  Info<< "Reached refinement limit." << endl;
858  }
859 
860  return returnReduce(nRefine-oldNRefine, sumOp<label>());
861 }
862 
863 
864 Foam::labelList Foam::meshRefinement::getRefineCandidateFaces
865 (
866  const labelList& refineCell
867 ) const
868 {
869  labelList testFaces(mesh_.nFaces());
870 
871  label nTest = 0;
872 
873  forAll(surfaceIndex_, facei)
874  {
875  if (surfaceIndex_[facei] != -1)
876  {
877  label own = mesh_.faceOwner()[facei];
878 
879  if (mesh_.isInternalFace(facei))
880  {
881  const label nei = mesh_.faceNeighbour()[facei];
882 
883  if (refineCell[own] == -1 || refineCell[nei] == -1)
884  {
885  testFaces[nTest++] = facei;
886  }
887  }
888  else
889  {
890  if (refineCell[own] == -1)
891  {
892  testFaces[nTest++] = facei;
893  }
894  }
895  }
896  }
897  testFaces.setSize(nTest);
898 
899  return testFaces;
900 }
901 
902 
903 Foam::label Foam::meshRefinement::markSurfaceRefinement
904 (
905  const label nAllowRefine,
906  const labelList& neiLevel,
907  const pointField& neiCc,
908 
909  labelList& refineCell,
910  label& nRefine
911 ) const
912 {
913  const labelList& cellLevel = meshCutter_.cellLevel();
914  const pointField& cellCentres = mesh_.cellCentres();
915 
916  const label oldNRefine = nRefine;
917 
918  // Use cached surfaceIndex_ to detect if any intersection. If so
919  // re-intersect to determine level wanted.
920 
921  // Collect candidate faces
922  // ~~~~~~~~~~~~~~~~~~~~~~~
923 
924  labelList testFaces(getRefineCandidateFaces(refineCell));
925 
926  // Collect segments
927  // ~~~~~~~~~~~~~~~~
928 
929  pointField start(testFaces.size());
930  pointField end(testFaces.size());
931  labelList minLevel(testFaces.size());
932 
933  forAll(testFaces, i)
934  {
935  const label facei = testFaces[i];
936 
937  const label own = mesh_.faceOwner()[facei];
938 
939  if (mesh_.isInternalFace(facei))
940  {
941  label nei = mesh_.faceNeighbour()[facei];
942 
943  start[i] = cellCentres[own];
944  end[i] = cellCentres[nei];
945  minLevel[i] = min(cellLevel[own], cellLevel[nei]);
946  }
947  else
948  {
949  const label bFacei = facei - mesh_.nInternalFaces();
950 
951  start[i] = cellCentres[own];
952  end[i] = neiCc[bFacei];
953  minLevel[i] = min(cellLevel[own], neiLevel[bFacei]);
954  }
955  }
956 
957  // Extend segments a bit
958  {
959  const vectorField smallVec(rootSmall*(end-start));
960  start -= smallVec;
961  end += smallVec;
962  }
963 
964 
965  // Do test for higher intersections
966  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
967 
968  labelList surfaceHit;
969  labelList surfaceMinLevel;
970  surfaces_.findHigherIntersection
971  (
972  start,
973  end,
974  minLevel,
975 
976  surfaceHit,
977  surfaceMinLevel
978  );
979 
980 
981  // Mark cells for refinement
982  // ~~~~~~~~~~~~~~~~~~~~~~~~~
983 
984  forAll(testFaces, i)
985  {
986  const label facei = testFaces[i];
987 
988  const label surfi = surfaceHit[i];
989 
990  if (surfi != -1)
991  {
992  // Found intersection with surface with higher wanted
993  // refinement. Check if the level field on that surface
994  // specifies an even higher level. Note:this is weird. Should
995  // do the check with the surfaceMinLevel whilst intersecting the
996  // surfaces?
997 
998  const label own = mesh_.faceOwner()[facei];
999 
1000  if (surfaceMinLevel[i] > cellLevel[own])
1001  {
1002  // Owner needs refining
1003  if
1004  (
1005  !markForRefine
1006  (
1007  surfi,
1008  nAllowRefine,
1009  refineCell[own],
1010  nRefine
1011  )
1012  )
1013  {
1014  break;
1015  }
1016  }
1017 
1018  if (mesh_.isInternalFace(facei))
1019  {
1020  const label nei = mesh_.faceNeighbour()[facei];
1021 
1022  if (surfaceMinLevel[i] > cellLevel[nei])
1023  {
1024  // Neighbour needs refining
1025  if
1026  (
1027  !markForRefine
1028  (
1029  surfi,
1030  nAllowRefine,
1031  refineCell[nei],
1032  nRefine
1033  )
1034  )
1035  {
1036  break;
1037  }
1038  }
1039  }
1040  }
1041  }
1042 
1043  if
1044  (
1045  returnReduce(nRefine, sumOp<label>())
1046  > returnReduce(nAllowRefine, sumOp<label>())
1047  )
1048  {
1049  Info<< "Reached refinement limit." << endl;
1050  }
1051 
1052  return returnReduce(nRefine-oldNRefine, sumOp<label>());
1053 }
1054 
1055 
1056 // Count number of matches of first argument in second argument
1057 Foam::label Foam::meshRefinement::countMatches
1058 (
1059  const List<point>& normals1,
1060  const List<point>& normals2,
1061  const scalar tol
1062 ) const
1063 {
1064  label nMatches = 0;
1065 
1066  forAll(normals1, i)
1067  {
1068  const vector& n1 = normals1[i];
1069 
1070  forAll(normals2, j)
1071  {
1072  const vector& n2 = normals2[j];
1073 
1074  if (magSqr(n1-n2) < tol)
1075  {
1076  nMatches++;
1077  break;
1078  }
1079  }
1080  }
1081 
1082  return nMatches;
1083 }
1084 
1085 
1086 // Mark cells for surface curvature based refinement.
1087 Foam::label Foam::meshRefinement::markSurfaceCurvatureRefinement
1088 (
1089  const scalar curvature,
1090  const label nAllowRefine,
1091  const labelList& neiLevel,
1092  const pointField& neiCc,
1093 
1094  labelList& refineCell,
1095  label& nRefine
1096 ) const
1097 {
1098  const labelList& cellLevel = meshCutter_.cellLevel();
1099  const pointField& cellCentres = mesh_.cellCentres();
1100 
1101  const label oldNRefine = nRefine;
1102 
1103  // 1. local test: any cell on more than one surface gets refined
1104  // (if its current level is < max of the surface max level)
1105 
1106  // 2. 'global' test: any cell on only one surface with a neighbour
1107  // on a different surface gets refined (if its current level etc.)
1108 
1109 
1110  const PackedBoolList isMasterFace(syncTools::getMasterFaces(mesh_));
1111 
1112 
1113  // Collect candidate faces (i.e. intersecting any surface and
1114  // owner/neighbour not yet refined.
1115  labelList testFaces(getRefineCandidateFaces(refineCell));
1116 
1117  // Collect segments
1118  pointField start(testFaces.size());
1119  pointField end(testFaces.size());
1120  labelList minLevel(testFaces.size());
1121 
1122  forAll(testFaces, i)
1123  {
1124  const label facei = testFaces[i];
1125  const label own = mesh_.faceOwner()[facei];
1126 
1127  if (mesh_.isInternalFace(facei))
1128  {
1129  const label nei = mesh_.faceNeighbour()[facei];
1130 
1131  start[i] = cellCentres[own];
1132  end[i] = cellCentres[nei];
1133  minLevel[i] = min(cellLevel[own], cellLevel[nei]);
1134  }
1135  else
1136  {
1137  const label bFacei = facei - mesh_.nInternalFaces();
1138 
1139  start[i] = cellCentres[own];
1140  end[i] = neiCc[bFacei];
1141 
1142  if (!isMasterFace[facei])
1143  {
1144  Swap(start[i], end[i]);
1145  }
1146 
1147  minLevel[i] = min(cellLevel[own], neiLevel[bFacei]);
1148  }
1149  }
1150 
1151  // Extend segments a bit
1152  {
1153  const vectorField smallVec(rootSmall*(end-start));
1154  start -= smallVec;
1155  end += smallVec;
1156  }
1157 
1158 
1159  // Test for all intersections (with surfaces of higher max level than
1160  // minLevel) and cache per cell the interesting inter
1161  labelListList cellSurfLevels(mesh_.nCells());
1162  List<vectorList> cellSurfNormals(mesh_.nCells());
1163 
1164  {
1165  // Per segment the normals of the surfaces
1166  List<vectorList> surfaceNormal;
1167 
1168  // Per segment the list of levels of the surfaces
1169  labelListList surfaceLevel;
1170 
1171  surfaces_.findAllHigherIntersections
1172  (
1173  start,
1174  end,
1175  minLevel, // max level of surface has to be bigger
1176  // than min level of neighbouring cells
1177 
1178  surfaces_.maxLevel(),
1179 
1180  surfaceNormal,
1181  surfaceLevel
1182  );
1183 
1184 
1185  // Sort the data according to surface location. This will guarantee
1186  // that on coupled faces both sides visit the intersections in
1187  // the same order so will decide the same
1188  labelList visitOrder;
1189  forAll(surfaceNormal, pointi)
1190  {
1191  vectorList& pNormals = surfaceNormal[pointi];
1192  labelList& pLevel = surfaceLevel[pointi];
1193 
1194  sortedOrder(pNormals, visitOrder, normalLess(pNormals));
1195 
1196  pNormals = List<point>(pNormals, visitOrder);
1197  pLevel = UIndirectList<label>(pLevel, visitOrder);
1198  }
1199 
1200  // Clear out unnecessary data
1201  start.clear();
1202  end.clear();
1203  minLevel.clear();
1204 
1205  // Convert face-wise data to cell.
1206  forAll(testFaces, i)
1207  {
1208  const label facei = testFaces[i];
1209  const label own = mesh_.faceOwner()[facei];
1210 
1211  const vectorList& fNormals = surfaceNormal[i];
1212  const labelList& fLevels = surfaceLevel[i];
1213 
1214  forAll(fNormals, hiti)
1215  {
1216  if (fLevels[hiti] > cellLevel[own])
1217  {
1218  cellSurfLevels[own].append(fLevels[hiti]);
1219  cellSurfNormals[own].append(fNormals[hiti]);
1220  }
1221 
1222  if (mesh_.isInternalFace(facei))
1223  {
1224  label nei = mesh_.faceNeighbour()[facei];
1225  if (fLevels[hiti] > cellLevel[nei])
1226  {
1227  cellSurfLevels[nei].append(fLevels[hiti]);
1228  cellSurfNormals[nei].append(fNormals[hiti]);
1229  }
1230  }
1231  }
1232  }
1233  }
1234 
1235 
1236 
1237  // Bit of statistics
1238  if (debug)
1239  {
1240  label nSet = 0;
1241  label nNormals = 9;
1242  forAll(cellSurfNormals, celli)
1243  {
1244  const vectorList& normals = cellSurfNormals[celli];
1245  if (normals.size())
1246  {
1247  nSet++;
1248  nNormals += normals.size();
1249  }
1250  }
1251  reduce(nSet, sumOp<label>());
1252  reduce(nNormals, sumOp<label>());
1253  Info<< "markSurfaceCurvatureRefinement :"
1254  << " cells:" << mesh_.globalData().nTotalCells()
1255  << " of which with normals:" << nSet
1256  << " ; total normals stored:" << nNormals
1257  << endl;
1258  }
1259 
1260 
1261 
1262  bool reachedLimit = false;
1263 
1264 
1265  // 1. Check normals per cell
1266  // ~~~~~~~~~~~~~~~~~~~~~~~~~
1267 
1268  for
1269  (
1270  label celli = 0;
1271  !reachedLimit && celli < cellSurfNormals.size();
1272  celli++
1273  )
1274  {
1275  const vectorList& normals = cellSurfNormals[celli];
1276  const labelList& levels = cellSurfLevels[celli];
1277 
1278  // n^2 comparison of all normals in a cell
1279  for (label i = 0; !reachedLimit && i < normals.size(); i++)
1280  {
1281  for (label j = i+1; !reachedLimit && j < normals.size(); j++)
1282  {
1283  if ((normals[i] & normals[j]) < curvature)
1284  {
1285  const label maxLevel = max(levels[i], levels[j]);
1286 
1287  if (cellLevel[celli] < maxLevel)
1288  {
1289  if
1290  (
1291  !markForRefine
1292  (
1293  maxLevel,
1294  nAllowRefine,
1295  refineCell[celli],
1296  nRefine
1297  )
1298  )
1299  {
1300  if (debug)
1301  {
1302  Pout<< "Stopped refining since reaching my cell"
1303  << " limit of " << mesh_.nCells()+7*nRefine
1304  << endl;
1305  }
1306  reachedLimit = true;
1307  break;
1308  }
1309  }
1310  }
1311  }
1312  }
1313  }
1314 
1315 
1316 
1317  // 2. Find out a measure of surface curvature
1318  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1319  // Look at normals between neighbouring surfaces
1320  // Loop over all faces. Could only be checkFaces, except if they're coupled
1321 
1322  // Internal faces
1323  for
1324  (
1325  label facei = 0;
1326  !reachedLimit && facei < mesh_.nInternalFaces();
1327  facei++
1328  )
1329  {
1330  const label own = mesh_.faceOwner()[facei];
1331  const label nei = mesh_.faceNeighbour()[facei];
1332 
1333  const vectorList& ownNormals = cellSurfNormals[own];
1334  const labelList& ownLevels = cellSurfLevels[own];
1335  const vectorList& neiNormals = cellSurfNormals[nei];
1336  const labelList& neiLevels = cellSurfLevels[nei];
1337 
1338 
1339  // Special case: owner normals set is a subset of the neighbour
1340  // normals. Do not do curvature refinement since other cell's normals
1341  // do not add any information. Avoids weird corner extensions of
1342  // refinement regions.
1343  bool ownisSubset =
1344  countMatches(ownNormals, neiNormals)
1345  == ownNormals.size();
1346 
1347  bool neiisSubset =
1348  countMatches(neiNormals, ownNormals)
1349  == neiNormals.size();
1350 
1351 
1352  if (!ownisSubset && !neiisSubset)
1353  {
1354  // n^2 comparison of between ownNormals and neiNormals
1355  for (label i = 0; !reachedLimit && i < ownNormals.size(); i++)
1356  {
1357  for (label j = 0; !reachedLimit && j < neiNormals.size(); j++)
1358  {
1359  // Have valid data on both sides. Check curvature.
1360  if ((ownNormals[i] & neiNormals[j]) < curvature)
1361  {
1362  // See which side to refine.
1363  if (cellLevel[own] < ownLevels[i])
1364  {
1365  if
1366  (
1367  !markForRefine
1368  (
1369  ownLevels[i],
1370  nAllowRefine,
1371  refineCell[own],
1372  nRefine
1373  )
1374  )
1375  {
1376  if (debug)
1377  {
1378  Pout<< "Stopped refining since reaching"
1379  << " my cell limit of "
1380  << mesh_.nCells()+7*nRefine << endl;
1381  }
1382  reachedLimit = true;
1383  break;
1384  }
1385  }
1386  if (cellLevel[nei] < neiLevels[j])
1387  {
1388  if
1389  (
1390  !markForRefine
1391  (
1392  neiLevels[j],
1393  nAllowRefine,
1394  refineCell[nei],
1395  nRefine
1396  )
1397  )
1398  {
1399  if (debug)
1400  {
1401  Pout<< "Stopped refining since reaching"
1402  << " my cell limit of "
1403  << mesh_.nCells()+7*nRefine << endl;
1404  }
1405  reachedLimit = true;
1406  break;
1407  }
1408  }
1409  }
1410  }
1411  }
1412  }
1413  }
1414 
1415 
1416  // Send over surface normal to neighbour cell.
1417  List<vectorList> neiSurfaceNormals;
1418  syncTools::swapBoundaryCellList(mesh_, cellSurfNormals, neiSurfaceNormals);
1419 
1420  // Boundary faces
1421  for
1422  (
1423  label facei = mesh_.nInternalFaces();
1424  !reachedLimit && facei < mesh_.nFaces();
1425  facei++
1426  )
1427  {
1428  label own = mesh_.faceOwner()[facei];
1429  label bFacei = facei - mesh_.nInternalFaces();
1430 
1431  const vectorList& ownNormals = cellSurfNormals[own];
1432  const labelList& ownLevels = cellSurfLevels[own];
1433  const vectorList& neiNormals = neiSurfaceNormals[bFacei];
1434 
1435  // Special case: owner normals set is a subset of the neighbour
1436  // normals. Do not do curvature refinement since other cell's normals
1437  // do not add any information. Avoids weird corner extensions of
1438  // refinement regions.
1439  bool ownisSubset =
1440  countMatches(ownNormals, neiNormals)
1441  == ownNormals.size();
1442 
1443  bool neiisSubset =
1444  countMatches(neiNormals, ownNormals)
1445  == neiNormals.size();
1446 
1447 
1448  if (!ownisSubset && !neiisSubset)
1449  {
1450  // n^2 comparison of between ownNormals and neiNormals
1451  for (label i = 0; !reachedLimit && i < ownNormals.size(); i++)
1452  {
1453  for (label j = 0; !reachedLimit && j < neiNormals.size(); j++)
1454  {
1455  // Have valid data on both sides. Check curvature.
1456  if ((ownNormals[i] & neiNormals[j]) < curvature)
1457  {
1458  if (cellLevel[own] < ownLevels[i])
1459  {
1460  if
1461  (
1462  !markForRefine
1463  (
1464  ownLevels[i],
1465  nAllowRefine,
1466  refineCell[own],
1467  nRefine
1468  )
1469  )
1470  {
1471  if (debug)
1472  {
1473  Pout<< "Stopped refining since reaching"
1474  << " my cell limit of "
1475  << mesh_.nCells()+7*nRefine
1476  << endl;
1477  }
1478  reachedLimit = true;
1479  break;
1480  }
1481  }
1482  }
1483  }
1484  }
1485  }
1486  }
1487 
1488 
1489  if
1490  (
1491  returnReduce(nRefine, sumOp<label>())
1492  > returnReduce(nAllowRefine, sumOp<label>())
1493  )
1494  {
1495  Info<< "Reached refinement limit." << endl;
1496  }
1497 
1498  return returnReduce(nRefine-oldNRefine, sumOp<label>());
1499 }
1500 
1501 
1503 (
1504  const scalar planarCos,
1505  const vector& point0,
1506  const vector& normal0,
1507 
1508  const vector& point1,
1509  const vector& normal1
1510 ) const
1511 {
1512  //- Hits differ and angles are oppositeish and
1513  // hits have a normal distance
1514  const vector d = point1 - point0;
1515  const scalar magD = mag(d);
1516 
1517  if (magD > mergeDistance())
1518  {
1519  scalar cosAngle = (normal0 & normal1);
1520 
1521  vector avg = Zero;
1522  if (cosAngle < (-1+planarCos))
1523  {
1524  // Opposite normals
1525  avg = 0.5*(normal0-normal1);
1526  }
1527  else if (cosAngle > (1-planarCos))
1528  {
1529  avg = 0.5*(normal0+normal1);
1530  }
1531 
1532  if (avg != vector::zero)
1533  {
1534  avg /= mag(avg);
1535 
1536  // Check normal distance of intersection locations
1537  if (mag(avg&d) > mergeDistance())
1538  {
1539  return true;
1540  }
1541  else
1542  {
1543  return false;
1544  }
1545  }
1546  else
1547  {
1548  return false;
1549  }
1550  }
1551  else
1552  {
1553  return false;
1554  }
1555 }
1556 
1557 
1558 // Mark small gaps
1560 (
1561  const scalar planarCos,
1562  const vector& point0,
1563  const vector& normal0,
1564 
1565  const vector& point1,
1566  const vector& normal1
1567 ) const
1568 {
1569  //- Hits differ and angles are oppositeish and
1570  // hits have a normal distance
1571  vector d = point1 - point0;
1572  const scalar magD = mag(d);
1573 
1574  if (magD > mergeDistance())
1575  {
1576  scalar cosAngle = (normal0 & normal1);
1577 
1578  vector avg = Zero;
1579  if (cosAngle < (-1+planarCos))
1580  {
1581  // Opposite normals
1582  avg = 0.5*(normal0-normal1);
1583  }
1584  else if (cosAngle > (1-planarCos))
1585  {
1586  avg = 0.5*(normal0+normal1);
1587  }
1588 
1589  if (avg != vector::zero)
1590  {
1591  avg /= mag(avg);
1592  d /= magD;
1593 
1594  // Check average normal with respect to intersection locations
1595  if (mag(avg&d) > Foam::cos(degToRad(45.0)))
1596  {
1597  return true;
1598  }
1599  else
1600  {
1601  return false;
1602  }
1603  }
1604  else
1605  {
1606  return false;
1607  }
1608  }
1609  else
1610  {
1611  return false;
1612  }
1613 }
1614 
1615 
1616 bool Foam::meshRefinement::checkProximity
1617 (
1618  const scalar planarCos,
1619  const label nAllowRefine,
1620 
1621  const label surfaceLevel, // current intersection max level
1622  const vector& surfaceLocation, // current intersection location
1623  const vector& surfaceNormal, // current intersection normal
1624 
1625  const label celli,
1626 
1627  label& cellMaxLevel, // cached max surface level for this cell
1628  vector& cellMaxLocation, // cached surface normal for this cell
1629  vector& cellMaxNormal, // cached surface normal for this cell
1630 
1632  label& nRefine
1633 ) const
1634 {
1635  const labelList& cellLevel = meshCutter_.cellLevel();
1636 
1637  // Test if surface applicable
1638  if (surfaceLevel > cellLevel[celli])
1639  {
1640  if (cellMaxLevel == -1)
1641  {
1642  // First visit of cell. Store
1643  cellMaxLevel = surfaceLevel;
1644  cellMaxLocation = surfaceLocation;
1645  cellMaxNormal = surfaceNormal;
1646  }
1647  else
1648  {
1649  // Second or more visit.
1650  // Check if
1651  // - different location
1652  // - opposite surface
1653 
1654  bool closeSurfaces = isNormalGap
1655  (
1656  planarCos,
1657  cellMaxLocation,
1658  cellMaxNormal,
1659  surfaceLocation,
1660  surfaceNormal
1661  );
1662 
1663  // Set normal to that of highest surface. Not really necessary
1664  // over here but we reuse cellMax info when doing coupled faces.
1665  if (surfaceLevel > cellMaxLevel)
1666  {
1667  cellMaxLevel = surfaceLevel;
1668  cellMaxLocation = surfaceLocation;
1669  cellMaxNormal = surfaceNormal;
1670  }
1671 
1672 
1673  if (closeSurfaces)
1674  {
1675  // Pout<< "Found gap:" << nl
1676  // << " location:" << surfaceLocation
1677  // << "\tnormal :" << surfaceNormal << nl
1679  // << "\tnormal :" << cellMaxNormal << nl
1680  // << "\tcos :" << (surfaceNormal&cellMaxNormal) << nl
1681  // << endl;
1682 
1683  return markForRefine
1684  (
1685  surfaceLevel, // mark with any non-neg number.
1686  nAllowRefine,
1687  refineCell[celli],
1688  nRefine
1689  );
1690  }
1691  }
1692  }
1693 
1694  // Did not reach refinement limit.
1695  return true;
1696 }
1697 
1698 
1699 Foam::label Foam::meshRefinement::markProximityRefinement
1700 (
1701  const scalar planarCos,
1702  const label nAllowRefine,
1703  const labelList& neiLevel,
1704  const pointField& neiCc,
1705 
1706  labelList& refineCell,
1707  label& nRefine
1708 ) const
1709 {
1710  const labelList& cellLevel = meshCutter_.cellLevel();
1711  const pointField& cellCentres = mesh_.cellCentres();
1712 
1713  const label oldNRefine = nRefine;
1714 
1715  // 1. local test: any cell on more than one surface gets refined
1716  // (if its current level is < max of the surface max level)
1717 
1718  // 2. 'global' test: any cell on only one surface with a neighbour
1719  // on a different surface gets refined (if its current level etc.)
1720 
1721 
1722  // Collect candidate faces (i.e. intersecting any surface and
1723  // owner/neighbour not yet refined.
1724  labelList testFaces(getRefineCandidateFaces(refineCell));
1725 
1726  // Collect segments
1727  pointField start(testFaces.size());
1728  pointField end(testFaces.size());
1729  labelList minLevel(testFaces.size());
1730 
1731  forAll(testFaces, i)
1732  {
1733  const label facei = testFaces[i];
1734  const label own = mesh_.faceOwner()[facei];
1735 
1736  if (mesh_.isInternalFace(facei))
1737  {
1738  const label nei = mesh_.faceNeighbour()[facei];
1739 
1740  start[i] = cellCentres[own];
1741  end[i] = cellCentres[nei];
1742  minLevel[i] = min(cellLevel[own], cellLevel[nei]);
1743  }
1744  else
1745  {
1746  const label bFacei = facei - mesh_.nInternalFaces();
1747 
1748  start[i] = cellCentres[own];
1749  end[i] = neiCc[bFacei];
1750  minLevel[i] = min(cellLevel[own], neiLevel[bFacei]);
1751  }
1752  }
1753 
1754  // Extend segments a bit
1755  {
1756  const vectorField smallVec(rootSmall*(end-start));
1757  start -= smallVec;
1758  end += smallVec;
1759  }
1760 
1761 
1762  // Test for all intersections (with surfaces of higher gap level than
1763  // minLevel) and cache per cell the max surface level and the local normal
1764  // on that surface.
1765  labelList cellMaxLevel(mesh_.nCells(), -1);
1766  vectorField cellMaxNormal(mesh_.nCells(), Zero);
1767  pointField cellMaxLocation(mesh_.nCells(), Zero);
1768 
1769  {
1770  // Per segment the normals of the surfaces
1771  List<vectorList> surfaceLocation;
1772  List<vectorList> surfaceNormal;
1773  // Per segment the list of levels of the surfaces
1774  labelListList surfaceLevel;
1775 
1776  surfaces_.findAllHigherIntersections
1777  (
1778  start,
1779  end,
1780  minLevel, // gap level of surface has to be bigger
1781  // than min level of neighbouring cells
1782 
1783  surfaces_.gapLevel(),
1784 
1785  surfaceLocation,
1786  surfaceNormal,
1787  surfaceLevel
1788  );
1789  // Clear out unnecessary data
1790  start.clear();
1791  end.clear();
1792  minLevel.clear();
1793 
1795  // OBJstream str
1796  //(
1797  // mesh_.time().path()
1798  // / "findAllHigherIntersections_"
1799  // + mesh_.time().name()
1800  // + ".obj"
1801  //);
1803  // OBJstream str2
1804  //(
1805  // mesh_.time().path()
1806  // / "findAllHigherIntersections2_"
1807  // + mesh_.time().name()
1808  // + ".obj"
1809  //);
1810 
1811  forAll(testFaces, i)
1812  {
1813  label facei = testFaces[i];
1814  label own = mesh_.faceOwner()[facei];
1815 
1816  const labelList& fLevels = surfaceLevel[i];
1817  const vectorList& fPoints = surfaceLocation[i];
1818  const vectorList& fNormals = surfaceNormal[i];
1819 
1820  forAll(fLevels, hiti)
1821  {
1822  checkProximity
1823  (
1824  planarCos,
1825  nAllowRefine,
1826 
1827  fLevels[hiti],
1828  fPoints[hiti],
1829  fNormals[hiti],
1830 
1831  own,
1832  cellMaxLevel[own],
1833  cellMaxLocation[own],
1834  cellMaxNormal[own],
1835 
1836  refineCell,
1837  nRefine
1838  );
1839  }
1840 
1841  if (mesh_.isInternalFace(facei))
1842  {
1843  label nei = mesh_.faceNeighbour()[facei];
1844 
1845  forAll(fLevels, hiti)
1846  {
1847  checkProximity
1848  (
1849  planarCos,
1850  nAllowRefine,
1851 
1852  fLevels[hiti],
1853  fPoints[hiti],
1854  fNormals[hiti],
1855 
1856  nei,
1857  cellMaxLevel[nei],
1858  cellMaxLocation[nei],
1859  cellMaxNormal[nei],
1860 
1861  refineCell,
1862  nRefine
1863  );
1864  }
1865  }
1866  }
1867  }
1868 
1869  // 2. Find out a measure of surface curvature and region edges.
1870  // Send over surface region and surface normal to neighbour cell.
1871 
1872  labelList neiBndMaxLevel(mesh_.nFaces()-mesh_.nInternalFaces());
1873  pointField neiBndMaxLocation(mesh_.nFaces()-mesh_.nInternalFaces());
1874  vectorField neiBndMaxNormal(mesh_.nFaces()-mesh_.nInternalFaces());
1875 
1876  for (label facei = mesh_.nInternalFaces(); facei < mesh_.nFaces(); facei++)
1877  {
1878  const label bFacei = facei-mesh_.nInternalFaces();
1879  const label own = mesh_.faceOwner()[facei];
1880 
1881  neiBndMaxLevel[bFacei] = cellMaxLevel[own];
1882  neiBndMaxLocation[bFacei] = cellMaxLocation[own];
1883  neiBndMaxNormal[bFacei] = cellMaxNormal[own];
1884  }
1885  syncTools::swapBoundaryFaceList(mesh_, neiBndMaxLevel);
1886  syncTools::swapBoundaryFaceList(mesh_, neiBndMaxLocation);
1887  syncTools::swapBoundaryFaceList(mesh_, neiBndMaxNormal);
1888 
1889  // Loop over all faces. Could only be checkFaces.. except if they're coupled
1890 
1891  // Internal faces
1892  for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
1893  {
1894  const label own = mesh_.faceOwner()[facei];
1895  const label nei = mesh_.faceNeighbour()[facei];
1896 
1897  if (cellMaxLevel[own] != -1 && cellMaxLevel[nei] != -1)
1898  {
1899  // Have valid data on both sides. Check planarCos.
1900  if
1901  (
1902  isNormalGap
1903  (
1904  planarCos,
1905  cellMaxLocation[own],
1906  cellMaxNormal[own],
1907  cellMaxLocation[nei],
1908  cellMaxNormal[nei]
1909  )
1910  )
1911  {
1912  // See which side to refine
1913  if (cellLevel[own] < cellMaxLevel[own])
1914  {
1915  if
1916  (
1917  !markForRefine
1918  (
1919  cellMaxLevel[own],
1920  nAllowRefine,
1921  refineCell[own],
1922  nRefine
1923  )
1924  )
1925  {
1926  if (debug)
1927  {
1928  Pout<< "Stopped refining since reaching my cell"
1929  << " limit of " << mesh_.nCells()+7*nRefine
1930  << endl;
1931  }
1932  break;
1933  }
1934  }
1935 
1936  if (cellLevel[nei] < cellMaxLevel[nei])
1937  {
1938  if
1939  (
1940  !markForRefine
1941  (
1942  cellMaxLevel[nei],
1943  nAllowRefine,
1944  refineCell[nei],
1945  nRefine
1946  )
1947  )
1948  {
1949  if (debug)
1950  {
1951  Pout<< "Stopped refining since reaching my cell"
1952  << " limit of " << mesh_.nCells()+7*nRefine
1953  << endl;
1954  }
1955  break;
1956  }
1957  }
1958  }
1959  }
1960  }
1961  // Boundary faces
1962  for (label facei = mesh_.nInternalFaces(); facei < mesh_.nFaces(); facei++)
1963  {
1964  const label own = mesh_.faceOwner()[facei];
1965  const label bFacei = facei - mesh_.nInternalFaces();
1966 
1967  if (cellLevel[own] < cellMaxLevel[own] && neiBndMaxLevel[bFacei] != -1)
1968  {
1969  // Have valid data on both sides. Check planarCos.
1970  if
1971  (
1972  isNormalGap
1973  (
1974  planarCos,
1975  cellMaxLocation[own],
1976  cellMaxNormal[own],
1977  neiBndMaxLocation[bFacei],
1978  neiBndMaxNormal[bFacei]
1979  )
1980  )
1981  {
1982  if
1983  (
1984  !markForRefine
1985  (
1986  cellMaxLevel[own],
1987  nAllowRefine,
1988  refineCell[own],
1989  nRefine
1990  )
1991  )
1992  {
1993  if (debug)
1994  {
1995  Pout<< "Stopped refining since reaching my cell"
1996  << " limit of " << mesh_.nCells()+7*nRefine
1997  << endl;
1998  }
1999  break;
2000  }
2001  }
2002  }
2003  }
2004 
2005  if
2006  (
2007  returnReduce(nRefine, sumOp<label>())
2008  > returnReduce(nAllowRefine, sumOp<label>())
2009  )
2010  {
2011  Info<< "Reached refinement limit." << endl;
2012  }
2013 
2014  return returnReduce(nRefine-oldNRefine, sumOp<label>());
2015 }
2016 
2017 
2018 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
2019 
2020 // Calculate list of cells to refine. Gets for any edge (start - end)
2021 // whether it intersects the surface. Does more accurate test and checks
2022 // the wanted level on the surface intersected.
2023 // Does approximate precalculation of how many cells can be refined before
2024 // hitting overall limit maxGlobalCells.
2026 (
2027  const List<point>& insidePoints,
2028  const scalar curvature,
2029  const scalar planarAngle,
2030 
2031  const bool featureRefinement,
2032  const bool featureDistanceRefinement,
2033  const bool internalRefinement,
2034  const bool surfaceRefinement,
2035  const bool curvatureRefinement,
2036  const bool gapRefinement,
2037  const label maxGlobalCells,
2038  const label maxLocalCells
2039 ) const
2040 {
2041  const label totNCells = mesh_.globalData().nTotalCells();
2042 
2043  labelList cellsToRefine;
2044 
2045  if (totNCells >= maxGlobalCells)
2046  {
2047  Info<< "No cells marked for refinement since reached limit "
2048  << maxGlobalCells << '.' << endl;
2049  }
2050  else
2051  {
2052  // Every cell I refine adds 7 cells. Estimate the number of cells
2053  // I am allowed to refine.
2054  // Assume perfect distribution so can only refine as much the fraction
2055  // of the mesh I hold. This prediction step prevents us having to do
2056  // lots of reduces to keep count of the total number of cells selected
2057  // for refinement.
2058 
2059  // scalar fraction = scalar(mesh_.nCells())/totNCells;
2060  // label myMaxCells = label(maxGlobalCells*fraction);
2061  // label nAllowRefine = (myMaxCells - mesh_.nCells())/7;
2063  //
2064  // Pout<< "refineCandidates:" << nl
2065  // << " total cells:" << totNCells << nl
2066  // << " local cells:" << mesh_.nCells() << nl
2067  // << " local fraction:" << fraction << nl
2068  // << " local allowable cells:" << myMaxCells << nl
2069  // << " local allowable refinement:" << nAllowRefine << nl
2070  // << endl;
2071 
2072  //- Disable refinement shortcut. nAllowRefine is per processor limit.
2073  const label nAllowRefine = labelMax / Pstream::nProcs();
2074 
2075  // Marked for refinement (>= 0) or not (-1). Actual value is the
2076  // index of the surface it intersects.
2077  labelList refineCell(mesh_.nCells(), -1);
2078  label nRefine = 0;
2079 
2080 
2081  // Swap neighbouring cell centres and cell level
2082  labelList neiLevel(mesh_.nFaces()-mesh_.nInternalFaces());
2083  pointField neiCc(mesh_.nFaces()-mesh_.nInternalFaces());
2084  calcNeighbourData(neiLevel, neiCc);
2085 
2086 
2087 
2088  // Cells pierced by feature edges
2089  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2090 
2091  if (featureRefinement)
2092  {
2093  const label nFeatures = markFeatureRefinement
2094  (
2095  insidePoints,
2096  nAllowRefine,
2097 
2098  refineCell,
2099  nRefine
2100  );
2101 
2102  Info<< "Marked for refinement due to explicit features "
2103  << ": " << nFeatures << " cells." << endl;
2104  }
2105 
2106  // Inside distance-to-feature
2107  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
2108 
2109  if (featureDistanceRefinement)
2110  {
2111  const label nCellsFeat = markInternalDistanceToFeatureRefinement
2112  (
2113  nAllowRefine,
2114 
2115  refineCell,
2116  nRefine
2117  );
2118  Info<< "Marked for refinement due to distance to explicit features "
2119  ": " << nCellsFeat << " cells." << endl;
2120  }
2121 
2122  // Inside refinement shells
2123  // ~~~~~~~~~~~~~~~~~~~~~~~~
2124 
2125  if (internalRefinement)
2126  {
2127  const label nShell = markInternalRefinement
2128  (
2129  nAllowRefine,
2130 
2131  refineCell,
2132  nRefine
2133  );
2134  Info<< "Marked for refinement due to refinement shells "
2135  << ": " << nShell << " cells." << endl;
2136  }
2137 
2138  // Refinement based on intersection of surface
2139  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2140 
2141  if (surfaceRefinement)
2142  {
2143  const label nSurf = markSurfaceRefinement
2144  (
2145  nAllowRefine,
2146  neiLevel,
2147  neiCc,
2148 
2149  refineCell,
2150  nRefine
2151  );
2152  Info<< "Marked for refinement due to surface intersection "
2153  << ": " << nSurf << " cells." << endl;
2154  }
2155 
2156  // Refinement based on curvature of surface
2157  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2158 
2159  if
2160  (
2161  curvatureRefinement
2162  && (curvature >= -1 && curvature <= 1)
2163  && (surfaces_.minLevel() != surfaces_.maxLevel())
2164  )
2165  {
2166  const label nCurv = markSurfaceCurvatureRefinement
2167  (
2168  curvature,
2169  nAllowRefine,
2170  neiLevel,
2171  neiCc,
2172 
2173  refineCell,
2174  nRefine
2175  );
2176  Info<< "Marked for refinement due to curvature/regions "
2177  << ": " << nCurv << " cells." << endl;
2178  }
2179 
2180 
2181  const scalar planarCos = Foam::cos(degToRad(planarAngle));
2182 
2183  if
2184  (
2185  gapRefinement
2186  && (planarCos >= -1 && planarCos <= 1)
2187  && (max(surfaces_.gapLevel()) > -1)
2188  )
2189  {
2190  Info<< "Specified gap level : " << max(surfaces_.gapLevel())
2191  << ", planar angle " << planarAngle << endl;
2192 
2193  const label nGap = markProximityRefinement
2194  (
2195  planarCos,
2196  nAllowRefine,
2197  neiLevel,
2198  neiCc,
2199 
2200  refineCell,
2201  nRefine
2202  );
2203  Info<< "Marked for refinement due to close opposite surfaces "
2204  << ": " << nGap << " cells." << endl;
2205  }
2206 
2207 
2208 
2209  // Pack cells-to-refine
2210  // ~~~~~~~~~~~~~~~~~~~~
2211 
2212  cellsToRefine.setSize(nRefine);
2213  nRefine = 0;
2214 
2215  forAll(refineCell, celli)
2216  {
2217  if (refineCell[celli] != -1)
2218  {
2219  cellsToRefine[nRefine++] = celli;
2220  }
2221  }
2222  }
2223 
2224  return cellsToRefine;
2225 }
2226 
2227 
2229 (
2230  const labelList& cellsToRefine
2231 )
2232 {
2233  // Mesh changing engine.
2234  polyTopoChange meshMod(mesh_);
2235 
2236  // Play refinement commands into mesh changer.
2237  meshCutter_.setRefinement(cellsToRefine, meshMod);
2238 
2239  // Create mesh (no inflation), return map from old to new mesh.
2240  autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(mesh_, false);
2241 
2242  // Update fields
2243  mesh_.topoChange(map);
2244 
2245  // Optionally inflate mesh
2246  if (map().hasMotionPoints())
2247  {
2248  mesh_.movePoints(map().preMotionPoints());
2249  }
2250  else
2251  {
2252  // Delete mesh volumes.
2253  mesh_.clearOut();
2254  }
2255 
2256  // Reset the instance for if in overwrite mode
2257  mesh_.setInstance(name());
2258 
2259  // Update intersection info
2260  topoChange(map, getChangedFaces(map, cellsToRefine));
2261 
2262  return map;
2263 }
2264 
2265 
2268 (
2269  const string& msg,
2270  decompositionMethod& decomposer,
2271  fvMeshDistribute& distributor,
2272  const labelList& cellsToRefine,
2273  const scalar maxLoadUnbalance
2274 )
2275 {
2276  // Do all refinement
2277  refine(cellsToRefine);
2278 
2279  if (debug&meshRefinement::MESH)
2280  {
2281  Pout<< "Writing refined but unbalanced " << msg
2282  << " mesh to time " << name() << endl;
2283  write
2284  (
2285  debugType(debug),
2286  writeType(writeLevel() | WRITEMESH),
2287  mesh_.time().path()/name()
2288  );
2289  Pout<< "Dumped debug data in = "
2290  << mesh_.time().cpuTimeIncrement() << " s" << endl;
2291 
2292  // test all is still synced across proc patches
2293  checkData();
2294  }
2295 
2296  Info<< "Refined mesh in = "
2297  << mesh_.time().cpuTimeIncrement() << " s" << endl;
2298  printMeshInfo(debug, "After refinement " + msg);
2299 
2300 
2301  // Load balancing
2302  // ~~~~~~~~~~~~~~
2303 
2305 
2306  if (Pstream::nProcs() > 1)
2307  {
2308  const scalar nIdealCells =
2309  mesh_.globalData().nTotalCells()
2310  / Pstream::nProcs();
2311 
2312  const scalar unbalance = returnReduce
2313  (
2314  mag(1.0-mesh_.nCells()/nIdealCells),
2315  maxOp<scalar>()
2316  );
2317 
2318  if (unbalance <= maxLoadUnbalance)
2319  {
2320  Info<< "Skipping balancing since max unbalance " << unbalance
2321  << " is less than allowable " << maxLoadUnbalance
2322  << endl;
2323  }
2324  else
2325  {
2326  scalarField cellWeights(mesh_.nCells(), 1);
2327 
2328  distMap = balance
2329  (
2330  false, // keepZoneFaces
2331  false, // keepBaffles
2332  cellWeights,
2333  decomposer,
2334  distributor
2335  );
2336 
2337  Info<< "Balanced mesh in = "
2338  << mesh_.time().cpuTimeIncrement() << " s" << endl;
2339 
2340  printMeshInfo(debug, "After balancing " + msg);
2341 
2342 
2343  if (debug&meshRefinement::MESH)
2344  {
2345  Pout<< "Writing balanced " << msg
2346  << " mesh to time " << name() << endl;
2347  write
2348  (
2349  debugType(debug),
2350  writeType(writeLevel() | WRITEMESH),
2351  mesh_.time().path()/name()
2352  );
2353  Pout<< "Dumped debug data in = "
2354  << mesh_.time().cpuTimeIncrement() << " s" << endl;
2355 
2356  // test all is still synced across proc patches
2357  checkData();
2358  }
2359  }
2360  }
2361 
2362  return distMap;
2363 }
2364 
2365 
2366 // Do load balancing followed by refinement of consistent set of cells.
2369 (
2370  const string& msg,
2371  decompositionMethod& decomposer,
2372  fvMeshDistribute& distributor,
2373  const labelList& initCellsToRefine,
2374  const scalar maxLoadUnbalance
2375 )
2376 {
2377  labelList cellsToRefine(initCellsToRefine);
2378 
2379  //{
2380  // globalIndex globalCells(mesh_.nCells());
2381  //
2382  // Info<< "** Distribution before balancing/refining:" << endl;
2383  // for (label proci = 0; proci < Pstream::nProcs(); proci++)
2384  // {
2385  // Info<< " " << proci << '\t'
2386  // << globalCells.localSize(proci) << endl;
2387  // }
2388  // Info<< endl;
2389  //}
2390  //{
2391  // globalIndex globalCells(cellsToRefine.size());
2392  //
2393  // Info<< "** Cells to be refined:" << endl;
2394  // for (label proci = 0; proci < Pstream::nProcs(); proci++)
2395  // {
2396  // Info<< " " << proci << '\t'
2397  // << globalCells.localSize(proci) << endl;
2398  // }
2399  // Info<< endl;
2400  //}
2401 
2402 
2403  // Load balancing
2404  // ~~~~~~~~~~~~~~
2405 
2407 
2408  if (Pstream::nProcs() > 1)
2409  {
2410  // First check if we need to balance at all. Precalculate number of
2411  // cells after refinement and see what maximum difference is.
2412  const scalar nNewCells =
2413  scalar(mesh_.nCells() + 7*cellsToRefine.size());
2414  const scalar nIdealNewCells =
2415  returnReduce(nNewCells, sumOp<scalar>())/Pstream::nProcs();
2416  const scalar unbalance = returnReduce
2417  (
2418  mag(1.0-nNewCells/nIdealNewCells),
2419  maxOp<scalar>()
2420  );
2421 
2422  if (unbalance <= maxLoadUnbalance)
2423  {
2424  Info<< "Skipping balancing since max unbalance " << unbalance
2425  << " is less than allowable " << maxLoadUnbalance
2426  << endl;
2427  }
2428  else
2429  {
2430  scalarField cellWeights(mesh_.nCells(), 1);
2431  forAll(cellsToRefine, i)
2432  {
2433  cellWeights[cellsToRefine[i]] += 7;
2434  }
2435 
2436  distMap = balance
2437  (
2438  false, // keepZoneFaces
2439  false, // keepBaffles
2440  cellWeights,
2441  decomposer,
2442  distributor
2443  );
2444 
2445  // Update cells to refine
2446  distMap().distributeCellIndices(cellsToRefine);
2447 
2448  Info<< "Balanced mesh in = "
2449  << mesh_.time().cpuTimeIncrement() << " s" << endl;
2450  }
2451 
2452  //{
2453  // globalIndex globalCells(mesh_.nCells());
2454  //
2455  // Info<< "** Distribution after balancing:" << endl;
2456  // for (label proci = 0; proci < Pstream::nProcs(); proci++)
2457  // {
2458  // Info<< " " << proci << '\t'
2459  // << globalCells.localSize(proci) << endl;
2460  // }
2461  // Info<< endl;
2462  //}
2463 
2464  printMeshInfo(debug, "After balancing " + msg);
2465 
2466  if (debug&meshRefinement::MESH)
2467  {
2468  Pout<< "Writing balanced " << msg
2469  << " mesh to time " << name() << endl;
2470  write
2471  (
2472  debugType(debug),
2473  writeType(writeLevel() | WRITEMESH),
2474  mesh_.time().path()/name()
2475  );
2476  Pout<< "Dumped debug data in = "
2477  << mesh_.time().cpuTimeIncrement() << " s" << endl;
2478 
2479  // test all is still synced across proc patches
2480  checkData();
2481  }
2482  }
2483 
2484 
2485  // Refinement
2486  // ~~~~~~~~~~
2487 
2488  refine(cellsToRefine);
2489 
2490  if (debug&meshRefinement::MESH)
2491  {
2492  Pout<< "Writing refined " << msg
2493  << " mesh to time " << name() << endl;
2494  write
2495  (
2496  debugType(debug),
2497  writeType(writeLevel() | WRITEMESH),
2498  mesh_.time().path()/name()
2499  );
2500  Pout<< "Dumped debug data in = "
2501  << mesh_.time().cpuTimeIncrement() << " s" << endl;
2502 
2503  // test all is still synced across proc patches
2504  checkData();
2505  }
2506 
2507  Info<< "Refined mesh in = "
2508  << mesh_.time().cpuTimeIncrement() << " s" << endl;
2509 
2510  //{
2511  // globalIndex globalCells(mesh_.nCells());
2512  //
2513  // Info<< "** After refinement distribution:" << endl;
2514  // for (label proci = 0; proci < Pstream::nProcs(); proci++)
2515  // {
2516  // Info<< " " << proci << '\t'
2517  // << globalCells.localSize(proci) << endl;
2518  // }
2519  // Info<< endl;
2520  //}
2521 
2522  printMeshInfo(debug, "After refinement " + msg);
2523 
2524  return distMap;
2525 }
2526 
2527 
2528 // ************************************************************************* //
scalar y
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:125
void setSize(const label)
Reset size of List.
Definition: List.C:281
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:85
virtual Ostream & write(const char)=0
Write character.
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static const Form zero
Definition: VectorSpace.H:113
static const direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:99
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Abstract base class for decomposition.
Sends/receives parts of mesh+fvfields to neighbouring processors. Used in load balancing.
label nTotalFaces() const
Return total number of faces in decomposed mesh. Not.
autoPtr< polyDistributionMap > refineAndBalance(const string &msg, decompositionMethod &decomposer, fvMeshDistribute &distributor, const labelList &cellsToRefine, const scalar maxLoadUnbalance)
Refine some cells and rebalance.
labelList refineCandidates(const List< point > &insidePoints, const scalar curvature, const scalar planarAngle, const bool featureRefinement, const bool featureDistanceRefinement, const bool internalRefinement, const bool surfaceRefinement, const bool curvatureRefinement, const bool gapRefinement, const label maxGlobalCells, const label maxLocalCells) const
Calculate list of cells to refine.
const fvMesh & mesh() const
Reference to mesh.
autoPtr< polyTopoChangeMap > refine(const labelList &cellsToRefine)
Refine some cells.
bool isGap(const scalar, const vector &, const vector &, const vector &, const vector &) const
Is local topology a small gap?
autoPtr< polyDistributionMap > balanceAndRefine(const string &msg, decompositionMethod &decomposer, fvMeshDistribute &distributor, const labelList &cellsToRefine, const scalar maxLoadUnbalance)
Balance before refining some cells.
bool isNormalGap(const scalar, const vector &, const vector &, const vector &, const vector &) const
Is local topology a small gap normal to the test vector.
To compare normals.
normalLess(const vectorList &values)
bool operator()(const label a, const label b) const
labelList cmptType
Component type.
vectorList cmptType
Component type.
Traits class for primitives.
Definition: pTraits.H:53
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1387
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1563
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:403
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1393
Direct mesh changes based on v1.3 polyTopoChange syntax.
autoPtr< polyTopoChangeMap > changeMesh(polyMesh &mesh, const bool inflate, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Inplace changes mesh without change of patches.
label nInternalFaces() const
label nFaces() const
const cellList & cells() const
Container with cells to refine. Refinement given as single direction.
Definition: refineCell.H:57
Contains information about location on a triSurface.
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &l)
Swap coupled boundary face values.
Definition: syncTools.H:436
static PackedBoolList getMasterFaces(const polyMesh &)
Get per face whether it is uncoupled or a master of a.
Definition: syncTools.C:153
static void syncFaceList(const polyMesh &mesh, UList< T > &l, const CombineOp &cop)
Synchronise values on all mesh faces.
Definition: syncTools.H:387
static void syncBoundaryFaceList(const polyMesh &, UList< T > &, const CombineOp &cop, const TransformOp &top, const bool parRun=Pstream::parRun())
Synchronise values on boundary faces only.
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData)
Swap to obtain neighbour cell values for all boundary faces.
label patchi
const cellShapeList & cells
volScalarField & b
Definition: createFields.H:27
labelList getChangedFaces(const polyMesh &mesh, const labelHashSet &patchIDs)
Get initial set of changed faces.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
const doubleScalar e
Definition: doubleScalar.H:105
List< label > labelList
A List of labels.
Definition: labelList.H:56
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
List< cell > cellList
list of cells
Definition: cellList.H:42
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
labelList findIndices(const ListType &, typename ListType::const_reference, const label start=0)
Find all occurrences of given element. Linear search.
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
vector point
Point is a vector.
Definition: point.H:41
List< vector > vectorList
A List of vectors.
Definition: vectorList.H:50
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
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)
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Field< vector > vectorField
Specialisation of Field<T> for vector.
static bool less(const vector &x, const vector &y)
To compare normals.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
void Swap(T &a, T &b)
Definition: Swap.H:43
static const label labelMax
Definition: label.H:62
void sortedOrder(const UList< T > &, labelList &order)
Generate the (stable) sort order for the list.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dimensionedScalar cos(const dimensionedScalar &ds)
uint8_t direction
Definition: direction.H:45
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
insidePoints((1 2 3))