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