Go to the documentation of this file.
55 #ifndef backgroundMeshDecomposition_H
56 #define backgroundMeshDecomposition_H
103 const Time& runTime_;
140 scalar minCellSizeLimit_;
151 scalar maxCellWeightCoeff_;
156 void initialRefinement();
168 scalar& weightEstimate
180 void buildPatchAndTree();
195 ClassName(
"backgroundMeshDecomposition");
228 template<
class Po
intType>
246 const scalar radiusSqr
266 template<
class Po
intType>
282 bool includeOwnProcessor =
false
288 const scalar& radiusSqr
294 const scalar radiusSqr
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
A class for handling file names.
Standard boundBox with 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.
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...
An enumeration wrapper for classification of a location as being inside/outside of a volume.
const fvMesh & mesh() const
Return access to the underlying mesh.
~backgroundMeshDecomposition()=default
Destructor.
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?
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
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.
labelList processorPosition(const List< PointType > &pts) const
What processor is the given position on?
labelList overlapProcessors(const point ¢re, const scalar radiusSqr) const
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Container with cells to refine. Refinement given as single direction.
PrimitivePatch<::Foam::List< face >, const pointField > bPatch
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
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,.
autoPtr< mapDistribute > distributePoints(List< PointType > &points) const
Distribute supplied the points to the appropriate processor.
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.
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.