VTK  9.3.0
vtkProbeFilter.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
59#ifndef vtkProbeFilter_h
60#define vtkProbeFilter_h
61
62#include "vtkDataSetAlgorithm.h"
63#include "vtkDataSetAttributes.h" // needed for vtkDataSetAttributes::FieldList
64#include "vtkFiltersCoreModule.h" // For export macro
65
66#include <vector> // For std::vector
67
68VTK_ABI_NAMESPACE_BEGIN
70class vtkCharArray;
71class vtkGenericCell;
72class vtkIdTypeArray;
73class vtkImageData;
74class vtkPointData;
76
77class VTKFILTERSCORE_EXPORT vtkProbeFilter : public vtkDataSetAlgorithm
78{
79public:
82 void PrintSelf(ostream& os, vtkIndent indent) override;
83
85
94
102
104
109 vtkSetMacro(CategoricalData, vtkTypeBool);
110 vtkGetMacro(CategoricalData, vtkTypeBool);
111 vtkBooleanMacro(CategoricalData, vtkTypeBool);
113
115
125 vtkSetMacro(SpatialMatch, vtkTypeBool);
126 vtkGetMacro(SpatialMatch, vtkTypeBool);
127 vtkBooleanMacro(SpatialMatch, vtkTypeBool);
129
131
137
139
144 vtkSetStringMacro(ValidPointMaskArrayName);
145 vtkGetStringMacro(ValidPointMaskArrayName);
147
149
153 vtkSetMacro(PassCellArrays, vtkTypeBool);
154 vtkBooleanMacro(PassCellArrays, vtkTypeBool);
155 vtkGetMacro(PassCellArrays, vtkTypeBool);
158
162 vtkSetMacro(PassPointArrays, vtkTypeBool);
163 vtkBooleanMacro(PassPointArrays, vtkTypeBool);
164 vtkGetMacro(PassPointArrays, vtkTypeBool);
166
168
172 vtkSetMacro(PassFieldArrays, vtkTypeBool);
173 vtkBooleanMacro(PassFieldArrays, vtkTypeBool);
174 vtkGetMacro(PassFieldArrays, vtkTypeBool);
176
178
183 vtkSetMacro(Tolerance, double);
184 vtkGetMacro(Tolerance, double);
186
188
196 vtkSetMacro(SnapToCellWithClosestPoint, bool);
197 vtkBooleanMacro(SnapToCellWithClosestPoint, bool);
198 vtkGetMacro(SnapToCellWithClosestPoint, bool);
200
202
207 vtkSetMacro(ComputeTolerance, bool);
208 vtkBooleanMacro(ComputeTolerance, bool);
209 vtkGetMacro(ComputeTolerance, bool);
211
213
220 vtkGetObjectMacro(FindCellStrategy, vtkFindCellStrategy);
222
224
233 vtkGetObjectMacro(CellLocatorPrototype, vtkAbstractCellLocator);
235
236protected:
238 ~vtkProbeFilter() override;
239
243
249
253 void Probe(vtkDataSet* input, vtkDataSet* source, vtkDataSet* output);
254
260
264 virtual void InitializeForProbing(vtkDataSet* input, vtkDataSet* output);
266 virtual void InitializeOutputArrays(vtkPointData* outPD, vtkIdType numPts);
267
272 void DoProbing(vtkDataSet* input, int srcIdx, vtkDataSet* source, vtkDataSet* output);
273
275
279
281
282 double Tolerance;
285
289
290 // Support various methods to support the FindCell() operation
293
296
297private:
298 vtkProbeFilter(const vtkProbeFilter&) = delete;
299 void operator=(const vtkProbeFilter&) = delete;
300
301 // Probe only those points that are marked as not-probed by the MaskPoints
302 // array.
303 void ProbeEmptyPoints(vtkDataSet* input, int srcIdx, vtkDataSet* source, vtkDataSet* output);
304
305 // A faster implementation for vtkImageData input.
306 void ProbePointsImageData(
307 vtkImageData* input, int srcIdx, vtkDataSet* source, vtkImageData* output);
308 void ProbeImagePointsInCell(vtkGenericCell* cell, vtkIdType cellId, vtkDataSet* source,
309 int srcBlockId, const double start[3], const double spacing[3], const int dim[3],
310 vtkPointData* outPD, char* maskArray, double* wtsBuff);
311
312 class ProbeImageDataWorklet;
313
314 // A faster implementation for vtkImageData source.
315 void ProbeImageDataPoints(
316 vtkDataSet* input, int srcIdx, vtkImageData* sourceImage, vtkDataSet* output);
317 void ProbeImageDataPointsSMP(vtkDataSet* input, vtkImageData* source, int srcIdx,
318 vtkPointData* outPD, char* maskArray, vtkIdList* pointIds, vtkIdType startId, vtkIdType endId,
319 bool baseThread);
320
321 class ProbeImageDataPointsWorklet;
322
323 class ProbeEmptyPointsWorklet;
324
325 std::vector<vtkDataArray*> InputCellArrays;
326 std::vector<vtkDataArray*> SourceCellArrays;
327};
328
329VTK_ABI_NAMESPACE_END
330#endif
an abstract base class for locators which find cells
Proxy object to connect input/output ports.
dynamic, self-adjusting array of char
general representation of visualization data
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:53
helper class to manage the vtkPointSet::FindCell() METHOD
provides thread-safe access to cells
list of point or cell ids
Definition vtkIdList.h:23
dynamic, self-adjusting array of vtkIdType
topologically and geometrically regular array of data
a simple class to control print indentation
Definition vtkIndent.h:29
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
represent and manipulate point attribute data
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.
bool SnapToCellWithClosestPoint
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
virtual void InitializeSourceArrays(vtkDataSet *source)
int vtkTypeBool
Definition vtkABI.h:64
boost::graph_traits< vtkGraph * >::vertex_descriptor source(boost::graph_traits< vtkGraph * >::edge_descriptor e, vtkGraph *)
int vtkIdType
Definition vtkType.h:315