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=========================================================================*/
68#ifndef vtkClipClosedSurface_h
69#define vtkClipClosedSurface_h
70
71#include "vtkFiltersGeneralModule.h" // For export macro
73
76class vtkDoubleArray;
77class vtkIdTypeArray;
78class vtkCellArray;
79class vtkPointData;
80class vtkCellData;
81class vtkPolygon;
82class vtkIdList;
83class vtkCCSEdgeLocator;
84
85enum
86{
90};
91
92class VTKFILTERSGENERAL_EXPORT vtkClipClosedSurface : public vtkPolyDataAlgorithm
93{
94public:
97 void PrintSelf(ostream& os, vtkIndent indent) override;
98
100
104 vtkGetObjectMacro(ClippingPlanes, vtkPlaneCollection);
106
108
113 vtkSetMacro(Tolerance, double);
114 vtkGetMacro(Tolerance, double);
116
118
122 vtkSetMacro(PassPointData, vtkTypeBool);
123 vtkBooleanMacro(PassPointData, vtkTypeBool);
124 vtkGetMacro(PassPointData, vtkTypeBool);
126
128
132 vtkSetMacro(GenerateOutline, vtkTypeBool);
133 vtkBooleanMacro(GenerateOutline, vtkTypeBool);
134 vtkGetMacro(GenerateOutline, vtkTypeBool);
136
138
142 vtkSetMacro(GenerateFaces, vtkTypeBool);
143 vtkBooleanMacro(GenerateFaces, vtkTypeBool);
144 vtkGetMacro(GenerateFaces, vtkTypeBool);
146
148
157 vtkSetClampMacro(ScalarMode, int, VTK_CCS_SCALAR_MODE_NONE, VTK_CCS_SCALAR_MODE_LABELS);
158 void SetScalarModeToNone() { this->SetScalarMode(VTK_CCS_SCALAR_MODE_NONE); }
159 void SetScalarModeToColors() { this->SetScalarMode(VTK_CCS_SCALAR_MODE_COLORS); }
160 void SetScalarModeToLabels() { this->SetScalarMode(VTK_CCS_SCALAR_MODE_LABELS); }
161 vtkGetMacro(ScalarMode, int);
164
166
172 vtkSetVector3Macro(BaseColor, double);
173 vtkGetVector3Macro(BaseColor, double);
175
177
182 vtkSetVector3Macro(ClipColor, double);
183 vtkGetVector3Macro(ClipColor, double);
185
187
192 vtkSetMacro(ActivePlaneId, int);
193 vtkGetMacro(ActivePlaneId, int);
195
197
202 vtkSetVector3Macro(ActivePlaneColor, double);
203 vtkGetVector3Macro(ActivePlaneColor, double);
205
207
213 vtkSetMacro(TriangulationErrorDisplay, vtkTypeBool);
214 vtkBooleanMacro(TriangulationErrorDisplay, vtkTypeBool);
215 vtkGetMacro(TriangulationErrorDisplay, vtkTypeBool);
217
218protected:
221
223
224 double Tolerance;
225
231 double BaseColor[3];
232 double ClipColor[3];
233 double ActivePlaneColor[3];
234
236
238
240 vtkInformationVector* outputVector, int requestFromOutputPort, vtkMTimeType* mtime) override;
241
243 vtkInformationVector* outputVector) override;
244
248 void ClipLines(vtkPoints* points, vtkDoubleArray* pointScalars, vtkPointData* pointData,
249 vtkCCSEdgeLocator* edgeLocator, vtkCellArray* inputCells, vtkCellArray* outputLines,
250 vtkCellData* inCellData, vtkCellData* outLineData);
251
259 vtkCCSEdgeLocator* edgeLocator, int triangulate, vtkCellArray* inputCells,
260 vtkCellArray* outputPolys, vtkCellArray* outputLines, vtkCellData* inCellData,
261 vtkCellData* outPolyData, vtkCellData* outLineData);
262
270 vtkCCSEdgeLocator* edgeLocator, double tol, vtkIdType i0, vtkIdType i1, double v0, double v1,
271 vtkIdType& i);
272
279
290 vtkCellArray* outputPolys, const double normal[3]);
291
298 static void BreakPolylines(vtkCellArray* inputLines, vtkCellArray* outputLines,
299 vtkUnsignedCharArray* inputScalars, vtkIdType firstLineScalar,
300 vtkUnsignedCharArray* outputScalars, const unsigned char color[3]);
301
307 static void CopyPolygons(vtkCellArray* inputPolys, vtkCellArray* outputPolys,
308 vtkUnsignedCharArray* inputScalars, vtkIdType firstPolyScalar,
309 vtkUnsignedCharArray* outputScalars, const unsigned char color[3]);
310
315 static void BreakTriangleStrips(vtkCellArray* inputStrips, vtkCellArray* outputPolys,
316 vtkUnsignedCharArray* inputScalars, vtkIdType firstStripScalar,
317 vtkUnsignedCharArray* outputScalars, const unsigned char color[3]);
318
325 vtkPolyData* output, vtkPoints* points, vtkPointData* pointData, int outputPointDataType);
326
330 static void CreateColorValues(const double color1[3], const double color2[3],
331 const double color3[3], unsigned char colors[3][3]);
332
333private:
335 void operator=(const vtkClipClosedSurface&) = delete;
336};
337
338#endif
object to represent cell connectivity
Definition: vtkCellArray.h:290
represent and manipulate cell attribute data
Definition: vtkCellData.h:142
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:140
dynamic, self-adjusting array of vtkIdType
a simple class to control print indentation
Definition: vtkIndent.h:113
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:142
represent and manipulate 3D points
Definition: vtkPoints.h:143
Superclass for algorithms that produce only polydata as output.
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:195
a cell that represents an n-sided polygon
Definition: vtkPolygon.h:149
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