meshRefinement.H
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-2018 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 Class
25  Foam::meshRefinement
26 
27 Description
28  Helper class which maintains intersections of (changing) mesh with
29  (static) surfaces.
30 
31  Maintains
32  - per face any intersections of the cc-cc segment with any of the surfaces
33 
34 SourceFiles
35  meshRefinement.C
36  meshRefinementBaffles.C
37  meshRefinementMerge.C
38  meshRefinementProblemCells.C
39  meshRefinementRefine.C
40 
41 \*---------------------------------------------------------------------------*/
42 
43 #ifndef meshRefinement_H
44 #define meshRefinement_H
45 
46 #include "hexRef8.H"
47 #include "mapPolyMesh.H"
48 #include "autoPtr.H"
49 #include "labelPair.H"
50 #include "indirectPrimitivePatch.H"
51 #include "pointFieldsFwd.H"
52 #include "Tuple2.H"
53 #include "pointIndexHit.H"
54 
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 
57 namespace Foam
58 {
59 
60 // Class forward declarations
61 class fvMesh;
62 class mapDistributePolyMesh;
63 class decompositionMethod;
64 class refinementSurfaces;
65 class refinementFeatures;
66 class shellSurfaces;
67 class removeCells;
68 class fvMeshDistribute;
69 class searchableSurface;
70 class regionSplit;
71 class globalIndex;
72 class removePoints;
73 class localPointRegion;
74 
75 class snapParameters;
76 
77 /*---------------------------------------------------------------------------*\
78  Class meshRefinement Declaration
79 \*---------------------------------------------------------------------------*/
80 
81 class meshRefinement
82 {
83 
84 public:
85 
86  // Public data types
87 
88  //- Enumeration for what to debug
89  enum IOdebugType
90  {
96  };
97 
99  enum debugType
100  {
101  MESH = 1<<IOMESH,
105  LAYERINFO = 1<<IOLAYERINFO
106  };
107 
108  //- Enumeration for what to output
109  enum IOoutputType
110  {
112  };
115  enum outputType
116  {
118  };
119 
120  //- Enumeration for what to write
121  enum IOwriteType
122  {
128  };
131  enum writeType
132  {
137  WRITELAYERFIELDS = 1<<IOWRITELAYERFIELDS
138  };
139 
140  //- Enumeration for how the userdata is to be mapped upon refinement.
141  enum mapType
142  {
144  KEEPALL = 2,
145  REMOVE = 4
146  };
147 
148 
149 private:
150 
151  // Static data members
152 
153  //- Control of writing level
154  static writeType writeLevel_;
155 
156  //- Control of output/log level
157  static outputType outputLevel_;
158 
159 
160  // Private data
161 
162  //- Reference to mesh
163  fvMesh& mesh_;
164 
165  //- Tolerance used for sorting coordinates (used in 'less' routine)
166  const scalar mergeDistance_;
167 
168  //- Overwrite the mesh?
169  const bool overwrite_;
170 
171  //- Instance of mesh upon construction. Used when in overwrite_ mode.
172  const word oldInstance_;
173 
174  //- All surface-intersection interaction
175  const refinementSurfaces& surfaces_;
176 
177  //- All feature-edge interaction
178  const refinementFeatures& features_;
179 
180  //- All shell-refinement interaction
181  const shellSurfaces& shells_;
182 
183  //- Refinement engine
184  hexRef8 meshCutter_;
185 
186  //- Per cc-cc vector the index of the surface hit
187  labelIOList surfaceIndex_;
188 
189  //- User supplied face based data.
190  List<Tuple2<mapType, labelList>> userFaceData_;
191 
192  //- Meshed patches - are treated differently. Stored as wordList since
193  // order changes.
194  wordList meshedPatches_;
195 
196 
197  // Private Member Functions
198 
199  //- Add patchfield of given type to all fields on mesh
200  template<class GeoField>
201  static void addPatchFields(fvMesh&, const word& patchFieldType);
202 
203  //- Reorder patchfields of all fields on mesh
204  template<class GeoField>
205  static void reorderPatchFields(fvMesh&, const labelList& oldToNew);
206 
207  //- Find out which faces have changed given cells (old mesh labels)
208  // that were marked for refinement.
209  static labelList getChangedFaces
210  (
211  const mapPolyMesh&,
212  const labelList& oldCellsToRefine
213  );
214 
215  //- Calculate coupled boundary end vector and refinement level
216  void calcNeighbourData
217  (
218  labelList& neiLevel,
219  pointField& neiCc
220  ) const;
221 
222  //- Find any intersection of surface. Store in surfaceIndex_.
223  void updateIntersections(const labelList& changedFaces);
224 
225  //- Remove cells. Put exposedFaces into exposedPatchIDs.
226  autoPtr<mapPolyMesh> doRemoveCells
227  (
228  const labelList& cellsToRemove,
229  const labelList& exposedFaces,
230  const labelList& exposedPatchIDs,
231  removeCells& cellRemover
232  );
233 
234  // Get cells which are inside any closed surface. Note that
235  // all closed surfaces
236  // will have already been oriented to have keepPoint outside.
237  labelList getInsideCells(const word&) const;
238 
239  // Do all to remove inside cells
240  autoPtr<mapPolyMesh> removeInsideCells
241  (
242  const string& msg,
243  const label exposedPatchi
244  );
245 
246 
247  // Refinement candidate selection
248 
249  //- Mark cell for refinement (if not already marked). Return false
250  // if refinelimit hit. Keeps running count (in nRefine) of cells
251  // marked for refinement
252  static bool markForRefine
253  (
254  const label markValue,
255  const label nAllowRefine,
256  label& cellValue,
257  label& nRefine
258  );
259 
260  //- Mark every cell with level of feature passing through it
261  // (or -1 if not passed through). Uses tracking.
262  void markFeatureCellLevel
263  (
264  const pointField& keepPoints,
265  labelList& maxFeatureLevel
266  ) const;
267 
268  //- Calculate list of cells to refine based on intersection of
269  // features.
270  label markFeatureRefinement
271  (
272  const pointField& keepPoints,
273  const label nAllowRefine,
274 
276  label& nRefine
277  ) const;
278 
279  //- Mark cells for distance-to-feature based refinement.
280  label markInternalDistanceToFeatureRefinement
281  (
282  const label nAllowRefine,
284  label& nRefine
285  ) const;
286 
287  //- Mark cells for refinement-shells based refinement.
288  label markInternalRefinement
289  (
290  const label nAllowRefine,
292  label& nRefine
293  ) const;
294 
295  //- Collect faces that are intersected and whose neighbours aren't
296  // yet marked for refinement.
297  labelList getRefineCandidateFaces
298  (
299  const labelList& refineCell
300  ) const;
301 
302  //- Mark cells for surface intersection based refinement.
303  label markSurfaceRefinement
304  (
305  const label nAllowRefine,
306  const labelList& neiLevel,
307  const pointField& neiCc,
309  label& nRefine
310  ) const;
311 
312  //- Helper: count number of normals1 that are in normals2
313  label countMatches
314  (
315  const List<point>& normals1,
316  const List<point>& normals2,
317  const scalar tol = 1e-6
318  ) const;
319 
320  //- Mark cells for surface curvature based refinement. Marks if
321  // local curvature > curvature or if on different regions
322  // (markDifferingRegions)
323  label markSurfaceCurvatureRefinement
324  (
325  const scalar curvature,
326  const label nAllowRefine,
327  const labelList& neiLevel,
328  const pointField& neiCc,
330  label& nRefine
331  ) const;
332 
333  //- Mark cell if local a gap topology or
334  bool checkProximity
335  (
336  const scalar planarCos,
337  const label nAllowRefine,
338 
339  const label surfaceLevel,
340  const vector& surfaceLocation,
341  const vector& surfaceNormal,
342 
343  const label celli,
344 
345  label& cellMaxLevel,
346  vector& cellMaxLocation,
347  vector& cellMaxNormal,
348 
350  label& nRefine
351  ) const;
352 
353  //- Mark cells for surface proximity based refinement.
354  label markProximityRefinement
355  (
356  const scalar curvature,
357  const label nAllowRefine,
358  const labelList& neiLevel,
359  const pointField& neiCc,
360 
362  label& nRefine
363  ) const;
364 
365  // Baffle handling
366 
367  //- Get faces to repatch. Returns map from face to patch.
368  Map<labelPair> getZoneBafflePatches
369  (
370  const bool allowBoundary,
371  const labelList& globalToMasterPatch,
372  const labelList& globalToSlavePatch
373  ) const;
374 
375  //- Determine patches for baffles on all intersected unnamed faces
376  void getBafflePatches
377  (
378  const labelList& globalToMasterPatch,
379  const labelList& neiLevel,
380  const pointField& neiCc,
381  labelList& ownPatch,
382  labelList& neiPatch
383  ) const;
384 
385  //- Repatches external face or creates baffle for internal face
386  // with user specified patches (might be different for both sides).
387  // Returns label of added face.
388  label createBaffle
389  (
390  const label facei,
391  const label ownPatch,
392  const label neiPatch,
393  polyTopoChange& meshMod
394  ) const;
395 
396  // Problem cell handling
397 
398  //- Helper function to mark face as being on 'boundary'. Used by
399  // markFacesOnProblemCells
400  void markBoundaryFace
401  (
402  const label facei,
403  boolList& isBoundaryFace,
404  boolList& isBoundaryEdge,
405  boolList& isBoundaryPoint
406  ) const;
407 
408  void findNearest
409  (
410  const labelList& meshFaces,
411  List<pointIndexHit>& nearestInfo,
412  labelList& nearestSurface,
413  labelList& nearestRegion,
414  vectorField& nearestNormal
415  ) const;
416 
417  Map<label> findEdgeConnectedProblemCells
418  (
419  const scalarField& perpendicularAngle,
420  const labelList&
421  ) const;
422 
423  bool isCollapsedFace
424  (
425  const pointField&,
426  const pointField& neiCc,
427  const scalar minFaceArea,
428  const scalar maxNonOrtho,
429  const label facei
430  ) const;
431 
432  bool isCollapsedCell
433  (
434  const pointField&,
435  const scalar volFraction,
436  const label celli
437  ) const;
438 
439  //- Returns list with for every internal face -1 or the patch
440  // they should be baffled into. If removeEdgeConnectedCells is set
441  // removes cells based on perpendicularAngle.
442  labelList markFacesOnProblemCells
443  (
444  const dictionary& motionDict,
445  const bool removeEdgeConnectedCells,
446  const scalarField& perpendicularAngle,
447  const labelList& globalToMasterPatch
448  ) const;
449 
450  //- Returns list with for every face the label of the nearest
451  // patch. Any unreached face (disconnected mesh?) becomes
452  // adaptPatchIDs[0]
453  labelList nearestPatch(const labelList& adaptPatchIDs) const;
454 
455  //- Returns list with for every internal face -1 or the patch
456  // they should be baffled into.
457  labelList markFacesOnProblemCellsGeometric
458  (
459  const snapParameters& snapParams,
460  const dictionary& motionDict
461  ) const;
462 
463 
464  // Baffle merging
465 
466  //- Extract those baffles (duplicate) faces that are on the edge
467  // of a baffle region. These are candidates for merging.
468  List<labelPair> freeStandingBaffles
469  (
470  const List<labelPair>&,
471  const scalar freeStandingAngle
472  ) const;
473 
474 
475  // Zone handling
476 
477  //- Finds zone per cell for cells inside closed named surfaces.
478  // (uses geometric test for insideness)
479  // Adapts namedSurfaceIndex so all faces on boundary of cellZone
480  // have corresponding faceZone.
481  void findCellZoneGeometric
482  (
483  const pointField& neiCc,
484  const labelList& closedNamedSurfaces,
485  labelList& namedSurfaceIndex,
486  const labelList& surfaceToCellZone,
487  labelList& cellToZone
488  ) const;
489 
490  //- Finds zone per cell for cells inside named surfaces that have
491  // an inside point specified.
492  void findCellZoneInsideWalk
493  (
494  const labelList& locationSurfaces,
495  const labelList& namedSurfaceIndex,
496  const labelList& surfaceToCellZone,
497  labelList& cellToZone
498  ) const;
499 
500  //- Determines cell zone from cell region information.
501  bool calcRegionToZone
502  (
503  const label surfZoneI,
504  const label ownRegion,
505  const label neiRegion,
506 
507  labelList& regionToCellZone
508  ) const;
509 
510  //- Finds zone per cell. Uses topological walk with all faces
511  // marked in namedSurfaceIndex regarded as blocked.
512  void findCellZoneTopo
513  (
514  const point& keepPoint,
515  const labelList& namedSurfaceIndex,
516  const labelList& surfaceToCellZone,
517  labelList& cellToZone
518  ) const;
519 
520  //- Make namedSurfaceIndex consistent with cellToZone
521  // - clear out any blocked faces in between same cell zone.
522  void makeConsistentFaceIndex
523  (
524  const labelList& cellToZone,
525  labelList& namedSurfaceIndex
526  ) const;
527 
528  //- Remove any loose standing cells
529  void handleSnapProblems
530  (
531  const snapParameters& snapParams,
532  const bool useTopologicalSnapDetection,
533  const bool removeEdgeConnectedCells,
534  const scalarField& perpendicularAngle,
535  const dictionary& motionDict,
536  Time& runTime,
537  const labelList& globalToMasterPatch,
538  const labelList& globalToSlavePatch
539  );
540 
541 
542  // Some patch utilities
543 
544  //- Get all faces in faceToZone that have no cellZone on
545  // either side.
546  labelList freeStandingBaffleFaces
547  (
548  const labelList& faceToZone,
549  const labelList& cellToZone,
550  const labelList& neiCellZone
551  ) const;
552 
553  //- Determine per patch edge the number of master faces. Used
554  // to detect non-manifold situations.
555  void calcPatchNumMasterFaces
556  (
557  const PackedBoolList& isMasterFace,
558  const indirectPrimitivePatch& patch,
559  labelList& nMasterFaces
560  ) const;
561 
562  //- Determine per patch face the (singly-) connected zone it
563  // is in. Return overall number of zones.
564  label markPatchZones
565  (
566  const indirectPrimitivePatch& patch,
567  const labelList& nMasterFaces,
568  labelList& faceToZone
569  ) const;
570 
571  //- Make faces consistent.
572  void consistentOrientation
573  (
574  const PackedBoolList& isMasterFace,
575  const indirectPrimitivePatch& patch,
576  const labelList& nMasterFaces,
577  const labelList& faceToZone,
578  const Map<label>& zoneToOrientation,
579  boolList& meshFlipMap
580  ) const;
581 
582 
583  //- Disallow default bitwise copy construct
585 
586  //- Disallow default bitwise assignment
587  void operator=(const meshRefinement&);
588 
589 public:
590 
591  //- Runtime type information
592  ClassName("meshRefinement");
593 
594 
595  // Constructors
596 
597  //- Construct from components
599  (
600  fvMesh& mesh,
601  const scalar mergeDistance,
602  const bool overwrite,
603  const refinementSurfaces&,
604  const refinementFeatures&,
605  const shellSurfaces&
606  );
607 
608 
609  // Member Functions
610 
611  // Access
612 
613  //- Reference to mesh
614  const fvMesh& mesh() const
615  {
616  return mesh_;
617  }
618  fvMesh& mesh()
619  {
620  return mesh_;
621  }
623  scalar mergeDistance() const
624  {
625  return mergeDistance_;
626  }
627 
628  //- Overwrite the mesh?
629  bool overwrite() const
630  {
631  return overwrite_;
632  }
633 
634  //- (points)instance of mesh upon construction
635  const word& oldInstance() const
636  {
637  return oldInstance_;
638  }
639 
640  //- Reference to surface search engines
641  const refinementSurfaces& surfaces() const
642  {
643  return surfaces_;
644  }
645 
646  //- Reference to feature edge mesh
647  const refinementFeatures& features() const
648  {
649  return features_;
650  }
651 
652  //- Reference to refinement shells (regions)
653  const shellSurfaces& shells() const
654  {
655  return shells_;
656  }
657 
658  //- Reference to meshcutting engine
659  const hexRef8& meshCutter() const
660  {
661  return meshCutter_;
662  }
663 
664  //- Per start-end edge the index of the surface hit
665  const labelList& surfaceIndex() const
666  {
667  return surfaceIndex_;
668  }
671  {
672  return surfaceIndex_;
673  }
674 
675  //- Additional face data that is maintained across
676  // topo changes. Every entry is a list over all faces.
677  // Bit of a hack. Additional flag to say whether to maintain master
678  // only (false) or increase set to account for face-from-face.
680  {
681  return userFaceData_;
682  }
685  {
686  return userFaceData_;
687  }
688 
689 
690  // Other
691 
692  //- Count number of intersections (local)
693  label countHits() const;
694 
695  //- Redecompose according to cell count
696  // keepZoneFaces : find all faceZones from zoned surfaces and keep
697  // owner and neighbour together
698  // keepBaffles : find all baffles and keep them together
700  (
701  const bool keepZoneFaces,
702  const bool keepBaffles,
703  const scalarField& cellWeights,
704  decompositionMethod& decomposer,
705  fvMeshDistribute& distributor
706  );
707 
708  //- Get faces with intersection.
709  labelList intersectedFaces() const;
710 
711  //- Get points on surfaces with intersection and boundary faces.
713 
714  //- Create patch from set of patches
716  (
717  const polyMesh&,
718  const labelList&
719  );
720 
721  //- Helper function to make a pointVectorField with correct
722  // bcs for mesh movement:
723  // - adaptPatchIDs : fixedValue
724  // - processor : calculated (so free to move)
725  // - cyclic/wedge/symmetry : slip
726  // - other : slip
728  (
729  const pointMesh& pMesh,
730  const labelList& adaptPatchIDs
731  );
732 
733  //- Helper function: check that face zones are synced
734  static void checkCoupledFaceZones(const polyMesh&);
735 
736  //- Helper: calculate edge weights (1/length)
737  static void calculateEdgeWeights
738  (
739  const polyMesh& mesh,
740  const PackedBoolList& isMasterEdge,
741  const labelList& meshPoints,
742  const edgeList& edges,
743  scalarField& edgeWeights,
744  scalarField& invSumWeight
745  );
746 
747  //- Helper: weighted sum (over all subset of mesh points) by
748  // summing contribution from (master) edges
749  template<class Type>
750  static void weightedSum
751  (
752  const polyMesh& mesh,
753  const PackedBoolList& isMasterEdge,
754  const labelList& meshPoints,
755  const edgeList& edges,
756  const scalarField& edgeWeights,
757  const Field<Type>& data,
759  );
760 
761 
762  // Refinement
763 
764  //- Is local topology a small gap?
765  bool isGap
766  (
767  const scalar,
768  const vector&,
769  const vector&,
770  const vector&,
771  const vector&
772  ) const;
773 
774  //- Is local topology a small gap normal to the test vector
775  bool isNormalGap
776  (
777  const scalar,
778  const vector&,
779  const vector&,
780  const vector&,
781  const vector&
782  ) const;
783 
784  //- Calculate list of cells to refine.
786  (
787  const pointField& keepPoints,
788  const scalar curvature,
789  const scalar planarAngle,
790 
791  const bool featureRefinement,
792  const bool featureDistanceRefinement,
793  const bool internalRefinement,
794  const bool surfaceRefinement,
795  const bool curvatureRefinement,
796  const bool gapRefinement,
797  const label maxGlobalCells,
798  const label maxLocalCells
799  ) const;
800 
801  //- Refine some cells
802  autoPtr<mapPolyMesh> refine(const labelList& cellsToRefine);
803 
804  //- Refine some cells and rebalance
806  (
807  const string& msg,
808  decompositionMethod& decomposer,
809  fvMeshDistribute& distributor,
810  const labelList& cellsToRefine,
811  const scalar maxLoadUnbalance
812  );
813 
814  //- Balance before refining some cells
816  (
817  const string& msg,
818  decompositionMethod& decomposer,
819  fvMeshDistribute& distributor,
820  const labelList& cellsToRefine,
821  const scalar maxLoadUnbalance
822  );
823 
824 
825  // Baffle handling
826 
827  //- Split off unreachable areas of mesh.
828  void baffleAndSplitMesh
829  (
830  const bool handleSnapProblems,
831 
832  // How to remove problem snaps
833  const snapParameters& snapParams,
834  const bool useTopologicalSnapDetection,
835  const bool removeEdgeConnectedCells,
836  const scalarField& perpendicularAngle,
837 
838  // How to handle free-standing baffles
839  const bool mergeFreeStanding,
840  const scalar freeStandingAngle,
841 
842  const dictionary& motionDict,
843  Time& runTime,
844  const labelList& globalToMasterPatch,
845  const labelList& globalToSlavePatch,
846  const point& keepPoint
847  );
848 
849  //- Split off (with optional buffer layers) unreachable areas
850  // of mesh. Does not introduce baffles.
852  (
853  const label nBufferLayers,
854  const labelList& globalToMasterPatch,
855  const labelList& globalToSlavePatch,
856  const point& keepPoint
857  );
858 
859  //- Find boundary points that connect to more than one cell
860  // region and split them.
862 
863  //- Find boundary points that connect to more than one cell
864  // region and split them.
866 
867  //- Create baffle for every internal face where ownPatch != -1.
868  // External faces get repatched according to ownPatch (neiPatch
869  // should be -1 for these)
871  (
872  const labelList& ownPatch,
873  const labelList& neiPatch
874  );
875 
876  //- Debug helper: check faceZones are not on processor patches
877  void checkZoneFaces() const;
878 
879  //- Create baffles for faces straddling zoned surfaces. Return
880  // baffles.
882  (
883  const labelList& globalToMasterPatch,
884  const labelList& globalToSlavePatch,
886  );
887 
888  //- Merge baffles. Gets pairs of faces.
890 
891  //- Put faces/cells into zones according to surface specification.
892  // Returns null if no zone surfaces present. Region containing
893  // the keepPoint will not be put into a cellZone.
895  (
896  const point& keepPoint,
897  const bool allowFreeStandingZoneFaces
898  );
899 
900 
901  // Other topo changes
902 
903  //- Helper:append patch to end of mesh.
904  static label appendPatch
905  (
906  fvMesh&,
907  const label insertPatchi,
908  const word&,
909  const dictionary&
910  );
911 
912  //- Helper:add patch to mesh. Update all registered fields.
913  // Used by addMeshedPatch to add patches originating from surfaces.
914  static label addPatch(fvMesh&, const word& name, const dictionary&);
915 
916  //- Add patch originating from meshing. Update meshedPatches_.
917  label addMeshedPatch(const word& name, const dictionary&);
918 
919  //- Get patchIDs for patches added in addMeshedPatch.
920  labelList meshedPatches() const;
921 
922  //- Select coupled faces that are not collocated
924 
925  //- Find region point is in. Uses optional perturbation to re-test.
926  static label findRegion
927  (
928  const polyMesh&,
929  const labelList& cellRegion,
930  const vector& perturbVec,
931  const point& p
932  );
933 
934  //- Split mesh. Keep part containing point.
936  (
937  const labelList& globalToMasterPatch,
938  const labelList& globalToSlavePatch,
939  const point& keepPoint
940  );
941 
942  //- Split faces into two
944  (
945  const labelList& splitFaces,
946  const labelPairList& splits
947  );
948 
949  //- Update local numbering for mesh redistribution
950  void distribute(const mapDistributePolyMesh&);
951 
952  //- Update for external change to mesh. changedFaces are in new mesh
953  // face labels.
954  void updateMesh
955  (
956  const mapPolyMesh&,
957  const labelList& changedFaces
958  );
959 
960  //- Helper: reorder list according to map.
961  template<class T>
962  static void updateList
963  (
964  const labelList& newToOld,
965  const T& nullValue,
966  List<T>& elems
967  );
968 
969 
970  // Restoring : is where other processes delete and reinsert data.
971 
972  //- Signal points/face/cells for which to store data
973  void storeData
974  (
975  const labelList& pointsToStore,
976  const labelList& facesToStore,
977  const labelList& cellsToStore
978  );
979 
980  //- Update local numbering + undo
981  // Data to restore given as new pointlabel + stored pointlabel
982  // (i.e. what was in pointsToStore)
983  void updateMesh
984  (
985  const mapPolyMesh&,
986  const labelList& changedFaces,
987  const Map<label>& pointsToRestore,
988  const Map<label>& facesToRestore,
989  const Map<label>& cellsToRestore
990  );
991 
992  // Merging coplanar faces and edges
993 
994  //- Merge coplanar faces. preserveFaces is != -1 for faces
995  // to be preserved
997  (
998  const scalar minCos,
999  const scalar concaveCos,
1000  const labelList& patchIDs,
1001  const dictionary& motionDict,
1002  const labelList& preserveFaces
1003  );
1004 
1006  (
1007  removePoints& pointRemover,
1008  const boolList& pointCanBeDeleted
1009  );
1010 
1012  (
1013  removePoints& pointRemover,
1014  const labelList& facesToRestore
1015  );
1016 
1018  (
1019  const labelList& candidateFaces,
1020  const labelHashSet& set
1021  ) const;
1022 
1023  // Pick up faces of cells of faces in set.
1025  (
1026  const labelHashSet& set
1027  ) const;
1028 
1029  //- Merge edges, maintain mesh quality. Return global number
1030  // of edges merged
1032  (
1033  const scalar minCos,
1034  const dictionary& motionDict
1035  );
1036 
1037 
1038  // Debug/IO
1039 
1040  //- Debugging: check that all faces still obey start()>end()
1041  void checkData();
1042 
1043  static void testSyncPointList
1044  (
1045  const string& msg,
1046  const polyMesh& mesh,
1047  const List<scalar>& fld
1048  );
1049 
1050  static void testSyncPointList
1051  (
1052  const string& msg,
1053  const polyMesh& mesh,
1054  const List<point>& fld
1055  );
1056 
1057  //- Compare two lists over all boundary faces
1058  template<class T>
1060  (
1061  const scalar mergeDistance,
1062  const string&,
1063  const UList<T>&,
1064  const UList<T>&
1065  ) const;
1066 
1067  //- Print list according to (collected and) sorted coordinate
1068  template<class T>
1069  static void collectAndPrint
1070  (
1071  const UList<point>& points,
1072  const UList<T>& data
1073  );
1074 
1075  //- Determine master point for subset of points. If coupled
1076  // chooses only one
1078  (
1079  const polyMesh& mesh,
1080  const labelList& meshPoints
1081  );
1082 
1083  //- Determine master edge for subset of edges. If coupled
1084  // chooses only one
1086  (
1087  const polyMesh& mesh,
1088  const labelList& meshEdges
1089  );
1090 
1091  //- Print some mesh stats.
1092  void printMeshInfo(const bool, const string&) const;
1093 
1094  //- Replacement for Time::timeName() : return oldInstance (if
1095  // overwrite_)
1096  word timeName() const;
1097 
1098  //- Set instance of all local IOobjects
1099  void setInstance(const fileName&);
1100 
1101  //- Write mesh and all data
1102  bool write() const;
1103 
1104  //- Write refinement level as volScalarFields for postprocessing
1105  void dumpRefinementLevel() const;
1106 
1107  //- Debug: Write intersection information to OBJ format
1108  void dumpIntersections(const fileName& prefix) const;
1109 
1110  //- Do any one of above IO functions
1111  void write
1112  (
1113  const debugType debugFlags,
1114  const writeType writeFlags,
1115  const fileName&
1116  ) const;
1117 
1118  //- Helper: calculate average
1119  template<class T>
1120  static T gAverage
1121  (
1122  const PackedBoolList& isMasterElem,
1123  const UList<T>& values
1124  );
1125 
1126  //- Get/set write level
1127  static writeType writeLevel();
1128  static void writeLevel(const writeType);
1129 
1130  //- Get/set output level
1131  static outputType outputLevel();
1132  static void outputLevel(const outputType);
1133 
1134 
1135  //- Helper: convert wordList into bit pattern using provided
1136  // NamedEnum
1137  template<class Enum>
1138  static int readFlags(const Enum& namedEnum, const wordList&);
1139 
1140 };
1141 
1142 
1143 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1144 
1145 } // End namespace Foam
1146 
1147 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1148 
1149 #ifdef NoRepository
1150  #include "meshRefinementTemplates.C"
1151 #endif
1152 
1153 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1154 
1155 #endif
1156 
1157 // ************************************************************************* //
const labelList & surfaceIndex() const
Per start-end edge the index of the surface hit.
void printMeshInfo(const bool, const string &) const
Print some mesh stats.
const refinementSurfaces & surfaces() const
Reference to surface search engines.
void checkZoneFaces() const
Debug helper: check faceZones are not on processor patches.
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
IOoutputType
Enumeration for what to output.
A class for handling file names.
Definition: fileName.H:69
const List< Tuple2< mapType, labelList > > & userFaceData() const
Additional face data that is maintained across.
labelList refineCandidates(const pointField &keepPoints, const scalar curvature, const scalar planarAngle, const bool featureRefinement, const bool featureDistanceRefinement, const bool internalRefinement, const bool surfaceRefinement, const bool curvatureRefinement, const bool gapRefinement, const label maxGlobalCells, const label maxLocalCells) const
Calculate list of cells to refine.
autoPtr< mapDistributePolyMesh > balanceAndRefine(const string &msg, decompositionMethod &decomposer, fvMeshDistribute &distributor, const labelList &cellsToRefine, const scalar maxLoadUnbalance)
Balance before refining some cells.
const word & oldInstance() const
(points)instance of mesh upon construction
label countHits() const
Count number of intersections (local)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
IOdebugType
Enumeration for what to debug.
bool overwrite() const
Overwrite the mesh?
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:60
static label findRegion(const polyMesh &, const labelList &cellRegion, const vector &perturbVec, const point &p)
Find region point is in. Uses optional perturbation to re-test.
static const NamedEnum< IOoutputType, 1 > IOoutputTypeNames
const shellSurfaces & shells() const
Reference to refinement shells (regions)
Given list of cells to remove insert all the topology changes.
Definition: removeCells.H:59
engineTime & runTime
autoPtr< mapPolyMesh > createZoneBaffles(const labelList &globalToMasterPatch, const labelList &globalToSlavePatch, List< labelPair > &)
Create baffles for faces straddling zoned surfaces. Return.
autoPtr< mapPolyMesh > createBaffles(const labelList &ownPatch, const labelList &neiPatch)
Create baffle for every internal face where ownPatch != -1.
labelList intersectedFaces() const
Get faces with intersection.
static writeType writeLevel()
Get/set write level.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
void baffleAndSplitMesh(const bool handleSnapProblems, const snapParameters &snapParams, const bool useTopologicalSnapDetection, const bool removeEdgeConnectedCells, const scalarField &perpendicularAngle, const bool mergeFreeStanding, const scalar freeStandingAngle, const dictionary &motionDict, Time &runTime, const labelList &globalToMasterPatch, const labelList &globalToSlavePatch, const point &keepPoint)
Split off unreachable areas of mesh.
autoPtr< mapPolyMesh > splitMesh(const label nBufferLayers, const labelList &globalToMasterPatch, const labelList &globalToSlavePatch, const point &keepPoint)
Split off (with optional buffer layers) unreachable areas.
set value to -1 any face that was refined
labelList intersectedPoints() const
Get points on surfaces with intersection and boundary faces.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
labelList collectFaces(const labelList &candidateFaces, const labelHashSet &set) const
autoPtr< mapPolyMesh > splitMeshRegions(const labelList &globalToMasterPatch, const labelList &globalToSlavePatch, const point &keepPoint)
Split mesh. Keep part containing point.
Takes mesh with &#39;baffles&#39; (= boundary faces sharing points). Determines for selected points on bounda...
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
const fvMesh & mesh() const
Reference to mesh.
IOwriteType
Enumeration for what to write.
Container for data on surfaces used for surface-driven refinement. Contains all the data about the le...
static label addPatch(fvMesh &, const word &name, const dictionary &)
Helper:add patch to mesh. Update all registered fields.
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
label mergeEdgesUndo(const scalar minCos, const dictionary &motionDict)
Merge edges, maintain mesh quality. Return global number.
autoPtr< mapPolyMesh > mergeBaffles(const List< labelPair > &)
Merge baffles. Gets pairs of faces.
Encapsulates queries for features.
static label appendPatch(fvMesh &, const label insertPatchi, const word &, const dictionary &)
Helper:append patch to end of mesh.
void distribute(const mapDistributePolyMesh &)
Update local numbering for mesh redistribution.
static void updateList(const labelList &newToOld, const T &nullValue, List< T > &elems)
Helper: reorder list according to map.
A list of faces which address into the list of points.
Sends/receives parts of mesh+fvfields to neighbouring processors. Used in load balancing.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
ClassName("meshRefinement")
Runtime type information.
const pointField & points
static void weightedSum(const polyMesh &mesh, const PackedBoolList &isMasterEdge, const labelList &meshPoints, const edgeList &edges, const scalarField &edgeWeights, const Field< Type > &data, Field< Type > &sum)
Helper: weighted sum (over all subset of mesh points) by.
static void calculateEdgeWeights(const polyMesh &mesh, const PackedBoolList &isMasterEdge, const labelList &meshPoints, const edgeList &edges, scalarField &edgeWeights, scalarField &invSumWeight)
Helper: calculate edge weights (1/length)
autoPtr< mapDistributePolyMesh > refineAndBalance(const string &msg, decompositionMethod &decomposer, fvMeshDistribute &distributor, const labelList &cellsToRefine, const scalar maxLoadUnbalance)
Refine some cells and rebalance.
autoPtr< mapPolyMesh > splitFaces(const labelList &splitFaces, const labelPairList &splits)
Split faces into two.
Refinement of (split) hexes using polyTopoChange.
Definition: hexRef8.H:64
A class for handling words, derived from string.
Definition: word.H:59
static PackedBoolList getMasterEdges(const polyMesh &mesh, const labelList &meshEdges)
Determine master edge for subset of edges. If coupled.
Simple container to keep together snap specific information.
static const NamedEnum< IOdebugType, 5 > IOdebugTypeNames
static void checkCoupledFaceZones(const polyMesh &)
Helper function: check that face zones are synced.
void testSyncBoundaryFaceList(const scalar mergeDistance, const string &, const UList< T > &, const UList< T > &) const
Compare two lists over all boundary faces.
Abstract base class for decomposition.
Encapsulates queries for volume refinement (&#39;refine all cells within shell&#39;).
Definition: shellSurfaces.H:52
autoPtr< mapDistributePolyMesh > balance(const bool keepZoneFaces, const bool keepBaffles, const scalarField &cellWeights, decompositionMethod &decomposer, fvMeshDistribute &distributor)
Redecompose according to cell count.
Container with cells to refine. Refinement given as single direction.
Definition: refineCell.H:56
autoPtr< mapPolyMesh > doRemovePoints(removePoints &pointRemover, const boolList &pointCanBeDeleted)
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:61
const refinementFeatures & features() const
Reference to feature edge mesh.
word timeName() const
Replacement for Time::timeName() : return oldInstance (if.
mapType
Enumeration for how the userdata is to be mapped upon refinement.
const hexRef8 & meshCutter() const
Reference to meshcutting engine.
void checkData()
Debugging: check that all faces still obey start()>end()
static int readFlags(const Enum &namedEnum, const wordList &)
Helper: convert wordList into bit pattern using provided.
Database for solution data, solver performance and other reduced data.
Definition: data.H:52
label addMeshedPatch(const word &name, const dictionary &)
Add patch originating from meshing. Update meshedPatches_.
labelList meshedPatches() const
Get patchIDs for patches added in addMeshedPatch.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
autoPtr< mapPolyMesh > doRestorePoints(removePoints &pointRemover, const labelList &facesToRestore)
static autoPtr< indirectPrimitivePatch > makePatch(const polyMesh &, const labelList &)
Create patch from set of patches.
static outputType outputLevel()
Get/set output level.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void updateMesh(const mapPolyMesh &, const labelList &changedFaces)
Update for external change to mesh. changedFaces are in new mesh.
label mergePatchFacesUndo(const scalar minCos, const scalar concaveCos, const labelList &patchIDs, const dictionary &motionDict, const labelList &preserveFaces)
Merge coplanar faces. preserveFaces is != -1 for faces.
scalar mergeDistance() const
Helper class which maintains intersections of (changing) mesh with (static) surfaces.
autoPtr< mapPolyMesh > refine(const labelList &cellsToRefine)
Refine some cells.
void setInstance(const fileName &)
Set instance of all local IOobjects.
static void collectAndPrint(const UList< point > &points, const UList< T > &data)
Print list according to (collected and) sorted coordinate.
autoPtr< mapPolyMesh > dupNonManifoldPoints()
Find boundary points that connect to more than one cell.
A bit-packed bool list.
static T gAverage(const PackedBoolList &isMasterElem, const UList< T > &values)
Helper: calculate average.
void storeData(const labelList &pointsToStore, const labelList &facesToStore, const labelList &cellsToStore)
Signal points/face/cells for which to store data.
autoPtr< mapPolyMesh > zonify(const point &keepPoint, const bool allowFreeStandingZoneFaces)
Put faces/cells into zones according to surface specification.
static void testSyncPointList(const string &msg, const polyMesh &mesh, const List< scalar > &fld)
have slaves (upon refinement) from master
labelList growFaceCellFace(const labelHashSet &set) const
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Direct mesh changes based on v1.3 polyTopoChange syntax.
void selectSeparatedCoupledFaces(boolList &) const
Select coupled faces that are not collocated.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:98
Removes selected points from mesh and updates faces using these points.
Definition: removePoints.H:58
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
bool isGap(const scalar, const vector &, const vector &, const vector &, const vector &) const
Is local topology a small gap?
Contains information about location on a triSurface.
bool write() const
Write mesh and all data.
static tmp< pointVectorField > makeDisplacementField(const pointMesh &pMesh, const labelList &adaptPatchIDs)
Helper function to make a pointVectorField with correct.
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
static PackedBoolList getMasterPoints(const polyMesh &mesh, const labelList &meshPoints)
Determine master point for subset of points. If coupled.
void dumpIntersections(const fileName &prefix) const
Debug: Write intersection information to OBJ format.
void dumpRefinementLevel() const
Write refinement level as volScalarFields for postprocessing.
bool isNormalGap(const scalar, const vector &, const vector &, const vector &, const vector &) const
Is local topology a small gap normal to the test vector.
Namespace for OpenFOAM.
static const NamedEnum< IOwriteType, 5 > IOwriteTypeNames
A HashTable to objects of type <T> with a label key.
Definition: Map.H:49