VTK  9.2.5
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 =========================================================================*/
180 #ifndef vtkProbeFilter_h
181 #define vtkProbeFilter_h
182 
183 #include "vtkDataSetAlgorithm.h"
184 #include "vtkDataSetAttributes.h" // needed for vtkDataSetAttributes::FieldList
185 #include "vtkFiltersCoreModule.h" // For export macro
186 
188 class vtkCell;
189 class vtkCharArray;
190 class vtkIdTypeArray;
191 class vtkImageData;
192 class vtkPointData;
193 class vtkFindCellStrategy;
194 
195 class VTKFILTERSCORE_EXPORT vtkProbeFilter : public vtkDataSetAlgorithm
196 {
197 public:
198  static vtkProbeFilter* New();
200  void PrintSelf(ostream& os, vtkIndent indent) override;
201 
203 
212 
220 
222 
227  vtkSetMacro(CategoricalData, vtkTypeBool);
228  vtkGetMacro(CategoricalData, vtkTypeBool);
229  vtkBooleanMacro(CategoricalData, vtkTypeBool);
231 
233 
243  vtkSetMacro(SpatialMatch, vtkTypeBool);
244  vtkGetMacro(SpatialMatch, vtkTypeBool);
245  vtkBooleanMacro(SpatialMatch, vtkTypeBool);
247 
249 
255 
257 
262  vtkSetStringMacro(ValidPointMaskArrayName);
263  vtkGetStringMacro(ValidPointMaskArrayName);
265 
267 
271  vtkSetMacro(PassCellArrays, vtkTypeBool);
272  vtkBooleanMacro(PassCellArrays, vtkTypeBool);
273  vtkGetMacro(PassCellArrays, vtkTypeBool);
276 
280  vtkSetMacro(PassPointArrays, vtkTypeBool);
281  vtkBooleanMacro(PassPointArrays, vtkTypeBool);
282  vtkGetMacro(PassPointArrays, vtkTypeBool);
284 
286 
290  vtkSetMacro(PassFieldArrays, vtkTypeBool);
291  vtkBooleanMacro(PassFieldArrays, vtkTypeBool);
292  vtkGetMacro(PassFieldArrays, vtkTypeBool);
294 
296 
301  vtkSetMacro(Tolerance, double);
302  vtkGetMacro(Tolerance, double);
304 
306 
311  vtkSetMacro(ComputeTolerance, bool);
312  vtkBooleanMacro(ComputeTolerance, bool);
313  vtkGetMacro(ComputeTolerance, bool);
315 
317 
324  vtkGetObjectMacro(FindCellStrategy, vtkFindCellStrategy);
326 
328 
337  vtkGetObjectMacro(CellLocatorPrototype, vtkAbstractCellLocator);
339 
340 protected:
342  ~vtkProbeFilter() override;
343 
347 
353 
357  void Probe(vtkDataSet* input, vtkDataSet* source, vtkDataSet* output);
358 
364 
368  virtual void InitializeForProbing(vtkDataSet* input, vtkDataSet* output);
369  virtual void InitializeOutputArrays(vtkPointData* outPD, vtkIdType numPts);
370 
375  void DoProbing(vtkDataSet* input, int srcIdx, vtkDataSet* source, vtkDataSet* output);
376 
378 
382 
384 
385  double Tolerance;
387 
391 
392  // Support various methods to support the FindCell() operation
395 
398 
399 private:
400  vtkProbeFilter(const vtkProbeFilter&) = delete;
401  void operator=(const vtkProbeFilter&) = delete;
402 
403  // Probe only those points that are marked as not-probed by the MaskPoints
404  // array.
405  void ProbeEmptyPoints(vtkDataSet* input, int srcIdx, vtkDataSet* source, vtkDataSet* output);
406 
407  // A faster implementation for vtkImageData input.
408  void ProbePointsImageData(
409  vtkImageData* input, int srcIdx, vtkDataSet* source, vtkImageData* output);
410  void ProbeImagePointsInCell(vtkCell* cell, vtkIdType cellId, vtkDataSet* source, int srcBlockId,
411  const double start[3], const double spacing[3], const int dim[3], vtkPointData* outPD,
412  char* maskArray, double* wtsBuff);
413 
414  class ProbeImageDataWorklet;
415 
416  // A faster implementation for vtkImageData source.
417  void ProbeImageDataPoints(
418  vtkDataSet* input, int srcIdx, vtkImageData* sourceImage, vtkDataSet* output);
419  void ProbeImageDataPointsSMP(vtkDataSet* input, vtkImageData* source, int srcIdx,
420  vtkPointData* outPD, char* maskArray, vtkIdList* pointIds, vtkIdType startId, vtkIdType endId,
421  bool baseThread);
422 
423  class ProbeImageDataPointsWorklet;
424 
425  class vtkVectorOfArrays;
426  vtkVectorOfArrays* CellArrays;
427 };
428 
429 #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:150
dynamic, self-adjusting array of char
Definition: vtkCharArray.h:71
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:172
helper class to manage the vtkPointSet::FindCell() METHOD
list of point or cell ids
Definition: vtkIdList.h:143
dynamic, self-adjusting array of vtkIdType
topologically and geometrically regular array of data
Definition: vtkImageData.h:163
a simple class to control print indentation
Definition: vtkIndent.h:119
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
represent and manipulate point attribute data
Definition: vtkPointData.h:151
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.
vtkDataObject * GetSource()
Specify the data set that will be probed at the input points.
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
vtkIdTypeArray * GetValidPoints()
Get the list of point ids in the output that contain attribute data interpolated from the source.
vtkDataSetAttributes::FieldList * CellList
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.
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.
vtkTypeBool PassPointArrays
vtkTypeBool PassFieldArrays
static vtkProbeFilter * New()
@ 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