OpenFOAM logo
The Open Source CFD Toolbox
  Source Guide OpenCFD Solutions Contact OpenFOAM

meshRefinement Class Reference

Helper class which maintains intersections of (changing) mesh with (static) surfaces. More...

Collaboration diagram for meshRefinement:

List of all members.


Public Types

enum  writeFlag { MESH = 1, SCALARLEVELS = 2, OBJINTERSECTIONS = 4 }
 Enumeration for debug dumping. More...
enum  mapType { MASTERONLY = 1, KEEPALL = 2, REMOVE = 4 }
 Enumeration for how the userdata is to be mapped upon refinement. More...

Public Member Functions

 ClassName ("meshRefinement")
 Runtime type information.
 meshRefinement (fvMesh &mesh, const scalar mergeDistance, const bool overwrite, const refinementSurfaces &, const shellSurfaces &)
 Construct from components.
const fvMeshmesh () const
 reference to mesh
fvMeshmesh ()
scalar mergeDistance () const
bool overwrite () const
 Overwrite the mesh?
const wordoldInstance () const
 (points)instance of mesh upon construction
const refinementSurfacessurfaces () const
 reference to surface search engines
const shellSurfacesshells () const
 reference to refinement shells (regions)
const hexRef8meshCutter () const
 reference to meshcutting engine
const labelListsurfaceIndex () const
 per start-end edge the index of the surface hit
labelListsurfaceIndex ()
const List< Tuple2< mapType,
labelList > > & 
userFaceData () const
 Additional face data that is maintained across.
List< Tuple2< mapType,
labelList > > & 
userFaceData ()
label countHits () const
 Count number of intersections (local).
labelList decomposeCombineRegions (const boolList &blockedFace, const List< labelPair > &explicitConnections, decompositionMethod &) const
 Helper function to get decomposition such that all connected.
autoPtr< mapDistributePolyMeshbalance (const bool keepZoneFaces, const bool keepBaffles, decompositionMethod &decomposer, fvMeshDistribute &distributor)
 Redecompose according to cell count.
labelList intersectedFaces () const
 Get faces with intersection.
labelList intersectedPoints () const
 Get points on surfaces with intersection and boundary faces.
labelList refineCandidates (const point &keepPoint, const scalar curvature, const PtrList< featureEdgeMesh > &featureMeshes, const labelList &featureLevels, const bool featureRefinement, const bool internalRefinement, const bool surfaceRefinement, const bool curvatureRefinement, const label maxGlobalCells, const label maxLocalCells) const
 Calculate list of cells to refine.
autoPtr< mapPolyMeshrefine (const labelList &cellsToRefine)
 Refine some cells.
autoPtr< mapDistributePolyMeshrefineAndBalance (const string &msg, decompositionMethod &decomposer, fvMeshDistribute &distributor, const labelList &cellsToRefine)
 Refine some cells and rebalance.
void baffleAndSplitMesh (const bool handleSnapProblems, const bool removeEdgeConnectedCells, const scalarField &perpendicularAngle, const bool mergeFreeStanding, const dictionary &motionDict, Time &runTime, const labelList &globalToPatch, const point &keepPoint)
 Split off unreachable areas of mesh.
autoPtr< mapPolyMeshsplitMesh (const label nBufferLayers, const labelList &globalToPatch, const point &keepPoint)
 Split off (with optional buffer layers) unreachable areas.
autoPtr< mapPolyMeshdupNonManifoldPoints ()
 Find boundary points that connect to more than one cell.
autoPtr< mapPolyMeshcreateBaffles (const labelList &ownPatch, const labelList &neiPatch)
 Create baffle for every internal face where ownPatch != -1.
List< labelPairgetDuplicateFaces (const labelList &testFaces) const
 Return a list of coupled face pairs, i.e. faces that.
autoPtr< mapPolyMeshmergeBaffles (const List< labelPair > &)
 Merge baffles. Gets pairs of faces.
autoPtr< mapPolyMeshzonify (const point &keepPoint)
 Put faces/cells into zones according to surface specification.
label addMeshedPatch (const word &name, const word &type)
 Add patch originating from meshing. Update meshedPatches_.
labelList meshedPatches () const
 Get patchIDs for patches added in addMeshedPatch.
autoPtr< mapPolyMeshsplitMeshRegions (const point &keepPoint)
 Split mesh. Keep part containing point.
void distribute (const mapDistributePolyMesh &)
 Update local numbering for mesh redistribution.
void updateMesh (const mapPolyMesh &, const labelList &changedFaces)
 Update for external change to mesh. changedFaces are in new mesh.
void storeData (const labelList &pointsToStore, const labelList &facesToStore, const labelList &cellsToStore)
 Signal points/face/cells for which to store data.
