ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / mesh / autoMesh / autoHexMesh / meshRefinement / meshRefinement.H
blobbe9994523a6e7fb0e57559c675b6d4c7188cf946
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
19     for more details.
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/>.
24 Class
25     Foam::meshRefinement
27 Description
28     Helper class which maintains intersections of (changing) mesh with
29     (static) surfaces.
31     Maintains
32     - per face any intersections of the cc-cc segment with any of the surfaces
34 SourceFiles
35     meshRefinement.C
36     meshRefinementBaffles.C
37     meshRefinementMerge.C
38     meshRefinementProblemCells.C
39     meshRefinementRefine.C
41 \*---------------------------------------------------------------------------*/
43 #ifndef meshRefinement_H
44 #define meshRefinement_H
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"
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 namespace Foam
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;
74 /*---------------------------------------------------------------------------*\
75                            Class meshRefinement Declaration
76 \*---------------------------------------------------------------------------*/
78 class meshRefinement
81 public:
83     // Public data types
85         //- Enumeration for debug dumping
86         enum writeFlag
87         {
88             MESH = 1,
89             SCALARLEVELS = 2,
90             OBJINTERSECTIONS = 4
91         };
94         //- Enumeration for how the userdata is to be mapped upon refinement.
95         enum mapType
96         {
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 */
100         };
103 private:
105     // Private data
107         //- Reference to mesh
108         fvMesh& 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
129         hexRef8 meshCutter_;
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
138         //  order changes.
139         wordList meshedPatches_;
142     // Private Member Functions
144         //- Reorder list according to map.
145         template<class T>
146         static void updateList
147         (
148             const labelList& newToOld,
149             const T& nullValue,
150             List<T>& elems
151         );
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
164         (
165             const mapPolyMesh&,
166             const labelList& oldCellsToRefine
167         );
169         //- Calculate coupled boundary end vector and refinement level
170         void calcNeighbourData
171         (
172             labelList& neiLevel,
173             pointField& neiCc
174         ) const;
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
181         (
182             const labelList& cellsToRemove,
183             const labelList& exposedFaces,
184             const labelList& exposedPatchIDs,
185             removeCells& cellRemover
186         );
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
195         (
196             const string& msg,
197             const label exposedPatchI
198         );
200         // For decomposeCombineRegions
202             //- Used in decomposeCombineRegions. Given global region per cell
203             //  determines master processor/cell for regions straddling
204             //  procboundaries.
205             void getCoupledRegionMaster
206             (
207                 const globalIndex& globalCells,
208                 const boolList& blockedFace,
209                 const regionSplit& globalRegion,
210                 Map<label>& regionToMaster
211             ) const;
213             //- Determine regions that are local to me or coupled ones that
214             //  are owned by me. Determine representative location.
215             void calcLocalRegions
216             (
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
225             ) const;
227             //- Convert region into global index.
228             static label getShiftedRegion
229             (
230                 const globalIndex& indexer,
231                 const Map<label>& globalToLocalRegion,
232                 const Map<label>& coupledRegionToShifted,
233                 const label globalRegion
234             );
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
241             (
242                 const labelList& globalRegion,
243                 const Map<label>& globalToLocalRegion,
244                 const Map<label>& coupledRegionToMaster,
245                 labelListList& regionRegions
246             ) const;
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
255             (
256                 const label markValue,
257                 const label nAllowRefine,
258                 label& cellValue,
259                 label& nRefine
260             );
262             //- Calculate list of cells to refine based on intersection of
263             //  features.
264             label markFeatureRefinement
265             (
266                 const point& keepPoint,
267                 const label nAllowRefine,
269                 labelList& refineCell,
270                 label& nRefine
271             ) const;
273             //- Mark cells for refinement-shells based refinement.
274             label markInternalRefinement
275             (
276                 const label nAllowRefine,
277                 labelList& refineCell,
278                 label& nRefine
279             ) const;
281             //- Collect faces that are intersected and whose neighbours aren't
282             //  yet marked  for refinement.
283             labelList getRefineCandidateFaces
284             (
285                 const labelList& refineCell
286             ) const;
288             //- Mark cells for surface intersection based refinement.
289             label markSurfaceRefinement
290             (
291                 const label nAllowRefine,
292                 const labelList& neiLevel,
293                 const pointField& neiCc,
294                 labelList& refineCell,
295                 label& nRefine
296             ) const;
298             //- Mark cell if local curvature > curvature or
299             //  markDifferingRegions = true and intersections with different
300             //  regions.
301             bool checkCurvature
302             (
303                 const scalar curvature,
304                 const label nAllowRefine,
305                 const label surfaceLevel,
306                 const vector& surfaceNormal,
307                 const label cellI,
308                 label& cellMaxLevel,
309                 vector& cellMaxNormal,
310                 labelList& refineCell,
311                 label& nRefine
312             ) const;
315             //- Mark cells for surface curvature based refinement. Marks if
316             //  local curvature > curvature or if on different regions
317             //  (markDifferingRegions)
318             label markSurfaceCurvatureRefinement
319             (
320                 const scalar curvature,
321                 const label nAllowRefine,
322                 const labelList& neiLevel,
323                 const pointField& neiCc,
324                 labelList& refineCell,
325                 label& nRefine
326             ) const;
328         // Baffle handling
330             //- Get faces to repatch. Returns map from face to patch.
331             Map<label> getZoneBafflePatches
332             (
333                 const bool allowBoundary,
334                 const labelList& globalToPatch
335             ) const;
337             //- Determine patches for baffles
338             void getBafflePatches
339             (
340                 const labelList& globalToPatch,
341                 const labelList& neiLevel,
342                 const pointField& neiCc,
343                 labelList& ownPatch,
344                 labelList& neiPatch
345             ) const;
347             //- Determine patch for baffle using some heuristic (and not
348             //  surface)
349             label getBafflePatch
350             (
351                 const labelList& facePatch,
352                 const label faceI
353             ) const;
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.
358             label createBaffle
359             (
360                 const label faceI,
361                 const label ownPatch,
362                 const label neiPatch,
363                 polyTopoChange& meshMod
364             ) const;
366         // Problem cell handling
368             //- Helper function to mark face as being on 'boundary'. Used by
369             //  markFacesOnProblemCells
370             void markBoundaryFace
371             (
372                 const label faceI,
373                 boolList& isBoundaryFace,
374                 boolList& isBoundaryEdge,
375                 boolList& isBoundaryPoint
376             ) const;
378             void findNearest
379             (
380                 const labelList& meshFaces,
381                 List<pointIndexHit>& nearestInfo,
382                 labelList& nearestSurface,
383                 labelList& nearestRegion,
384                 vectorField& nearestNormal
385             ) const;
387             Map<label> findEdgeConnectedProblemCells
388             (
389                 const scalarField& perpendicularAngle,
390                 const labelList&
391             ) const;
393             bool isCollapsedFace
394             (
395                 const pointField&,
396                 const pointField& neiCc,
397                 const scalar minFaceArea,
398                 const scalar maxNonOrtho,
399                 const label faceI
400             ) const;
402             bool isCollapsedCell
403             (
404                 const pointField&,
405                 const scalar volFraction,
406                 const label cellI
407             ) const;
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
413             (
414                 const dictionary& motionDict,
415                 const bool removeEdgeConnectedCells,
416                 const scalarField& perpendicularAngle,
417                 const labelList& globalToPatch
418             ) const;
421         // Baffle merging
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;
428         // Zone handling
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
435             (
436                 const labelList& closedNamedSurfaces,
437                 labelList& namedSurfaceIndex,
438                 const labelList& surfaceToCellZone,
439                 labelList& cellToZone
440             ) const;
442             //- Finds zone per cell for cells inside named surfaces that have
443             //  an inside point specified.
444             void findCellZoneInsideWalk
445             (
446                 const labelList& locationSurfaces,
447                 const labelList& namedSurfaceIndex,
448                 const labelList& surfaceToCellZone,
449                 labelList& cellToZone
450             ) const;
452             //- Determines cell zone from cell region information.
453             bool calcRegionToZone
454             (
455                 const label surfZoneI,
456                 const label ownRegion,
457                 const label neiRegion,
459                 labelList& regionToCellZone
460             ) const;
462             //- Finds zone per cell. Uses topological walk with all faces
463             //  marked in namedSurfaceIndex regarded as blocked.
464             void findCellZoneTopo
465             (
466                 const point& keepPoint,
467                 const labelList& namedSurfaceIndex,
468                 const labelList& surfaceToCellZone,
469                 labelList& cellToZone
470             ) const;
472             void makeConsistentFaceIndex
473             (
474                 const labelList& cellToZone,
475                 labelList& namedSurfaceIndex
476             ) const;
479         //- Disallow default bitwise copy construct
480         meshRefinement(const meshRefinement&);
482         //- Disallow default bitwise assignment
483         void operator=(const meshRefinement&);
485 public:
487     //- Runtime type information
488     ClassName("meshRefinement");
491     // Constructors
493         //- Construct from components
494         meshRefinement
495         (
496             fvMesh& mesh,
497             const scalar mergeDistance,
498             const bool overwrite,
499             const refinementSurfaces&,
500             const refinementFeatures&,
501             const shellSurfaces&
502         );
505     // Member Functions
507         // Access
509             //- reference to mesh
510             const fvMesh& mesh() const
511             {
512                 return mesh_;
513             }
514             fvMesh& mesh()
515             {
516                 return mesh_;
517             }
519             scalar mergeDistance() const
520             {
521                 return mergeDistance_;
522             }
524             //- Overwrite the mesh?
525             bool overwrite() const
526             {
527                 return overwrite_;
528             }
530             //- (points)instance of mesh upon construction
531             const word& oldInstance() const
532             {
533                 return oldInstance_;
534             }
536             //- reference to surface search engines
537             const refinementSurfaces& surfaces() const
538             {
539                 return surfaces_;
540             }
542             //- reference to feature edge mesh
543             const refinementFeatures& features() const
544             {
545                 return features_;
546             }
548             //- reference to refinement shells (regions)
549             const shellSurfaces& shells() const
550             {
551                 return shells_;
552             }
554             //- reference to meshcutting engine
555             const hexRef8& meshCutter() const
556             {
557                 return meshCutter_;
558             }
560             //- per start-end edge the index of the surface hit
561             const labelList& surfaceIndex() const
562             {
563                 return surfaceIndex_;
564             }
566             labelList& surfaceIndex()
567             {
568                 return surfaceIndex_;
569             }
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
576             {
577                 return userFaceData_;
578             }
580             List<Tuple2<mapType, labelList> >& userFaceData()
581             {
582                 return userFaceData_;
583             }
586         // Other
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
597             (
598                 const scalarField& cellWeights,
599                 const boolList& blockedFace,
600                 const List<labelPair>& explicitConnections,
601                 decompositionMethod&
602             ) const;
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
609             (
610                 const bool keepZoneFaces,
611                 const bool keepBaffles,
612                 const scalarField& cellWeights,
613                 decompositionMethod& decomposer,
614                 fvMeshDistribute& distributor
615             );
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
625             (
626                 const polyMesh&,
627                 const labelList&
628             );
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
635             //  - other                 : slip
636             static tmp<pointVectorField> makeDisplacementField
637             (
638                 const pointMesh& pMesh,
639                 const labelList& adaptPatchIDs
640             );
642             //- Helper function: check that face zones are synced
643             static void checkCoupledFaceZones(const polyMesh&);
646         // Refinement
648             //- Calculate list of cells to refine.
649             labelList refineCandidates
650             (
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
660             ) const;
662             //- Refine some cells
663             autoPtr<mapPolyMesh> refine(const labelList& cellsToRefine);
665             //- Refine some cells and rebalance
666             autoPtr<mapDistributePolyMesh> refineAndBalance
667             (
668                 const string& msg,
669                 decompositionMethod& decomposer,
670                 fvMeshDistribute& distributor,
671                 const labelList& cellsToRefine,
672                 const scalar maxLoadUnbalance
673             );
675             //- Balance before refining some cells
676             autoPtr<mapDistributePolyMesh> balanceAndRefine
677             (
678                 const string& msg,
679                 decompositionMethod& decomposer,
680                 fvMeshDistribute& distributor,
681                 const labelList& cellsToRefine,
682                 const scalar maxLoadUnbalance
683             );
686         // Baffle handling
688             //- Split off unreachable areas of mesh.
689             void baffleAndSplitMesh
690             (
691                 const bool handleSnapProblems,
692                 const bool removeEdgeConnectedCells,
693                 const scalarField& perpendicularAngle,
694                 const bool mergeFreeStanding,
695                 const dictionary& motionDict,
696                 Time& runTime,
697                 const labelList& globalToPatch,
698                 const point& keepPoint
699             );
701             //- Split off (with optional buffer layers) unreachable areas
702             //  of mesh. Does not introduce baffles.
703             autoPtr<mapPolyMesh> splitMesh
704             (
705                 const label nBufferLayers,
706                 const labelList& globalToPatch,
707                 const point& keepPoint
708             );
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
718             (
719                 const labelList& ownPatch,
720                 const labelList& neiPatch
721             );
723             //- Create baffles for faces straddling zoned surfaces. Return
724             //  baffles.
725             autoPtr<mapPolyMesh> createZoneBaffles
726             (
727                 const labelList& globalToPatch,
728                 List<labelPair>&
729             );
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
742             (
743                 const point& keepPoint,
744                 const bool allowFreeStandingZoneFaces
745             );
748         // Other topo changes
750             //- Helper:append patch to end of mesh.
751             static label appendPatch
752             (
753                 fvMesh&,
754                 const label insertPatchI,
755                 const word&,
756                 const dictionary&
757             );
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
776             //  face labels.
777             void updateMesh
778             (
779                 const mapPolyMesh&,
780                 const labelList& changedFaces
781             );
784             // Restoring : is where other processes delete and reinsert data.
786                 //- Signal points/face/cells for which to store data
787                 void storeData
788                 (
789                     const labelList& pointsToStore,
790                     const labelList& facesToStore,
791                     const labelList& cellsToStore
792                 );
794                 //- Update local numbering + undo
795                 //  Data to restore given as new pointlabel + stored pointlabel
796                 //  (i.e. what was in pointsToStore)
797                 void updateMesh
798                 (
799                     const mapPolyMesh&,
800                     const labelList& changedFaces,
801                     const Map<label>& pointsToRestore,
802                     const Map<label>& facesToRestore,
803                     const Map<label>& cellsToRestore
804                 );
806             // Merging coplanar faces and edges
808                 label mergePatchFacesUndo
809                 (
810                     const scalar minCos,
811                     const scalar concaveCos,
812                     const labelList& patchIDs,
813                     const dictionary& motionDict
814                 );
816                 autoPtr<mapPolyMesh> doRemovePoints
817                 (
818                     removePoints& pointRemover,
819                     const boolList& pointCanBeDeleted
820                 );
822                 autoPtr<mapPolyMesh> doRestorePoints
823                 (
824                     removePoints& pointRemover,
825                     const labelList& facesToRestore
826                 );
828                 labelList collectFaces
829                 (
830                     const labelList& candidateFaces,
831                     const labelHashSet& set
832                 ) const;
834                 // Pick up faces of cells of faces in set.
835                 labelList growFaceCellFace
836                 (
837                     const labelHashSet& set
838                 ) const;
840                 label mergeEdgesUndo
841                 (
842                     const scalar minCos,
843                     const dictionary& motionDict
844                 );
847         // Debug/IO
849             //- Debugging: check that all faces still obey start()>end()
850             void checkData();
852             //- Compare two lists over all boundary faces
853             template<class T>
854             void testSyncBoundaryFaceList
855             (
856                 const scalar mergeDistance,
857                 const string&,
858                 const UList<T>&,
859                 const UList<T>&
860             ) const;
862             //- Print some mesh stats.
863             void printMeshInfo(const bool, const string&) const;
865             //- Replacement for Time::timeName() : return oldInstance (if
866             //  overwrite_)
867             word timeName() const;
869             //- Set instance of all local IOobjects
870             void setInstance(const fileName&);
872             //- Write mesh and all data
873             bool write() const;
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
882             //  writeFlag values.
883             void write(const label flag, const fileName&) const;
887 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
889 } // End namespace Foam
891 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
893 #ifdef NoRepository
894 #   include "meshRefinementTemplates.C"
895 #endif
897 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
899 #endif
901 // ************************************************************************* //