FreeFOAM The Cross-Platform CFD Toolkit
meshRefinement Class Reference

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

#include <autoMesh/meshRefinement.H>


Detailed Description

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 scalarField &cellWeights, 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, const scalarField &cellWeights, 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, const scalar maxLoadUnbalance)
 Refine some cells and rebalance.
autoPtr< mapDistributePolyMeshbalanceAndRefine (const string &msg, decompositionMethod &decomposer, fvMeshDistribute &distributor, const labelList &cellsToRefine, const scalar maxLoadUnbalance)
 Balance before refining some cells.
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, const bool allowFreeStandingZoneFaces)
 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.

Member Enumeration Documentation

enum writeFlag

Enumeration for debug dumping.

Enumerator:
MESH 
SCALARLEVELS 
OBJINTERSECTIONS 

Definition at line 85 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 94 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 846 of file meshRefinement.C.

References Foam::identity().


Member Function Documentation

ClassName ( "meshRefinement"  )

Runtime type information.

const fvMesh& mesh ( ) const
inline

reference to mesh

Definition at line 496 of file meshRefinement.H.

fvMesh& mesh ( )
inline

Definition at line 500 of file meshRefinement.H.

scalar mergeDistance ( ) const
inline

Definition at line 505 of file meshRefinement.H.

bool overwrite ( ) const
inline

Overwrite the mesh?

Definition at line 511 of file meshRefinement.H.

const word& oldInstance ( ) const
inline

(points)instance of mesh upon construction

Definition at line 517 of file meshRefinement.H.

const refinementSurfaces& surfaces ( ) const
inline

reference to surface search engines

Definition at line 523 of file meshRefinement.H.

const shellSurfaces& shells ( ) const
inline

reference to refinement shells (regions)

Definition at line 529 of file meshRefinement.H.

const hexRef8& meshCutter ( ) const
inline

reference to meshcutting engine

Definition at line 535 of file meshRefinement.H.

const labelList& surfaceIndex ( ) const
inline

per start-end edge the index of the surface hit

Definition at line 541 of file meshRefinement.H.

labelList& surfaceIndex ( )
inline

Definition at line 546 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 555 of file meshRefinement.H.

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

Definition at line 560 of file meshRefinement.H.

Foam::label countHits ( ) const

Count number of intersections (local)

Definition at line 930 of file meshRefinement.C.

References forAll, PackedList< nBits >::get(), and syncTools::getMasterFaces().

