VTK
vtkContourFilter.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkContourFilter.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=========================================================================*/
53#ifndef vtkContourFilter_h
54#define vtkContourFilter_h
55
56#include "vtkFiltersCoreModule.h" // For export macro
58
59#include "vtkContourValues.h" // Needed for inline methods
60
62class vtkScalarTree;
68
69class VTKFILTERSCORE_EXPORT vtkContourFilter : public vtkPolyDataAlgorithm
70{
71public:
73 void PrintSelf(ostream& os, vtkIndent indent) VTK_OVERRIDE;
74
80
82
85 void SetValue(int i, double value);
86 double GetValue(int i);
87 double *GetValues();
88 void GetValues(double *contourValues);
89 void SetNumberOfContours(int number);
90 int GetNumberOfContours();
91 void GenerateValues(int numContours, double range[2]);
92 void GenerateValues(int numContours, double rangeStart, double rangeEnd);
94
98 vtkMTimeType GetMTime() VTK_OVERRIDE;
99
101
111 vtkSetMacro(ComputeNormals,int);
112 vtkGetMacro(ComputeNormals,int);
113 vtkBooleanMacro(ComputeNormals,int);
115
117
125 vtkSetMacro(ComputeGradients,int);
126 vtkGetMacro(ComputeGradients,int);
127 vtkBooleanMacro(ComputeGradients,int);
129
131
134 vtkSetMacro(ComputeScalars,int);
135 vtkGetMacro(ComputeScalars,int);
136 vtkBooleanMacro(ComputeScalars,int);
138
140
143 vtkSetMacro(UseScalarTree,int);
144 vtkGetMacro(UseScalarTree,int);
145 vtkBooleanMacro(UseScalarTree,int);
147
149
152 virtual void SetScalarTree(vtkScalarTree*);
153 vtkGetObjectMacro(ScalarTree,vtkScalarTree);
155
157
161 void SetLocator(vtkIncrementalPointLocator *locator);
162 vtkGetObjectMacro(Locator,vtkIncrementalPointLocator);
164
169 void CreateDefaultLocator();
170
172
176 void SetArrayComponent( int );
177 int GetArrayComponent();
179
180
182
189 vtkSetMacro(GenerateTriangles,int);
190 vtkGetMacro(GenerateTriangles,int);
191 vtkBooleanMacro(GenerateTriangles,int);
193
195
200 void SetOutputPointsPrecision(int precision);
201 int GetOutputPointsPrecision() const;
203
204protected:
206 ~vtkContourFilter() VTK_OVERRIDE;
207
208 void ReportReferences(vtkGarbageCollector*) VTK_OVERRIDE;
209
210 int RequestData(vtkInformation* request,
211 vtkInformationVector** inputVector,
212 vtkInformationVector* outputVector) VTK_OVERRIDE;
213 int RequestUpdateExtent(vtkInformation*,
215 vtkInformationVector*) VTK_OVERRIDE;
216 int FillInputPortInformation(int port, vtkInformation *info) VTK_OVERRIDE;
217
218 vtkContourValues *ContourValues;
219 int ComputeNormals;
220 int ComputeGradients;
221 int ComputeScalars;
223 int UseScalarTree;
224 vtkScalarTree *ScalarTree;
225 int OutputPointsPrecision;
226 int GenerateTriangles;
227
228 vtkSynchronizedTemplates2D *SynchronizedTemplates2D;
229 vtkSynchronizedTemplates3D *SynchronizedTemplates3D;
230 vtkGridSynchronizedTemplates3D *GridSynchronizedTemplates;
231 vtkRectilinearSynchronizedTemplates *RectilinearSynchronizedTemplates;
232 vtkCallbackCommand *InternalProgressCallbackCommand;
233
234 static void InternalProgressCallbackFunction(vtkObject *caller,
235 unsigned long eid,
236 void *clientData,
237 void *callData);
238
239private:
240 vtkContourFilter(const vtkContourFilter&) VTK_DELETE_FUNCTION;
241 void operator=(const vtkContourFilter&) VTK_DELETE_FUNCTION;
242};
243
248inline void vtkContourFilter::SetValue(int i, double value)
249{this->ContourValues->SetValue(i,value);}
250
254inline double vtkContourFilter::GetValue(int i)
255{return this->ContourValues->GetValue(i);}
256
262{return this->ContourValues->GetValues();}
263
269inline void vtkContourFilter::GetValues(double *contourValues)
270{this->ContourValues->GetValues(contourValues);}
271
278{this->ContourValues->SetNumberOfContours(number);}
279
284{return this->ContourValues->GetNumberOfContours();}
285
290inline void vtkContourFilter::GenerateValues(int numContours, double range[2])
291{this->ContourValues->GenerateValues(numContours, range);}
292
297inline void vtkContourFilter::GenerateValues(int numContours, double
298 rangeStart, double rangeEnd)
299{this->ContourValues->GenerateValues(numContours, rangeStart, rangeEnd);}
300
301
302#endif
303
304
supports function callbacks
generate isosurfaces/isolines from scalar values
double GetValue(int i)
Get the ith contour value.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkMTimeType GetMTime() override
Modified GetMTime Because we delegate to vtkContourValues.
void SetNumberOfContours(int number)
Set the number of contours to place into the list.
int GetNumberOfContours()
Get the number of contours in the list of contour values.
void GenerateValues(int numContours, double range[2])
Generate numContours equally spaced contour values between specified range.
static vtkContourFilter * New()
Construct object with initial range (0,1) and single contour value of 0.0.
double * GetValues()
Get a pointer to an array of contour values.
helper object to manage setting and generating contour values
void SetValue(int i, double value)
Set the ith contour value.
Detect and break reference loops.
generate isosurface from structured grids
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:40
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
abstract base class for most VTK objects
Definition: vtkObject.h:60
Superclass for algorithms that produce only polydata as output.
generate isosurface from rectilinear grid
organize data according to scalar values (used to accelerate contouring operations)
Definition: vtkScalarTree.h:55
generate isoline(s) from a structured points set
generate isosurface from structured points
@ info
Definition: vtkX3D.h:376
@ value
Definition: vtkX3D.h:220
@ port
Definition: vtkX3D.h:447
@ range
Definition: vtkX3D.h:238
vtkSetMacro(IgnoreDriverBugs, bool)
Updates the extensions string.
vtkBooleanMacro(IgnoreDriverBugs, bool)
Updates the extensions string.
vtkTypeUInt64 vtkMTimeType
Definition: vtkType.h:248