44#ifndef cellCellStencils_inverseDistance_H
45#define cellCellStencils_inverseDistance_H
70 public cellCellStencil
75 inverseDistance(
const inverseDistance&) =
delete;
78 void operator=(
const inverseDistance&) =
delete;
156 const unsigned int val
166 const unsigned int val
204 const label destMesh,
205 const label currentDonorMesh,
206 const label newDonorMesh
258 const scalar wantedFraction,
269 const bool allowHoleDonors
Minimal example by using system/controlDict.functions:
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
A dynamic list of packed unsigned integers, with the number of bits per item specified by the <Width>...
Buffers for inter-processor communications streams (UOPstream, UIPstream).
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
A bounding box defined in terms of min/max extrema points.
static const labelIOList & zoneID(const fvMesh &)
Helper: get reference to registered zoneID. Loads volScalarField.
Inverse-distance-weighted interpolation stencil.
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.
vector smallVec_
Small perturbation vector for geometric tests.
labelList cellTypes_
Per cell the cell type.
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.
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.
virtual void stencilWeights(const point &sample, const pointList &donorCcs, scalarList &weights) const
Calculate inverse distance weights for a single acceptor.
static point position(const boundBox &bb, const labelVector &nDivs, const label boxI)
Convert index back into coordinate.
static label index(const labelVector &nDivs, const labelVector &)
Convert ijk indices into single index.
virtual ~inverseDistance()
Destructor.
virtual const labelUList & interpolationCells() const
Indices of interpolated cells.
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.
void seedCell(const label cellI, const scalar wantedFraction, bitSet &isFront, scalarField &fraction) const
Seed faces of cell with wantedFraction (if higher than current).
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.
static labelVector index3(const labelVector &nDivs, const label)
Convert single index into ijk.
const bool allowHoleDonors_
Allow holes as donors.
bool betterDonor(const label destMesh, const label currentDonorMesh, const label newDonorMesh) const
If multiple donors meshes: decide which is best.
virtual const mapDistribute & cellInterpolationMap() const
Return a communication schedule.
labelListList cellStencil_
Interpolation stencil.
virtual const scalarList & cellInterpolationWeight() const
Per interpolated cell the interpolation factor. (0 = use.
autoPtr< globalIndex > compactedRegionSplit(const fvMesh &mesh, const globalIndex &globalFaces, labelList &cellRegion) const
virtual const labelUList & cellTypes() const
Return the cell type list.
const bool allowInterpolatedDonors_
Allow interpolared as donors.
autoPtr< mapDistribute > cellInterpolationMap_
Fetch interpolated cells.
virtual bool update()
Update stencils. Return false if nothing changed.
void holeExtrapolationStencil(const globalIndex &globalCells)
Make holes next to live ones type SPECIAL and include in interpolation stencil.
const dictionary dict_
Dictionary of motion control parameters.
TypeName("inverseDistance")
Runtime type information.
volScalarField cellInterpolationWeight_
Amount of interpolation.
virtual const labelListList & cellStencil() const
Per interpolated cell the neighbour cells (in terms of slots as.
scalarListList cellInterpolationWeights_
Interpolation weights.
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.
labelList interpolationCells_
Indices of interpolated cells.
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.
virtual const scalarListList & cellInterpolationWeights() const
Weights for cellStencil.
virtual void createStencil(const globalIndex &, const bool allowHoleDonors)
Create stencil starting from the donor containing the acceptor allowHoleDonors : allow donors of type...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Holds a reference to the original mesh (the baseMesh) and optionally to a subset of that mesh (the su...
Mesh data needed to do the Finite Volume discretisation.
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Class containing processor-to-processor mapping information.
Standard boundBox with extra functionality for use in octree.
const bitSet isBlockedFace(intersectedFaces())
List< scalarList > scalarListList
List of scalarList.
List< labelList > labelListList
List of labelList.
List< label > labelList
A List of labels.
Vector< label > labelVector
Vector of labels.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
List< point > pointList
List of point.
vector point
Point is a vector.
List< bool > boolList
A List of bools.
UList< label > labelUList
A UList of labels.
List< scalar > scalarList
List of scalar.
wordList patchTypes(nPatches)
List< treeBoundBox > meshBb(1, treeBoundBox(coarseMesh.points()).extend(rndGen, 1e-3))
Data for a set of interpolated/donor set.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.