backgroundMeshDecomposition Class Reference

Store a background polyMesh to use for the decomposition of space and queries for parallel conformalVoronoiMesh. More...

Public Member Functions

 ClassName ("backgroundMeshDecomposition")
 Runtime type information. More...
 
 backgroundMeshDecomposition (const Time &runTime, Random &rndGen, const conformationSurfaces &geometryToConformTo, const dictionary &coeffsDict)
 Construct from components in foamyHexMesh operation. More...
 
 backgroundMeshDecomposition (const backgroundMeshDecomposition &)=delete
 Disallow default bitwise copy construction. More...
 
 ~backgroundMeshDecomposition ()
 Destructor. More...
 
autoPtr< polyDistributionMapdistribute (volScalarField &cellWeights)
 Redistribute the background mesh based on a supplied weight field,. More...
 
template<class PointType >
autoPtr< distributionMapdistributePoints (List< PointType > &points) const
 Distribute supplied the points to the appropriate processor. More...
 
bool positionOnThisProcessor (const point &pt) const
 Is the given position inside the domain of this decomposition. More...
 
boolList positionOnThisProcessor (const List< point > &pts) const
 Are the given positions inside the domain of this decomposition. More...
 
bool overlapsThisProcessor (const treeBoundBox &box) const
 Does the given box overlap the faces of the boundary of this. More...
 
bool overlapsThisProcessor (const point &centre, const scalar radiusSqr) const
 Does the given sphere overlap the faces of the boundary of this. More...
 
