VTK
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 *
19 * Copyright (c) 2000 - 2009, Lawrence Livermore National Security, LLC
20 * Produced at the Lawrence Livermore National Laboratory
21 * LLNL-CODE-400124
22 * All rights reserved.
23 *
24 * This file was adapted from the VisIt clipper (vtkVisItClipper). For details,
25 * see https://visit.llnl.gov/. The full copyright notice is contained in the
26 * file COPYRIGHT located at the root of the VisIt distribution or at
27 * http://www.llnl.gov/visit/copyright.html.
28 *
29 *****************************************************************************/
30 
31 
96 #ifndef vtkTableBasedClipDataSet_h
97 #define vtkTableBasedClipDataSet_h
98 
99 #include "vtkFiltersGeneralModule.h" // For export macro
101 
102 class vtkCallbackCommand;
103 class vtkImplicitFunction;
105 
106 class VTKFILTERSGENERAL_EXPORT vtkTableBasedClipDataSet : public vtkUnstructuredGridAlgorithm
107 {
108 public:
110  void PrintSelf( ostream & os, vtkIndent indent ) override;
111 
117 
121  vtkMTimeType GetMTime() override;
122 
124 
131  vtkSetMacro( InsideOut, vtkTypeBool );
132  vtkGetMacro( InsideOut, vtkTypeBool );
133  vtkBooleanMacro( InsideOut, vtkTypeBool );
135 
137 
143  vtkSetMacro( Value, double );
144  vtkGetMacro( Value, double );
146 
148 
153  vtkSetMacro( UseValueAsOffset, bool );
154  vtkGetMacro( UseValueAsOffset, bool );
155  vtkBooleanMacro( UseValueAsOffset, bool );
157 
159 
165  vtkGetObjectMacro( ClipFunction, vtkImplicitFunction );
167 
169 
175  vtkSetMacro( GenerateClipScalars, vtkTypeBool );
176  vtkGetMacro( GenerateClipScalars, vtkTypeBool );
177  vtkBooleanMacro( GenerateClipScalars, vtkTypeBool );
179 
181 
190  vtkGetObjectMacro( Locator, vtkIncrementalPointLocator );
192 
194 
199  vtkSetClampMacro( MergeTolerance, double, 0.0001, 0.25 );
200  vtkGetMacro( MergeTolerance, double );
202 
208 
210 
214  vtkSetMacro( GenerateClippedOutput, vtkTypeBool );
215  vtkGetMacro( GenerateClippedOutput, vtkTypeBool );
216  vtkBooleanMacro( GenerateClippedOutput, vtkTypeBool );
218 
223 
225 
230  vtkSetClampMacro(OutputPointsPrecision, int, SINGLE_PRECISION, DEFAULT_PRECISION);
231  vtkGetMacro(OutputPointsPrecision, int);
233 
234 protected:
237 
241 
247  void ClipDataSet( vtkDataSet * pDataSet,
248  vtkDataArray * clipAray, vtkUnstructuredGrid * unstruct );
249 
254  void ClipImageData( vtkDataSet * inputGrd, vtkDataArray * clipAray,
255  double isoValue, vtkUnstructuredGrid * outputUG );
256 
263  void ClipPolyData( vtkDataSet * inputGrd, vtkDataArray * clipAray,
264  double isoValue, vtkUnstructuredGrid * outputUG );
265 
272  void ClipRectilinearGridData( vtkDataSet * inputGrd, vtkDataArray * clipAray,
273  double isoValue, vtkUnstructuredGrid * outputUG );
274 
281  void ClipStructuredGridData( vtkDataSet * inputGrd, vtkDataArray * clipAray,
282  double isoValue, vtkUnstructuredGrid * outputUG );
283 
290  void ClipUnstructuredGridData( vtkDataSet * inputGrd, vtkDataArray * clipAray,
291  double isoValue, vtkUnstructuredGrid * outputUG );
292 
293 
297  static void InternalProgressCallbackFunction( vtkObject *, unsigned long,
298  void * clientdata, void * );
299 
304 
305 
310  double Value;
315 
317 
318 private:
320  void operator= ( const vtkTableBasedClipDataSet & ) = delete;
321 };
322 
323 #endif
vtkTableBasedClipDataSet::New
static vtkTableBasedClipDataSet * New()
Create an instance with a user-specified implicit function, turning off IVARs InsideOut and GenerateC...
vtkTableBasedClipDataSet::InsideOut
vtkTypeBool InsideOut
Definition: vtkTableBasedClipDataSet.h:306
vtkTableBasedClipDataSet::GetClippedOutput
vtkUnstructuredGrid * GetClippedOutput()
Return the clipped output.
vtkTableBasedClipDataSet::InternalProgressCallbackFunction
static void InternalProgressCallbackFunction(vtkObject *, unsigned long, void *clientdata, void *)
Register a callback function with the InternalProgressObserver.
vtkTableBasedClipDataSet::GenerateClipScalars
vtkTypeBool GenerateClipScalars
Definition: vtkTableBasedClipDataSet.h:307
vtkTableBasedClipDataSet::GetMTime
vtkMTimeType GetMTime() override
Get the MTime for which the point locator and clip function are considered.
vtkAlgorithm
Superclass for all sources, filters, and sinks in VTK.
Definition: vtkAlgorithm.h:60
vtkTableBasedClipDataSet::FillInputPortInformation
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
vtkTableBasedClipDataSet::ClipPolyData
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...
vtkInformationVector
Store zero or more vtkInformation instances.
Definition: vtkInformationVector.h:42
vtkTableBasedClipDataSet::ClipDataSet
void ClipDataSet(vtkDataSet *pDataSet, vtkDataArray *clipAray, vtkUnstructuredGrid *unstruct)
This function resorts to the sibling class vtkClipDataSet to handle special grids (such as cylinders ...
vtkObject
abstract base class for most VTK objects
Definition: vtkObject.h:60
vtkTableBasedClipDataSet::Value
double Value
Definition: vtkTableBasedClipDataSet.h:310
vtkDataArray
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:55
vtkTableBasedClipDataSet::OutputPointsPrecision
int OutputPointsPrecision
Definition: vtkTableBasedClipDataSet.h:316
vtkTableBasedClipDataSet::InternalProgressObserver
vtkCallbackCommand * InternalProgressObserver
Definition: vtkTableBasedClipDataSet.h:312
vtkTableBasedClipDataSet::~vtkTableBasedClipDataSet
~vtkTableBasedClipDataSet() override
vtkTableBasedClipDataSet::InternalProgressCallback
void InternalProgressCallback(vtkAlgorithm *algorithm)
The actual operation executed by the callback function.
vtkImplicitFunction
abstract interface for implicit functions
Definition: vtkImplicitFunction.h:61
vtkTableBasedClipDataSet::ClipRectilinearGridData
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...
vtkTableBasedClipDataSet::ClipFunction
vtkImplicitFunction * ClipFunction
Definition: vtkTableBasedClipDataSet.h:313
vtkX3D::port
@ port
Definition: vtkX3D.h:447
vtkTableBasedClipDataSet::SetClipFunction
virtual void SetClipFunction(vtkImplicitFunction *)
Set/Get the implicit function with which to perform the clipping operation.
vtkTableBasedClipDataSet::ClipStructuredGridData
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::SetLocator
void SetLocator(vtkIncrementalPointLocator *locator)
Set/Get a point locator locator for merging duplicate points.
vtkTableBasedClipDataSet::vtkTableBasedClipDataSet
vtkTableBasedClipDataSet(vtkImplicitFunction *cf=nullptr)
vtkIndent
a simple class to control print indentation
Definition: vtkIndent.h:40
vtkTableBasedClipDataSet::CreateDefaultLocator
void CreateDefaultLocator()
Create a default point locator when none is specified.
vtkTableBasedClipDataSet::GenerateClippedOutput
vtkTypeBool GenerateClippedOutput
Definition: vtkTableBasedClipDataSet.h:308
vtkIncrementalPointLocator
Abstract class in support of both point location and point insertion.
Definition: vtkIncrementalPointLocator.h:52
vtkTableBasedClipDataSet::Locator
vtkIncrementalPointLocator * Locator
Definition: vtkTableBasedClipDataSet.h:314
vtkTableBasedClipDataSet::RequestData
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
vtkDataSet
abstract class to specify dataset behavior
Definition: vtkDataSet.h:63
vtkInformation
Store vtkAlgorithm input/output information.
Definition: vtkInformation.h:87
vtkX3D::info
@ info
Definition: vtkX3D.h:376
vtkTableBasedClipDataSet
Clip any dataset with a user-specified implicit function or an input scalar point data array.
Definition: vtkTableBasedClipDataSet.h:107
vtkTableBasedClipDataSet::PrintSelf
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkCallbackCommand
supports function callbacks
Definition: vtkCallbackCommand.h:51
vtkUnstructuredGridAlgorithm.h
vtkUnstructuredGridAlgorithm
Superclass for algorithms that produce only unstructured grid as output.
Definition: vtkUnstructuredGridAlgorithm.h:41
vtkUnstructuredGrid
dataset represents arbitrary combinations of all possible cell types
Definition: vtkUnstructuredGrid.h:89
vtkTableBasedClipDataSet::ClipImageData
void ClipImageData(vtkDataSet *inputGrd, vtkDataArray *clipAray, double isoValue, vtkUnstructuredGrid *outputUG)
This function takes a vtkImageData as a vtkRectilinearGrid, which is then clipped by ClipRectilinearG...
vtkTableBasedClipDataSet::UseValueAsOffset
bool UseValueAsOffset
Definition: vtkTableBasedClipDataSet.h:309
vtkTableBasedClipDataSet::ClipUnstructuredGridData
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...
vtkTableBasedClipDataSet::MergeTolerance
double MergeTolerance
Definition: vtkTableBasedClipDataSet.h:311
vtkTypeBool
int vtkTypeBool
Definition: vtkABI.h:69
vtkMTimeType
vtkTypeUInt32 vtkMTimeType
Definition: vtkType.h:302