VTK  9.3.20240425
vtkTemporalPathLineFilter.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
27#ifndef vtkTemporalPathLineFilter_h
28#define vtkTemporalPathLineFilter_h
29
30#include "vtkDeprecation.h" // For VTK_DEPRECATED_IN_9_3_0
31#include "vtkFiltersGeneralModule.h" // For export macro
33#include "vtkTemporalAlgorithm.h" // For temporal algorithm
34
35#ifndef __VTK_WRAP__
36#define vtkPolyDataAlgorithm vtkTemporalAlgorithm<vtkPolyDataAlgorithm>
37#endif
38
39VTK_ABI_NAMESPACE_BEGIN
40class vtkPoints;
41class vtkCellArray;
42class vtkMergePoints;
43class vtkFloatArray;
44
45VTK_ABI_NAMESPACE_END
46#include "vtkSmartPointer.h" // for memory safety
47#include <set> // Because we want to use it
48VTK_ABI_NAMESPACE_BEGIN
49class ParticleTrail;
50class vtkTemporalPathLineFilterInternals;
52
53class VTKFILTERSGENERAL_EXPORT vtkTemporalPathLineFilter : public vtkPolyDataAlgorithm
54{
55public:
57
61#ifndef __VTK_WRAP__
62#undef vtkPassInputTypeAlgorithm
63#endif
65 void PrintSelf(ostream& os, vtkIndent indent) override;
67
68#if defined(__VTK_WRAP__) || defined(__WRAP_GCCXML)
70#endif
71
73
77 vtkSetMacro(MaskPoints, int);
78 vtkGetMacro(MaskPoints, int);
80
82
90 vtkSetMacro(MaxTrackLength, unsigned int);
91 vtkGetMacro(MaxTrackLength, unsigned int);
93
95
103 vtkSetStringMacro(IdChannelArray);
104 vtkGetStringMacro(IdChannelArray);
106
108
116 vtkSetVector3Macro(MaxStepDistance, double);
117 vtkGetVector3Macro(MaxStepDistance, double);
119
121
127 vtkSetMacro(KeepDeadTrails, bool);
128 vtkGetMacro(KeepDeadTrails, bool);
130
132
141 VTK_DEPRECATED_IN_9_3_0("Running backward will not be supported anymore.")
142 virtual void SetBackwardTime(bool backward);
143 VTK_DEPRECATED_IN_9_3_0("Running backward will not be supported anymore.")
144 vtkGetMacro(BackwardTime, bool);
146
151 void Flush();
152
158 void SetSelectionConnection(vtkAlgorithmOutput* algOutput);
159
165 void SetSelectionData(vtkDataSet* input);
166
167protected:
170
171 //
172 // Make sure the pipeline knows what type we expect as input
173 //
174 int FillInputPortInformation(int port, vtkInformation* info) override;
175 int FillOutputPortInformation(int port, vtkInformation* info) override;
176
177 int Initialize(vtkInformation* request, vtkInformationVector** inputVector,
178 vtkInformationVector* outputVector) override;
179 int Execute(vtkInformation* request, vtkInformationVector** inputVector,
180 vtkInformationVector* outputVector) override;
181 int Finalize(vtkInformation* request, vtkInformationVector** inputVector,
182 vtkInformationVector* outputVector) override;
183 void IncrementTrail(TrailPointer trail, vtkDataSet* input, vtkIdType i);
184
186
187 // void Initialize(vtkDataSet* input, vtkDataSet* selection,
188 // vtkPolyData* pathLines, vtkPolyData* particles);
189
190 // internal data variables
191 int NumberOfTimeSteps = 0;
192 int MaskPoints = 200;
193 unsigned int MaxTrackLength = 10;
194 unsigned int LastTrackLength = 10;
195 char* IdChannelArray = nullptr;
196 double MaxStepDistance[3] = { 1, 1, 1 };
198 bool KeepDeadTrails = false;
199 bool BackwardTime = false;
200 //
201
208
209 int CurrentTimeIndex = 0;
210
211 //
212private:
213 void AccumulateTrails(vtkDataSet* input, vtkDataSet* selection);
214 void PostExecute(vtkDataSet* input, vtkPolyData* pathLines, vtkPolyData* particles);
215 void InitializeExecute(vtkDataSet* input, vtkPolyData* pathLines);
216
218 void operator=(const vtkTemporalPathLineFilter&) = delete;
219};
220
221VTK_ABI_NAMESPACE_END
222#endif
Proxy object to connect input/output ports.
object to represent cell connectivity
abstract class to specify dataset behavior
Definition vtkDataSet.h:166
dynamic, self-adjusting array of float
a simple class to control print indentation
Definition vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
merge exactly coincident points
represent and manipulate 3D points
Definition vtkPoints.h:139
Superclass for algorithms that produce only polydata as output.
concrete dataset represents vertices, lines, polygons, and triangle strips
Hold a reference to a vtkObjectBase instance.
Generate a Polydata Pointset from any Dataset.
vtkSmartPointer< vtkTemporalPathLineFilterInternals > Internals
vtkSmartPointer< vtkFloatArray > TrailId
vtkSmartPointer< vtkPoints > LineCoordinates
vtkSmartPointer< vtkCellArray > Vertices
void PrintSelf(ostream &os, vtkIndent indent) override
Standard Type-Macro.
vtkSmartPointer< vtkCellArray > PolyLines
vtkSmartPointer< vtkPoints > VertexCoordinates
static vtkTemporalPathLineFilter * New()
Standard Type-Macro.
#define VTK_DEPRECATED_IN_9_3_0(reason)
#define vtkCreateWrappedTemporalAlgorithmInterface()
vtkSmartPointer< ParticleTrail > TrailPointer
int vtkIdType
Definition vtkType.h:315