VTK  9.4.20241222
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
41#include "vtkFiltersParallelDIY2Module.h" // For export macro
42#include "vtkSmartPointer.h" // For sampling line
43
44#include <memory> // for unique_ptr
45
46VTK_ABI_NAMESPACE_BEGIN
47class vtkDataSet;
49class vtkIdList;
51class vtkPoints;
52class vtkPolyData;
53class vtkVector3d;
54
55class VTKFILTERSPARALLELDIY2_EXPORT vtkProbeLineFilter : public vtkDataObjectAlgorithm
56{
57public:
59 void PrintSelf(ostream& os, vtkIndent indent) override;
60
62
64
68 vtkGetObjectMacro(Controller, vtkMultiProcessController);
70
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
187protected:
190
193 int FillInputPortInformation(int port, vtkInformation* info) override;
194
201 vtkPoints* points, vtkIdList* pointIds, vtkDataObject* input, double tol) const;
202
204
210 const vtkVector3d& p1, const vtkVector3d& p2, vtkDataObject* input, double tolerance) const;
212 const vtkVector3d& p1, const vtkVector3d& p2, vtkDataObject* input, double tolerance) const;
214
216
226 const vtkVector3d& p1, const vtkVector3d& p2, vtkDataSet* dataset, double tolerance) const;
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
244private:
245 vtkProbeLineFilter(const vtkProbeLineFilter&) = delete;
246 void operator=(const vtkProbeLineFilter&) = delete;
247
248 struct vtkInternals;
249 std::unique_ptr<vtkInternals> Internal;
250};
251
252VTK_ABI_NAMESPACE_END
253#endif
Proxy object to connect input/output ports.
Superclass for algorithms that produce only data object as output.
general representation of visualization data
abstract class to specify dataset behavior
Definition vtkDataSet.h:165
A dataset containing a grid of vtkHyperTree instances arranged as a rectilinear grid.
list of point or cell ids
Definition vtkIdList.h:133
a simple class to control print indentation
Definition vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Multiprocessing communication superclass.
represent and manipulate 3D points
Definition vtkPoints.h:139
concrete dataset represents vertices, lines, polygons, and triangle strips
probe dataset along a line in parallel
vtkSmartPointer< vtkPolyData > SampleLineAtEachCell(const vtkVector3d &p1, const vtkVector3d &p2, vtkDataObject *input, double tolerance) const
Generate sampling points and their probed data between p1 and p2 according to SamplingPattern.
vtkSmartPointer< vtkPolyData > IntersectCells(const vtkVector3d &p1, const vtkVector3d &p2, vtkHyperTreeGrid *dataset, double tolerance) const
Compute all intersections between a segment and a given dataset / htg.
virtual void SetSourceConnection(vtkAlgorithmOutput *input)
Set the source for creating the lines to probe from.
SamplingPatternEnum
Sampling pattern enumeration.
virtual void SetController(vtkMultiProcessController *)
Set and get the controller.
int RequestDataObject(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
~vtkProbeLineFilter() override
vtkSmartPointer< vtkPolyData > CreateSamplingPolyLine(vtkPoints *points, vtkIdList *pointIds, vtkDataObject *input, double tol) const
Given a line / polyline cell defined by points and pointIds, return the probing of the input dataset ...
vtkSmartPointer< vtkPolyData > IntersectCells(const vtkVector3d &p1, const vtkVector3d &p2, vtkDataSet *dataset, double tolerance) const
Compute all intersections between a segment and a given dataset / htg.
static vtkProbeLineFilter * New()
vtkSmartPointer< vtkPolyData > SampleLineUniformly(const vtkVector3d &p1, const vtkVector3d &p2, vtkDataObject *input, double tolerance) const
Generate sampling points and their probed data between p1 and p2 according to SamplingPattern.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
Hold a reference to a vtkObjectBase instance.