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
68#ifndef vtkProbeFilter_h
69#define vtkProbeFilter_h
70
71#include "vtkDataSetAlgorithm.h"
72#include "vtkDataSetAttributes.h" // needed for vtkDataSetAttributes::FieldList
73#include "vtkFiltersCoreModule.h" // For export macro
74
75#include <vector> // For std::vector
76
77VTK_ABI_NAMESPACE_BEGIN
79class vtkCharArray;
80class vtkGenericCell;
81class vtkIdTypeArray;
82class vtkImageData;
83class vtkPointData;
85
86class VTKFILTERSCORE_EXPORT vtkProbeFilter : public vtkDataSetAlgorithm
87{
88public:
91 void PrintSelf(ostream& os, vtkIndent indent) override;
92
94
103
111
113
118 vtkSetMacro(CategoricalData, vtkTypeBool);
119 vtkGetMacro(CategoricalData, vtkTypeBool);
120 vtkBooleanMacro(CategoricalData, vtkTypeBool);
122
124
134 vtkSetMacro(SpatialMatch, vtkTypeBool);
135 vtkGetMacro(SpatialMatch, vtkTypeBool);
136 vtkBooleanMacro(SpatialMatch, vtkTypeBool);
138
140
146
148
153 vtkSetStringMacro(ValidPointMaskArrayName);
154 vtkGetStringMacro(ValidPointMaskArrayName);
156
158
162 vtkSetMacro(PassCellArrays, vtkTypeBool);
163 vtkBooleanMacro(PassCellArrays, vtkTypeBool);
164 vtkGetMacro(PassCellArrays, vtkTypeBool);
167
171 vtkSetMacro(PassPointArrays, vtkTypeBool);
172 vtkBooleanMacro(PassPointArrays, vtkTypeBool);
173 vtkGetMacro(PassPointArrays, vtkTypeBool);
175
177
181 vtkSetMacro(PassFieldArrays, vtkTypeBool);
182 vtkBooleanMacro(PassFieldArrays, vtkTypeBool);
183 vtkGetMacro(PassFieldArrays, vtkTypeBool);
185
187
192 vtkSetMacro(Tolerance, double);
193 vtkGetMacro(Tolerance, double);
195
197
205 vtkSetMacro(SnapToCellWithClosestPoint, bool);
206 vtkBooleanMacro(SnapToCellWithClosestPoint, bool);
207 vtkGetMacro(SnapToCellWithClosestPoint, bool);
209
211
216 vtkSetMacro(ComputeTolerance, bool);
217 vtkBooleanMacro(ComputeTolerance, bool);
218 vtkGetMacro(ComputeTolerance, bool);
220
222
229 vtkGetObjectMacro(FindCellStrategy, vtkFindCellStrategy);
231
233
242 vtkGetObjectMacro(CellLocatorPrototype, vtkAbstractCellLocator);
244
245protected:
247 ~vtkProbeFilter() override;
248
252
258
262 void Probe(vtkDataSet* input, vtkDataSet* source, vtkDataSet* output);
263
269
273 virtual void InitializeForProbing(vtkDataSet* input, vtkDataSet* output);
275 virtual void InitializeOutputArrays(vtkPointData* outPD, vtkIdType numPts);
276
281 void DoProbing(vtkDataSet* input, int srcIdx, vtkDataSet* source, vtkDataSet* output);
282
284
288
290
291 double Tolerance;
294
298
299 // Support various methods to support the FindCell() operation
302
305
306private:
307 vtkProbeFilter(const vtkProbeFilter&) = delete;
308 void operator=(const vtkProbeFilter&) = delete;
309
310 // Probe only those points that are marked as not-probed by the MaskPoints
311 // array.
312 void ProbeEmptyPoints(vtkDataSet* input, int srcIdx, vtkDataSet* source, vtkDataSet* output);
313
314 // A faster implementation for vtkImageData input.
315 void ProbePointsImageData(
316 vtkImageData* input, int srcIdx, vtkDataSet* source, vtkImageData* output);
317 void ProbeImagePointsInCell(vtkGenericCell* cell, vtkIdType cellId, vtkDataSet* source,
318 int srcBlockId, const double start[3], const double spacing[3], const int dim[3],
319 vtkPointData* outPD, char* maskArray, double* wtsBuff);
320
321 class ProbeImageDataWorklet;
322
323 // A faster implementation for vtkImageData source.
324 void ProbeImageDataPoints(
325 vtkDataSet* input, int srcIdx, vtkImageData* sourceImage, vtkDataSet* output);
326 void ProbeImageDataPointsSMP(vtkDataSet* input, vtkImageData* source, int srcIdx,
327 vtkPointData* outPD, char* maskArray, vtkIdList* pointIds, vtkIdType startId, vtkIdType endId,
328 bool baseThread);
329
330 class ProbeImageDataPointsWorklet;
331
332 class ProbeEmptyPointsWorklet;
333
334 std::vector<vtkDataArray*> InputCellArrays;
335 std::vector<vtkDataArray*> SourceCellArrays;
336};
337
338VTK_ABI_NAMESPACE_END
339#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:62
helper class to manage the vtkPointSet::FindCell() METHOD
provides thread-safe access to cells
list of point or cell ids
Definition vtkIdList.h:32
dynamic, self-adjusting array of vtkIdType
topologically and geometrically regular array of data
a simple class to control print indentation
Definition vtkIndent.h:38
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