VTK  9.2.6
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=========================================================================*/
54
55#ifndef vtkClipClosedSurface_h
56#define vtkClipClosedSurface_h
57
58#include "vtkFiltersGeneralModule.h" // For export macro
60
63class vtkDoubleArray;
64class vtkIdTypeArray;
65class vtkCellArray;
66class vtkPointData;
67class vtkCellData;
68class vtkPolygon;
69class vtkIdList;
70class vtkCCSEdgeLocator;
71
72enum
73{
77};
78
79class VTKFILTERSGENERAL_EXPORT vtkClipClosedSurface : public vtkPolyDataAlgorithm
80{
81public:
83
88 void PrintSelf(ostream& os, vtkIndent indent) override;
90
92
98
100
105 vtkSetMacro(Tolerance, double);
106 vtkGetMacro(Tolerance, double);
108
110
115 vtkBooleanMacro(PassPointData, vtkTypeBool);
118
120
125 vtkBooleanMacro(GenerateOutline, vtkTypeBool);
128
130
135 vtkBooleanMacro(GenerateFaces, vtkTypeBool);
138
140
153 vtkGetMacro(ScalarMode, int);
156
158
164 vtkSetVector3Macro(BaseColor, double);
165 vtkGetVector3Macro(BaseColor, double);
167
169
174 vtkSetVector3Macro(ClipColor, double);
175 vtkGetVector3Macro(ClipColor, double);
177
179
184 vtkSetMacro(ActivePlaneId, int);
185 vtkGetMacro(ActivePlaneId, int);
187
189
194 vtkSetVector3Macro(ActivePlaneColor, double);
195 vtkGetVector3Macro(ActivePlaneColor, double);
197
199
209
210protected:
213
215
216 double Tolerance;
217
223 double BaseColor[3];
224 double ClipColor[3];
226
228
230
232 vtkInformationVector* outputVector, int requestFromOutputPort, vtkMTimeType* mtime) override;
233
235 vtkInformationVector* outputVector) override;
236
240 void ClipLines(vtkPoints* points, vtkDoubleArray* pointScalars, vtkPointData* pointData,
241 vtkCCSEdgeLocator* edgeLocator, vtkCellArray* inputCells, vtkCellArray* outputLines,
242 vtkCellData* inCellData, vtkCellData* outLineData);
243
250 void ClipAndContourPolys(vtkPoints* points, vtkDoubleArray* pointScalars, vtkPointData* pointData,
251 vtkCCSEdgeLocator* edgeLocator, int triangulate, vtkCellArray* inputCells,
252 vtkCellArray* outputPolys, vtkCellArray* outputLines, vtkCellData* inCellData,
253 vtkCellData* outPolyData, vtkCellData* outLineData);
254
261 static int InterpolateEdge(vtkPoints* points, vtkPointData* pointData,
262 vtkCCSEdgeLocator* edgeLocator, double tol, vtkIdType i0, vtkIdType i1, double v0, double v1,
263 vtkIdType& i);
264
270 int TriangulatePolygon(vtkIdList* polygon, vtkPoints* points, vtkCellArray* triangles);
271
281 void TriangulateContours(vtkPolyData* data, vtkIdType firstLine, vtkIdType numLines,
282 vtkCellArray* outputPolys, const double normal[3]);
283
290 static void BreakPolylines(vtkCellArray* inputLines, vtkCellArray* outputLines,
291 vtkUnsignedCharArray* inputScalars, vtkIdType firstLineScalar,
292 vtkUnsignedCharArray* outputScalars, const unsigned char color[3]);
293
299 static void CopyPolygons(vtkCellArray* inputPolys, vtkCellArray* outputPolys,
300 vtkUnsignedCharArray* inputScalars, vtkIdType firstPolyScalar,
301 vtkUnsignedCharArray* outputScalars, const unsigned char color[3]);
302
307 static void BreakTriangleStrips(vtkCellArray* inputStrips, vtkCellArray* outputPolys,
308 vtkUnsignedCharArray* inputScalars, vtkIdType firstStripScalar,
309 vtkUnsignedCharArray* outputScalars, const unsigned char color[3]);
310
317 vtkPolyData* output, vtkPoints* points, vtkPointData* pointData, int outputPointDataType);
318
322 static void CreateColorValues(const double color1[3], const double color2[3],
323 const double color3[3], unsigned char colors[3][3]);
324
325private:
327 void operator=(const vtkClipClosedSurface&) = delete;
328};
329
330#endif
object to represent cell connectivity
represent and manipulate cell attribute data
Definition vtkCellData.h:42
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
Standard methods for instantiation, obtaining type information, and printing.
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()
Standard methods for instantiation, obtaining type information, and printing.
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.
virtual void SetScalarMode(int)
Set whether to add cell scalars, so that new faces and outlines can be distinguished from original fa...
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...
dynamic, self-adjusting array of double
list of point or cell ids
Definition vtkIdList.h:34
dynamic, self-adjusting array of vtkIdType
a simple class to control print indentation
Definition vtkIndent.h:40
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
maintain a list of planes
represent and manipulate point attribute data
represent and manipulate 3D points
Definition vtkPoints.h:40
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition vtkPolyData.h:91
a cell that represents an n-sided polygon
Definition vtkPolygon.h:43
dynamic, self-adjusting array of unsigned char
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