VTK  9.6.1
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
58
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#include "vtkWrappingHints.h" // For VTK_MARSHALAUTO
66
67#include <vector> // For std::vector
68
69VTK_ABI_NAMESPACE_BEGIN
71class vtkCharArray;
72class vtkGenericCell;
73class vtkIdTypeArray;
74class vtkImageData;
75class vtkPointData;
77
78class VTKFILTERSCORE_EXPORT VTK_MARSHALAUTO vtkProbeFilter : public vtkDataSetAlgorithm
79{
80public:
83 void PrintSelf(ostream& os, vtkIndent indent) override;
84
86
95
103
105
112 vtkBooleanMacro(CategoricalData, vtkTypeBool);
114
116
128 vtkBooleanMacro(SpatialMatch, vtkTypeBool);
130
132
138
140
145 vtkSetStringMacro(ValidPointMaskArrayName);
146 vtkGetStringMacro(ValidPointMaskArrayName);
148
150
155 vtkBooleanMacro(PassCellArrays, vtkTypeBool);
159
164 vtkBooleanMacro(PassPointArrays, vtkTypeBool);
167
169
174 vtkBooleanMacro(PassFieldArrays, vtkTypeBool);
177
179
184 vtkSetMacro(Tolerance, double);
185 vtkGetMacro(Tolerance, double);
187
189
197 vtkSetMacro(SnapToCellWithClosestPoint, bool);
198 vtkBooleanMacro(SnapToCellWithClosestPoint, bool);
199 vtkGetMacro(SnapToCellWithClosestPoint, bool);
201
203
208 vtkSetMacro(ComputeTolerance, bool);
209 vtkBooleanMacro(ComputeTolerance, bool);
210 vtkGetMacro(ComputeTolerance, bool);
212
214
223
225
236
237protected:
239 ~vtkProbeFilter() override;
240
244
250
254 void Probe(vtkDataSet* input, vtkDataSet* source, vtkDataSet* output);
255
261
265 virtual void InitializeForProbing(vtkDataSet* input, vtkDataSet* output);
267 virtual void InitializeOutputArrays(vtkPointData* outPD, vtkIdType numPts);
268
273 void DoProbing(vtkDataSet* input, int srcIdx, vtkDataSet* source, vtkDataSet* output);
274
276
280
282
283 double Tolerance;
286
290
291 // Support various methods to support the FindCell() operation
294
297
298private:
299 vtkProbeFilter(const vtkProbeFilter&) = delete;
300 void operator=(const vtkProbeFilter&) = delete;
301
302 // Probe only those points that are marked as not-probed by the MaskPoints
303 // array.
304 void ProbeEmptyPoints(vtkDataSet* input, int srcIdx, vtkDataSet* source, vtkDataSet* output);
305
306 // A faster implementation for vtkImageData input.
307 void ProbePointsImageData(
308 vtkImageData* input, int srcIdx, vtkDataSet* source, vtkImageData* output);
309 void ProbeImagePointsInCell(vtkGenericCell* cell, vtkIdType cellId, vtkDataSet* source,
310 int srcBlockId, vtkImageData* input, 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
vtkDataSetAttributesFieldList FieldList
abstract class to specify dataset behavior
Definition vtkDataSet.h:57
helper class to manage the vtkPointSet::FindCell() method
provides thread-safe access to cells
list of point or cell ids
Definition vtkIdList.h:24
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
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:368
#define VTK_MARSHALAUTO