VTK  9.1.0
vtkDIYKdTreeUtilities.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkDIYKdTreeUtilities.h
5
6 Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7 All rights reserved.
8 See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9
10 This software is distributed WITHOUT ANY WARRANTY; without even
11 the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12 PURPOSE. See the above copyright notice for more information.
13
14=========================================================================*/
24#ifndef vtkDIYKdTreeUtilities_h
25#define vtkDIYKdTreeUtilities_h
26
27#include "vtkBoundingBox.h" // for vtkBoundingBox
28#include "vtkDIYExplicitAssigner.h" // for vtkDIYExplicitAssigner
29#include "vtkFiltersParallelDIY2Module.h" // for export macros
30#include "vtkObject.h"
31#include "vtkSmartPointer.h" // for vtkSmartPointer
32
33#include <memory> // for std::shared_ptr
34#include <vector> // for std::vector
35
36class vtkDataObject;
37class vtkDataSet;
38class vtkIntArray;
41class vtkPoints;
43
44class VTKFILTERSPARALLELDIY2_EXPORT vtkDIYKdTreeUtilities : public vtkObject
45{
46public:
48 void PrintSelf(ostream& os, vtkIndent indent) override;
49
63 static std::vector<vtkBoundingBox> GenerateCuts(vtkDataObject* dobj, int number_of_partitions,
64 bool use_cell_centers, vtkMultiProcessController* controller = nullptr,
65 const double* local_bounds = nullptr);
66
71 static std::vector<vtkBoundingBox> GenerateCuts(const std::vector<vtkDataObject*>& dobjs,
72 int number_of_partitions, bool use_cell_centers,
73 vtkMultiProcessController* controller = nullptr, const double* local_bounds = nullptr);
74
88 static std::vector<vtkBoundingBox> GenerateCuts(
89 const std::vector<vtkSmartPointer<vtkPoints>>& points, int number_of_partitions,
90 vtkMultiProcessController* controller = nullptr, const double* local_bounds = nullptr);
91
112 vtkMultiProcessController* controller, std::shared_ptr<diy::Assigner> block_assigner = nullptr);
113
123 vtkMultiProcessController* controller, vtkIdType* mb_offset = nullptr);
124
133 static std::vector<int> ComputeAssignments(int num_blocks, int num_ranks);
134
140 static vtkDIYExplicitAssigner CreateAssigner(diy::mpi::communicator& comm, int num_blocks);
141
150 static void ResizeCuts(std::vector<vtkBoundingBox>& cuts, int size);
151
152protected:
155
156private:
158 void operator=(const vtkDIYKdTreeUtilities&) = delete;
159};
160
161#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:113
dynamic, self-adjusting array of int
Definition: vtkIntArray.h:149
Multiprocessing communication superclass.
abstract base class for most VTK objects
Definition: vtkObject.h:82
composite dataset to encapsulates a dataset consisting of partitions.
represent and manipulate 3D points
Definition: vtkPoints.h:143
Hold a reference to a vtkObjectBase instance.
dataset represents arbitrary combinations of all possible cell types
int vtkIdType
Definition: vtkType.h:332