VTK
vtkBandedPolyDataContourFilter.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkBandedPolyDataContourFilter.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=========================================================================*/
47#ifndef vtkBandedPolyDataContourFilter_h
48#define vtkBandedPolyDataContourFilter_h
49
50#include "vtkFiltersModelingModule.h" // For export macro
52
53#include "vtkContourValues.h" // Needed for inline methods
54
55class vtkPoints;
56class vtkCellArray;
57class vtkPointData;
58class vtkDataArray;
59class vtkFloatArray;
60class vtkDoubleArray;
61
62#define VTK_SCALAR_MODE_INDEX 0
63#define VTK_SCALAR_MODE_VALUE 1
64
65class VTKFILTERSMODELING_EXPORT vtkBandedPolyDataContourFilter : public vtkPolyDataAlgorithm
66{
67public:
69 void PrintSelf(ostream& os, vtkIndent indent);
70
75
77
83 void SetValue(int i, double value);
84 double GetValue(int i);
85 double *GetValues();
86 void GetValues(double *contourValues);
87 void SetNumberOfContours(int number);
88 int GetNumberOfContours();
89 void GenerateValues(int numContours, double range[2]);
90 void GenerateValues(int numContours, double rangeStart, double rangeEnd);
92
94
100 vtkSetMacro(Clipping,int);
101 vtkGetMacro(Clipping,int);
102 vtkBooleanMacro(Clipping,int);
104
106
112 vtkSetClampMacro(ScalarMode,int,VTK_SCALAR_MODE_INDEX,VTK_SCALAR_MODE_VALUE);
113 vtkGetMacro(ScalarMode,int);
115 {this->SetScalarMode(VTK_SCALAR_MODE_INDEX);}
117 {this->SetScalarMode(VTK_SCALAR_MODE_VALUE);}
119
121
127 vtkSetMacro(GenerateContourEdges,int);
128 vtkGetMacro(GenerateContourEdges,int);
129 vtkBooleanMacro(GenerateContourEdges,int);
131
133
139 vtkSetMacro(ClipTolerance,double);
140 vtkGetMacro(ClipTolerance,double);
142
144
148 vtkSetMacro(Component,int);
149 vtkGetMacro(Component,int);
151
157
163
164protected:
167
169
171 int IsContourValue(double val);
172 int ClipEdge(int v1, int v2, vtkPoints *pts, vtkDataArray *inScalars,
173 vtkDoubleArray *outScalars,
174 vtkPointData *inPD, vtkPointData *outPD, vtkIdType edgePts[]);
175 int InsertCell(vtkCellArray *cells, int npts, vtkIdType *pts,
176 int cellId, double s, vtkFloatArray *newS);
178 int cellId, double s, vtkFloatArray *newS);
179 int ComputeClippedIndex(double s);
180 int InsertNextScalar(vtkFloatArray* scalars, int cellId, int idx);
181 // data members
183
187
188 // sorted and cleaned contour values
189 double *ClipValues;
191 int ClipIndex[2]; //indices outside of this range (inclusive) are clipped
192 double ClipTolerance; //specify numerical accuracy during clipping
193 double InternalClipTolerance; //used to clean up numerical problems
194
195 //the second output
197
198private:
200 void operator=(const vtkBandedPolyDataContourFilter&) VTK_DELETE_FUNCTION;
201};
202
208 {this->ContourValues->SetValue(i,value);}
209
214 {return this->ContourValues->GetValue(i);}
215
221 {return this->ContourValues->GetValues();}
222
228inline void vtkBandedPolyDataContourFilter::GetValues(double *contourValues)
229 {this->ContourValues->GetValues(contourValues);}
230
237 {this->ContourValues->SetNumberOfContours(number);}
238
243 {return this->ContourValues->GetNumberOfContours();}
244
250 double range[2])
251 {this->ContourValues->GenerateValues(numContours, range);}
252
258 double rangeStart,
259 double rangeEnd)
260 {this->ContourValues->GenerateValues(numContours, rangeStart, rangeEnd);}
261
262
263#endif
generate filled contours for vtkPolyData
double * GetValues()
Get a pointer to an array of contour values.
double GetValue(int i)
Get the ith contour value.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *)
This is called by the superclass.
int InsertCell(vtkCellArray *cells, int npts, vtkIdType *pts, int cellId, double s, vtkFloatArray *newS)
int InsertLine(vtkCellArray *cells, vtkIdType pt1, vtkIdType pt2, 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.
void GenerateValues(int numContours, double range[2])
Generate numContours equally spaced contour values between specified range.
int GetNumberOfContours()
Get the number of contours in the list of contour values.
vtkMTimeType GetMTime()
Overload GetMTime because we delegate to vtkContourValues so its modified time must be taken into acc...
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.
void PrintSelf(ostream &os, vtkIndent indent)
Methods invoked by print to print information about the object including superclasses.
object to represent cell connectivity
Definition: vtkCellArray.h:51
helper object to manage setting and generating contour values
double * GetValues()
Return a pointer to a list of contour values.
int GetNumberOfContours()
Return the number of contours in the.
void GenerateValues(int numContours, double range[2])
Generate numContours equally spaced contour values between specified range.
void SetNumberOfContours(const int number)
Set the number of contours to place into the list.
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
Definition: vtkDataArray.h:55
dynamic, self-adjusting array of double
dynamic, self-adjusting array of float
Definition: vtkFloatArray.h:42
a simple class to control print indentation
Definition: vtkIndent.h:40
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
represent and manipulate point attribute data
Definition: vtkPointData.h:38
represent and manipulate 3D points
Definition: vtkPoints.h:40
Superclass for algorithms that produce only polydata as output.
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:86
@ value
Definition: vtkX3D.h:220
@ range
Definition: vtkX3D.h:238
#define VTK_SCALAR_MODE_VALUE
#define VTK_SCALAR_MODE_INDEX
vtkSetMacro(IgnoreDriverBugs, bool)
Updates the extensions string.
vtkBooleanMacro(IgnoreDriverBugs, bool)
Updates the extensions string.
int vtkIdType
Definition: vtkType.h:287
vtkTypeUInt64 vtkMTimeType
Definition: vtkType.h:248