VTK  9.4.20250310
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
79#ifndef vtkClipClosedSurface_h
80#define vtkClipClosedSurface_h
81
82#include "vtkFiltersGeneralModule.h" // For export macro
84
85VTK_ABI_NAMESPACE_BEGIN
88class vtkDoubleArray;
89class vtkIdTypeArray;
90class vtkCellArray;
91class vtkPointData;
92class vtkCellData;
93class vtkPolygon;
94class vtkIdList;
95class vtkCCSEdgeLocator;
96
97enum
98{
103
104class VTKFILTERSGENERAL_EXPORT vtkClipClosedSurface : public vtkPolyDataAlgorithm
105{
106public:
108
113 void PrintSelf(ostream& os, vtkIndent indent) override;
115
117
121 vtkGetObjectMacro(ClippingPlanes, vtkPlaneCollection);
123
125
130 vtkSetMacro(Tolerance, double);
131 vtkGetMacro(Tolerance, double);
133
135
139 vtkSetMacro(PassPointData, vtkTypeBool);
140 vtkBooleanMacro(PassPointData, vtkTypeBool);
141 vtkGetMacro(PassPointData, vtkTypeBool);
143
145
149 vtkSetMacro(GenerateOutline, vtkTypeBool);
150 vtkBooleanMacro(GenerateOutline, vtkTypeBool);
151 vtkGetMacro(GenerateOutline, vtkTypeBool);
153
155
159 vtkSetMacro(GenerateFaces, vtkTypeBool);
160 vtkBooleanMacro(GenerateFaces, vtkTypeBool);
161 vtkGetMacro(GenerateFaces, vtkTypeBool);
163
165
174 vtkSetClampMacro(ScalarMode, int, VTK_CCS_SCALAR_MODE_NONE, VTK_CCS_SCALAR_MODE_LABELS);
175 void SetScalarModeToNone() { this->SetScalarMode(VTK_CCS_SCALAR_MODE_NONE); }
176 void SetScalarModeToColors() { this->SetScalarMode(VTK_CCS_SCALAR_MODE_COLORS); }
177 void SetScalarModeToLabels() { this->SetScalarMode(VTK_CCS_SCALAR_MODE_LABELS); }
178 vtkGetMacro(ScalarMode, int);
181
183
189 vtkSetVector3Macro(BaseColor, double);
190 vtkGetVector3Macro(BaseColor, double);
192
194
199 vtkSetVector3Macro(ClipColor, double);
200 vtkGetVector3Macro(ClipColor, double);
202
204
209 vtkSetMacro(ActivePlaneId, int);
210 vtkGetMacro(ActivePlaneId, int);
212
214
219 vtkSetVector3Macro(ActivePlaneColor, double);
220 vtkGetVector3Macro(ActivePlaneColor, double);
222
224
230 vtkSetMacro(TriangulationErrorDisplay, vtkTypeBool);
231 vtkBooleanMacro(TriangulationErrorDisplay, vtkTypeBool);
232 vtkGetMacro(TriangulationErrorDisplay, vtkTypeBool);
234
236
245 vtkSetMacro(InsideOut, vtkTypeBool);
246 vtkGetMacro(InsideOut, vtkTypeBool);
247 vtkBooleanMacro(InsideOut, vtkTypeBool);
249
251
256 vtkSetMacro(GenerateClipFaceOutput, vtkTypeBool);
257 vtkGetMacro(GenerateClipFaceOutput, vtkTypeBool);
258 vtkBooleanMacro(GenerateClipFaceOutput, vtkTypeBool);
260
265
266protected:
269
271
272 double Tolerance;
273
279 double BaseColor[3];
280 double ClipColor[3];
281 double ActivePlaneColor[3];
282 vtkTypeBool InsideOut = false;
283 vtkTypeBool GenerateClipFaceOutput = false;
284
286
288
290 vtkInformationVector* outputVector, int requestFromOutputPort, vtkMTimeType* mtime) override;
291
293 vtkInformationVector* outputVector) override;
294
298 void ClipLines(vtkPoints* points, vtkDoubleArray* pointScalars, vtkPointData* pointData,
299 vtkCCSEdgeLocator* edgeLocator, vtkCellArray* inputCells, vtkCellArray* outputLines,
300 vtkCellData* inCellData, vtkCellData* outLineData);
301
308 void ClipAndContourPolys(vtkPoints* points, vtkDoubleArray* pointScalars, vtkPointData* pointData,
309 vtkCCSEdgeLocator* edgeLocator, int triangulate, vtkCellArray* inputCells,
310 vtkCellArray* outputPolys, vtkCellArray* outputLines, vtkCellData* inCellData,
311 vtkCellData* outPolyData, vtkCellData* outLineData);
312
319 static int InterpolateEdge(vtkPoints* points, vtkPointData* pointData,
320 vtkCCSEdgeLocator* edgeLocator, double tol, vtkIdType i0, vtkIdType i1, double v0, double v1,
321 vtkIdType& i);
322
328 int TriangulatePolygon(vtkIdList* polygon, vtkPoints* points, vtkCellArray* triangles);
329
339 void TriangulateContours(vtkPolyData* data, vtkIdType firstLine, vtkIdType numLines,
340 vtkCellArray* outputPolys, const double normal[3]);
341
348 static void BreakPolylines(vtkCellArray* inputLines, vtkCellArray* outputLines,
349 vtkUnsignedCharArray* inputScalars, vtkIdType firstLineScalar,
350 vtkUnsignedCharArray* outputScalars, const unsigned char color[3]);
351
357 static void CopyPolygons(vtkCellArray* inputPolys, vtkCellArray* outputPolys,
358 vtkUnsignedCharArray* inputScalars, vtkIdType firstPolyScalar,
359 vtkUnsignedCharArray* outputScalars, const unsigned char color[3]);
360
365 static void BreakTriangleStrips(vtkCellArray* inputStrips, vtkCellArray* outputPolys,
366 vtkUnsignedCharArray* inputScalars, vtkIdType firstStripScalar,
367 vtkUnsignedCharArray* outputScalars, const unsigned char color[3]);
368
375 vtkPolyData* output, vtkPoints* points, vtkPointData* pointData, int outputPointDataType);
376
380 static void CreateColorValues(const double color1[3], const double color2[3],
381 const double color3[3], unsigned char colors[3][3]);
382
383private:
385 void operator=(const vtkClipClosedSurface&) = delete;
386};
387
388VTK_ABI_NAMESPACE_END
389#endif
object to represent cell connectivity
represent and manipulate cell attribute data
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:133
dynamic, self-adjusting array of vtkIdType
a simple class to control print indentation
Definition vtkIndent.h:108
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:139
Superclass for algorithms that produce only polydata as output.
concrete dataset represents vertices, lines, polygons, and triangle strips
a cell that represents an n-sided polygon
Definition vtkPolygon.h:131
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:332
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:287