pointIndexHit findLine (const point &start, const point &end) const
 Find nearest intersection of line between start and end, (exposing. More...
 
pointIndexHit findLineAny (const point &start, const point &end) const
 Find any intersection of line between start and end, (exposing. More...
 
template<class PointType >
labelList processorPosition (const List< PointType > &pts) const
 What processor is the given position on? More...
 
labelList processorNearestPosition (const List< point > &pts) const
 What is the nearest processor to the given position? More...
 
List< List< pointIndexHit > > intersectsProcessors (const List< point > &starts, const List< point > &ends, bool includeOwnProcessor=false) const
 Which processors are intersected by the line segment, returns all. More...
 
bool overlapsOtherProcessors (const point &centre, const scalar &radiusSqr) const
 
labelList overlapProcessors (const point &centre, const scalar radiusSqr) const
 
const fvMeshmesh () const
 Return access to the underlying mesh. More...
 
const indexedOctree< treeDataBPatch > & tree () const
 Return access to the underlying tree. More...
 
const treeBoundBoxprocBounds () const
 Return the boundBox of this processor. More...
 
const labelListcellLevel () const
 Return the cell level of the underlying mesh. More...
 
const labelListpointLevel () const
 Return the point level of the underlying mesh. More...
 
const decompositionMethoddecomposer () const
 Return the current decomposition method. More...
 
void operator= (const backgroundMeshDecomposition &)=delete
 Disallow default bitwise assignment. More...
 

Static Public Member Functions

static autoPtr< distributionMapbuildMap (const List< label > &toProc)
 Build a distributionMap for the supplied destination processor data. More...
 

Detailed Description

Store a background polyMesh to use for the decomposition of space and queries for parallel conformalVoronoiMesh.

The requirements are:

  • To have a decomposition of space which can quickly interrogate an arbitrary location from any processor to reliably and unambiguously determine which processor owns the space that the point is in, i.e. as the vertices move, or need inserted as part of the surface conformation, send them to the correct proc.
  • To be able to be dynamically built, refined and redistributed to other procs the partitioning as the meshing progresses to balance the load.
  • To be able to query whether a sphere (the circumsphere of a Delaunay tet) overlaps any part of the space defined by the structure, and whether a ray (Voronoi edge) penetrates any part of the space defined by the structure, this is what determines if points get referred to a processor.
Source files

Definition at line 92 of file backgroundMeshDecomposition.H.

Constructor & Destructor Documentation

◆ backgroundMeshDecomposition() [1/2]

backgroundMeshDecomposition ( const Time runTime,
Random rndGen,
const conformationSurfaces geometryToConformTo,
const dictionary coeffsDict 
)

Construct from components in foamyHexMesh operation.

◆ backgroundMeshDecomposition() [2/2]

Disallow default bitwise copy construction.

◆ ~backgroundMeshDecomposition()

Destructor.

Member Function Documentation

◆ ClassName()

ClassName ( "backgroundMeshDecomposition"  )

Runtime type information.

◆ buildMap()

static autoPtr<distributionMap> buildMap ( const List< label > &  toProc)
static

Build a distributionMap for the supplied destination processor data.

◆ distribute()

autoPtr<polyDistributionMap> distribute ( volScalarField cellWeights)

Redistribute the background mesh based on a supplied weight field,.

returning a map to use to redistribute vertices.

◆ distributePoints()

autoPtr<distributionMap> distributePoints ( List< PointType > &  points) const

Distribute supplied the points to the appropriate processor.

◆ positionOnThisProcessor() [1/2]

bool positionOnThisProcessor ( const point pt) const

Is the given position inside the domain of this decomposition.

◆ positionOnThisProcessor() [2/2]

boolList positionOnThisProcessor ( const List< point > &  pts) const

Are the given positions inside the domain of this decomposition.

◆ overlapsThisProcessor() [1/2]

bool overlapsThisProcessor ( const treeBoundBox box) const

Does the given box overlap the faces of the boundary of this.

processor

◆ overlapsThisProcessor() [2/2]

bool overlapsThisProcessor ( const point centre,
const scalar  radiusSqr 
) const

Does the given sphere overlap the faces of the boundary of this.

processor

◆ findLine()

pointIndexHit findLine ( const point start,
const point end 
) const

Find nearest intersection of line between start and end, (exposing.

underlying indexedOctree)

◆ findLineAny()

pointIndexHit findLineAny ( const point start,
const point end 
) const

Find any intersection of line between start and end, (exposing.

underlying indexedOctree)

◆ processorPosition()

labelList processorPosition ( const List< PointType > &  pts) const

What processor is the given position on?

◆ processorNearestPosition()

labelList processorNearestPosition ( const List< point > &  pts) const

What is the nearest processor to the given position?

◆ intersectsProcessors()

List<List<pointIndexHit> > intersectsProcessors ( const List< point > &  starts,
const List< point > &  ends,
bool  includeOwnProcessor = false 
) const

Which processors are intersected by the line segment, returns all.

processors whose boundary patch is intersected by the sphere. By default this does not return the processor that the query is launched from, it is assumed that the point is on that processor. The index data member of the pointIndexHit is replaced with the processor index.

◆ overlapsOtherProcessors()

bool overlapsOtherProcessors ( const point centre,
const scalar &  radiusSqr 
) const

◆ overlapProcessors()

labelList overlapProcessors ( const point centre,
const scalar  radiusSqr 
) const

◆ mesh()

const Foam::fvMesh & mesh ( ) const
inline

Return access to the underlying mesh.

Definition at line 28 of file backgroundMeshDecompositionI.H.

◆ tree()

const Foam::indexedOctree< Foam::treeDataBPatch > & tree ( ) const
inline

Return access to the underlying tree.

Definition at line 35 of file backgroundMeshDecompositionI.H.

◆ procBounds()

const Foam::treeBoundBox & procBounds ( ) const
inline

Return the boundBox of this processor.

Definition at line 42 of file backgroundMeshDecompositionI.H.

References UPstream::myProcNo().

Here is the call graph for this function:

◆ cellLevel()

const Foam::labelList & cellLevel ( ) const
inline

Return the cell level of the underlying mesh.

Definition at line 48 of file backgroundMeshDecompositionI.H.

References hexRef8::cellLevel().

Here is the call graph for this function:

◆ pointLevel()

const Foam::labelList & pointLevel ( ) const
inline

Return the point level of the underlying mesh.

Definition at line 54 of file backgroundMeshDecompositionI.H.

References hexRef8::pointLevel().

Here is the call graph for this function:

◆ decomposer()

const Foam::decompositionMethod & decomposer ( ) const
inline

Return the current decomposition method.

Definition at line 61 of file backgroundMeshDecompositionI.H.

◆ operator=()

void operator= ( const backgroundMeshDecomposition )
delete

Disallow default bitwise assignment.


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