Foam::labelList decomposeCombineRegions ( const scalarField cellWeights,
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 950 of file meshRefinement.C.

References decompositionMethod::decompose(), HashTable< T, Key, Hash >::end(), HashTable< T, label, Hash< label > >::end(), HashTable< T, Key, Hash >::find(), forAll, forAllConstIter, Pstream::mapCombineGather(), and Pstream::mapCombineScatter().

Foam::autoPtr< Foam::mapDistributePolyMesh > balance ( const bool  keepZoneFaces,
const bool  keepBaffles,
const scalarField cellWeights,
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 1091 of file meshRefinement.C.

References fvMeshDistribute::countCells(), decompositionMethod::decompose(), fvMeshDistribute::distribute(), Foam::endl(), refinementSurfaces::faceZoneNames(), ZoneMesh< ZoneType, MeshType >::findZoneID(), forAll, Foam::identity(), Foam::Info, Pstream::listCombineGather(), Pstream::listCombineScatter(), Pstream::parRun(), Foam::Pout, Foam::reduce(), Foam::returnReduce(), List< T >::size(), and syncTools::syncFaceList().

Foam::labelList intersectedFaces ( ) const

Get faces with intersection.

Definition at line 1263 of file meshRefinement.C.

References forAll.

Foam::labelList intersectedPoints ( ) const

Get points on surfaces with intersection and boundary faces.

Definition at line 1290 of file meshRefinement.C.

References f(), and forAll.

Foam::autoPtr< Foam::indirectPrimitivePatch > makePatch ( const polyMesh mesh,
const labelList patchIDs 
)
static
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 1405 of file meshRefinement.C.

References IOobject::AUTO_WRITE, pointMesh::boundary(), Foam::dimLength, forAll, mesh, IOobject::NO_READ, PtrList< T >::size(), objectRegistry::time(), Time::timeName(), and Vector< scalar >::zero.

Referenced by autoLayerDriver::addLayers(), and autoSnapDriver::doSnap().

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 1048 of file meshRefinementRefine.C.

References Foam::endl(), forAll, Foam::Info, Pstream::nProcs(), and List< T >::setSize().

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

Refine some cells.

Definition at line 1208 of file meshRefinementRefine.C.

References polyTopoChange::changeMesh().

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

Refine some cells and rebalance.

Definition at line 1251 of file meshRefinementRefine.C.

References Foam::endl(), Foam::Info, Foam::mag(), Pstream::nProcs(), Foam::Pout, Foam::returnReduce(), and timeName.

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

Balance before refining some cells.

Definition at line 1351 of file meshRefinementRefine.C.

References Foam::endl(), forAll, Foam::Info, Foam::mag(), Pstream::nProcs(), Foam::Pout, Foam::returnReduce(), List< T >::size(), and timeName.

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 
)
Foam::autoPtr< Foam::mapPolyMesh > splitMesh ( const label  nBufferLayers,
const labelList globalToPatch,
const point keepPoint 
)
Foam::autoPtr< Foam::mapPolyMesh > dupNonManifoldPoints ( )

Find boundary points that connect to more than one cell.

region and split them.

Definition at line 1975 of file meshRefinementBaffles.C.

References polyTopoChange::changeMesh(), Foam::endl(), Foam::Info, localPointRegion::meshPointMap(), Foam::returnReduce(), duplicatePoints::setRefinement(), and HashTable< T, Key, Hash >::size().

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 383 of file meshRefinementBaffles.C.

References Foam::abort(), polyTopoChange::changeMesh(), Foam::endl(), Foam::FatalError, FatalErrorIn, forAll, HashSet< Key, Hash >::insert(), List< T >::size(), faceSet::sync(), syncTools::syncFaceList(), and HashTable< T, Key, Hash >::toc().

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

Add patch originating from meshing. Update meshedPatches_.

Definition at line 1729 of file meshRefinement.C.

References Foam::findIndex(), and Foam::name().

Foam::labelList meshedPatches ( ) const

Get patchIDs for patches added in addMeshedPatch.

Definition at line 1756 of file meshRefinement.C.

References Foam::abort(), Foam::FatalError, FatalErrorIn, and forAll.

Foam::autoPtr< Foam::mapPolyMesh > splitMeshRegions ( const point keepPoint)
void distribute ( const mapDistributePolyMesh map)

Update local numbering for mesh redistribution.

Definition at line 1855 of file meshRefinement.C.

References autoPtr< T >::clear(), mapDistributePolyMesh::distributeFaceData(), E(), treeBoundBox::extend(), forAll, and autoPtr< T >::valid().

void updateMesh ( const mapPolyMesh map,
const labelList changedFaces 
)

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

face labels.

Definition at line 1913 of file meshRefinement.C.

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

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

Definition at line 1925 of file meshRefinement.C.

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 1942 of file meshRefinement.C.

References mapPolyMesh::faceMap(), forAll, mapPolyMesh::reverseFaceMap(), List< T >::size(), and List< T >::transfer().

Foam::label mergePatchFaces ( const scalar  minCos,
const scalar  concaveCos,
const labelList patchIDs 
)

Merge faces on the same patch (usually from exposing refinement)

Returns global number of faces merged.

Definition at line 35 of file meshRefinementMerge.C.

References polyTopoChange::changeMesh(), polyPatch::coupled(), Foam::endl(), forAll, combineFaces::getMergeSets(), Foam::Info, patches, Foam::returnReduce(), combineFaces::setRefinement(), polyPatch::start(), and combineFaces::updateMesh().

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

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

Definition at line 251 of file meshRefinement.C.

References Foam::endl(), localPointRegion::findDuplicateFaces(), forAll, Foam::identity(), Foam::Pout, syncTools::swapBoundaryFaceList(), and WarningIn.

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

Compare two lists over all boundary faces.

Definition at line 63 of file meshRefinementTemplates.C.

References Foam::abort(), Foam::FatalError, FatalErrorIn, forAll, Foam::mag(), patchIdentifier::name(), patches, UList< T >::size(), polyPatch::start(), and Foam::T().

void printMeshInfo ( const bool  debug,
const string msg 
) const
Foam::word timeName ( ) const

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

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

overwrite_)

Definition at line 2123 of file meshRefinement.C.

Referenced by autoHexMeshDriver::writeMesh().

void setInstance ( const fileName inst)

Set instance of all local IOobjects.

Definition at line 441 of file meshRefinement.C.

void dumpRefinementLevel ( ) const

Write refinement level as volScalarFields for postprocessing.

Definition at line 2136 of file meshRefinement.C.

References IOobject::AUTO_WRITE, Foam::dimless, forAll, MeshObject< polyMesh, pointMesh >::New(), timeName, and regIOobject::write().

void dumpIntersections ( const fileName prefix) const

Debug: Write intersection information to OBJ format.

Definition at line 2192 of file meshRefinement.C.

References Foam::endl(), forAll, OFstream::name(), Foam::nl, Foam::Pout, List< T >::size(), Foam::system(), and Foam::meshTools::writeOBJ().

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 2272 of file meshRefinement.C.


The documentation for this class was generated from the following files: