VTK  9.3.0
vtkClipClosedSurface.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
51#ifndef vtkClipClosedSurface_h
52#define vtkClipClosedSurface_h
53
54#include "vtkFiltersGeneralModule.h" // For export macro
56
57VTK_ABI_NAMESPACE_BEGIN
60class vtkDoubleArray;
61class vtkIdTypeArray;
62class vtkCellArray;
63class vtkPointData;
64class vtkCellData;
65class vtkPolygon;
66class vtkIdList;
67class vtkCCSEdgeLocator;
68
69enum
70{
74};
75
76class VTKFILTERSGENERAL_EXPORT vtkClipClosedSurface : public vtkPolyDataAlgorithm
77{
78public:
80
85 void PrintSelf(ostream& os, vtkIndent indent) override;
87
89
93 vtkGetObjectMacro(ClippingPlanes, vtkPlaneCollection);
95
97
102 vtkSetMacro(Tolerance, double);
103 vtkGetMacro(Tolerance, double);
105
107
111 vtkSetMacro(PassPointData, vtkTypeBool);
112 vtkBooleanMacro(PassPointData, vtkTypeBool);
113 vtkGetMacro(PassPointData, vtkTypeBool);
115
117
121 vtkSetMacro(GenerateOutline, vtkTypeBool);
122 vtkBooleanMacro(GenerateOutline, vtkTypeBool);
123 vtkGetMacro(GenerateOutline, vtkTypeBool);
125
127
131 vtkSetMacro(GenerateFaces, vtkTypeBool);
132 vtkBooleanMacro(GenerateFaces, vtkTypeBool);
133 vtkGetMacro(GenerateFaces, vtkTypeBool);
135
137
146 vtkSetClampMacro(ScalarMode, int, VTK_CCS_SCALAR_MODE_NONE, VTK_CCS_SCALAR_MODE_LABELS);
147 void SetScalarModeToNone() { this->SetScalarMode(VTK_CCS_SCALAR_MODE_NONE); }
148 void SetScalarModeToColors() { this->SetScalarMode(VTK_CCS_SCALAR_MODE_COLORS); }
149 void SetScalarModeToLabels() { this->SetScalarMode(VTK_CCS_SCALAR_MODE_LABELS); }
150 vtkGetMacro(ScalarMode, int);
153
155
161 vtkSetVector3Macro(BaseColor, double);
162 vtkGetVector3Macro(BaseColor, double);
164
166
171 vtkSetVector3Macro(ClipColor, double);
172 vtkGetVector3Macro(ClipColor, double);
174
176
181 vtkSetMacro(ActivePlaneId, int);
182 vtkGetMacro(ActivePlaneId, int);
184
186
191 vtkSetVector3Macro(ActivePlaneColor, double);
192 vtkGetVector3Macro(ActivePlaneColor, double);
194
196
202 vtkSetMacro(TriangulationErrorDisplay, vtkTypeBool);
203 vtkBooleanMacro(TriangulationErrorDisplay, vtkTypeBool);
204 vtkGetMacro(TriangulationErrorDisplay, vtkTypeBool);
206
208
217 vtkSetMacro(InsideOut, vtkTypeBool);
218 vtkGetMacro(InsideOut, vtkTypeBool);
219 vtkBooleanMacro(InsideOut, vtkTypeBool);
221
223
228 vtkSetMacro(GenerateClipFaceOutput, vtkTypeBool);
229 vtkGetMacro(GenerateClipFaceOutput, vtkTypeBool);
230 vtkBooleanMacro(GenerateClipFaceOutput, vtkTypeBool);
232
237
238protected:
241
243
244 double Tolerance;
245
251 double BaseColor[3];
252 double ClipColor[3];
253 double ActivePlaneColor[3];
254 vtkTypeBool InsideOut = false;
255 vtkTypeBool GenerateClipFaceOutput = false;
256
258
260
262 vtkInformationVector* outputVector, int requestFromOutputPort, vtkMTimeType* mtime) override;
263
265 vtkInformationVector* outputVector) override;
266
270 void ClipLines(vtkPoints* points, vtkDoubleArray* pointScalars, vtkPointData* pointData,
271 vtkCCSEdgeLocator* edgeLocator, vtkCellArray* inputCells, vtkCellArray* outputLines,
272 vtkCellData* inCellData, vtkCellData* outLineData);
273
280 void ClipAndContourPolys(vtkPoints* points, vtkDoubleArray* pointScalars, vtkPointData* pointData,
281 vtkCCSEdgeLocator* edgeLocator, int triangulate, vtkCellArray* inputCells,
282 vtkCellArray* outputPolys, vtkCellArray* outputLines, vtkCellData* inCellData,
283 vtkCellData* outPolyData, vtkCellData* outLineData);
284
291 static int InterpolateEdge(vtkPoints* points, vtkPointData* pointData,
292 vtkCCSEdgeLocator* edgeLocator, double tol, vtkIdType i0, vtkIdType i1, double v0, double v1,
293 vtkIdType& i);
294
300 int TriangulatePolygon(vtkIdList* polygon, vtkPoints* points, vtkCellArray* triangles);
301
311 void TriangulateContours(vtkPolyData* data, vtkIdType firstLine, vtkIdType numLines,
312 vtkCellArray* outputPolys, const double normal[3]);
313
320 static void BreakPolylines(vtkCellArray* inputLines, vtkCellArray* outputLines,
321 vtkUnsignedCharArray* inputScalars, vtkIdType firstLineScalar,
322 vtkUnsignedCharArray* outputScalars, const unsigned char color[3]);
323
329 static void CopyPolygons(vtkCellArray* inputPolys, vtkCellArray* outputPolys,
330 vtkUnsignedCharArray* inputScalars, vtkIdType firstPolyScalar,
331 vtkUnsignedCharArray* outputScalars, const unsigned char color[3]);
332
337 static void BreakTriangleStrips(vtkCellArray* inputStrips, vtkCellArray* outputPolys,
338 vtkUnsignedCharArray* inputScalars, vtkIdType firstStripScalar,
339 vtkUnsignedCharArray* outputScalars, const unsigned char color[3]);
340
347 vtkPolyData* output, vtkPoints* points, vtkPointData* pointData, int outputPointDataType);
348
352 static void CreateColorValues(const double color1[3], const double color2[3],
353 const double color3[3], unsigned char colors[3][3]);
354
355private:
357 void operator=(const vtkClipClosedSurface&) = delete;
358};
359
360VTK_ABI_NAMESPACE_END
361#endif
object to represent cell connectivity
represent and manipulate cell attribute data
Definition vtkCellData.h:31
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
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.
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...
vtkPolyData * GetClipFaceOutput()
Return the clip face triangulated output.
dynamic, self-adjusting array of double
list of point or cell ids
Definition vtkIdList.h:23
dynamic, self-adjusting array of vtkIdType
a simple class to control print indentation
Definition vtkIndent.h:29
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:29
Superclass for algorithms that produce only polydata as output.
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition vtkPolyData.h:80
a cell that represents an n-sided polygon
Definition vtkPolygon.h:32
dynamic, self-adjusting array of unsigned char
int vtkTypeBool
Definition vtkABI.h:64
@ VTK_CCS_SCALAR_MODE_NONE
@ VTK_CCS_SCALAR_MODE_LABELS
@ VTK_CCS_SCALAR_MODE_COLORS
int vtkIdType
Definition vtkType.h:315
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:270