FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

meshRefinement.H

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 Class
00025     Foam::meshRefinement
00026 
00027 Description
00028     Helper class which maintains intersections of (changing) mesh with
00029     (static) surfaces.
00030 
00031     Maintains
00032     - per face any intersections of the cc-cc segment with any of the surfaces
00033 
00034 SourceFiles
00035     meshRefinement.C
00036     meshRefinementBaffles.C
00037     meshRefinementMerge.C
00038     meshRefinementProblemCells.C
00039     meshRefinementRefine.C
00040 
00041 \*---------------------------------------------------------------------------*/
00042 
00043 #ifndef meshRefinement_H
00044 #define meshRefinement_H
00045 
00046 #include <dynamicMesh/hexRef8.H>
00047 #include <OpenFOAM/mapPolyMesh.H>
00048 #include <OpenFOAM/autoPtr.H>
00049 #include <OpenFOAM/labelPair.H>
00050 #include <OpenFOAM/indirectPrimitivePatch.H>
00051 #include <OpenFOAM/pointFieldsFwd.H>
00052 #include <OpenFOAM/Tuple2.H>
00053 #include <meshTools/pointIndexHit.H>
00054 
00055 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00056 
00057 namespace Foam
00058 {
00059 
00060 // Class forward declarations
00061 class fvMesh;
00062 class mapDistributePolyMesh;
00063 class decompositionMethod;
00064 class refinementSurfaces;
00065 class shellSurfaces;
00066 class removeCells;
00067 class featureEdgeMesh;
00068 class fvMeshDistribute;
00069 class searchableSurface;
00070 class regionSplit;
00071 class globalIndex;
00072 
00073 /*---------------------------------------------------------------------------*\
00074                            Class meshRefinement Declaration
00075 \*---------------------------------------------------------------------------*/
00076 
00077 class meshRefinement
00078 {
00079 
00080 public:
00081 
00082     // Public data types
00083 
00084         //- Enumeration for debug dumping
00085         enum writeFlag
00086         {
00087             MESH = 1,
00088             SCALARLEVELS = 2,
00089             OBJINTERSECTIONS = 4
00090         };
00091 
00092 
00093         //- Enumeration for how the userdata is to be mapped upon refinement.
00094         enum mapType
00095         {
00096             MASTERONLY = 1, 
00097             KEEPALL = 2,    
00098             REMOVE = 4      
00099         };
00100 
00101 
00102 private:
00103 
00104     // Private data
00105 
00106         //- Reference to mesh
00107         fvMesh& mesh_;
00108 
00109         //- tolerance used for sorting coordinates (used in 'less' routine)
00110         const scalar mergeDistance_;
00111 
00112         //- overwrite the mesh?
00113         const bool overwrite_;
00114 
00115         //- Instance of mesh upon construction. Used when in overwrite_ mode.
00116         const word oldInstance_;
00117 
00118         //- All surface-intersection interaction
00119         const refinementSurfaces& surfaces_;
00120 
00121         //- All shell-refinement interaction
00122         const shellSurfaces& shells_;
00123 
00124         //- refinement engine
00125         hexRef8 meshCutter_;
00126 
00127         //- per cc-cc vector the index of the surface hit
00128         labelIOList surfaceIndex_;
00129 
00130         //- user supplied face based data.
00131         List<Tuple2<mapType, labelList> > userFaceData_;
00132 
00133         //- Meshed patches - are treated differently. Stored as wordList since
00134         //  order changes.
00135         wordList meshedPatches_;
00136 
00137 
00138     // Private Member Functions
00139 
00140         //- Reorder list according to map.
00141         template<class T>
00142         static void updateList
00143         (
00144             const labelList& newToOld,
00145             const T& nullValue,
00146             List<T>& elems
00147         );
00148 
00149         //- Add patchfield of given type to all fields on mesh
00150         template<class GeoField>
00151         static void addPatchFields(fvMesh&, const word& patchFieldType);
00152 
00153         //- Reorder patchfields of all fields on mesh
00154         template<class GeoField>
00155         static void reorderPatchFields(fvMesh&, const labelList& oldToNew);
00156 
00157         //- Find out which faces have changed given cells (old mesh labels)
00158         //  that were marked for refinement.
00159         static labelList getChangedFaces
00160         (
00161             const mapPolyMesh&,
00162             const labelList& oldCellsToRefine
00163         );
00164 
00165         //- Calculate coupled boundary end vector and refinement level
00166         void calcNeighbourData
00167         (
00168             labelList& neiLevel,
00169             pointField& neiCc
00170         ) const;
00171 
00172         //- Find any intersection of surface. Store in surfaceIndex_.
00173         void updateIntersections(const labelList& changedFaces);
00174 
00175         //- Remove cells. Put exposedFaces into exposedPatchIDs.
00176         autoPtr<mapPolyMesh> doRemoveCells
00177         (
00178             const labelList& cellsToRemove,
00179             const labelList& exposedFaces,
00180             const labelList& exposedPatchIDs,
00181             removeCells& cellRemover
00182         );
00183 
00184         // Get cells which are inside any closed surface. Note that
00185         // all closed surfaces
00186         // will have already been oriented to have keepPoint outside.
00187         labelList getInsideCells(const word&) const;
00188 
00189         // Do all to remove inside cells
00190         autoPtr<mapPolyMesh> removeInsideCells
00191         (
00192             const string& msg,
00193             const label exposedPatchI
00194         );
00195 
00196         // For decomposeCombineRegions
00197 
00198             //- Used in decomposeCombineRegions. Given global region per cell
00199             //  determines master processor/cell for regions straddling
00200             //  procboundaries.
00201             void getCoupledRegionMaster
00202             (
00203                 const globalIndex& globalCells,
00204                 const boolList& blockedFace,
00205                 const regionSplit& globalRegion,
00206                 Map<label>& regionToMaster
00207             ) const;
00208 
00209             //- Determine regions that are local to me or coupled ones that
00210             //  are owned by me. Determine representative location.
00211             void calcLocalRegions
00212             (
00213                 const globalIndex& globalCells,
00214                 const labelList& globalRegion,
00215                 const Map<label>& coupledRegionToMaster,
00216                 const scalarField& cellWeights,
00217 
00218                 Map<label>& globalToLocalRegion,
00219                 pointField& localPoints,
00220                 scalarField& localWeights
00221             ) const;
00222 
00223             //- Convert region into global index.
00224             static label getShiftedRegion
00225             (
00226                 const globalIndex& indexer,
00227                 const Map<label>& globalToLocalRegion,
00228                 const Map<label>& coupledRegionToShifted,
00229                 const label globalRegion
00230             );
00231 
00232             //- helper: add element if not in list. Linear search.
00233             static void addUnique(const label, labelList&);
00234 
00235             //- Calculate region connectivity. Major communication.
00236             void calcRegionRegions
00237             (
00238                 const labelList& globalRegion,
00239                 const Map<label>& globalToLocalRegion,
00240                 const Map<label>& coupledRegionToMaster,
00241                 labelListList& regionRegions
00242             ) const;
00243 
00244 
00245         // Refinement candidate selection
00246 
00247             //- Mark cell for refinement (if not already marked). Return false
00248             // if refinelimit hit. Keeps running count (in nRefine) of cells
00249             // marked for refinement
00250             static bool markForRefine
00251             (
00252                 const label markValue,
00253                 const label nAllowRefine,
00254                 label& cellValue,
00255                 label& nRefine
00256             );
00257 
00258             //- Calculate list of cells to refine based on intersection of
00259             //  features.
00260             label markFeatureRefinement
00261             (
00262                 const point& keepPoint,
00263                 const PtrList<featureEdgeMesh>& featureMeshes,
00264                 const labelList& featureLevels,
00265                 const label nAllowRefine,
00266 
00267                 labelList& refineCell,
00268                 label& nRefine
00269             ) const;
00270 
00271             //- Mark cells for refinement-shells based refinement.
00272             label markInternalRefinement
00273             (
00274                 const label nAllowRefine,
00275                 labelList& refineCell,
00276                 label& nRefine
00277             ) const;
00278 
00279             //- Collect faces that are intersected and whose neighbours aren't
00280             //  yet marked  for refinement.
00281             labelList getRefineCandidateFaces
00282             (
00283                 const labelList& refineCell
00284             ) const;
00285 
00286             //- Mark cells for surface intersection based refinement.
00287             label markSurfaceRefinement
00288             (
00289                 const label nAllowRefine,
00290                 const labelList& neiLevel,
00291                 const pointField& neiCc,
00292                 labelList& refineCell,
00293                 label& nRefine
00294             ) const;
00295 
00296             //- Mark cell if local curvature > curvature or
00297             //  markDifferingRegions = true and intersections with different
00298             //  regions.
00299             bool checkCurvature
00300             (
00301                 const scalar curvature,
00302                 const label nAllowRefine,
00303                 const label surfaceLevel,
00304                 const vector& surfaceNormal,
00305                 const label cellI,
00306                 label& cellMaxLevel,
00307                 vector& cellMaxNormal,
00308                 labelList& refineCell,
00309                 label& nRefine
00310             ) const;
00311 
00312 
00313             //- Mark cells for surface curvature based refinement. Marks if
00314             //  local curvature > curvature or if on different regions
00315             //  (markDifferingRegions)
00316             label markSurfaceCurvatureRefinement
00317             (
00318                 const scalar curvature,
00319                 const label nAllowRefine,
00320                 const labelList& neiLevel,
00321                 const pointField& neiCc,
00322                 labelList& refineCell,
00323                 label& nRefine
00324             ) const;
00325 
00326         // Baffle handling
00327 
00328             //- Determine patches for baffles
00329             void getBafflePatches
00330             (
00331                 const labelList& globalToPatch,
00332                 const labelList& neiLevel,
00333                 const pointField& neiCc,
00334                 labelList& ownPatch,
00335                 labelList& neiPatch
00336             ) const;
00337 
00338             //- Determine patch for baffle using some heuristic (and not
00339             //  surface)
00340             label getBafflePatch
00341             (
00342                 const labelList& facePatch,
00343                 const label faceI
00344             ) const;
00345 
00346             //- Repatches external face or creates baffle for internal face
00347             //  with user specified patches (might be different for both sides).
00348             //  Returns label of added face.
00349             label createBaffle
00350             (
00351                 const label faceI,
00352                 const label ownPatch,
00353                 const label neiPatch,
00354                 polyTopoChange& meshMod
00355             ) const;
00356 
00357         // Problem cell handling
00358 
00359             //- Helper function to mark face as being on 'boundary'. Used by
00360             //  markFacesOnProblemCells
00361             void markBoundaryFace
00362             (
00363                 const label faceI,
00364                 boolList& isBoundaryFace,
00365                 boolList& isBoundaryEdge,
00366                 boolList& isBoundaryPoint
00367             ) const;
00368 
00369             void findNearest
00370             (
00371                 const labelList& meshFaces,
00372                 List<pointIndexHit>& nearestInfo,
00373                 labelList& nearestSurface,
00374                 labelList& nearestRegion,
00375                 vectorField& nearestNormal
00376             ) const;
00377 
00378             Map<label> findEdgeConnectedProblemCells
00379             (
00380                 const scalarField& perpendicularAngle,
00381                 const labelList&
00382             ) const;
00383 
00384             bool isCollapsedFace
00385             (
00386                 const pointField&,
00387                 const pointField& neiCc,
00388                 const scalar minFaceArea,
00389                 const scalar maxNonOrtho,
00390                 const label faceI
00391             ) const;
00392 
00393             bool isCollapsedCell
00394             (
00395                 const pointField&,
00396                 const scalar volFraction,
00397                 const label cellI
00398             ) const;
00399 
00400             //- Returns list with for every internal face -1 or the patch
00401             //  they should be baffled into. If removeEdgeConnectedCells is set
00402             //  removes cells based on perpendicularAngle.
00403             labelList markFacesOnProblemCells
00404             (
00405                 const dictionary& motionDict,
00406                 const bool removeEdgeConnectedCells,
00407                 const scalarField& perpendicularAngle,
00408                 const labelList& globalToPatch
00409             ) const;
00410 
00412             //labelList markFacesOnProblemCellsGeometric
00413             //(
00414             //    const dictionary& motionDict
00415             //) const;
00416 
00417 
00418         // Baffle merging
00419 
00420             //- Extract those baffles (duplicate) faces that are on the edge
00421             //  of a baffle region. These are candidates for merging.
00422             List<labelPair> filterDuplicateFaces(const List<labelPair>&) const;
00423 
00424 
00425         // Zone handling
00426 
00427             //- Finds zone per cell for cells inside closed named surfaces.
00428             //  (uses geometric test for insideness)
00429             //  Adapts namedSurfaceIndex so all faces on boundary of cellZone
00430             //  have corresponding faceZone.
00431             void findCellZoneGeometric
00432             (
00433                 const labelList& closedNamedSurfaces,
00434                 labelList& namedSurfaceIndex,
00435                 const labelList& surfaceToCellZone,
00436                 labelList& cellToZone
00437             ) const;
00438 
00439             //- Determines cell zone from cell region information.
00440             bool calcRegionToZone
00441             (
00442                 const label surfZoneI,
00443                 const label ownRegion,
00444                 const label neiRegion,
00445 
00446                 labelList& regionToCellZone
00447             ) const;
00448 
00449             //- Finds zone per cell. Uses topological walk with all faces
00450             //  marked in namedSurfaceIndex regarded as blocked.
00451             void findCellZoneTopo
00452             (
00453                 const point& keepPoint,
00454                 const labelList& namedSurfaceIndex,
00455                 const labelList& surfaceToCellZone,
00456                 labelList& cellToZone
00457             ) const;
00458 
00459             void makeConsistentFaceIndex
00460             (
00461                 const labelList& cellToZone,
00462                 labelList& namedSurfaceIndex
00463             ) const;
00464 
00465 
00466         //- Disallow default bitwise copy construct
00467         meshRefinement(const meshRefinement&);
00468 
00469         //- Disallow default bitwise assignment
00470         void operator=(const meshRefinement&);
00471 
00472 public:
00473 
00474     //- Runtime type information
00475     ClassName("meshRefinement");
00476 
00477 
00478     // Constructors
00479 
00480         //- Construct from components
00481         meshRefinement
00482         (
00483             fvMesh& mesh,
00484             const scalar mergeDistance,
00485             const bool overwrite,
00486             const refinementSurfaces&,
00487             const shellSurfaces&
00488         );
00489 
00490 
00491     // Member Functions
00492 
00493         // Access
00494 
00495             //- reference to mesh
00496             const fvMesh& mesh() const
00497             {
00498                 return mesh_;
00499             }
00500             fvMesh& mesh()
00501             {
00502                 return mesh_;
00503             }
00504 
00505             scalar mergeDistance() const
00506             {
00507                 return mergeDistance_;
00508             }
00509 
00510             //- Overwrite the mesh?
00511             bool overwrite() const
00512             {
00513                 return overwrite_;
00514             }
00515 
00516             //- (points)instance of mesh upon construction
00517             const word& oldInstance() const
00518             {
00519                 return oldInstance_;
00520             }
00521 
00522             //- reference to surface search engines
00523             const refinementSurfaces& surfaces() const
00524             {
00525                 return surfaces_;
00526             }
00527 
00528             //- reference to refinement shells (regions)
00529             const shellSurfaces& shells() const
00530             {
00531                 return shells_;
00532             }
00533 
00534             //- reference to meshcutting engine
00535             const hexRef8& meshCutter() const
00536             {
00537                 return meshCutter_;
00538             }
00539 
00540             //- per start-end edge the index of the surface hit
00541             const labelList& surfaceIndex() const
00542             {
00543                 return surfaceIndex_;
00544             }
00545 
00546             labelList& surfaceIndex()
00547             {
00548                 return surfaceIndex_;
00549             }
00550 
00551             //- Additional face data that is maintained across
00552             //  topo changes. Every entry is a list over all faces.
00553             //  Bit of a hack. Additional flag to say whether to maintain master
00554             //  only (false) or increase set to account for face-from-face.
00555             const List<Tuple2<mapType, labelList> >& userFaceData() const
00556             {
00557                 return userFaceData_;
00558             }
00559 
00560             List<Tuple2<mapType, labelList> >& userFaceData()
00561             {
00562                 return userFaceData_;
00563             }
00564 
00565 
00566         // Other
00567 
00568             //- Count number of intersections (local)
00569             label countHits() const;
00570 
00571             //- Helper function to get decomposition such that all connected
00572             //  regions get moved onto one processor. Used to prevent baffles
00573             //  straddling processor boundaries. explicitConnections is to
00574             //  keep pairs of non-coupled boundary faces together
00575             //  (e.g. to keep baffles together)
00576             labelList decomposeCombineRegions
00577             (
00578                 const scalarField& cellWeights,
00579                 const boolList& blockedFace,
00580                 const List<labelPair>& explicitConnections,
00581                 decompositionMethod&
00582             ) const;
00583 
00584             //- Redecompose according to cell count
00585             //  keepZoneFaces : find all faceZones from zoned surfaces and keep
00586             //                  owner and neighbour together
00587             //  keepBaffles   : find all baffles and keep them together
00588             autoPtr<mapDistributePolyMesh> balance
00589             (
00590                 const bool keepZoneFaces,
00591                 const bool keepBaffles,
00592                 const scalarField& cellWeights,
00593                 decompositionMethod& decomposer,
00594                 fvMeshDistribute& distributor
00595             );
00596 
00597             //- Get faces with intersection.
00598             labelList intersectedFaces() const;
00599 
00600             //- Get points on surfaces with intersection and boundary faces.
00601             labelList intersectedPoints() const;
00602 
00603             //- Create patch from set of patches
00604             static autoPtr<indirectPrimitivePatch> makePatch
00605             (
00606                 const polyMesh&,
00607                 const labelList&
00608             );
00609 
00610             //- Helper function to make a pointVectorField with correct
00611             //  bcs for mesh movement:
00612             //  - adaptPatchIDs         : fixedValue
00613             //  - processor             : calculated (so free to move)
00614             //  - cyclic/wedge/symmetry : slip
00615             //  - other                 : slip
00616             static tmp<pointVectorField> makeDisplacementField
00617             (
00618                 const pointMesh& pMesh,
00619                 const labelList& adaptPatchIDs
00620             );
00621 
00622             //- Helper function: check that face zones are synced
00623             static void checkCoupledFaceZones(const polyMesh&);
00624 
00625 
00626         // Refinement
00627 
00628             //- Calculate list of cells to refine.
00629             labelList refineCandidates
00630             (
00631                 const point& keepPoint,
00632                 const scalar curvature,
00633 
00634                 const PtrList<featureEdgeMesh>& featureMeshes,
00635                 const labelList& featureLevels,
00636 
00637                 const bool featureRefinement,
00638                 const bool internalRefinement,
00639                 const bool surfaceRefinement,
00640                 const bool curvatureRefinement,
00641                 const label maxGlobalCells,
00642                 const label maxLocalCells
00643             ) const;
00644 
00645             //- Refine some cells
00646             autoPtr<mapPolyMesh> refine(const labelList& cellsToRefine);
00647 
00648             //- Refine some cells and rebalance
00649             autoPtr<mapDistributePolyMesh> refineAndBalance
00650             (
00651                 const string& msg,
00652                 decompositionMethod& decomposer,
00653                 fvMeshDistribute& distributor,
00654                 const labelList& cellsToRefine,
00655                 const scalar maxLoadUnbalance
00656             );
00657 
00658             //- Balance before refining some cells
00659             autoPtr<mapDistributePolyMesh> balanceAndRefine
00660             (
00661                 const string& msg,
00662                 decompositionMethod& decomposer,
00663                 fvMeshDistribute& distributor,
00664                 const labelList& cellsToRefine,
00665                 const scalar maxLoadUnbalance
00666             );
00667 
00668 
00669         // Baffle handling
00670 
00671             //- Split off unreachable areas of mesh.
00672             void baffleAndSplitMesh
00673             (
00674                 const bool handleSnapProblems,
00675                 const bool removeEdgeConnectedCells,
00676                 const scalarField& perpendicularAngle,
00677                 const bool mergeFreeStanding,
00678                 const dictionary& motionDict,
00679                 Time& runTime,
00680                 const labelList& globalToPatch,
00681                 const point& keepPoint
00682             );
00683 
00684             //- Split off (with optional buffer layers) unreachable areas
00685             //  of mesh. Does not introduce baffles.
00686             autoPtr<mapPolyMesh> splitMesh
00687             (
00688                 const label nBufferLayers,
00689                 const labelList& globalToPatch,
00690                 const point& keepPoint
00691             );
00692 
00693             //- Find boundary points that connect to more than one cell
00694             //  region and split them.
00695             autoPtr<mapPolyMesh> dupNonManifoldPoints();
00696 
00697             //- Create baffle for every internal face where ownPatch != -1.
00698             //  External faces get repatched according to ownPatch (neiPatch
00699             //  should be -1 for these)
00700             autoPtr<mapPolyMesh> createBaffles
00701             (
00702                 const labelList& ownPatch,
00703                 const labelList& neiPatch
00704             );
00705 
00706             //- Return a list of coupled face pairs, i.e. faces that
00707             //  use the same vertices.
00708             List<labelPair> getDuplicateFaces(const labelList& testFaces) const;
00709 
00710             //- Merge baffles. Gets pairs of faces.
00711             autoPtr<mapPolyMesh> mergeBaffles(const List<labelPair>&);
00712 
00713             //- Put faces/cells into zones according to surface specification.
00714             //  Returns null if no zone surfaces present. Region containing
00715             //  the keepPoint will not be put into a cellZone.
00716             autoPtr<mapPolyMesh> zonify
00717             (
00718                 const point& keepPoint,
00719                 const bool allowFreeStandingZoneFaces
00720             );
00721 
00722 
00723         // Other topo changes
00724 
00725             //- Helper:add patch to mesh. Update all registered fields.
00726             //  Use addMeshedPatch to add patches originating from surfaces.
00727             static label addPatch(fvMesh&, const word& name, const word& type);
00728 
00729             //- Add patch originating from meshing. Update meshedPatches_.
00730             label addMeshedPatch(const word& name, const word& type);
00731 
00732             //- Get patchIDs for patches added in addMeshedPatch.
00733             labelList meshedPatches() const;
00734 
00735             //- Split mesh. Keep part containing point.
00736             autoPtr<mapPolyMesh> splitMeshRegions(const point& keepPoint);
00737 
00738             //- Update local numbering for mesh redistribution
00739             void distribute(const mapDistributePolyMesh&);
00740 
00741             //- Update for external change to mesh. changedFaces are in new mesh
00742             //  face labels.
00743             void updateMesh
00744             (
00745                 const mapPolyMesh&,
00746                 const labelList& changedFaces
00747             );
00748 
00749 
00750             // Restoring : is where other processes delete and reinsert data.
00751 
00752                 //- Signal points/face/cells for which to store data
00753                 void storeData
00754                 (
00755                     const labelList& pointsToStore,
00756                     const labelList& facesToStore,
00757                     const labelList& cellsToStore
00758                 );
00759 
00760                 //- Update local numbering + undo
00761                 //  Data to restore given as new pointlabel + stored pointlabel
00762                 //  (i.e. what was in pointsToStore)
00763                 void updateMesh
00764                 (
00765                     const mapPolyMesh&,
00766                     const labelList& changedFaces,
00767                     const Map<label>& pointsToRestore,
00768                     const Map<label>& facesToRestore,
00769                     const Map<label>& cellsToRestore
00770                 );
00771 
00772 
00773             //- Merge faces on the same patch (usually from exposing refinement)
00774             //  Returns global number of faces merged.
00775             label mergePatchFaces
00776             (
00777                 const scalar minCos,
00778                 const scalar concaveCos,
00779                 const labelList& patchIDs
00780             );
00781 
00782             //- Remove points not used by any face or points used
00783             // by only two faces where the edges are in line
00784             autoPtr<mapPolyMesh> mergeEdges(const scalar minCos);
00785 
00786 
00787         // Debug/IO
00788 
00789             //- Debugging: check that all faces still obey start()>end()
00790             void checkData();
00791 
00792             //- Compare two lists over all boundary faces
00793             template<class T>
00794             void testSyncBoundaryFaceList
00795             (
00796                 const scalar mergeDistance,
00797                 const string&,
00798                 const UList<T>&,
00799                 const UList<T>&
00800             ) const;
00801 
00802             //- Print some mesh stats.
00803             void printMeshInfo(const bool, const string&) const;
00804 
00805             //- Replacement for Time::timeName() : return oldInstance (if
00806             //  overwrite_)
00807             word timeName() const;
00808 
00809             //- Set instance of all local IOobjects
00810             void setInstance(const fileName&);
00811 
00812             //- Write mesh and all data
00813             bool write() const;
00814 
00815             //- Write refinement level as volScalarFields for postprocessing
00816             void dumpRefinementLevel() const;
00817 
00818             //- Debug: Write intersection information to OBJ format
00819             void dumpIntersections(const fileName& prefix) const;
00820 
00821             //- Do any one of above IO functions. flag is combination of
00822             //  writeFlag values.
00823             void write(const label flag, const fileName&) const;
00824 };
00825 
00826 
00827 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00828 
00829 } // End namespace Foam
00830 
00831 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00832 
00833 #ifdef NoRepository
00834 #   include <autoMesh/meshRefinementTemplates.C>
00835 #endif
00836 
00837 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00838 
00839 #endif
00840 
00841 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines