VTK  9.1.0
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=========================================================================*/
113#ifndef vtkContour3DLinearGrid_h
114#define vtkContour3DLinearGrid_h
115
116#include "vtkContourValues.h" // Needed for inline methods
118#include "vtkFiltersCoreModule.h" // For export macro
119
120class vtkPolyData;
122class vtkScalarTree;
123struct vtkScalarTreeMap;
124
125class VTKFILTERSCORE_EXPORT vtkContour3DLinearGrid : public vtkDataObjectAlgorithm
126{
127public:
129
134 void PrintSelf(ostream& os, vtkIndent indent) override;
136
138
141 void SetValue(int i, double value);
142 double GetValue(int i);
143 double* GetValues();
144 void GetValues(double* contourValues);
145 void SetNumberOfContours(int number);
146 vtkIdType GetNumberOfContours();
147 void GenerateValues(int numContours, double range[2]);
148 void GenerateValues(int numContours, double rangeStart, double rangeEnd);
150
152
157 vtkSetMacro(MergePoints, vtkTypeBool);
158 vtkGetMacro(MergePoints, vtkTypeBool);
159 vtkBooleanMacro(MergePoints, vtkTypeBool);
161
163
167 vtkSetMacro(InterpolateAttributes, vtkTypeBool);
168 vtkGetMacro(InterpolateAttributes, vtkTypeBool);
169 vtkBooleanMacro(InterpolateAttributes, vtkTypeBool);
171
173
178 vtkSetMacro(ComputeNormals, vtkTypeBool);
179 vtkGetMacro(ComputeNormals, vtkTypeBool);
180 vtkBooleanMacro(ComputeNormals, vtkTypeBool);
182
184
190 vtkSetMacro(ComputeScalars, vtkTypeBool);
191 vtkGetMacro(ComputeScalars, vtkTypeBool);
192 vtkBooleanMacro(ComputeScalars, vtkTypeBool);
194
196
201 void SetOutputPointsPrecision(int precision);
204
210
212
217 vtkSetMacro(UseScalarTree, vtkTypeBool);
218 vtkGetMacro(UseScalarTree, vtkTypeBool);
219 vtkBooleanMacro(UseScalarTree, vtkTypeBool);
221
223
228 vtkGetObjectMacro(ScalarTree, vtkScalarTree);
230
232
240 vtkSetMacro(SequentialProcessing, vtkTypeBool);
241 vtkGetMacro(SequentialProcessing, vtkTypeBool);
242 vtkBooleanMacro(SequentialProcessing, vtkTypeBool);
244
249 int GetNumberOfThreadsUsed() { return this->NumberOfThreadsUsed; }
250
259 bool GetLargeIds() { return this->LargeIds; }
260
267 static bool CanFullyProcessDataObject(vtkDataObject* object, const char* scalarArrayName);
268
269protected:
272
281 bool LargeIds; // indicate whether integral ids are large(==true) or not
282
283 // Manage scalar trees, including mapping scalar tree to input dataset
286 struct vtkScalarTreeMap* ScalarTreeMap;
287
288 // Process the data: input unstructured grid and output polydata
289 void ProcessPiece(vtkUnstructuredGrid* input, vtkDataArray* inScalars, vtkPolyData* output);
290
292 vtkInformationVector* outputVector) override;
294 vtkInformationVector* outputVector) override;
296
297private:
299 void operator=(const vtkContour3DLinearGrid&) = delete;
300};
301
306inline void vtkContour3DLinearGrid::SetValue(int i, double value)
307{
308 this->ContourValues->SetValue(i, value);
309}
310
315{
316 return this->ContourValues->GetValue(i);
317}
318
324{
325 return this->ContourValues->GetValues();
326}
327
333inline void vtkContour3DLinearGrid::GetValues(double* contourValues)
334{
335 this->ContourValues->GetValues(contourValues);
336}
337
344{
345 this->ContourValues->SetNumberOfContours(number);
346}
347
352{
353 return this->ContourValues->GetNumberOfContours();
354}
355
360inline void vtkContour3DLinearGrid::GenerateValues(int numContours, double range[2])
361{
362 this->ContourValues->GenerateValues(numContours, range);
363}
364
370 int numContours, double rangeStart, double rangeEnd)
371{
372 this->ContourValues->GenerateValues(numContours, rangeStart, rangeEnd);
373}
374
375#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
vtkContourValues * ContourValues
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
Definition: vtkDataArray.h:159
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:113
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:195
organize data according to scalar values (used to accelerate contouring operations)
Definition: vtkScalarTree.h:55
dataset represents arbitrary combinations of all possible cell types
@ info
Definition: vtkX3D.h:382
@ value
Definition: vtkX3D.h:226
@ port
Definition: vtkX3D.h:453
@ range
Definition: vtkX3D.h:244
int vtkTypeBool
Definition: vtkABI.h:69
int vtkIdType
Definition: vtkType.h:332
vtkTypeUInt32 vtkMTimeType
Definition: vtkType.h:287