VTK  9.0.1
vtkPlaneCutter.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkPlaneCutter.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 =========================================================================*/
62 #ifndef vtkPlaneCutter_h
63 #define vtkPlaneCutter_h
64 
65 #include "vtkDataSetAlgorithm.h"
66 #include "vtkFiltersCoreModule.h" // For export macro
67 #include "vtkSmartPointer.h" // For SmartPointer
68 #include <vector> // For vector
69 
70 class vtkCellArray;
71 class vtkCellData;
72 class vtkImageData;
74 class vtkPlane;
75 class vtkPointData;
76 class vtkPoints;
77 class vtkSphereTree;
78 class vtkStructuredGrid;
80 
81 class VTKFILTERSCORE_EXPORT vtkPlaneCutter : public vtkDataSetAlgorithm
82 {
83 public:
85 
88  static vtkPlaneCutter* New();
90  void PrintSelf(ostream& os, vtkIndent indent) override;
92 
96  vtkMTimeType GetMTime() override;
97 
99 
104  virtual void SetPlane(vtkPlane*);
105  vtkGetObjectMacro(Plane, vtkPlane);
107 
109 
115  vtkSetMacro(ComputeNormals, bool);
116  vtkGetMacro(ComputeNormals, bool);
117  vtkBooleanMacro(ComputeNormals, bool);
119 
121 
126  vtkSetMacro(InterpolateAttributes, bool);
127  vtkGetMacro(InterpolateAttributes, bool);
128  vtkBooleanMacro(InterpolateAttributes, bool);
130 
132 
137  vtkSetMacro(GeneratePolygons, bool);
138  vtkGetMacro(GeneratePolygons, bool);
139  vtkBooleanMacro(GeneratePolygons, bool);
141 
143 
149  vtkSetMacro(BuildTree, bool);
150  vtkGetMacro(BuildTree, bool);
151  vtkBooleanMacro(BuildTree, bool);
153 
155 
161  vtkSetMacro(BuildHierarchy, bool);
162  vtkGetMacro(BuildHierarchy, bool);
163  vtkBooleanMacro(BuildHierarchy, bool);
165 
171 
172 protected:
174  ~vtkPlaneCutter() override;
175 
180  bool BuildTree;
182 
183  // Helpers
184  std::vector<vtkSmartPointer<vtkSphereTree> > SphereTrees;
185 
186  // Pipeline-related methods
192 
193  virtual int ExecuteDataSet(vtkDataSet* input, vtkSphereTree* tree, vtkMultiPieceDataSet* output);
194 
195  static void AddNormalArray(double* planeNormal, vtkDataSet* ds);
197 
198 private:
199  vtkPlaneCutter(const vtkPlaneCutter&) = delete;
200  void operator=(const vtkPlaneCutter&) = delete;
201 };
202 
203 #endif
object to represent cell connectivity
Definition: vtkCellArray.h:180
represent and manipulate cell attribute data
Definition: vtkCellData.h:33
Superclass for algorithms that produce output of the same type as input.
abstract class to specify dataset behavior
Definition: vtkDataSet.h:57
topologically and geometrically regular array of data
Definition: vtkImageData.h:42
a simple class to control print indentation
Definition: vtkIndent.h:34
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
composite dataset to encapsulates pieces of dataset.
cut any dataset with a plane and generate a polygonal cut surface
static void InitializeOutput(vtkMultiPieceDataSet *output)
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called within ProcessRequest when a request asks the algorithm to do its work.
static vtkPlaneCutter * New()
Standard construction and print methods.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int RequestUpdateExtent(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called within ProcessRequest when each filter in the pipeline decides what portion of its inp...
virtual int ExecuteDataSet(vtkDataSet *input, vtkSphereTree *tree, vtkMultiPieceDataSet *output)
int RequestDataObject(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called within ProcessRequest to when a request asks the algorithm to create empty output data...
bool InterpolateAttributes
static void AddNormalArray(double *planeNormal, vtkDataSet *ds)
int FillOutputPortInformation(int port, vtkInformation *info) override
Fill the output port information objects for this algorithm.
std::vector< vtkSmartPointer< vtkSphereTree > > SphereTrees
~vtkPlaneCutter() override
virtual void SetPlane(vtkPlane *)
Specify the plane (an implicit function) to perform the cutting.
vtkTypeBool ProcessRequest(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
See vtkAlgorithm for details.
vtkMTimeType GetMTime() override
The modified time depends on the delegated cut plane.
vtkPlane * Plane
perform various plane computations
Definition: vtkPlane.h:32
represent and manipulate point attribute data
Definition: vtkPointData.h:32
represent and manipulate 3D points
Definition: vtkPoints.h:34
class to build and traverse sphere trees
Definition: vtkSphereTree.h:70
topologically regular array of data
dataset represents arbitrary combinations of all possible cell types
@ info
Definition: vtkX3D.h:382
@ port
Definition: vtkX3D.h:453
int vtkTypeBool
Definition: vtkABI.h:69
vtkTypeUInt32 vtkMTimeType
Definition: vtkType.h:293