VTK  9.1.0
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 =========================================================================*/
174 #ifndef vtkProbeFilter_h
175 #define vtkProbeFilter_h
176 
177 #include "vtkDataSetAlgorithm.h"
178 #include "vtkDataSetAttributes.h" // needed for vtkDataSetAttributes::FieldList
179 #include "vtkFiltersCoreModule.h" // For export macro
180 
182 class vtkCell;
183 class vtkCharArray;
184 class vtkIdTypeArray;
185 class vtkImageData;
186 class vtkPointData;
187 class vtkFindCellStrategy;
188 
189 class VTKFILTERSCORE_EXPORT vtkProbeFilter : public vtkDataSetAlgorithm
190 {
191 public:
192  static vtkProbeFilter* New();
194  void PrintSelf(ostream& os, vtkIndent indent) override;
195 
197 
206 
214 
216 
221  vtkSetMacro(CategoricalData, vtkTypeBool);
222  vtkGetMacro(CategoricalData, vtkTypeBool);
223  vtkBooleanMacro(CategoricalData, vtkTypeBool);
225 
227 
237  vtkSetMacro(SpatialMatch, vtkTypeBool);
238  vtkGetMacro(SpatialMatch, vtkTypeBool);
239  vtkBooleanMacro(SpatialMatch, vtkTypeBool);
241 
243 
249 
251 
256  vtkSetStringMacro(ValidPointMaskArrayName);
257  vtkGetStringMacro(ValidPointMaskArrayName);
259 
261 
265  vtkSetMacro(PassCellArrays, vtkTypeBool);
266  vtkBooleanMacro(PassCellArrays, vtkTypeBool);
267  vtkGetMacro(PassCellArrays, vtkTypeBool);
270 
274  vtkSetMacro(PassPointArrays, vtkTypeBool);
275  vtkBooleanMacro(PassPointArrays, vtkTypeBool);
276  vtkGetMacro(PassPointArrays, vtkTypeBool);
278 
280 
284  vtkSetMacro(PassFieldArrays, vtkTypeBool);
285  vtkBooleanMacro(PassFieldArrays, vtkTypeBool);
286  vtkGetMacro(PassFieldArrays, vtkTypeBool);
288 
290 
295  vtkSetMacro(Tolerance, double);
296  vtkGetMacro(Tolerance, double);
298 
300 
305  vtkSetMacro(ComputeTolerance, bool);
306  vtkBooleanMacro(ComputeTolerance, bool);
307  vtkGetMacro(ComputeTolerance, bool);
309 
311 
318  vtkGetObjectMacro(FindCellStrategy, vtkFindCellStrategy);
320 
322 
331  vtkGetObjectMacro(CellLocatorPrototype, vtkAbstractCellLocator);
333 
334 protected:
336  ~vtkProbeFilter() override;
337 
341 
347 
351  void Probe(vtkDataSet* input, vtkDataSet* source, vtkDataSet* output);
352 
358 
362  virtual void InitializeForProbing(vtkDataSet* input, vtkDataSet* output);
363  virtual void InitializeOutputArrays(vtkPointData* outPD, vtkIdType numPts);
364 
369  void DoProbing(vtkDataSet* input, int srcIdx, vtkDataSet* source, vtkDataSet* output);
370 
372 
376 
378 
379  double Tolerance;
381 
385 
386  // Support various methods to support the FindCell() operation
389 
392 
393 private:
394  vtkProbeFilter(const vtkProbeFilter&) = delete;
395  void operator=(const vtkProbeFilter&) = delete;
396 
397  // Probe only those points that are marked as not-probed by the MaskPoints
398  // array.
399  void ProbeEmptyPoints(vtkDataSet* input, int srcIdx, vtkDataSet* source, vtkDataSet* output);
400 
401  // A faster implementation for vtkImageData input.
402  void ProbePointsImageData(
403  vtkImageData* input, int srcIdx, vtkDataSet* source, vtkImageData* output);
404  void ProbeImagePointsInCell(vtkCell* cell, vtkIdType cellId, vtkDataSet* source, int srcBlockId,
405  const double start[3], const double spacing[3], const int dim[3], vtkPointData* outPD,
406  char* maskArray, double* wtsBuff);
407 
408  class ProbeImageDataWorklet;
409 
410  // A faster implementation for vtkImageData source.
411  void ProbeImageDataPoints(
412  vtkDataSet* input, int srcIdx, vtkImageData* sourceImage, vtkDataSet* output);
413  void ProbeImageDataPointsSMP(vtkDataSet* input, vtkImageData* source, int srcIdx,
414  vtkPointData* outPD, char* maskArray, vtkIdList* pointIds, vtkIdType startId, vtkIdType endId,
415  bool baseThread);
416 
417  class ProbeImageDataPointsWorklet;
418 
419  class vtkVectorOfArrays;
420  vtkVectorOfArrays* CellArrays;
421 };
422 
423 #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:147
dynamic, self-adjusting array of char
Definition: vtkCharArray.h:68
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:166
helper class to manage the vtkPointSet::FindCell() METHOD
list of point or cell ids
Definition: vtkIdList.h:140
dynamic, self-adjusting array of vtkIdType
topologically and geometrically regular array of data
Definition: vtkImageData.h:157
a simple class to control print indentation
Definition: vtkIndent.h:113
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
represent and manipulate point attribute data
Definition: vtkPointData.h:142
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