VTK  9.3.20240423
vtkPlaneCutter.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
62#ifndef vtkPlaneCutter_h
63#define vtkPlaneCutter_h
64
66#include "vtkFiltersCoreModule.h" // For export macro
67#include "vtkSmartPointer.h" // For SmartPointer
68#include <map> // For std::map
69
70VTK_ABI_NAMESPACE_BEGIN
72class vtkPlane;
73class vtkPolyData;
74class vtkSphereTree;
75
76class VTKFILTERSCORE_EXPORT vtkPlaneCutter : public vtkDataObjectAlgorithm
77{
78public:
80
85 void PrintSelf(ostream& os, vtkIndent indent) override;
87
92
94
99 virtual void SetPlane(vtkPlane*);
100 vtkGetObjectMacro(Plane, vtkPlane);
102
104
110 vtkSetMacro(ComputeNormals, bool);
111 vtkGetMacro(ComputeNormals, bool);
112 vtkBooleanMacro(ComputeNormals, bool);
114
116
121 vtkSetMacro(InterpolateAttributes, bool);
122 vtkGetMacro(InterpolateAttributes, bool);
123 vtkBooleanMacro(InterpolateAttributes, bool);
125
127
132 vtkSetMacro(GeneratePolygons, bool);
133 vtkGetMacro(GeneratePolygons, bool);
134 vtkBooleanMacro(GeneratePolygons, bool);
136
138
144 vtkSetMacro(BuildTree, bool);
145 vtkGetMacro(BuildTree, bool);
146 vtkBooleanMacro(BuildTree, bool);
148
150
156 vtkSetMacro(BuildHierarchy, bool);
157 vtkGetMacro(BuildHierarchy, bool);
158 vtkBooleanMacro(BuildHierarchy, bool);
160
162
170 vtkSetMacro(MergePoints, bool);
171 vtkGetMacro(MergePoints, bool);
172 vtkBooleanMacro(MergePoints, bool);
174
176
181 vtkSetClampMacro(OutputPointsPrecision, int, SINGLE_PRECISION, DEFAULT_PRECISION);
182 vtkGetMacro(OutputPointsPrecision, int);
184
185protected:
187 ~vtkPlaneCutter() override;
188
197
198 // Support delegation to vtkPolyDataPlaneCutter/vtk3DLinearGridPlaneCutter.
200
201 std::map<vtkDataSet*, vtkSmartPointer<vtkSphereTree>> SphereTrees;
202 std::map<vtkDataSet*, bool> CanBeFullyProcessed;
204 {
207
209 : Input(nullptr)
210 , LastMTime(0)
211 {
212 }
214 : Input(input)
215 , LastMTime(mtime)
216 {
217 }
218 };
220
221 // Pipeline-related methods
225 int FillInputPortInformation(int port, vtkInformation* info) override;
226 int FillOutputPortInformation(int port, vtkInformation* info) override;
227
230
231 static void AddNormalArray(double* planeNormal, vtkPolyData* polyData);
232
233private:
234 vtkPlaneCutter(const vtkPlaneCutter&) = delete;
235 void operator=(const vtkPlaneCutter&) = delete;
236};
237
238VTK_ABI_NAMESPACE_END
239#endif
Superclass for algorithms that produce only data object as output.
provides implementation for most abstract methods in the superclass vtkCompositeDataSet
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
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
cut any dataset with a plane and generate a polygonal cut surface
int ExecuteDataSet(vtkDataSet *input, vtkPolyData *output)
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
static void AddNormalArray(double *planeNormal, vtkPolyData *polyData)
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
void PrintSelf(ostream &os, vtkIndent indent) override
Standard construction and print methods.
int RequestUpdateExtent(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
int RequestDataObject(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
vtkInputInfo InputInfo
int ExecuteDataObjectTree(vtkDataObjectTree *input, vtkDataObjectTree *output)
static vtkPlaneCutter * New()
Standard construction and print methods.
std::map< vtkDataSet *, bool > CanBeFullyProcessed
int FillOutputPortInformation(int port, vtkInformation *info) override
Fill the output port information objects for this algorithm.
std::map< vtkDataSet *, vtkSmartPointer< vtkSphereTree > > SphereTrees
~vtkPlaneCutter() override
virtual void SetPlane(vtkPlane *)
Specify the plane (an implicit function) to perform the cutting.
vtkMTimeType GetMTime() override
The modified time depends on the delegated cut plane.
perform various plane computations
Definition vtkPlane.h:135
concrete dataset represents vertices, lines, polygons, and triangle strips
class to build and traverse sphere trees
vtkInputInfo(vtkDataObject *input, vtkMTimeType mtime)
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:270