VTK  9.4.20250420
vtkHyperTreeGridFeatureEdges.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
28#ifndef vtkHyperTreeGridFeatureEdges_h
29#define vtkHyperTreeGridFeatureEdges_h
30
31#include "vtkFiltersHyperTreeModule.h" // For export macro
33#include "vtkMergePoints.h" // For vtkMergePoints
34#include "vtkSmartPointer.h" // For vtkNew
35
36VTK_ABI_NAMESPACE_BEGIN
37
38class vtkCellArray;
39class vtkDataObject;
44class vtkInformation;
45class vtkPoints;
46
47class VTKFILTERSHYPERTREE_EXPORT vtkHyperTreeGridFeatureEdges : public vtkHyperTreeGridAlgorithm
48{
49public:
52 void PrintSelf(ostream& os, vtkIndent indent) override;
53
55
60 vtkSetMacro(MergePoints, bool);
61 vtkGetMacro(MergePoints, bool);
63
64protected:
66 ~vtkHyperTreeGridFeatureEdges() override = default;
67
72
77
82 bool MergePoints = false;
83
84private:
86 void operator=(const vtkHyperTreeGridFeatureEdges&) = delete;
87
89
96 void Process1DHTG(vtkHyperTreeGrid* input, vtkPoints* outPoints, vtkCellArray* outCells);
97 void Process2DHTG(vtkHyperTreeGrid* input, vtkPoints* outPoints, vtkCellArray* outCells);
98 void Process3DHTG(vtkHyperTreeGrid* input, vtkPoints* outPoints, vtkCellArray* outCells);
100
102
107 void RecursivelyProcess1DHTGTree(vtkHyperTreeGrid* input, vtkPoints* outPoints,
109 void RecursivelyProcess2DHTGTree(vtkHyperTreeGrid* input, vtkPoints* outPoints,
111 void RecursivelyProcess3DHTGTree(vtkHyperTreeGrid* input, vtkPoints* outPoints,
114
116
121 vtkSmartPointer<vtkPoints> Build2DCellPoints(
125
136 bool ShouldAddEdge2D(
137 vtkHyperTreeGridNonOrientedVonNeumannSuperCursor* cursor, unsigned int edgeId);
138
149 bool ShouldAddEdge3D(vtkHyperTreeGridNonOrientedMooreSuperCursor* cursor, unsigned int edgeId);
150
157 void InsertNewEdge(double* edgePt1, double* edgePt2, vtkPoints* outPoints, vtkCellArray* outCells,
158 vtkIdType cellId);
159
160 unsigned int OrientationAxe1D = 0;
161 const unsigned int* OrientationAxes2D = nullptr;
162
167};
168
169VTK_ABI_NAMESPACE_END
170#endif /* vtkHyperTreeGridFeatureEdges_h */
object to represent cell connectivity
general representation of visualization data
Superclass for algorithms that produce a hyper tree grid as output.
Generates feature edges from an Hyper Tree Grid.
int FillOutputPortInformation(int, vtkInformation *) override
For this algorithm, the output is a vtkPolyData instance.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
static vtkHyperTreeGridFeatureEdges * New()
~vtkHyperTreeGridFeatureEdges() override=default
int ProcessTrees(vtkHyperTreeGrid *, vtkDataObject *) override
Main routine to generate feature edges.
A dataset containing a grid of vtkHyperTree instances arranged as a rectilinear grid.
a simple class to control print indentation
Definition vtkIndent.h:108
Store vtkAlgorithm input/output information.
represent and manipulate 3D points
Definition vtkPoints.h:139
Hold a reference to a vtkObjectBase instance.
int vtkIdType
Definition vtkType.h:332