VTK  9.3.20240424
vtkBandedPolyDataContourFilter.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
101#ifndef vtkBandedPolyDataContourFilter_h
102#define vtkBandedPolyDataContourFilter_h
103
104#include "vtkFiltersModelingModule.h" // For export macro
105#include "vtkPolyDataAlgorithm.h"
106
107#include "vtkContourValues.h" // Needed for inline methods
108
109VTK_ABI_NAMESPACE_BEGIN
110class vtkPoints;
111class vtkCellArray;
112class vtkPointData;
113class vtkDataArray;
114class vtkFloatArray;
115class vtkDoubleArray;
116struct vtkBandedPolyDataContourFilterInternals;
117
118#define VTK_SCALAR_MODE_INDEX 0
119#define VTK_SCALAR_MODE_VALUE 1
120
121class VTKFILTERSMODELING_EXPORT vtkBandedPolyDataContourFilter : public vtkPolyDataAlgorithm
122{
123public:
125 void PrintSelf(ostream& os, vtkIndent indent) override;
126
131
133
139 void SetValue(int i, double value);
140 double GetValue(int i);
141 double* GetValues();
142 void GetValues(double* contourValues);
143 void SetNumberOfContours(int number);
144 vtkIdType GetNumberOfContours();
145 void GenerateValues(int numContours, double range[2]);
146 void GenerateValues(int numContours, double rangeStart, double rangeEnd);
148
150
156 vtkSetMacro(Clipping, vtkTypeBool);
157 vtkGetMacro(Clipping, vtkTypeBool);
158 vtkBooleanMacro(Clipping, vtkTypeBool);
160
162
168 vtkSetClampMacro(ScalarMode, int, VTK_SCALAR_MODE_INDEX, VTK_SCALAR_MODE_VALUE);
169 vtkGetMacro(ScalarMode, int);
170 void SetScalarModeToIndex() { this->SetScalarMode(VTK_SCALAR_MODE_INDEX); }
171 void SetScalarModeToValue() { this->SetScalarMode(VTK_SCALAR_MODE_VALUE); }
173
175
181 vtkSetMacro(GenerateContourEdges, vtkTypeBool);
182 vtkGetMacro(GenerateContourEdges, vtkTypeBool);
183 vtkBooleanMacro(GenerateContourEdges, vtkTypeBool);
185
187
193 vtkSetMacro(ClipTolerance, double);
194 vtkGetMacro(ClipTolerance, double);
196
198
202 vtkSetMacro(Component, int);
203 vtkGetMacro(Component, int);
205
211
217
218protected:
221
223
224 int ClipEdge(int v1, int v2, vtkPoints* pts, vtkDataArray* inScalars, vtkDoubleArray* outScalars,
225 vtkPointData* inPD, vtkPointData* outPD, vtkIdType edgePts[]);
227 vtkCellArray* cells, int npts, const vtkIdType* pts, int cellId, double s, vtkFloatArray* newS);
229 vtkCellArray* cells, vtkIdType pt1, vtkIdType pt2, int cellId, double s, vtkFloatArray* newS);
230 int ComputeClippedIndex(double s);
231 int InsertNextScalar(vtkFloatArray* scalars, int cellId, int idx);
232 // data members
234
238 double ClipTolerance; // specify numerical accuracy during clipping
239 // the second output
241
242 vtkBandedPolyDataContourFilterInternals* Internal;
243
244private:
246 void operator=(const vtkBandedPolyDataContourFilter&) = delete;
247};
248
253inline void vtkBandedPolyDataContourFilter::SetValue(int i, double value)
254{
255 this->ContourValues->SetValue(i, value);
256}
257
262{
263 return this->ContourValues->GetValue(i);
264}
265
271{
272 return this->ContourValues->GetValues();
273}
274
280inline void vtkBandedPolyDataContourFilter::GetValues(double* contourValues)
281{
282 this->ContourValues->GetValues(contourValues);
283}
284
291{
292 this->ContourValues->SetNumberOfContours(number);
293}
294
299{
300 return this->ContourValues->GetNumberOfContours();
301}
302
307inline void vtkBandedPolyDataContourFilter::GenerateValues(int numContours, double range[2])
308{
309 this->ContourValues->GenerateValues(numContours, range);
310}
311
317 int numContours, double rangeStart, double rangeEnd)
318{
319 this->ContourValues->GenerateValues(numContours, rangeStart, rangeEnd);
320}
321
322VTK_ABI_NAMESPACE_END
323#endif
generate filled contours for vtkPolyData
vtkMTimeType GetMTime() override
Overload GetMTime because we delegate to vtkContourValues so its modified time must be taken into acc...
double * GetValues()
Get a pointer to an array of contour values.
double GetValue(int i)
Get the ith contour value.
int InsertLine(vtkCellArray *cells, vtkIdType pt1, vtkIdType pt2, int cellId, double s, vtkFloatArray *newS)
int InsertCell(vtkCellArray *cells, int npts, const vtkIdType *pts, int cellId, double s, vtkFloatArray *newS)
vtkSmartPointer< vtkContourValues > ContourValues
void SetScalarModeToIndex()
Control whether the cell scalars are output as an integer index or a scalar value.
~vtkBandedPolyDataContourFilter() override
vtkPolyData * GetContourEdgesOutput()
Get the second output which contains the edges dividing the contour bands.
static vtkBandedPolyDataContourFilter * New()
Construct object with no contours defined.
void GenerateValues(int numContours, double range[2])
Generate numContours equally spaced contour values between specified range.
void SetScalarModeToValue()
Control whether the cell scalars are output as an integer index or a scalar value.
int InsertNextScalar(vtkFloatArray *scalars, int cellId, int idx)
int ClipEdge(int v1, int v2, vtkPoints *pts, vtkDataArray *inScalars, vtkDoubleArray *outScalars, vtkPointData *inPD, vtkPointData *outPD, vtkIdType edgePts[])
void SetValue(int i, double value)
Methods to set / get contour values.
void SetNumberOfContours(int number)
Set the number of contours to place into the list.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
vtkBandedPolyDataContourFilterInternals * Internal
vtkIdType GetNumberOfContours()
Get the number of contours in the list of contour values.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
object to represent cell connectivity
double * GetValues()
Return a pointer to a list of contour values.
int GetNumberOfContours()
Return the number of contours in the.
void SetNumberOfContours(int number)
Set the number of contours to place into the list.
void GenerateValues(int numContours, double range[2])
Generate numContours equally spaced contour values between specified range.
void SetValue(int i, double value)
Set the ith contour value.
double GetValue(int i)
Get the ith contour value.
abstract superclass for arrays of numeric data
dynamic, self-adjusting array of double
dynamic, self-adjusting array of float
a simple class to control print indentation
Definition vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
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
Hold a reference to a vtkObjectBase instance.
int vtkTypeBool
Definition vtkABI.h:64
#define VTK_SCALAR_MODE_VALUE
#define VTK_SCALAR_MODE_INDEX
int vtkIdType
Definition vtkType.h:315
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:270