VTK  9.1.0
vtkTableBasedClipDataSet.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkTableBasedClipDataSet.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 =========================================================================*/
15 
16 /*****************************************************************************
17  *
18  * Copyright (c) 2000 - 2009, Lawrence Livermore National Security, LLC
19  * Produced at the Lawrence Livermore National Laboratory
20  * LLNL-CODE-400124
21  * All rights reserved.
22  *
23  * This file was adapted from the VisIt clipper (vtkVisItClipper). For details,
24  * see https://visit.llnl.gov/. The full copyright notice is contained in the
25  * file COPYRIGHT located at the root of the VisIt distribution or at
26  * http://www.llnl.gov/visit/copyright.html.
27  *
28  *****************************************************************************/
29 
91 #ifndef vtkTableBasedClipDataSet_h
92 #define vtkTableBasedClipDataSet_h
93 
94 #include "vtkFiltersGeneralModule.h" // For export macro
96 
97 class vtkCallbackCommand;
100 
101 class VTKFILTERSGENERAL_EXPORT vtkTableBasedClipDataSet : public vtkUnstructuredGridAlgorithm
102 {
103 public:
105  void PrintSelf(ostream& os, vtkIndent indent) override;
106 
112 
116  vtkMTimeType GetMTime() override;
117 
119 
126  vtkSetMacro(InsideOut, vtkTypeBool);
127  vtkGetMacro(InsideOut, vtkTypeBool);
128  vtkBooleanMacro(InsideOut, vtkTypeBool);
130 
132 
138  vtkSetMacro(Value, double);
139  vtkGetMacro(Value, double);
141 
143 
148  vtkSetMacro(UseValueAsOffset, bool);
149  vtkGetMacro(UseValueAsOffset, bool);
150  vtkBooleanMacro(UseValueAsOffset, bool);
152 
154 
160  vtkGetObjectMacro(ClipFunction, vtkImplicitFunction);
162 
164 
170  vtkSetMacro(GenerateClipScalars, vtkTypeBool);
171  vtkGetMacro(GenerateClipScalars, vtkTypeBool);
172  vtkBooleanMacro(GenerateClipScalars, vtkTypeBool);
174 
176 
185  vtkGetObjectMacro(Locator, vtkIncrementalPointLocator);
187 
189 
194  vtkSetClampMacro(MergeTolerance, double, 0.0001, 0.25);
195  vtkGetMacro(MergeTolerance, double);
197 
203 
205 
209  vtkSetMacro(GenerateClippedOutput, vtkTypeBool);
210  vtkGetMacro(GenerateClippedOutput, vtkTypeBool);
211  vtkBooleanMacro(GenerateClippedOutput, vtkTypeBool);
213 
218 
220 
225  vtkSetClampMacro(OutputPointsPrecision, int, SINGLE_PRECISION, DEFAULT_PRECISION);
226  vtkGetMacro(OutputPointsPrecision, int);
228 
229 protected:
232 
235 
241  void ClipDataSet(vtkDataSet* pDataSet, vtkDataArray* clipAray, vtkUnstructuredGrid* unstruct);
242 
248  vtkDataSet* inputGrd, vtkDataArray* clipAray, double isoValue, vtkUnstructuredGrid* outputUG);
249 
257  vtkDataSet* inputGrd, vtkDataArray* clipAray, double isoValue, vtkUnstructuredGrid* outputUG);
258 
266  vtkDataSet* inputGrd, vtkDataArray* clipAray, double isoValue, vtkUnstructuredGrid* outputUG);
267 
275  vtkDataSet* inputGrd, vtkDataArray* clipAray, double isoValue, vtkUnstructuredGrid* outputUG);
276 
284  vtkDataSet* inputGrd, vtkDataArray* clipAray, double isoValue, vtkUnstructuredGrid* outputUG);
285 
289  static void InternalProgressCallbackFunction(vtkObject*, unsigned long, void* clientdata, void*);
290 
295 
300  double Value;
305 
307 
308 private:
310  void operator=(const vtkTableBasedClipDataSet&) = delete;
311 };
312 
313 #endif
Superclass for all sources, filters, and sinks in VTK.
Definition: vtkAlgorithm.h:64
supports function callbacks
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:50
abstract class to specify dataset behavior
Definition: vtkDataSet.h:57
abstract interface for implicit functions
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:34
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
abstract base class for most VTK objects
Definition: vtkObject.h:63
Clip any dataset with a user-specified implicit function or an input scalar point data array.
vtkTableBasedClipDataSet(vtkImplicitFunction *cf=nullptr)
void ClipDataSet(vtkDataSet *pDataSet, vtkDataArray *clipAray, vtkUnstructuredGrid *unstruct)
This function resorts to the sibling class vtkClipDataSet to handle special grids (such as cylinders ...
static vtkTableBasedClipDataSet * New()
Create an instance with a user-specified implicit function, turning off IVARs InsideOut and GenerateC...
void ClipPolyData(vtkDataSet *inputGrd, vtkDataArray *clipAray, double isoValue, vtkUnstructuredGrid *outputUG)
This function clips a vtkPolyData object based on a specified iso-value (isoValue) using a scalar poi...
void ClipStructuredGridData(vtkDataSet *inputGrd, vtkDataArray *clipAray, double isoValue, vtkUnstructuredGrid *outputUG)
This function clips a vtkStructuredGrid based on a specified iso-value (isoValue) using a scalar poin...
~vtkTableBasedClipDataSet() override
void SetLocator(vtkIncrementalPointLocator *locator)
Set/Get a point locator locator for merging duplicate points.
void ClipUnstructuredGridData(vtkDataSet *inputGrd, vtkDataArray *clipAray, double isoValue, vtkUnstructuredGrid *outputUG)
This function clips a vtkUnstructuredGrid based on a specified iso-value (isoValue) using a scalar po...
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
static void InternalProgressCallbackFunction(vtkObject *, unsigned long, void *clientdata, void *)
Register a callback function with the InternalProgressObserver.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
vtkImplicitFunction * ClipFunction
virtual void SetClipFunction(vtkImplicitFunction *)
Set/Get the implicit function with which to perform the clipping operation.
vtkUnstructuredGrid * GetClippedOutput()
Return the clipped output.
vtkMTimeType GetMTime() override
Get the MTime for which the point locator and clip function are considered.
vtkIncrementalPointLocator * Locator
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void InternalProgressCallback(vtkAlgorithm *algorithm)
The actual operation executed by the callback function.
void CreateDefaultLocator()
Create a default point locator when none is specified.
void ClipRectilinearGridData(vtkDataSet *inputGrd, vtkDataArray *clipAray, double isoValue, vtkUnstructuredGrid *outputUG)
This function clips a vtkRectilinearGrid based on a specified iso-value (isoValue) using a scalar poi...
void ClipImageData(vtkDataSet *inputGrd, vtkDataArray *clipAray, double isoValue, vtkUnstructuredGrid *outputUG)
This function takes a vtkImageData as a vtkRectilinearGrid, which is then clipped by ClipRectilinearG...
vtkCallbackCommand * InternalProgressObserver
Superclass for algorithms that produce only unstructured grid as output.
dataset represents arbitrary combinations of all possible cell types
@ info
Definition: vtkX3D.h:382
@ port
Definition: vtkX3D.h:453
int vtkTypeBool
Definition: vtkABI.h:69
vtkTypeUInt32 vtkMTimeType
Definition: vtkType.h:287