VTK  9.1.0
vtkProbeFilter.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkProbeFilter.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=========================================================================*/
74#ifndef vtkProbeFilter_h
75#define vtkProbeFilter_h
76
77#include "vtkDataSetAlgorithm.h"
78#include "vtkDataSetAttributes.h" // needed for vtkDataSetAttributes::FieldList
79#include "vtkFiltersCoreModule.h" // For export macro
80
82class vtkCell;
83class vtkCharArray;
84class vtkIdTypeArray;
85class vtkImageData;
86class vtkPointData;
88
89class VTKFILTERSCORE_EXPORT vtkProbeFilter : public vtkDataSetAlgorithm
90{
91public:
94 void PrintSelf(ostream& os, vtkIndent indent) override;
95
97
106
114
116
121 vtkSetMacro(CategoricalData, vtkTypeBool);
122 vtkGetMacro(CategoricalData, vtkTypeBool);
123 vtkBooleanMacro(CategoricalData, vtkTypeBool);
125
127
137 vtkSetMacro(SpatialMatch, vtkTypeBool);
138 vtkGetMacro(SpatialMatch, vtkTypeBool);
139 vtkBooleanMacro(SpatialMatch, vtkTypeBool);
141
143
149
151
156 vtkSetStringMacro(ValidPointMaskArrayName);
157 vtkGetStringMacro(ValidPointMaskArrayName);
159
161
165 vtkSetMacro(PassCellArrays, vtkTypeBool);
166 vtkBooleanMacro(PassCellArrays, vtkTypeBool);
167 vtkGetMacro(PassCellArrays, vtkTypeBool);
170
174 vtkSetMacro(PassPointArrays, vtkTypeBool);
175 vtkBooleanMacro(PassPointArrays, vtkTypeBool);
176 vtkGetMacro(PassPointArrays, vtkTypeBool);
178
180
184 vtkSetMacro(PassFieldArrays, vtkTypeBool);
185 vtkBooleanMacro(PassFieldArrays, vtkTypeBool);
186 vtkGetMacro(PassFieldArrays, vtkTypeBool);
188
190
195 vtkSetMacro(Tolerance, double);
196 vtkGetMacro(Tolerance, double);
198
200
205 vtkSetMacro(ComputeTolerance, bool);
206 vtkBooleanMacro(ComputeTolerance, bool);
207 vtkGetMacro(ComputeTolerance, bool);
209
211
218 vtkGetObjectMacro(FindCellStrategy, vtkFindCellStrategy);
220
222
231 vtkGetObjectMacro(CellLocatorPrototype, vtkAbstractCellLocator);
233
234protected:
236 ~vtkProbeFilter() override;
237
241
247
251 void Probe(vtkDataSet* input, vtkDataSet* source, vtkDataSet* output);
252
258
262 virtual void InitializeForProbing(vtkDataSet* input, vtkDataSet* output);
263 virtual void InitializeOutputArrays(vtkPointData* outPD, vtkIdType numPts);
264
269 void DoProbing(vtkDataSet* input, int srcIdx, vtkDataSet* source, vtkDataSet* output);
270
272
276
278
279 double Tolerance;
281
285
286 // Support various methods to support the FindCell() operation
289
292
293private:
294 vtkProbeFilter(const vtkProbeFilter&) = delete;
295 void operator=(const vtkProbeFilter&) = delete;
296
297 // Probe only those points that are marked as not-probed by the MaskPoints
298 // array.
299 void ProbeEmptyPoints(vtkDataSet* input, int srcIdx, vtkDataSet* source, vtkDataSet* output);
300
301 // A faster implementation for vtkImageData input.
302 void ProbePointsImageData(
303 vtkImageData* input, int srcIdx, vtkDataSet* source, vtkImageData* output);
304 void ProbeImagePointsInCell(vtkCell* cell, vtkIdType cellId, vtkDataSet* source, int srcBlockId,
305 const double start[3], const double spacing[3], const int dim[3], vtkPointData* outPD,
306 char* maskArray, double* wtsBuff);
307
308 class ProbeImageDataWorklet;
309
310 // A faster implementation for vtkImageData source.
311 void ProbeImageDataPoints(
312 vtkDataSet* input, int srcIdx, vtkImageData* sourceImage, vtkDataSet* output);
313 void ProbeImageDataPointsSMP(vtkDataSet* input, vtkImageData* source, int srcIdx,
314 vtkPointData* outPD, char* maskArray, vtkIdList* pointIds, vtkIdType startId, vtkIdType endId,
315 bool baseThread);
316
317 class ProbeImageDataPointsWorklet;
318
319 class vtkVectorOfArrays;
320 vtkVectorOfArrays* CellArrays;
321};
322
323#endif
an abstract base class for locators which find cells
Proxy object to connect input/output ports.
abstract class to specify cell behavior
Definition: vtkCell.h:67
dynamic, self-adjusting array of char
Definition: vtkCharArray.h:58
general representation of visualization data
Definition: vtkDataObject.h:69
Superclass for algorithms that produce output of the same type as input.
helps manage arrays from multiple vtkDataSetAttributes.
abstract class to specify dataset behavior
Definition: vtkDataSet.h:66
helper class to manage the vtkPointSet::FindCell() METHOD
list of point or cell ids
Definition: vtkIdList.h:40
dynamic, self-adjusting array of vtkIdType
topologically and geometrically regular array of data
Definition: vtkImageData.h:57
a simple class to control print indentation
Definition: vtkIndent.h:43
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
represent and manipulate point attribute data
Definition: vtkPointData.h:42
sample data values at specified point locations
virtual void InitializeForProbing(vtkDataSet *input, vtkDataSet *output)
Initializes output and various arrays which keep track for probing status.
void SetSourceData(vtkDataObject *source)
Specify the data set that will be probed at the input points.
vtkFindCellStrategy * FindCellStrategy
vtkIdTypeArray * ValidPoints
virtual void SetCellLocatorPrototype(vtkAbstractCellLocator *)
Set/Get the prototype cell locator to perform the FindCell() operation.
virtual void SetFindCellStrategy(vtkFindCellStrategy *)
Set / get the strategy used to perform the FindCell() operation.
virtual void InitializeOutputArrays(vtkPointData *outPD, vtkIdType numPts)
vtkTypeBool CategoricalData
void Probe(vtkDataSet *input, vtkDataSet *source, vtkDataSet *output)
Equivalent to calling BuildFieldList(); InitializeForProbing(); DoProbing().
~vtkProbeFilter() override
int RequestUpdateExtent(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called within ProcessRequest when each filter in the pipeline decides what portion of its inp...
vtkCharArray * MaskPoints
char * ValidPointMaskArrayName
vtkTypeBool SpatialMatch
vtkDataSetAttributes::FieldList * CellList
vtkIdTypeArray * GetValidPoints()
Get the list of point ids in the output that contain attribute data interpolated from the source.
vtkTypeBool PassCellArrays
vtkDataSetAttributes::FieldList * PointList
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void SetSourceConnection(vtkAlgorithmOutput *algOutput)
Specify the data set that will be probed at the input points.
void PassAttributeData(vtkDataSet *input, vtkDataObject *source, vtkDataSet *output)
Call at end of RequestData() to pass attribute data respecting the PassCellArrays,...
void BuildFieldList(vtkDataSet *source)
Build the field lists.
vtkAbstractCellLocator * CellLocatorPrototype
int RequestInformation(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called within ProcessRequest when a request asks for Information.
static vtkProbeFilter * New()
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called within ProcessRequest when a request asks the algorithm to do its work.
void DoProbing(vtkDataSet *input, int srcIdx, vtkDataSet *source, vtkDataSet *output)
Probe appropriate points srcIdx is the index in the PointList for the given source.
vtkDataObject * GetSource()
Specify the data set that will be probed at the input points.
vtkTypeBool PassPointArrays
vtkTypeBool PassFieldArrays
@ spacing
Definition: vtkX3D.h:487
int vtkTypeBool
Definition: vtkABI.h:69
boost::graph_traits< vtkGraph * >::vertex_descriptor source(boost::graph_traits< vtkGraph * >::edge_descriptor e, vtkGraph *)
int vtkIdType
Definition: vtkType.h:332