Go to the documentation of this file.
44 #ifndef cellCellStencils_inverseDistance_H
45 #define cellCellStencils_inverseDistance_H
61 namespace cellCellStencils
134 const unsigned int val
144 const unsigned int val
189 const label destMesh,
190 const label currentDonorMesh,
191 const label newDonorMesh
242 const scalar wantedFraction,
250 const scalar layerRelax,
virtual ~inverseDistance()
Destructor.
virtual const labelUList & cellTypes() const
Return the cell type list.
volScalarField cellInterpolationWeight_
Amount of interpolation.
Inverse-distance-weighted interpolation stencil.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
static point position(const boundBox &bb, const labelVector &nDivs, const label boxI)
Convert index back into coordinate.
void uncompactedRegionSplit(const fvMesh &mesh, const globalIndex &globalFaces, const label nZones, const labelList &zoneID, const labelList &cellTypes, const boolList &isBlockedFace, labelList &cellRegion) const
Replacement of regionSplit.
virtual void stencilWeights(const point &sample, const pointList &donorCcs, scalarList &weights) const
Calculate inverse distance weights for a single acceptor.
Standard boundBox with extra functionality for use in octree.
Buffers for inter-processor communications streams (UOPstream, UIPstream).
List< treeBoundBox > meshBb(1, treeBoundBox(boundBox(coarseMesh.points(), false)).extend(rndGen, 1e-3))
static label index(const labelVector &nDivs, const labelVector &)
Convert ijk indices into single index.
TypeName("inverseDistance")
Runtime type information.
wordList patchTypes(nPatches)
PtrList< fvMeshSubset > meshParts(nZones)
static bool overlaps(const boundBox &bb, const labelVector &nDivs, const PackedList< 2 > &voxels, const treeBoundBox &subBb, const unsigned int val)
Is any voxel inside subBb set to val.
scalarListList cellInterpolationWeights_
Interpolation weights.
const dictionary dict_
Dictionary of motion control parameters.
static treeBoundBox cellBb(const primitiveMesh &mesh, const label celli)
Calculate bounding box of cell.
void walkFront(const scalar layerRelax, const labelListList &allStencil, labelList &allCellTypes, scalarField &allWeight) const
Surround holes with layer(s) of interpolated cells.
const labelIOList & zoneID() const
Helper: get reference to registered zoneID. Loads volScalarField.
Class containing processor-to-processor mapping information.
virtual const scalarList & cellInterpolationWeight() const
Per interpolated cell the interpolation factor. (0 = use.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
static void markBoundaries(const fvMesh &mesh, const vector &smallVec, const boundBox &bb, const labelVector &nDivs, PackedList< 2 > &patchTypes, const labelList &cellMap, labelList &patchCellTypes)
Mark voxels of patchTypes with type of patch face.
static void fill(PackedList< 2 > &elems, const boundBox &bb, const labelVector &nDivs, const boundBox &subBb, const unsigned int val)
Fill all elements overlapping subBb with value val.
virtual const labelListList & cellStencil() const
Per interpolated cell the neighbour cells (in terms of slots as.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
virtual void createStencil(const globalIndex &)
Create stencil starting from the donor containing the acceptor.
Mesh data needed to do the Finite Volume discretisation.
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
virtual const mapDistribute & cellInterpolationMap() const
Return a communication schedule.
void seedCell(const label cellI, const scalar wantedFraction, bitSet &isFront, scalarField &fraction) const
Seed faces of cell with wantedFraction (if higher than current)
virtual bool update()
Update stencils. Return false if nothing changed.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
labelList interpolationCells_
Indices of interpolated cells.
Calculation of interpolation stencils.
labelList cellTypes_
Per cell the cell type.
virtual const labelUList & interpolationCells() const
Indices of interpolated cells.
static labelVector index3(const labelVector &nDivs, const label)
Convert single index into ijk.
virtual const scalarListList & cellInterpolationWeights() const
Weights for cellStencil.
A bounding box defined in terms of min/max extrema points.
vector smallVec_
Small perturbation vector for geometric tests.
void findHoles(const globalIndex &globalCells, const fvMesh &mesh, const labelList &zoneID, const labelListList &stencil, labelList &cellTypes) const
Do flood filling to detect unreachable (from patches) sections.
autoPtr< globalIndex > compactedRegionSplit(const fvMesh &mesh, const globalIndex &globalFaces, labelList &cellRegion) const
void markDonors(const globalIndex &globalCells, PstreamBuffers &pBufs, const PtrList< fvMeshSubset > &meshParts, const List< treeBoundBoxList > &meshBb, const labelList &allCellTypes, const label srcI, const label tgtI, labelListList &allStencil, labelList &allDonor) const
Determine donors for all tgt cells.
autoPtr< mapDistribute > cellInterpolationMap_
Fetch interpolated cells.
bool betterDonor(const label destMesh, const label currentDonorMesh, const label newDonorMesh) const
If multiple donors meshes: decide which is best.
labelListList cellStencil_
Interpolation stencil.
void markPatchesAsHoles(PstreamBuffers &pBufs, const PtrList< fvMeshSubset > &meshParts, const List< treeBoundBoxList > &patchBb, const List< labelVector > &patchDivisions, const PtrList< PackedList< 2 >> &patchParts, const label srcI, const label tgtI, labelList &allCellTypes) const
Mark all cells overlapping (a voxel covered by) a src patch.
Minimal example by using system/controlDict.functions:
Cell-face mesh analysis engine.