VTK  9.3.0
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
32 #ifndef vtkBandedPolyDataContourFilter_h
33 #define vtkBandedPolyDataContourFilter_h
34 
35 #include "vtkFiltersModelingModule.h" // For export macro
36 #include "vtkPolyDataAlgorithm.h"
37 
38 #include "vtkContourValues.h" // Needed for inline methods
39 
40 VTK_ABI_NAMESPACE_BEGIN
41 class vtkPoints;
42 class vtkCellArray;
43 class vtkPointData;
44 class vtkDataArray;
45 class vtkFloatArray;
46 class vtkDoubleArray;
47 struct vtkBandedPolyDataContourFilterInternals;
48 
49 #define VTK_SCALAR_MODE_INDEX 0
50 #define VTK_SCALAR_MODE_VALUE 1
51 
52 class VTKFILTERSMODELING_EXPORT vtkBandedPolyDataContourFilter : public vtkPolyDataAlgorithm
53 {
54 public:
56  void PrintSelf(ostream& os, vtkIndent indent) override;
57 
62 
64 
70  void SetValue(int i, double value);
71  double GetValue(int i);
72  double* GetValues();
73  void GetValues(double* contourValues);
74  void SetNumberOfContours(int number);
75  vtkIdType GetNumberOfContours();
76  void GenerateValues(int numContours, double range[2]);
77  void GenerateValues(int numContours, double rangeStart, double rangeEnd);
79 
81 
87  vtkSetMacro(Clipping, vtkTypeBool);
88  vtkGetMacro(Clipping, vtkTypeBool);
89  vtkBooleanMacro(Clipping, vtkTypeBool);
91 
93 
99  vtkSetClampMacro(ScalarMode, int, VTK_SCALAR_MODE_INDEX, VTK_SCALAR_MODE_VALUE);
100  vtkGetMacro(ScalarMode, int);
101  void SetScalarModeToIndex() { this->SetScalarMode(VTK_SCALAR_MODE_INDEX); }
102  void SetScalarModeToValue() { this->SetScalarMode(VTK_SCALAR_MODE_VALUE); }
104 
106 
112  vtkSetMacro(GenerateContourEdges, vtkTypeBool);
113  vtkGetMacro(GenerateContourEdges, vtkTypeBool);
114  vtkBooleanMacro(GenerateContourEdges, vtkTypeBool);
116 
118 
124  vtkSetMacro(ClipTolerance, double);
125  vtkGetMacro(ClipTolerance, double);
127 
129 
133  vtkSetMacro(Component, int);
134  vtkGetMacro(Component, int);
136 
142 
147  vtkMTimeType GetMTime() override;
148 
149 protected:
152 
154 
155  int ClipEdge(int v1, int v2, vtkPoints* pts, vtkDataArray* inScalars, vtkDoubleArray* outScalars,
156  vtkPointData* inPD, vtkPointData* outPD, vtkIdType edgePts[]);
158  vtkCellArray* cells, int npts, const vtkIdType* pts, int cellId, double s, vtkFloatArray* newS);
160  vtkCellArray* cells, vtkIdType pt1, vtkIdType pt2, int cellId, double s, vtkFloatArray* newS);
161  int ComputeClippedIndex(double s);
162  int InsertNextScalar(vtkFloatArray* scalars, int cellId, int idx);
163  // data members
165 
169  double ClipTolerance; // specify numerical accuracy during clipping
170  // the second output
172 
173  vtkBandedPolyDataContourFilterInternals* Internal;
174 
175 private:
177  void operator=(const vtkBandedPolyDataContourFilter&) = delete;
178 };
179 
185 {
186  this->ContourValues->SetValue(i, value);
187 }
188 
193 {
194  return this->ContourValues->GetValue(i);
195 }
196 
202 {
203  return this->ContourValues->GetValues();
204 }
205 
211 inline void vtkBandedPolyDataContourFilter::GetValues(double* contourValues)
212 {
213  this->ContourValues->GetValues(contourValues);
214 }
215 
222 {
223  this->ContourValues->SetNumberOfContours(number);
224 }
225 
230 {
231  return this->ContourValues->GetNumberOfContours();
232 }
233 
238 inline void vtkBandedPolyDataContourFilter::GenerateValues(int numContours, double range[2])
239 {
240  this->ContourValues->GenerateValues(numContours, range);
241 }
242 
248  int numContours, double rangeStart, double rangeEnd)
249 {
250  this->ContourValues->GenerateValues(numContours, rangeStart, rangeEnd);
251 }
252 
253 VTK_ABI_NAMESPACE_END
254 #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)
vtkPolyData * GetContourEdgesOutput()
Get the second output which contains the edges dividing the contour bands.
static vtkBandedPolyDataContourFilter * New()
Construct object with no contours defined.
vtkSmartPointer< vtkContourValues > ContourValues
void SetScalarModeToIndex()
Control whether the cell scalars are output as an integer index or a scalar value.
~vtkBandedPolyDataContourFilter() override
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
Definition: vtkCellArray.h:176
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.
double * GetValues()
Return a pointer to a list of contour values.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:45
dynamic, self-adjusting array of double
dynamic, self-adjusting array of float
Definition: vtkFloatArray.h:31
a simple class to control print indentation
Definition: vtkIndent.h:29
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
represent and manipulate point attribute data
Definition: vtkPointData.h:30
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
@ value
Definition: vtkX3D.h:220
@ range
Definition: vtkX3D.h:238
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