1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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.
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
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/>.
28 Helper class which maintains intersections of (changing) mesh with
32 - per face any intersections of the cc-cc segment with any of the surfaces
36 meshRefinementBaffles.C
38 meshRefinementProblemCells.C
39 meshRefinementRefine.C
41 \*---------------------------------------------------------------------------*/
43 #ifndef meshRefinement_H
44 #define meshRefinement_H
47 #include "mapPolyMesh.H"
49 #include "labelPair.H"
50 #include "indirectPrimitivePatch.H"
51 #include "pointFieldsFwd.H"
53 #include "pointIndexHit.H"
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60 // Class forward declarations
62 class mapDistributePolyMesh;
63 class decompositionMethod;
64 class refinementSurfaces;
65 class refinementFeatures;
68 class fvMeshDistribute;
69 class searchableSurface;
74 /*---------------------------------------------------------------------------*\
75 Class meshRefinement Declaration
76 \*---------------------------------------------------------------------------*/
85 //- Enumeration for debug dumping
94 //- Enumeration for how the userdata is to be mapped upon refinement.
97 MASTERONLY = 1, /*!< maintain master only */
98 KEEPALL = 2, /*!< have slaves (upon refinement) from master */
99 REMOVE = 4 /*!< set value to -1 any face that was refined */
107 //- Reference to mesh
110 //- tolerance used for sorting coordinates (used in 'less' routine)
111 const scalar mergeDistance_;
113 //- overwrite the mesh?
114 const bool overwrite_;
116 //- Instance of mesh upon construction. Used when in overwrite_ mode.
117 const word oldInstance_;
119 //- All surface-intersection interaction
120 const refinementSurfaces& surfaces_;
122 //- All feature-edge interaction
123 const refinementFeatures& features_;
125 //- All shell-refinement interaction
126 const shellSurfaces& shells_;
128 //- refinement engine
131 //- per cc-cc vector the index of the surface hit
132 labelIOList surfaceIndex_;
134 //- user supplied face based data.
135 List<Tuple2<mapType, labelList> > userFaceData_;
137 //- Meshed patches - are treated differently. Stored as wordList since
139 wordList meshedPatches_;
142 // Private Member Functions
144 //- Reorder list according to map.
146 static void updateList
148 const labelList& newToOld,
153 //- Add patchfield of given type to all fields on mesh
154 template<class GeoField>
155 static void addPatchFields(fvMesh&, const word& patchFieldType);
157 //- Reorder patchfields of all fields on mesh
158 template<class GeoField>
159 static void reorderPatchFields(fvMesh&, const labelList& oldToNew);
161 //- Find out which faces have changed given cells (old mesh labels)
162 // that were marked for refinement.
163 static labelList getChangedFaces
166 const labelList& oldCellsToRefine
169 //- Calculate coupled boundary end vector and refinement level
170 void calcNeighbourData
176 //- Find any intersection of surface. Store in surfaceIndex_.
177 void updateIntersections(const labelList& changedFaces);
179 //- Remove cells. Put exposedFaces into exposedPatchIDs.
180 autoPtr<mapPolyMesh> doRemoveCells
182 const labelList& cellsToRemove,
183 const labelList& exposedFaces,
184 const labelList& exposedPatchIDs,
185 removeCells& cellRemover
188 // Get cells which are inside any closed surface. Note that
189 // all closed surfaces
190 // will have already been oriented to have keepPoint outside.
191 labelList getInsideCells(const word&) const;
193 // Do all to remove inside cells
194 autoPtr<mapPolyMesh> removeInsideCells
197 const label exposedPatchI
200 // For decomposeCombineRegions
202 //- Used in decomposeCombineRegions. Given global region per cell
203 // determines master processor/cell for regions straddling
205 void getCoupledRegionMaster
207 const globalIndex& globalCells,
208 const boolList& blockedFace,
209 const regionSplit& globalRegion,
210 Map<label>& regionToMaster
213 //- Determine regions that are local to me or coupled ones that
214 // are owned by me. Determine representative location.
215 void calcLocalRegions
217 const globalIndex& globalCells,
218 const labelList& globalRegion,
219 const Map<label>& coupledRegionToMaster,
220 const scalarField& cellWeights,
222 Map<label>& globalToLocalRegion,
223 pointField& localPoints,
224 scalarField& localWeights
227 //- Convert region into global index.
228 static label getShiftedRegion
230 const globalIndex& indexer,
231 const Map<label>& globalToLocalRegion,
232 const Map<label>& coupledRegionToShifted,
233 const label globalRegion
236 //- helper: add element if not in list. Linear search.
237 static void addUnique(const label, labelList&);
239 //- Calculate region connectivity. Major communication.
240 void calcRegionRegions
242 const labelList& globalRegion,
243 const Map<label>& globalToLocalRegion,
244 const Map<label>& coupledRegionToMaster,
245 labelListList& regionRegions
249 // Refinement candidate selection
251 //- Mark cell for refinement (if not already marked). Return false
252 // if refinelimit hit. Keeps running count (in nRefine) of cells
253 // marked for refinement
254 static bool markForRefine
256 const label markValue,
257 const label nAllowRefine,
262 //- Calculate list of cells to refine based on intersection of
264 label markFeatureRefinement
266 const point& keepPoint,
267 const label nAllowRefine,
269 labelList& refineCell,
273 //- Mark cells for refinement-shells based refinement.
274 label markInternalRefinement
276 const label nAllowRefine,
277 labelList& refineCell,
281 //- Collect faces that are intersected and whose neighbours aren't
282 // yet marked for refinement.
283 labelList getRefineCandidateFaces
285 const labelList& refineCell
288 //- Mark cells for surface intersection based refinement.
289 label markSurfaceRefinement
291 const label nAllowRefine,
292 const labelList& neiLevel,
293 const pointField& neiCc,
294 labelList& refineCell,
298 //- Mark cell if local curvature > curvature or
299 // markDifferingRegions = true and intersections with different
303 const scalar curvature,
304 const label nAllowRefine,
305 const label surfaceLevel,
306 const vector& surfaceNormal,
309 vector& cellMaxNormal,
310 labelList& refineCell,
315 //- Mark cells for surface curvature based refinement. Marks if
316 // local curvature > curvature or if on different regions
317 // (markDifferingRegions)
318 label markSurfaceCurvatureRefinement
320 const scalar curvature,
321 const label nAllowRefine,
322 const labelList& neiLevel,
323 const pointField& neiCc,
324 labelList& refineCell,
330 //- Get faces to repatch. Returns map from face to patch.
331 Map<label> getZoneBafflePatches
333 const bool allowBoundary,
334 const labelList& globalToPatch
337 //- Determine patches for baffles
338 void getBafflePatches
340 const labelList& globalToPatch,
341 const labelList& neiLevel,
342 const pointField& neiCc,
347 //- Determine patch for baffle using some heuristic (and not
351 const labelList& facePatch,
355 //- Repatches external face or creates baffle for internal face
356 // with user specified patches (might be different for both sides).
357 // Returns label of added face.
361 const label ownPatch,
362 const label neiPatch,
363 polyTopoChange& meshMod
366 // Problem cell handling
368 //- Helper function to mark face as being on 'boundary'. Used by
369 // markFacesOnProblemCells
370 void markBoundaryFace
373 boolList& isBoundaryFace,
374 boolList& isBoundaryEdge,
375 boolList& isBoundaryPoint
380 const labelList& meshFaces,
381 List<pointIndexHit>& nearestInfo,
382 labelList& nearestSurface,
383 labelList& nearestRegion,
384 vectorField& nearestNormal
387 Map<label> findEdgeConnectedProblemCells
389 const scalarField& perpendicularAngle,
396 const pointField& neiCc,
397 const scalar minFaceArea,
398 const scalar maxNonOrtho,
405 const scalar volFraction,
409 //- Returns list with for every internal face -1 or the patch
410 // they should be baffled into. If removeEdgeConnectedCells is set
411 // removes cells based on perpendicularAngle.
412 labelList markFacesOnProblemCells
414 const dictionary& motionDict,
415 const bool removeEdgeConnectedCells,
416 const scalarField& perpendicularAngle,
417 const labelList& globalToPatch
423 //- Extract those baffles (duplicate) faces that are on the edge
424 // of a baffle region. These are candidates for merging.
425 List<labelPair> filterDuplicateFaces(const List<labelPair>&) const;
430 //- Finds zone per cell for cells inside closed named surfaces.
431 // (uses geometric test for insideness)
432 // Adapts namedSurfaceIndex so all faces on boundary of cellZone
433 // have corresponding faceZone.
434 void findCellZoneGeometric
436 const labelList& closedNamedSurfaces,
437 labelList& namedSurfaceIndex,
438 const labelList& surfaceToCellZone,
439 labelList& cellToZone
442 //- Finds zone per cell for cells inside named surfaces that have
443 // an inside point specified.
444 void findCellZoneInsideWalk
446 const labelList& locationSurfaces,
447 const labelList& namedSurfaceIndex,
448 const labelList& surfaceToCellZone,
449 labelList& cellToZone
452 //- Determines cell zone from cell region information.
453 bool calcRegionToZone
455 const label surfZoneI,
456 const label ownRegion,
457 const label neiRegion,
459 labelList& regionToCellZone
462 //- Finds zone per cell. Uses topological walk with all faces
463 // marked in namedSurfaceIndex regarded as blocked.
464 void findCellZoneTopo
466 const point& keepPoint,
467 const labelList& namedSurfaceIndex,
468 const labelList& surfaceToCellZone,
469 labelList& cellToZone
472 void makeConsistentFaceIndex
474 const labelList& cellToZone,
475 labelList& namedSurfaceIndex
479 //- Disallow default bitwise copy construct
480 meshRefinement(const meshRefinement&);
482 //- Disallow default bitwise assignment
483 void operator=(const meshRefinement&);
487 //- Runtime type information
488 ClassName("meshRefinement");
493 //- Construct from components
497 const scalar mergeDistance,
498 const bool overwrite,
499 const refinementSurfaces&,
500 const refinementFeatures&,
509 //- reference to mesh
510 const fvMesh& mesh() const
519 scalar mergeDistance() const
521 return mergeDistance_;
524 //- Overwrite the mesh?
525 bool overwrite() const
530 //- (points)instance of mesh upon construction
531 const word& oldInstance() const
536 //- reference to surface search engines
537 const refinementSurfaces& surfaces() const
542 //- reference to feature edge mesh
543 const refinementFeatures& features() const
548 //- reference to refinement shells (regions)
549 const shellSurfaces& shells() const
554 //- reference to meshcutting engine
555 const hexRef8& meshCutter() const
560 //- per start-end edge the index of the surface hit
561 const labelList& surfaceIndex() const
563 return surfaceIndex_;
566 labelList& surfaceIndex()
568 return surfaceIndex_;
571 //- Additional face data that is maintained across
572 // topo changes. Every entry is a list over all faces.
573 // Bit of a hack. Additional flag to say whether to maintain master
574 // only (false) or increase set to account for face-from-face.
575 const List<Tuple2<mapType, labelList> >& userFaceData() const
577 return userFaceData_;
580 List<Tuple2<mapType, labelList> >& userFaceData()
582 return userFaceData_;
588 //- Count number of intersections (local)
589 label countHits() const;
591 //- Helper function to get decomposition such that all connected
592 // regions get moved onto one processor. Used to prevent baffles
593 // straddling processor boundaries. explicitConnections is to
594 // keep pairs of non-coupled boundary faces together
595 // (e.g. to keep baffles together)
596 labelList decomposeCombineRegions
598 const scalarField& cellWeights,
599 const boolList& blockedFace,
600 const List<labelPair>& explicitConnections,
604 //- Redecompose according to cell count
605 // keepZoneFaces : find all faceZones from zoned surfaces and keep
606 // owner and neighbour together
607 // keepBaffles : find all baffles and keep them together
608 autoPtr<mapDistributePolyMesh> balance
610 const bool keepZoneFaces,
611 const bool keepBaffles,
612 const scalarField& cellWeights,
613 decompositionMethod& decomposer,
614 fvMeshDistribute& distributor
617 //- Get faces with intersection.
618 labelList intersectedFaces() const;
620 //- Get points on surfaces with intersection and boundary faces.
621 labelList intersectedPoints() const;
623 //- Create patch from set of patches
624 static autoPtr<indirectPrimitivePatch> makePatch
630 //- Helper function to make a pointVectorField with correct
631 // bcs for mesh movement:
632 // - adaptPatchIDs : fixedValue
633 // - processor : calculated (so free to move)
634 // - cyclic/wedge/symmetry : slip
636 static tmp<pointVectorField> makeDisplacementField
638 const pointMesh& pMesh,
639 const labelList& adaptPatchIDs
642 //- Helper function: check that face zones are synced
643 static void checkCoupledFaceZones(const polyMesh&);
648 //- Calculate list of cells to refine.
649 labelList refineCandidates
651 const point& keepPoint,
652 const scalar curvature,
654 const bool featureRefinement,
655 const bool internalRefinement,
656 const bool surfaceRefinement,
657 const bool curvatureRefinement,
658 const label maxGlobalCells,
659 const label maxLocalCells
662 //- Refine some cells
663 autoPtr<mapPolyMesh> refine(const labelList& cellsToRefine);
665 //- Refine some cells and rebalance
666 autoPtr<mapDistributePolyMesh> refineAndBalance
669 decompositionMethod& decomposer,
670 fvMeshDistribute& distributor,
671 const labelList& cellsToRefine,
672 const scalar maxLoadUnbalance
675 //- Balance before refining some cells
676 autoPtr<mapDistributePolyMesh> balanceAndRefine
679 decompositionMethod& decomposer,
680 fvMeshDistribute& distributor,
681 const labelList& cellsToRefine,
682 const scalar maxLoadUnbalance
688 //- Split off unreachable areas of mesh.
689 void baffleAndSplitMesh
691 const bool handleSnapProblems,
692 const bool removeEdgeConnectedCells,
693 const scalarField& perpendicularAngle,
694 const bool mergeFreeStanding,
695 const dictionary& motionDict,
697 const labelList& globalToPatch,
698 const point& keepPoint
701 //- Split off (with optional buffer layers) unreachable areas
702 // of mesh. Does not introduce baffles.
703 autoPtr<mapPolyMesh> splitMesh
705 const label nBufferLayers,
706 const labelList& globalToPatch,
707 const point& keepPoint
710 //- Find boundary points that connect to more than one cell
711 // region and split them.
712 autoPtr<mapPolyMesh> dupNonManifoldPoints();
714 //- Create baffle for every internal face where ownPatch != -1.
715 // External faces get repatched according to ownPatch (neiPatch
716 // should be -1 for these)
717 autoPtr<mapPolyMesh> createBaffles
719 const labelList& ownPatch,
720 const labelList& neiPatch
723 //- Create baffles for faces straddling zoned surfaces. Return
725 autoPtr<mapPolyMesh> createZoneBaffles
727 const labelList& globalToPatch,
731 //- Return a list of coupled face pairs, i.e. faces that
732 // use the same vertices.
733 List<labelPair> getDuplicateFaces(const labelList& testFaces) const;
735 //- Merge baffles. Gets pairs of faces.
736 autoPtr<mapPolyMesh> mergeBaffles(const List<labelPair>&);
738 //- Put faces/cells into zones according to surface specification.
739 // Returns null if no zone surfaces present. Region containing
740 // the keepPoint will not be put into a cellZone.
741 autoPtr<mapPolyMesh> zonify
743 const point& keepPoint,
744 const bool allowFreeStandingZoneFaces
748 // Other topo changes
750 //- Helper:append patch to end of mesh.
751 static label appendPatch
754 const label insertPatchI,
759 //- Helper:add patch to mesh. Update all registered fields.
760 // Used by addMeshedPatch to add patches originating from surfaces.
761 static label addPatch(fvMesh&, const word& name, const dictionary&);
763 //- Add patch originating from meshing. Update meshedPatches_.
764 label addMeshedPatch(const word& name, const dictionary&);
766 //- Get patchIDs for patches added in addMeshedPatch.
767 labelList meshedPatches() const;
769 //- Split mesh. Keep part containing point.
770 autoPtr<mapPolyMesh> splitMeshRegions(const point& keepPoint);
772 //- Update local numbering for mesh redistribution
773 void distribute(const mapDistributePolyMesh&);
775 //- Update for external change to mesh. changedFaces are in new mesh
780 const labelList& changedFaces
784 // Restoring : is where other processes delete and reinsert data.
786 //- Signal points/face/cells for which to store data
789 const labelList& pointsToStore,
790 const labelList& facesToStore,
791 const labelList& cellsToStore
794 //- Update local numbering + undo
795 // Data to restore given as new pointlabel + stored pointlabel
796 // (i.e. what was in pointsToStore)
800 const labelList& changedFaces,
801 const Map<label>& pointsToRestore,
802 const Map<label>& facesToRestore,
803 const Map<label>& cellsToRestore
806 // Merging coplanar faces and edges
808 label mergePatchFacesUndo
811 const scalar concaveCos,
812 const labelList& patchIDs,
813 const dictionary& motionDict
816 autoPtr<mapPolyMesh> doRemovePoints
818 removePoints& pointRemover,
819 const boolList& pointCanBeDeleted
822 autoPtr<mapPolyMesh> doRestorePoints
824 removePoints& pointRemover,
825 const labelList& facesToRestore
828 labelList collectFaces
830 const labelList& candidateFaces,
831 const labelHashSet& set
834 // Pick up faces of cells of faces in set.
835 labelList growFaceCellFace
837 const labelHashSet& set
843 const dictionary& motionDict
849 //- Debugging: check that all faces still obey start()>end()
852 //- Compare two lists over all boundary faces
854 void testSyncBoundaryFaceList
856 const scalar mergeDistance,
862 //- Print some mesh stats.
863 void printMeshInfo(const bool, const string&) const;
865 //- Replacement for Time::timeName() : return oldInstance (if
867 word timeName() const;
869 //- Set instance of all local IOobjects
870 void setInstance(const fileName&);
872 //- Write mesh and all data
875 //- Write refinement level as volScalarFields for postprocessing
876 void dumpRefinementLevel() const;
878 //- Debug: Write intersection information to OBJ format
879 void dumpIntersections(const fileName& prefix) const;
881 //- Do any one of above IO functions. flag is combination of
883 void write(const label flag, const fileName&) const;
887 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
889 } // End namespace Foam
891 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
894 # include "meshRefinementTemplates.C"
897 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
901 // ************************************************************************* //