VTK  9.3.1
vtkProbeLineFilter.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
37 #ifndef vtkProbeLineFilter_h
38 #define vtkProbeLineFilter_h
39 
40 #include "vtkDataObjectAlgorithm.h"
41 #include "vtkFiltersParallelDIY2Module.h" // For export macro
42 #include "vtkSmartPointer.h" // For sampling line
43 
44 #include <memory> // for unique_ptr
45 
46 VTK_ABI_NAMESPACE_BEGIN
47 class vtkDataSet;
48 class vtkHyperTreeGrid;
49 class vtkIdList;
51 class vtkPoints;
52 class vtkPolyData;
53 class vtkVector3d;
54 
55 class VTKFILTERSPARALLELDIY2_EXPORT vtkProbeLineFilter : public vtkDataObjectAlgorithm
56 {
57 public:
59  void PrintSelf(ostream& os, vtkIndent indent) override;
60 
61  static vtkProbeLineFilter* New();
62 
64 
67  virtual void SetController(vtkMultiProcessController*);
68  vtkGetObjectMacro(Controller, vtkMultiProcessController);
70 
75  virtual void SetSourceConnection(vtkAlgorithmOutput* input);
76 
78 
82  vtkSetMacro(PassCellArrays, bool);
83  vtkBooleanMacro(PassCellArrays, bool);
84  vtkGetMacro(PassCellArrays, bool);
87 
91  vtkSetMacro(PassPointArrays, bool);
92  vtkBooleanMacro(PassPointArrays, bool);
93  vtkGetMacro(PassPointArrays, bool);
95 
97 
101  vtkSetMacro(PassFieldArrays, bool);
102  vtkBooleanMacro(PassFieldArrays, bool);
103  vtkGetMacro(PassFieldArrays, bool);
105 
107 
112  vtkSetMacro(Tolerance, double);
113  vtkGetMacro(Tolerance, double);
115 
117 
122  vtkSetMacro(ComputeTolerance, bool);
123  vtkBooleanMacro(ComputeTolerance, bool);
124  vtkGetMacro(ComputeTolerance, bool);
126 
128 
140  vtkSetMacro(PassPartialArrays, bool);
141  vtkGetMacro(PassPartialArrays, bool);
142  vtkBooleanMacro(PassPartialArrays, bool);
144 
149  {
150  SAMPLE_LINE_AT_CELL_BOUNDARIES = 0,
151  SAMPLE_LINE_AT_SEGMENT_CENTERS = 1,
152  SAMPLE_LINE_UNIFORMLY = 2
153  };
154 
156 
160  vtkGetMacro(SamplingPattern, int);
161  vtkSetClampMacro(SamplingPattern, int, 0, 2);
163 
165 
170  vtkGetMacro(LineResolution, int);
171  vtkSetMacro(LineResolution, int);
173 
175 
182  vtkGetMacro(AggregateAsPolyData, bool);
183  vtkSetMacro(AggregateAsPolyData, bool);
184  vtkBooleanMacro(AggregateAsPolyData, bool);
186 
187 protected:
189  ~vtkProbeLineFilter() override;
190 
193  int FillInputPortInformation(int port, vtkInformation* info) override;
194 
200  vtkSmartPointer<vtkPolyData> CreateSamplingPolyLine(
201  vtkPoints* points, vtkIdList* pointIds, vtkDataObject* input, double tol) const;
202 
204 
209  vtkSmartPointer<vtkPolyData> SampleLineAtEachCell(
210  const vtkVector3d& p1, const vtkVector3d& p2, vtkDataObject* input, double tolerance) const;
211  vtkSmartPointer<vtkPolyData> SampleLineUniformly(
212  const vtkVector3d& p1, const vtkVector3d& p2, vtkDataObject* input, double tolerance) const;
214 
216 
225  vtkSmartPointer<vtkPolyData> IntersectCells(
226  const vtkVector3d& p1, const vtkVector3d& p2, vtkDataSet* dataset, double tolerance) const;
227  vtkSmartPointer<vtkPolyData> IntersectCells(const vtkVector3d& p1, const vtkVector3d& p2,
228  vtkHyperTreeGrid* dataset, double tolerance) const;
230 
231  vtkMultiProcessController* Controller = nullptr;
232 
233  int SamplingPattern = SAMPLE_LINE_AT_CELL_BOUNDARIES;
234  int LineResolution = 1000;
235 
236  bool AggregateAsPolyData = true;
237  bool PassPartialArrays = false;
238  bool PassCellArrays = false;
239  bool PassPointArrays = false;
240  bool PassFieldArrays = false;
241  bool ComputeTolerance = true;
242  double Tolerance = 1.0;
243 
244 private:
245  vtkProbeLineFilter(const vtkProbeLineFilter&) = delete;
246  void operator=(const vtkProbeLineFilter&) = delete;
247 
248  struct vtkInternals;
249  std::unique_ptr<vtkInternals> Internal;
250 };
251 
252 VTK_ABI_NAMESPACE_END
253 #endif
virtual int RequestDataObject(vtkInformation *, vtkInformationVector **, vtkInformationVector *)
This is called by the superclass.
Store vtkAlgorithm input/output information.
abstract class to specify dataset behavior
Definition: vtkDataSet.h:52
SamplingPatternEnum
Sampling pattern enumeration.
static vtkDataObjectAlgorithm * New()
A dataset containing a grid of vtkHyperTree instances arranged as a rectilinear grid.
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:79
Proxy object to connect input/output ports.
a simple class to control print indentation
Definition: vtkIndent.h:28
list of point or cell ids
Definition: vtkIdList.h:22
probe dataset along a line in parallel
Superclass for algorithms that produce only data object as output.
virtual int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *)
Store zero or more vtkInformation instances.
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
general representation of visualization data
Definition: vtkDataObject.h:54
represent and manipulate 3D points
Definition: vtkPoints.h:28
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
Multiprocessing communication superclass.