Go to the documentation of this file.
53 #ifndef backgroundMeshDecomposition_H
54 #define backgroundMeshDecomposition_H
68 #include "primitivePatch.H"
166 scalar& weightEstimate
190 ClassName(
"backgroundMeshDecomposition");
223 template<
typename Po
intType>
241 const scalar radiusSqr
261 template<
typename Po
intType>
277 bool includeOwnProcessor =
false
283 const scalar& radiusSqr
289 const scalar radiusSqr
Simple random number generator.
hexRef8 meshCutter_
Refinement object.
labelList selectRefinementCells(List< volumeType > &volumeStatus, volScalarField &cellWeights) const
Select cells for refinement at the surface of the geometry to be.
Random & rndGen_
Random number generator.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
A class for handling file names.
~backgroundMeshDecomposition()
Destructor.
void operator=(const backgroundMeshDecomposition &)
Disallow default bitwise assignment.
label minLevels_
Minimum normal level of refinement.
autoPtr< indexedOctree< treeDataBPatch > > bFTreePtr_
Search tree for the boundaryFaces_ patch.
label volRes_
How fine should the initial sample of the volume a box be to.
Standard boundBox + extra functionality for use in octree.
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.
const Time & runTime_
Method details dictionary.
autoPtr< bPatch > boundaryFacesPtr_
Patch containing an independent representation of the surface of the.
backgroundMeshDecomposition(const backgroundMeshDecomposition &)
Disallow default bitwise copy construct.
void buildPatchAndTree()
Build the surface patch and search tree.
scalar spanScale_
Scale of a cell span vs cell size used to decide to refine a cell.
fvMesh mesh_
Mesh stored on for this processor, specifiying the domain that it.
const conformationSurfaces & geometryToConformTo_
Reference to surface.
Mesh consisting of general polyhedral cells.
pointIndexHit findLineAny(const point &start, const point &end) const
Find any intersection of line between start and end, (exposing.
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const fvMesh & mesh() const
Return access to the underlying mesh.
scalar mergeDist_
merge distance required by fvMeshDistribute
bool refineCell(label cellI, volumeType volType, scalar &weightEstimate) const
Estimate the number of vertices that will be in this cell, returns.
Non-pointer based hierarchical recursive searching.
bool overlapsOtherProcessors(const point ¢re, const scalar &radiusSqr) const
ClassName("backgroundMeshDecomposition")
Runtime type information.
labelList processorNearestPosition(const List< point > &pts) const
What is the nearest processor to the given position?
void printMeshData(const polyMesh &mesh) const
Print details of the decomposed mesh.
A list of keyword definitions, which are a keyword followed by any number of values (e....
static autoPtr< mapDistribute > buildMap(const List< label > &toProc)
Build a mapDistribute for the supplied destination processor data.
Mesh data needed to do the Finite Volume discretisation.
treeBoundBox globalBackgroundBounds_
The overall bounds of all of the background meshes, used to test if.
labelList processorPosition(const List< PointType > &pts) const
What processor is the given position on?
labelList overlapProcessors(const point ¢re, const scalar radiusSqr) const
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Refinement of (split) hexes using polyTopoChange.
pointIndexHit findLine(const point &start, const point &end) const
Find nearest intersection of line between start and end, (exposing.
treeDataPrimitivePatch< bPatch > treeDataBPatch
PrimitivePatch< face, List, const pointField > bPatch
const treeBoundBox & procBounds() const
Return the boundBox of this processor.
Encapsulation of data needed to search on PrimitivePatches.
const indexedOctree< treeDataBPatch > & tree() const
Return access to the underlying tree.
Store a background polyMesh to use for the decomposition of space and queries for parallel conformalV...
const labelList & pointLevel() const
Return the point level of the underlying mesh.
autoPtr< mapDistributePolyMesh > distribute(volScalarField &cellWeights)
Redistribute the background mesh based on a supplied weight field,.
scalar minCellSizeLimit_
Smallest minimum cell size allowed, i.e. to avoid high initial.
autoPtr< mapDistribute > distributePoints(List< PointType > &points) const
Distribute supplied the points to the appropriate processor.
Generic GeometricField class.
const labelList & cellLevel() const
Return the cell level of the underlying mesh.
CGAL data structures used for 3D Delaunay meshing.
bool overlapsThisProcessor(const treeBoundBox &box) const
Does the given box overlap the faces of the boundary of this.
treeBoundBoxList allBackgroundMeshBounds_
The bounds of all background meshes on all processors.
cachedRandom rndGen(label(0), -1)
scalar maxCellWeightCoeff_
Allowed factor above the average cell weight before a background.
A list of faces which address into the list of points.
bool positionOnThisProcessor(const point &pt) const
Is the given position inside the domain of this decomposition.