void updateMesh (const mapPolyMesh &, const labelList &changedFaces, const Map< label > &pointsToRestore, const Map< label > &facesToRestore, const Map< label > &cellsToRestore)
 Update local numbering + undo.
label mergePatchFaces (const scalar minCos, const scalar concaveCos, const labelList &patchIDs)
 Merge faces on the same patch (usually from exposing refinement).
autoPtr< mapPolyMeshmergeEdges (const scalar minCos)
 Remove points not used by any face or points used.
void checkData ()
 Debugging: check that all faces still obey start()>end().
template<class T >
void testSyncBoundaryFaceList (const scalar mergeDistance, const string &, const UList< T > &, const UList< T > &) const
 Compare two lists over all boundary faces.
void printMeshInfo (const bool, const string &) const
 Print some mesh stats.
word timeName () const
 Replacement for Time::timeName() : return oldInstance (if.
void setInstance (const fileName &)
 Set instance of all local IOobjects.
bool write () const
 Write mesh and all data.
void dumpRefinementLevel () const
 Write refinement level as volScalarFields for postprocessing.
void dumpIntersections (const fileName &prefix) const
 Debug: Write intersection information to OBJ format.
void write (const label flag, const fileName &) const
 Do any one of above IO functions. flag is combination of.

Static Public Member Functions

static autoPtr
< indirectPrimitivePatch
makePatch (const polyMesh &, const labelList &)
 Create patch from set of patches.
static tmp< pointVectorFieldmakeDisplacementField (const pointMesh &pMesh, const labelList &adaptPatchIDs)
 Helper function to make a pointVectorField with correct.
static void checkCoupledFaceZones (const polyMesh &)
 Helper function: check that face zones are synced.
static label addPatch (fvMesh &, const word &name, const word &type)
 Helper:add patch to mesh. Update all registered fields.

Detailed Description

Helper class which maintains intersections of (changing) mesh with (static) surfaces.

Maintains

  • per face any intersections of the cc-cc segment with any of the surfaces

Source files

Definition at line 72 of file meshRefinement.H.


Member Enumeration Documentation

enum writeFlag

Enumeration for debug dumping.

Enumerator:
MESH 
SCALARLEVELS 
OBJINTERSECTIONS 

Definition at line 80 of file meshRefinement.H.

enum mapType

Enumeration for how the userdata is to be mapped upon refinement.

Enumerator:
MASTERONLY  maintain master only
KEEPALL  have slaves (upon refinement) from master
REMOVE  set value to -1 any face that was refined

Definition at line 89 of file meshRefinement.H.


Constructor & Destructor Documentation

meshRefinement ( fvMesh mesh,
const scalar  mergeDistance,
const bool  overwrite,
const refinementSurfaces surfaces,
const shellSurfaces shells 
)

Construct from components.

Definition at line 813 of file meshRefinement.C.


Member Function Documentation

ClassName ( "meshRefinement"   ) 

Runtime type information.

const fvMesh& mesh (  )  const [inline]

reference to mesh

Definition at line 521 of file meshRefinement.H.

Referenced by autoSnapDriver::repatchToSurface().

Here is the caller graph for this function:

fvMesh& mesh (  )  [inline]

Definition at line 525 of file meshRefinement.H.

scalar mergeDistance (  )  const [inline]

Definition at line 530 of file meshRefinement.H.

bool overwrite (  )  const [inline]

Overwrite the mesh?

Definition at line 536 of file meshRefinement.H.

Referenced by meshRefinement::mergePatchFaces().

Here is the caller graph for this function:

const word& oldInstance (  )  const [inline]

(points)instance of mesh upon construction

Definition at line 542 of file meshRefinement.H.

Referenced by meshRefinement::mergePatchFaces().

Here is the caller graph for this function:

const refinementSurfaces& surfaces (  )  const [inline]

reference to surface search engines

Definition at line 548 of file meshRefinement.H.

Referenced by autoSnapDriver::repatchToSurface().

Here is the caller graph for this function:

const shellSurfaces& shells (  )  const [inline]

reference to refinement shells (regions)

Definition at line 554 of file meshRefinement.H.

const hexRef8& meshCutter (  )  const [inline]

reference to meshcutting engine

Definition at line 560 of file meshRefinement.H.

const labelList& surfaceIndex (  )  const [inline]

per start-end edge the index of the surface hit

Definition at line 566 of file meshRefinement.H.

labelList& surfaceIndex (  )  [inline]

Definition at line 571 of file meshRefinement.H.

const List<Tuple2<mapType, labelList> >& userFaceData (  )  const [inline]

Additional face data that is maintained across.

topo changes. Every entry is a list over all faces. Bit of a hack. Additional flag to say whether to maintain master only (false) or increase set to account for face-from-face.

Definition at line 582 of file meshRefinement.H.

List<Tuple2<mapType, labelList> >& userFaceData (  )  [inline]

Definition at line 587 of file meshRefinement.H.

Foam::label countHits (  )  const

Count number of intersections (local).

Definition at line 897 of file meshRefinement.C.

References Foam::identity().

Here is the call graph for this function:

Foam::labelList decomposeCombineRegions ( const boolList blockedFace,
const List< labelPair > &  explicitConnections,
decompositionMethod decomposer 
) const

Helper function to get decomposition such that all connected.

regions get moved onto one processor. Used to prevent baffles straddling processor boundaries. explicitConnections is to keep pairs of non-coupled boundary faces together (e.g. to keep baffles together)

Definition at line 917 of file meshRefinement.C.

Foam::autoPtr< Foam::mapDistributePolyMesh > balance ( const bool  keepZoneFaces,
const bool  keepBaffles,
decompositionMethod decomposer,
fvMeshDistribute distributor 
)

Redecompose according to cell count.

keepZoneFaces : find all faceZones from zoned surfaces and keep owner and neighbour together keepBaffles : find all baffles and keep them together

Definition at line 1049 of file meshRefinement.C.

Foam::labelList intersectedFaces (  )  const

Get faces with intersection.

Definition at line 1199 of file meshRefinement.C.

Foam::labelList intersectedPoints (  )  const

Get points on surfaces with intersection and boundary faces.

Definition at line 1226 of file meshRefinement.C.

Foam::autoPtr< Foam::indirectPrimitivePatch > makePatch ( const polyMesh mesh,
const labelList patchIDs 
) [static]

Create patch from set of patches.

Definition at line 1295 of file meshRefinement.C.

Referenced by autoSnapDriver::repatchToSurface().

Here is the caller graph for this function:

Foam::tmp< Foam::pointVectorField > makeDisplacementField ( const pointMesh pMesh,
const labelList adaptPatchIDs 
) [static]

Helper function to make a pointVectorField with correct.

bcs for mesh movement:

  • adaptPatchIDs : fixedValue
  • processor : calculated (so free to move)
  • cyclic/wedge/symmetry : slip
  • other : slip

Definition at line 1341 of file meshRefinement.C.

void checkCoupledFaceZones ( const polyMesh mesh  )  [static]

Helper function: check that face zones are synced.

Definition at line 1396 of file meshRefinement.C.

Referenced by autoHexMeshDriver::autoHexMeshDriver().

Here is the caller graph for this function:

Foam::labelList refineCandidates ( const point keepPoint,
const scalar  curvature,
const PtrList< featureEdgeMesh > &  featureMeshes,
const labelList featureLevels,
const bool  featureRefinement,
const bool  internalRefinement,
const bool  surfaceRefinement,
const bool  curvatureRefinement,
const label  maxGlobalCells,
const label  maxLocalCells 
) const

Calculate list of cells to refine.

Disable refinement shortcut. nAllowRefine is per processor limit.

Definition at line 1039 of file meshRefinementRefine.C.

References Foam::endl(), forAll, polyMesh::globalData(), Foam::Info, refinementSurfaces::maxLevel(), refinementSurfaces::minLevel(), primitiveMesh::nCells(), primitiveMesh::nFaces(), primitiveMesh::nInternalFaces(), Pstream::nProcs(), and globalMeshData::nTotalCells().

Here is the call graph for this function:

Foam::autoPtr< Foam::mapPolyMesh > refine ( const labelList cellsToRefine  ) 

Refine some cells.

Definition at line 1199 of file meshRefinementRefine.C.

Foam::autoPtr< Foam::mapDistributePolyMesh > refineAndBalance ( const string msg,
decompositionMethod decomposer,
fvMeshDistribute distributor,
const labelList cellsToRefine 
)

Refine some cells and rebalance.

Definition at line 1242 of file meshRefinementRefine.C.

void baffleAndSplitMesh ( const bool  handleSnapProblems,
const bool  removeEdgeConnectedCells,
const scalarField perpendicularAngle,
const bool  mergeFreeStanding,
const dictionary motionDict,
Time runTime,
const labelList globalToPatch,
const point keepPoint 
)

Split off unreachable areas of mesh.

Definition at line 1443 of file meshRefinementBaffles.C.

Foam::autoPtr< Foam::mapPolyMesh > splitMesh ( const label  nBufferLayers,
const labelList globalToPatch,
const point keepPoint 
)

Split off (with optional buffer layers) unreachable areas.

of mesh. Does not introduce baffles.

Definition at line 1680 of file meshRefinementBaffles.C.

Foam::autoPtr< Foam::mapPolyMesh > dupNonManifoldPoints (  ) 

Find boundary points that connect to more than one cell.

region and split them.

Definition at line 1969 of file meshRefinementBaffles.C.

Foam::autoPtr< Foam::mapPolyMesh > createBaffles ( const labelList ownPatch,
const labelList neiPatch 
)

Create baffle for every internal face where ownPatch != -1.

External faces get repatched according to ownPatch (neiPatch should be -1 for these)

Redo the intersections on the newly create baffle faces. Note that

this changes also the cell centre positions.

Definition at line 375 of file meshRefinementBaffles.C.

Foam::List< Foam::labelPair > getDuplicateFaces ( const labelList testFaces  )  const

Return a list of coupled face pairs, i.e. faces that.

use the same vertices.

Definition at line 523 of file meshRefinementBaffles.C.

Foam::autoPtr< Foam::mapPolyMesh > mergeBaffles ( const List< labelPair > &  couples  ) 

Merge baffles. Gets pairs of faces.

Definition at line 733 of file meshRefinementBaffles.C.

Foam::autoPtr< Foam::mapPolyMesh > zonify ( const point keepPoint  ) 

Put faces/cells into zones according to surface specification.

Returns null if no zone surfaces present. Region containing the keepPoint will not be put into a cellZone.

Definition at line 2027 of file meshRefinementBaffles.C.

Foam::label addPatch ( fvMesh mesh,
const word name,
const word type 
) [static]

Helper:add patch to mesh. Update all registered fields.

Use addMeshedPatch to add patches originating from surfaces.

Definition at line 1493 of file meshRefinement.C.

Foam::label addMeshedPatch ( const word name,
const word type 
)

Add patch originating from meshing. Update meshedPatches_.

Definition at line 1664 of file meshRefinement.C.

Foam::labelList meshedPatches (  )  const

Get patchIDs for patches added in addMeshedPatch.

Definition at line 1691 of file meshRefinement.C.

Foam::autoPtr< Foam::mapPolyMesh > splitMeshRegions ( const point keepPoint  ) 

Split mesh. Keep part containing point.

Definition at line 1711 of file meshRefinement.C.

void distribute ( const mapDistributePolyMesh map  ) 

Update local numbering for mesh redistribution.

Definition at line 1790 of file meshRefinement.C.

void updateMesh ( const mapPolyMesh map,
const labelList changedFaces 
)

Update for external change to mesh. changedFaces are in new mesh.

face labels.

Definition at line 1848 of file meshRefinement.C.

Referenced by meshRefinement::mergePatchFaces(), and meshRefinement::storeData().

Here is the caller graph for this function:

void storeData ( const labelList pointsToStore,
const labelList facesToStore,
const labelList cellsToStore 
)

Signal points/face/cells for which to store data.

Definition at line 1860 of file meshRefinement.C.

References meshRefinement::updateMesh().

Here is the call graph for this function:

void updateMesh ( const mapPolyMesh map,
const labelList changedFaces,
const Map< label > &  pointsToRestore,
const Map< label > &  facesToRestore,
const Map< label > &  cellsToRestore 
)

Update local numbering + undo.

Data to restore given as new pointlabel + stored pointlabel (i.e. what was in pointsToStore)

Definition at line 1877 of file meshRefinement.C.

Foam::autoPtr< Foam::mapPolyMesh > mergeEdges ( const scalar  minCos  ) 

Remove points not used by any face or points used.

by only two faces where the edges are in line

Definition at line 135 of file meshRefinementMerge.C.

void checkData (  ) 

Debugging: check that all faces still obey start()>end().

Definition at line 242 of file meshRefinement.C.

Referenced by autoRefineDriver::doRefine().

Here is the caller graph for this function:

void testSyncBoundaryFaceList ( const scalar  mergeDistance,
const string msg,
const UList< T > &  faceData,
const UList< T > &  syncedFaceData 
) const [inline]

Compare two lists over all boundary faces.

Definition at line 55 of file meshRefinementTemplates.C.

void printMeshInfo ( const bool  debug,
const string msg 
) const

Print some mesh stats.

Definition at line 2010 of file meshRefinement.C.

Foam::word timeName (  )  const

Replacement for Time::timeName() : return oldInstance (if.

Return either time().constant() or oldInstance.

overwrite_)

Definition at line 2058 of file meshRefinement.C.

References Foam::Info.

void setInstance ( const fileName inst  ) 

Set instance of all local IOobjects.

Definition at line 414 of file meshRefinement.C.

bool write (  )  const

Write mesh and all data.

Definition at line 1971 of file meshRefinement.C.

void dumpRefinementLevel (  )  const

Write refinement level as volScalarFields for postprocessing.

Definition at line 2071 of file meshRefinement.C.

void dumpIntersections ( const fileName prefix  )  const

Debug: Write intersection information to OBJ format.

Definition at line 2126 of file meshRefinement.C.

void write ( const label  flag,
const fileName prefix 
) const

Do any one of above IO functions. flag is combination of.

writeFlag values.

Definition at line 2206 of file meshRefinement.C.


The documentation for this class was generated from the following files:
Copyright © 2000-2009 OpenCFD Ltd