VTK  9.3.20240423
vtkDIYKdTreeUtilities.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
15#ifndef vtkDIYKdTreeUtilities_h
16#define vtkDIYKdTreeUtilities_h
17
18#include "vtkBoundingBox.h" // for vtkBoundingBox
19#include "vtkDIYExplicitAssigner.h" // for vtkDIYExplicitAssigner
20#include "vtkFiltersParallelDIY2Module.h" // for export macros
21#include "vtkObject.h"
22#include "vtkSmartPointer.h" // for vtkSmartPointer
23
24#include <memory> // for std::shared_ptr
25#include <vector> // for std::vector
26
27VTK_ABI_NAMESPACE_BEGIN
28class vtkDataObject;
29class vtkDataSet;
30class vtkIntArray;
33class vtkPoints;
35
36class VTKFILTERSPARALLELDIY2_EXPORT vtkDIYKdTreeUtilities : public vtkObject
37{
38public:
40 void PrintSelf(ostream& os, vtkIndent indent) override;
41
55 static std::vector<vtkBoundingBox> GenerateCuts(vtkDataObject* dobj, int number_of_partitions,
56 bool use_cell_centers, vtkMultiProcessController* controller = nullptr,
57 const double* local_bounds = nullptr);
58
63 static std::vector<vtkBoundingBox> GenerateCuts(const std::vector<vtkDataObject*>& dobjs,
64 int number_of_partitions, bool use_cell_centers,
65 vtkMultiProcessController* controller = nullptr, const double* local_bounds = nullptr);
66
80 static std::vector<vtkBoundingBox> GenerateCuts(
81 const std::vector<vtkSmartPointer<vtkPoints>>& points, int number_of_partitions,
82 vtkMultiProcessController* controller = nullptr, const double* local_bounds = nullptr);
83
104 vtkMultiProcessController* controller, std::shared_ptr<diy::Assigner> block_assigner = nullptr);
105
115 vtkMultiProcessController* controller, vtkIdType* mb_offset = nullptr);
116
125 static std::vector<int> ComputeAssignments(int num_blocks, int num_ranks);
126
132 static vtkDIYExplicitAssigner CreateAssigner(diy::mpi::communicator& comm, int num_blocks);
133
142 static void ResizeCuts(std::vector<vtkBoundingBox>& cuts, int size);
143
144protected:
147
148private:
150 void operator=(const vtkDIYKdTreeUtilities&) = delete;
151};
152
153VTK_ABI_NAMESPACE_END
154#endif
assigner for use with DIY
collection of utility functions for DIY-based KdTree algorithm
static std::vector< vtkBoundingBox > GenerateCuts(const std::vector< vtkSmartPointer< vtkPoints > > &points, int number_of_partitions, vtkMultiProcessController *controller=nullptr, const double *local_bounds=nullptr)
Given a collection of points, this method will generate box cuts in the domain to approximately load ...
static vtkDIYExplicitAssigner CreateAssigner(diy::mpi::communicator &comm, int num_blocks)
Returns an assigner that assigns power-of-two blocks to an arbitrary number of ranks such that each r...
static void ResizeCuts(std::vector< vtkBoundingBox > &cuts, int size)
GenerateCuts returns a kd-tree with power of 2 nodes.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
static vtkSmartPointer< vtkPartitionedDataSet > Exchange(vtkPartitionedDataSet *parts, vtkMultiProcessController *controller, std::shared_ptr< diy::Assigner > block_assigner=nullptr)
Exchange parts in the partitioned dataset among ranks in the parallel group defined by the controller...
static std::vector< int > ComputeAssignments(int num_blocks, int num_ranks)
GenerateCuts returns a kd-tree with power of 2 nodes.
static std::vector< vtkBoundingBox > GenerateCuts(vtkDataObject *dobj, int number_of_partitions, bool use_cell_centers, vtkMultiProcessController *controller=nullptr, const double *local_bounds=nullptr)
Given a dataset (or a composite dataset), this method will generate box cuts in the domain to approxi...
static bool GenerateGlobalCellIds(vtkPartitionedDataSet *parts, vtkMultiProcessController *controller, vtkIdType *mb_offset=nullptr)
Generates and adds global cell ids to datasets in parts.
~vtkDIYKdTreeUtilities() override
static std::vector< vtkBoundingBox > GenerateCuts(const std::vector< vtkDataObject * > &dobjs, int number_of_partitions, bool use_cell_centers, vtkMultiProcessController *controller=nullptr, const double *local_bounds=nullptr)
Another variant to GenerateCuts that simply takes in a vector of dataobjects, each can be a dataset o...
general representation of visualization data
abstract class to specify dataset behavior
Definition vtkDataSet.h:166
a simple class to control print indentation
Definition vtkIndent.h:108
dynamic, self-adjusting array of int
Multiprocessing communication superclass.
abstract base class for most VTK objects
Definition vtkObject.h:162
composite dataset to encapsulates a dataset consisting of partitions.
represent and manipulate 3D points
Definition vtkPoints.h:139
Hold a reference to a vtkObjectBase instance.
dataset represents arbitrary combinations of all possible cell types
int vtkIdType
Definition vtkType.h:315