VTK  9.2.6
vtkContour3DLinearGrid.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkContour3DLinearGrid.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=========================================================================*/
116#ifndef vtkContour3DLinearGrid_h
117#define vtkContour3DLinearGrid_h
118
119#include "vtkContourValues.h" // Needed for inline methods
121#include "vtkFiltersCoreModule.h" // For export macro
122
123class vtkPolyData;
125class vtkScalarTree;
126struct vtkScalarTreeMap;
127
128class VTKFILTERSCORE_EXPORT vtkContour3DLinearGrid : public vtkDataObjectAlgorithm
129{
130public:
132
137 void PrintSelf(ostream& os, vtkIndent indent) override;
139
141
144 void SetValue(int i, double value);
145 double GetValue(int i);
146 double* GetValues();
147 void GetValues(double* contourValues);
148 void SetNumberOfContours(int number);
149 vtkIdType GetNumberOfContours();
150 void GenerateValues(int numContours, double range[2]);
151 void GenerateValues(int numContours, double rangeStart, double rangeEnd);
153
155
160 vtkSetMacro(MergePoints, vtkTypeBool);
161 vtkGetMacro(MergePoints, vtkTypeBool);
162 vtkBooleanMacro(MergePoints, vtkTypeBool);
164
166
170 vtkSetMacro(InterpolateAttributes, vtkTypeBool);
171 vtkGetMacro(InterpolateAttributes, vtkTypeBool);
172 vtkBooleanMacro(InterpolateAttributes, vtkTypeBool);
174
176
181 vtkSetMacro(ComputeNormals, vtkTypeBool);
182 vtkGetMacro(ComputeNormals, vtkTypeBool);
183 vtkBooleanMacro(ComputeNormals, vtkTypeBool);
185
187
193 vtkSetMacro(ComputeScalars, vtkTypeBool);
194 vtkGetMacro(ComputeScalars, vtkTypeBool);
195 vtkBooleanMacro(ComputeScalars, vtkTypeBool);
197
199
204 void SetOutputPointsPrecision(int precision);
207
213
215
220 vtkSetMacro(UseScalarTree, vtkTypeBool);
221 vtkGetMacro(UseScalarTree, vtkTypeBool);
222 vtkBooleanMacro(UseScalarTree, vtkTypeBool);
224
226
231 vtkGetObjectMacro(ScalarTree, vtkScalarTree);
233
235
243 vtkSetMacro(SequentialProcessing, vtkTypeBool);
244 vtkGetMacro(SequentialProcessing, vtkTypeBool);
245 vtkBooleanMacro(SequentialProcessing, vtkTypeBool);
247
252 int GetNumberOfThreadsUsed() { return this->NumberOfThreadsUsed; }
253
262 bool GetLargeIds() { return this->LargeIds; }
263
270 static bool CanFullyProcessDataObject(vtkDataObject* object, const char* scalarArrayName);
271
272protected:
275
284 bool LargeIds; // indicate whether integral ids are large(==true) or not
285
286 // Manage scalar trees, including mapping scalar tree to input dataset
289 struct vtkScalarTreeMap* ScalarTreeMap;
290
291 // Process the data: input unstructured grid and output polydata
292 void ProcessPiece(vtkUnstructuredGrid* input, vtkDataArray* inScalars, vtkPolyData* output);
293
295 vtkInformationVector* outputVector) override;
297 vtkInformationVector* outputVector) override;
298 int FillInputPortInformation(int port, vtkInformation* info) override;
299
300private:
302 void operator=(const vtkContour3DLinearGrid&) = delete;
303};
304
309inline void vtkContour3DLinearGrid::SetValue(int i, double value)
310{
311 this->ContourValues->SetValue(i, value);
312}
313
318{
319 return this->ContourValues->GetValue(i);
320}
321
327{
328 return this->ContourValues->GetValues();
329}
330
336inline void vtkContour3DLinearGrid::GetValues(double* contourValues)
337{
338 this->ContourValues->GetValues(contourValues);
339}
340
347{
348 this->ContourValues->SetNumberOfContours(number);
349}
350
355{
356 return this->ContourValues->GetNumberOfContours();
357}
358
363inline void vtkContour3DLinearGrid::GenerateValues(int numContours, double range[2])
364{
365 this->ContourValues->GenerateValues(numContours, range);
366}
367
373 int numContours, double rangeStart, double rangeEnd)
374{
375 this->ContourValues->GenerateValues(numContours, rangeStart, rangeEnd);
376}
377
378#endif
fast generation of isosurface from 3D linear cells
int GetOutputPointsPrecision() const
Set/get the desired precision for the output types.
static bool CanFullyProcessDataObject(vtkDataObject *object, const char *scalarArrayName)
Returns true if the data object passed in is fully supported by this filter, i.e.,...
bool GetLargeIds()
Inform the user as to whether large ids were used during filter execution.
static vtkContour3DLinearGrid * New()
Standard methods for construction, type info, and printing.
vtkMTimeType GetMTime() override
Overloaded GetMTime() because of delegation to the internal vtkContourValues class.
void SetOutputPointsPrecision(int precision)
Set/get the desired precision for the output types.
double * GetValues()
Get a pointer to an array of contour values.
virtual void SetScalarTree(vtkScalarTree *)
Specify the scalar tree to use.
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
void ProcessPiece(vtkUnstructuredGrid *input, vtkDataArray *inScalars, vtkPolyData *output)
void SetNumberOfContours(int number)
Set the number of contours to place into the list.
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
vtkIdType GetNumberOfContours()
Get the number of contours in the list of contour values.
int RequestDataObject(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
This is called by the superclass.
void SetValue(int i, double value)
Methods to set / get contour values.
struct vtkScalarTreeMap * ScalarTreeMap
~vtkContour3DLinearGrid() override
double GetValue(int i)
Get the ith contour value.
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods for construction, type info, and printing.
int GetNumberOfThreadsUsed()
Return the number of threads actually used during execution.
void GenerateValues(int numContours, double range[2])
Generate numContours equally spaced contour values between specified range.
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
Superclass for algorithms that produce only data object as output.
general representation of visualization data
a simple class to control print indentation
Definition vtkIndent.h:40
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition vtkPolyData.h:91
organize data according to scalar values (used to accelerate contouring operations)
dataset represents arbitrary combinations of all possible cell types
int vtkTypeBool
Definition vtkABI.h:69
int vtkIdType
Definition vtkType.h:332
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:287