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