VTK  9.2.6
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=========================================================================*/
71#ifndef vtkProbeFilter_h
72#define vtkProbeFilter_h
73
74#include "vtkDataSetAlgorithm.h"
75#include "vtkDataSetAttributes.h" // needed for vtkDataSetAttributes::FieldList
76#include "vtkFiltersCoreModule.h" // For export macro
77
79class vtkCell;
80class vtkCharArray;
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
202 vtkSetMacro(ComputeTolerance, bool);
203 vtkBooleanMacro(ComputeTolerance, bool);
204 vtkGetMacro(ComputeTolerance, bool);
206
208
215 vtkGetObjectMacro(FindCellStrategy, vtkFindCellStrategy);
217
219
228 vtkGetObjectMacro(CellLocatorPrototype, vtkAbstractCellLocator);
230
231protected:
233 ~vtkProbeFilter() override;
234
238
244
248 void Probe(vtkDataSet* input, vtkDataSet* source, vtkDataSet* output);
249
255
259 virtual void InitializeForProbing(vtkDataSet* input, vtkDataSet* output);
260 virtual void InitializeOutputArrays(vtkPointData* outPD, vtkIdType numPts);
261
266 void DoProbing(vtkDataSet* input, int srcIdx, vtkDataSet* source, vtkDataSet* output);
267
269
273
275
276 double Tolerance;
278
282
283 // Support various methods to support the FindCell() operation
286
289
290private:
291 vtkProbeFilter(const vtkProbeFilter&) = delete;
292 void operator=(const vtkProbeFilter&) = delete;
293
294 // Probe only those points that are marked as not-probed by the MaskPoints
295 // array.
296 void ProbeEmptyPoints(vtkDataSet* input, int srcIdx, vtkDataSet* source, vtkDataSet* output);
297
298 // A faster implementation for vtkImageData input.
299 void ProbePointsImageData(
300 vtkImageData* input, int srcIdx, vtkDataSet* source, vtkImageData* output);
301 void ProbeImagePointsInCell(vtkCell* cell, vtkIdType cellId, vtkDataSet* source, int srcBlockId,
302 const double start[3], const double spacing[3], const int dim[3], vtkPointData* outPD,
303 char* maskArray, double* wtsBuff);
304
305 class ProbeImageDataWorklet;
306
307 // A faster implementation for vtkImageData source.
308 void ProbeImageDataPoints(
309 vtkDataSet* input, int srcIdx, vtkImageData* sourceImage, vtkDataSet* output);
310 void ProbeImageDataPointsSMP(vtkDataSet* input, vtkImageData* source, int srcIdx,
311 vtkPointData* outPD, char* maskArray, vtkIdList* pointIds, vtkIdType startId, vtkIdType endId,
312 bool baseThread);
313
314 class ProbeImageDataPointsWorklet;
315
316 class vtkVectorOfArrays;
317 vtkVectorOfArrays* CellArrays;
318};
319
320#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:61
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:63
helper class to manage the vtkPointSet::FindCell() METHOD
list of point or cell ids
Definition vtkIdList.h:34
dynamic, self-adjusting array of vtkIdType
topologically and geometrically regular array of data
a simple class to control print indentation
Definition vtkIndent.h:40
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.
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
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