1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377
|
// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
// SPDX-FileCopyrightText: Copyright (c) Sandia Corporation
// SPDX-License-Identifier: BSD-3-Clause
/**
* @class vtkPDistributedDataFilter
* @brief Distribute data among processors
*
*
* This filter redistributes data among processors in a parallel
* application into spatially contiguous vtkUnstructuredGrids.
* The execution model anticipated is that all processes read in
* part of a large vtkDataSet. Each process sets the input of
* filter to be that DataSet. When executed, this filter builds
* in parallel a k-d tree, decomposing the space occupied by the
* distributed DataSet into spatial regions. It assigns each
* spatial region to a processor. The data is then redistributed
* and the output is a single vtkUnstructuredGrid containing the
* cells in the process' assigned regions.
*
* This filter is sometimes called "D3" for "distributed data decomposition".
*
* Enhancement: You can set the k-d tree decomposition, rather than
* have D3 compute it. This allows you to divide a dataset using
* the decomposition computed for another dataset. Obtain a description
* of the k-d tree cuts this way:
*
* @code{cpp}
* vtkBSPCuts *cuts = D3Object1->GetCuts()
* @endcode
*
* And set it this way:
*
* @code{cpp}
* D3Object2->SetCuts(cuts)
* @endcode
*
*
* It is desirable to have a field array of global node IDs
* for two reasons:
*
* 1. When merging together sub grids that were distributed
* across processors, global node IDs can be used to remove
* duplicate points and significantly reduce the size of the
* resulting output grid. If no such array is available,
* D3 will use a tolerance to merge points, which is much
* slower.
*
* 2. If ghost cells have been requested, D3 requires a
* global node ID array in order to request and transfer
* ghost cells in parallel among the processors. If there
* is no global node ID array, D3 will in parallel create
* a global node ID array, and the time to do this can be
* significant.
*
* D3 uses `vtkPointData::GetGlobalIds` to access global
* node ids from the input. If none is found,
* and ghost cells have been requested, D3 will create a
* temporary global node ID array before acquiring ghost cells.
*
* It is also desirable to have global element IDs (vtkCellData::GetGlobalIds).
* However, if they don't exist D3 can create them relatively quickly.
*
* @warning
* The Execute() method must be called by all processes in the
* parallel application, or it will hang. If you are not certain
* that your pipeline will execute identically on all processors,
* you may want to use this filter in an explicit execution mode.
*
* @sa
* vtkKdTree vtkPKdTree vtkBSPCuts
*/
#ifndef vtkPDistributedDataFilter_h
#define vtkPDistributedDataFilter_h
#include "vtkDistributedDataFilter.h"
#include "vtkFiltersParallelGeometryModule.h" // For export macro
VTK_ABI_NAMESPACE_BEGIN
class vtkBSPCuts;
class vtkDataArray;
class vtkFloatArray;
class vtkIdList;
class vtkIdTypeArray;
class vtkIntArray;
class vtkMultiProcessController;
class vtkPDistributedDataFilterSTLCloak;
class vtkPKdTree;
class vtkUnstructuredGrid;
class VTKFILTERSPARALLELGEOMETRY_EXPORT vtkPDistributedDataFilter : public vtkDistributedDataFilter
{
public:
vtkTypeMacro(vtkPDistributedDataFilter, vtkDistributedDataFilter);
void PrintSelf(ostream& os, vtkIndent indent) override;
static vtkPDistributedDataFilter* New();
protected:
vtkPDistributedDataFilter();
~vtkPDistributedDataFilter() override;
/**
* Build a vtkUnstructuredGrid for a spatial region from the
* data distributed across processes. Execute() must be called
* by all processes, or it will hang.
*/
int RequestData(vtkInformation*, vtkInformationVector**, vtkInformationVector*) override;
void SingleProcessExecute(vtkDataSet* input, vtkUnstructuredGrid* output);
/**
* Implementation for request data.
*/
int RequestDataInternal(vtkDataSet* input, vtkUnstructuredGrid* output);
private:
enum
{
DeleteNo = 0,
DeleteYes = 1
};
enum
{
DuplicateCellsNo = 0,
DuplicateCellsYes = 1
};
enum
{
GhostCellsNo = 0,
GhostCellsYes = 1
};
enum
{
UnsetGhostLevel = 99
};
/**
* ?
*/
int PartitionDataAndAssignToProcesses(vtkDataSet* set);
/**
* ?
*/
vtkUnstructuredGrid* RedistributeDataSet(
vtkDataSet* set, vtkDataSet* input, int filterOutDuplicateCells);
/**
* ?
*/
int ClipGridCells(vtkUnstructuredGrid* grid);
/**
* ?
*/
vtkUnstructuredGrid* AcquireGhostCells(vtkUnstructuredGrid* grid);
/**
* ?
*/
void ComputeMyRegionBounds();
/**
* ?
*/
int CheckFieldArrayTypes(vtkDataSet* set);
/**
* If any processes have 0 cell input data sets, then
* spread the input data sets around (quickly) before formal
* redistribution.
*/
vtkDataSet* TestFixTooFewInputFiles(vtkDataSet* input, int& duplicateCells);
/**
* ?
*/
vtkUnstructuredGrid* MPIRedistribute(
vtkDataSet* in, vtkDataSet* input, int filterOutDuplicateCells);
/**
* ?
*/
vtkIdList** GetCellIdsForProcess(int proc, int* nlists);
/**
* Fills in the Source and Target arrays which contain a schedule to allow
* each processor to talk to every other.
*/
void SetUpPairWiseExchange();
///@{
/**
* ?
*/
void FreeIntArrays(vtkIdTypeArray** ar);
static void FreeIdLists(vtkIdList** lists, int nlists);
static vtkIdType GetIdListSize(vtkIdList** lists, int nlists);
///@}
///@{
/**
* This transfers counts (array sizes) between processes.
*/
vtkIdTypeArray* ExchangeCounts(vtkIdType myCount, int tag);
vtkIdTypeArray* ExchangeCountsLean(vtkIdType myCount, int tag);
vtkIdTypeArray* ExchangeCountsFast(vtkIdType myCount, int tag);
///@}
///@{
/**
* This transfers id valued data arrays between processes.
*/
vtkIdTypeArray** ExchangeIdArrays(vtkIdTypeArray** arIn, int deleteSendArrays, int tag);
vtkIdTypeArray** ExchangeIdArraysLean(vtkIdTypeArray** arIn, int deleteSendArrays, int tag);
vtkIdTypeArray** ExchangeIdArraysFast(vtkIdTypeArray** arIn, int deleteSendArrays, int tag);
///@}
///@{
/**
* This transfers float valued data arrays between processes.
*/
vtkFloatArray** ExchangeFloatArrays(vtkFloatArray** myArray, int deleteSendArrays, int tag);
vtkFloatArray** ExchangeFloatArraysLean(vtkFloatArray** myArray, int deleteSendArrays, int tag);
vtkFloatArray** ExchangeFloatArraysFast(vtkFloatArray** myArray, int deleteSendArrays, int tag);
///@}
///@{
/**
* ?
*/
vtkUnstructuredGrid* ExchangeMergeSubGrids(vtkIdList** cellIds, int deleteCellIds,
vtkDataSet* myGrid, int deleteMyGrid, int filterOutDuplicateCells, int ghostCellFlag, int tag);
vtkUnstructuredGrid* ExchangeMergeSubGrids(vtkIdList*** cellIds, int* numLists, int deleteCellIds,
vtkDataSet* myGrid, int deleteMyGrid, int filterOutDuplicateCells, int ghostCellFlag, int tag);
vtkUnstructuredGrid* ExchangeMergeSubGridsLean(vtkIdList*** cellIds, int* numLists,
int deleteCellIds, vtkDataSet* myGrid, int deleteMyGrid, int filterOutDuplicateCells,
int ghostCellFlag, int tag);
vtkUnstructuredGrid* ExchangeMergeSubGridsFast(vtkIdList*** cellIds, int* numLists,
int deleteCellIds, vtkDataSet* myGrid, int deleteMyGrid, int filterOutDuplicateCells,
int ghostCellFlag, int tag);
///@}
///@{
/**
* ?
*/
char* MarshallDataSet(vtkUnstructuredGrid* extractedGrid, vtkIdType& size);
vtkUnstructuredGrid* UnMarshallDataSet(char* buf, vtkIdType size);
///@}
///@{
/**
* ?
*/
void ClipCellsToSpatialRegion(vtkUnstructuredGrid* grid);
#if 0
void ClipWithVtkClipDataSet(vtkUnstructuredGrid *grid, double *bounds,
vtkUnstructuredGrid **outside, vtkUnstructuredGrid **inside);
#endif
///@}
void ClipWithBoxClipDataSet(vtkUnstructuredGrid* grid, double* bounds,
vtkUnstructuredGrid** outside, vtkUnstructuredGrid** inside);
///@{
/**
* Accessors to the "GLOBALID" point and cell arrays of the dataset.
* Global ids are used by D3 to uniquely name all points and cells
* so that after shuffling data between processors, redundant information
* can be quickly eliminated.
*/
vtkIdTypeArray* GetGlobalNodeIdArray(vtkDataSet* set);
vtkIdType* GetGlobalNodeIds(vtkDataSet* set);
vtkIdTypeArray* GetGlobalElementIdArray(vtkDataSet* set);
vtkIdType* GetGlobalElementIds(vtkDataSet* set);
int AssignGlobalNodeIds(vtkUnstructuredGrid* grid);
int AssignGlobalElementIds(vtkDataSet* in);
vtkIdTypeArray** FindGlobalPointIds(vtkFloatArray** ptarray, vtkIdTypeArray* ids,
vtkUnstructuredGrid* grid, vtkIdType& numUniqueMissingPoints);
///@}
/**
* ?
*/
vtkIdTypeArray** MakeProcessLists(
vtkIdTypeArray** pointIds, vtkPDistributedDataFilterSTLCloak* procs);
/**
* ?
*/
vtkIdList** BuildRequestedGrids(vtkIdTypeArray** globalPtIds, vtkUnstructuredGrid* grid,
vtkPDistributedDataFilterSTLCloak* ptIdMap);
///@{
/**
* ?
*/
int InMySpatialRegion(float x, float y, float z);
int InMySpatialRegion(double x, double y, double z);
int StrictlyInsideMyBounds(float x, float y, float z);
int StrictlyInsideMyBounds(double x, double y, double z);
///@}
///@{
/**
* ?
*/
vtkIdTypeArray** GetGhostPointIds(
int ghostLevel, vtkUnstructuredGrid* grid, int AddCellsIAlreadyHave);
vtkUnstructuredGrid* AddGhostCellsUniqueCellAssignment(
vtkUnstructuredGrid* myGrid, vtkPDistributedDataFilterSTLCloak* globalToLocalMap);
vtkUnstructuredGrid* AddGhostCellsDuplicateCellAssignment(
vtkUnstructuredGrid* myGrid, vtkPDistributedDataFilterSTLCloak* globalToLocalMap);
vtkUnstructuredGrid* SetMergeGhostGrid(vtkUnstructuredGrid* ghostCellGrid,
vtkUnstructuredGrid* incomingGhostCells, int ghostLevel,
vtkPDistributedDataFilterSTLCloak* idMap);
///@}
///@{
/**
* ?
*/
vtkUnstructuredGrid* ExtractCells(vtkIdList* list, int deleteCellLists, vtkDataSet* in);
vtkUnstructuredGrid* ExtractCells(
vtkIdList** lists, int nlists, int deleteCellLists, vtkDataSet* in);
vtkUnstructuredGrid* ExtractZeroCellGrid(vtkDataSet* in);
///@}
///@{
/**
* ?
*/
static int GlobalPointIdIsUsed(
vtkUnstructuredGrid* grid, int ptId, vtkPDistributedDataFilterSTLCloak* globalToLocal);
static int LocalPointIdIsUsed(vtkUnstructuredGrid* grid, int ptId);
static vtkIdType FindId(vtkIdTypeArray* ids, vtkIdType gid, vtkIdType startLoc);
///@}
/**
* ?
*/
static vtkIdTypeArray* AddPointAndCells(vtkIdType gid, vtkIdType localId,
vtkUnstructuredGrid* grid, vtkIdType* gidCells, vtkIdTypeArray* ids);
///@{
/**
* ?
*/
static void AddConstantUnsignedCharPointArray(
vtkUnstructuredGrid* grid, const char* arrayName, unsigned char val);
static void AddConstantUnsignedCharCellArray(
vtkUnstructuredGrid* grid, const char* arrayName, unsigned char val);
///@}
/**
* ?
*/
static void RemoveRemoteCellsFromList(
vtkIdList* cellList, vtkIdType* gidCells, vtkIdType* remoteCells, vtkIdType nRemoteCells);
/**
* ?
*/
static vtkUnstructuredGrid* MergeGrids(vtkDataSet** sets, int nsets, int deleteDataSets,
int useGlobalNodeIds, float pointMergeTolerance, int useGlobalCellIds);
vtkPDistributedDataFilter(const vtkPDistributedDataFilter&) = delete;
void operator=(const vtkPDistributedDataFilter&) = delete;
};
VTK_ABI_NAMESPACE_END
#endif
|