VTK  9.1.0
vtkClipClosedSurface.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkClipClosedSurface.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=========================================================================*/
49#ifndef vtkClipClosedSurface_h
50#define vtkClipClosedSurface_h
51
52#include "vtkFiltersGeneralModule.h" // For export macro
54
57class vtkDoubleArray;
58class vtkIdTypeArray;
59class vtkCellArray;
60class vtkPointData;
61class vtkCellData;
62class vtkPolygon;
63class vtkIdList;
64class vtkCCSEdgeLocator;
65
66enum
67{
71};
72
73class VTKFILTERSGENERAL_EXPORT vtkClipClosedSurface : public vtkPolyDataAlgorithm
74{
75public:
78 void PrintSelf(ostream& os, vtkIndent indent) override;
79
81
85 vtkGetObjectMacro(ClippingPlanes, vtkPlaneCollection);
87
89
94 vtkSetMacro(Tolerance, double);
95 vtkGetMacro(Tolerance, double);
97
99
103 vtkSetMacro(PassPointData, vtkTypeBool);
104 vtkBooleanMacro(PassPointData, vtkTypeBool);
105 vtkGetMacro(PassPointData, vtkTypeBool);
107
109
113 vtkSetMacro(GenerateOutline, vtkTypeBool);
114 vtkBooleanMacro(GenerateOutline, vtkTypeBool);
115 vtkGetMacro(GenerateOutline, vtkTypeBool);
117
119
123 vtkSetMacro(GenerateFaces, vtkTypeBool);
124 vtkBooleanMacro(GenerateFaces, vtkTypeBool);
125 vtkGetMacro(GenerateFaces, vtkTypeBool);
127
129
138 vtkSetClampMacro(ScalarMode, int, VTK_CCS_SCALAR_MODE_NONE, VTK_CCS_SCALAR_MODE_LABELS);
139 void SetScalarModeToNone() { this->SetScalarMode(VTK_CCS_SCALAR_MODE_NONE); }
140 void SetScalarModeToColors() { this->SetScalarMode(VTK_CCS_SCALAR_MODE_COLORS); }
141 void SetScalarModeToLabels() { this->SetScalarMode(VTK_CCS_SCALAR_MODE_LABELS); }
142 vtkGetMacro(ScalarMode, int);
145
147
153 vtkSetVector3Macro(BaseColor, double);
154 vtkGetVector3Macro(BaseColor, double);
156
158
163 vtkSetVector3Macro(ClipColor, double);
164 vtkGetVector3Macro(ClipColor, double);
166
168
173 vtkSetMacro(ActivePlaneId, int);
174 vtkGetMacro(ActivePlaneId, int);
176
178
183 vtkSetVector3Macro(ActivePlaneColor, double);
184 vtkGetVector3Macro(ActivePlaneColor, double);
186
188
194 vtkSetMacro(TriangulationErrorDisplay, vtkTypeBool);
195 vtkBooleanMacro(TriangulationErrorDisplay, vtkTypeBool);
196 vtkGetMacro(TriangulationErrorDisplay, vtkTypeBool);
198
199protected:
202
204
205 double Tolerance;
206
212 double BaseColor[3];
213 double ClipColor[3];
214 double ActivePlaneColor[3];
215
217
219
221 vtkInformationVector* outputVector, int requestFromOutputPort, vtkMTimeType* mtime) override;
222
224 vtkInformationVector* outputVector) override;
225
229 void ClipLines(vtkPoints* points, vtkDoubleArray* pointScalars, vtkPointData* pointData,
230 vtkCCSEdgeLocator* edgeLocator, vtkCellArray* inputCells, vtkCellArray* outputLines,
231 vtkCellData* inCellData, vtkCellData* outLineData);
232
240 vtkCCSEdgeLocator* edgeLocator, int triangulate, vtkCellArray* inputCells,
241 vtkCellArray* outputPolys, vtkCellArray* outputLines, vtkCellData* inCellData,
242 vtkCellData* outPolyData, vtkCellData* outLineData);
243
251 vtkCCSEdgeLocator* edgeLocator, double tol, vtkIdType i0, vtkIdType i1, double v0, double v1,
252 vtkIdType& i);
253
260
271 vtkCellArray* outputPolys, const double normal[3]);
272
279 static void BreakPolylines(vtkCellArray* inputLines, vtkCellArray* outputLines,
280 vtkUnsignedCharArray* inputScalars, vtkIdType firstLineScalar,
281 vtkUnsignedCharArray* outputScalars, const unsigned char color[3]);
282
288 static void CopyPolygons(vtkCellArray* inputPolys, vtkCellArray* outputPolys,
289 vtkUnsignedCharArray* inputScalars, vtkIdType firstPolyScalar,
290 vtkUnsignedCharArray* outputScalars, const unsigned char color[3]);
291
296 static void BreakTriangleStrips(vtkCellArray* inputStrips, vtkCellArray* outputPolys,
297 vtkUnsignedCharArray* inputScalars, vtkIdType firstStripScalar,
298 vtkUnsignedCharArray* outputScalars, const unsigned char color[3]);
299
306 vtkPolyData* output, vtkPoints* points, vtkPointData* pointData, int outputPointDataType);
307
311 static void CreateColorValues(const double color1[3], const double color2[3],
312 const double color3[3], unsigned char colors[3][3]);
313
314private:
316 void operator=(const vtkClipClosedSurface&) = delete;
317};
318
319#endif
object to represent cell connectivity
Definition: vtkCellArray.h:181
represent and manipulate cell attribute data
Definition: vtkCellData.h:33
Clip a closed surface with a plane collection.
static void BreakTriangleStrips(vtkCellArray *inputStrips, vtkCellArray *outputPolys, vtkUnsignedCharArray *inputScalars, vtkIdType firstStripScalar, vtkUnsignedCharArray *outputScalars, const unsigned char color[3])
Break triangle strips and add the triangles to the output.
static void CopyPolygons(vtkCellArray *inputPolys, vtkCellArray *outputPolys, vtkUnsignedCharArray *inputScalars, vtkIdType firstPolyScalar, vtkUnsignedCharArray *outputScalars, const unsigned char color[3])
Copy polygons and their associated scalars to a new array.
void ClipLines(vtkPoints *points, vtkDoubleArray *pointScalars, vtkPointData *pointData, vtkCCSEdgeLocator *edgeLocator, vtkCellArray *inputCells, vtkCellArray *outputLines, vtkCellData *inCellData, vtkCellData *outLineData)
Method for clipping lines and copying the scalar data.
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
This is called by the superclass.
static void SqueezeOutputPoints(vtkPolyData *output, vtkPoints *points, vtkPointData *pointData, int outputPointDataType)
Squeeze the points and store them in the output.
vtkPlaneCollection * ClippingPlanes
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void SetScalarModeToLabels()
Set whether to add cell scalars, so that new faces and outlines can be distinguished from original fa...
void SetScalarModeToColors()
Set whether to add cell scalars, so that new faces and outlines can be distinguished from original fa...
virtual void SetClippingPlanes(vtkPlaneCollection *planes)
Set the vtkPlaneCollection that holds the clipping planes.
const char * GetScalarModeAsString()
Set whether to add cell scalars, so that new faces and outlines can be distinguished from original fa...
void TriangulateContours(vtkPolyData *data, vtkIdType firstLine, vtkIdType numLines, vtkCellArray *outputPolys, const double normal[3])
Given some closed contour lines, create a triangle mesh that fills those lines.
static void CreateColorValues(const double color1[3], const double color2[3], const double color3[3], unsigned char colors[3][3])
Take three colors as doubles, and convert to unsigned char.
void SetScalarModeToNone()
Set whether to add cell scalars, so that new faces and outlines can be distinguished from original fa...
int TriangulatePolygon(vtkIdList *polygon, vtkPoints *points, vtkCellArray *triangles)
A robust method for triangulating a polygon.
int ComputePipelineMTime(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector, int requestFromOutputPort, vtkMTimeType *mtime) override
A special version of ProcessRequest meant specifically for the pipeline modified time request.
~vtkClipClosedSurface() override
static vtkClipClosedSurface * New()
static int InterpolateEdge(vtkPoints *points, vtkPointData *pointData, vtkCCSEdgeLocator *edgeLocator, double tol, vtkIdType i0, vtkIdType i1, double v0, double v1, vtkIdType &i)
A helper function for interpolating a new point along an edge.
static void BreakPolylines(vtkCellArray *inputLines, vtkCellArray *outputLines, vtkUnsignedCharArray *inputScalars, vtkIdType firstLineScalar, vtkUnsignedCharArray *outputScalars, const unsigned char color[3])
Break polylines into individual lines, copying scalar values from inputScalars starting at firstLineS...
void ClipAndContourPolys(vtkPoints *points, vtkDoubleArray *pointScalars, vtkPointData *pointData, vtkCCSEdgeLocator *edgeLocator, int triangulate, vtkCellArray *inputCells, vtkCellArray *outputPolys, vtkCellArray *outputLines, vtkCellData *inCellData, vtkCellData *outPolyData, vtkCellData *outLineData)
Clip and contour polys in one step, in order to guarantee that the contour lines exactly match the ne...
vtkTypeBool TriangulationErrorDisplay
dynamic, self-adjusting array of double
list of point or cell ids
Definition: vtkIdList.h:31
dynamic, self-adjusting array of vtkIdType
a simple class to control print indentation
Definition: vtkIndent.h:34
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
maintain a list of planes
represent and manipulate point attribute data
Definition: vtkPointData.h:33
represent and manipulate 3D points
Definition: vtkPoints.h:34
Superclass for algorithms that produce only polydata as output.
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:86
a cell that represents an n-sided polygon
Definition: vtkPolygon.h:40
dynamic, self-adjusting array of unsigned char
@ points
Definition: vtkX3D.h:452
@ color
Definition: vtkX3D.h:227
@ data
Definition: vtkX3D.h:321
int vtkTypeBool
Definition: vtkABI.h:69
@ VTK_CCS_SCALAR_MODE_NONE
@ VTK_CCS_SCALAR_MODE_LABELS
@ VTK_CCS_SCALAR_MODE_COLORS
int vtkIdType
Definition: vtkType.h:332
vtkTypeUInt32 vtkMTimeType
Definition: vtkType.h:287