VTK  9.4.20241117
vtkVectorFieldTopology.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
18#ifndef vtkVectorFieldTopology_h
19#define vtkVectorFieldTopology_h
20
21#include "vtkFiltersFlowPathsModule.h" // For export macro
23#include "vtkStreamTracer.h" // for vtkStreamSurface::CELL_LENGTH_UNIT
24
25VTK_ABI_NAMESPACE_BEGIN
27class vtkImageData;
28class vtkPolyData;
31
32class VTKFILTERSFLOWPATHS_EXPORT vtkVectorFieldTopology : public vtkPolyDataAlgorithm
33{
34public:
37 void PrintSelf(ostream& os, vtkIndent indent) override;
38
40
46 vtkSetMacro(IntegrationStepUnit, int);
47 vtkGetMacro(IntegrationStepUnit, int);
49
51
54 vtkSetMacro(MaxNumSteps, int);
55 vtkGetMacro(MaxNumSteps, int);
57
59
63 vtkSetMacro(IntegrationStepSize, double);
64 vtkGetMacro(IntegrationStepSize, double);
66
68
72 vtkSetMacro(SeparatrixDistance, double);
73 vtkGetMacro(SeparatrixDistance, double);
75
77
80 vtkSetMacro(UseIterativeSeeding, bool);
81 vtkGetMacro(UseIterativeSeeding, bool);
83
85
88 vtkSetMacro(ComputeSurfaces, bool);
89 vtkGetMacro(ComputeSurfaces, bool);
91
93
96 vtkSetMacro(ExcludeBoundary, bool);
97 vtkGetMacro(ExcludeBoundary, bool);
99
101
104 vtkSetMacro(UseBoundarySwitchPoints, bool);
105 vtkGetMacro(UseBoundarySwitchPoints, bool);
107
109
118 vtkSetMacro(VectorAngleThreshold, double);
119 vtkGetMacro(VectorAngleThreshold, double);
121
123
126 vtkSetMacro(OffsetAwayFromBoundary, double);
127 vtkGetMacro(OffsetAwayFromBoundary, double);
129
131
134 vtkSetMacro(EpsilonCriticalPoint, double);
135 vtkGetMacro(EpsilonCriticalPoint, double);
137
143 void SetInterpolatorType(int interpType);
144
149
154
155protected:
158
159 int FillInputPortInformation(int port, vtkInformation* info) override;
160 int FillOutputPortInformation(int port, vtkInformation* info) override;
162
163private:
165 void operator=(const vtkVectorFieldTopology&) = delete;
166
170 int Validate();
171
179 int ImageDataPrepare(vtkDataSet* dataSetInput, vtkUnstructuredGrid* tridataset);
180
188 int UnstructuredGridPrepare(vtkDataSet* dataSetInput, vtkUnstructuredGrid* tridataset);
189
195 int RemoveBoundary(vtkSmartPointer<vtkUnstructuredGrid> tridataset);
196
204 int ComputeCriticalPoints2D(
206
214 int ComputeCriticalPoints3D(
216
227 static void InterpolateVector(
228 double x0, double x1, double x, const double v0[3], const double v1[3], double v[3]);
229
236 int ComputeBoundarySwitchPoints(
237 vtkPolyData* boundarySwitchPoints, vtkUnstructuredGrid* tridataset);
238
254 int ComputeSeparatricesBoundarySwitchPoints(vtkPolyData* boundarySwitchPoints,
255 vtkPolyData* separatrices, vtkDataSet* dataset, vtkPoints* interestPoints,
256 int integrationStepUnit, double dist, int maxNumSteps);
257
275 int ComputeSeparatricesBoundarySwitchLines(vtkPolyData* boundarySwitchLines,
276 vtkPolyData* separatrices, vtkDataSet* dataset, int integrationStepUnit, double dist,
277 int maxNumSteps, bool computeSurfaces, bool useIterativeSeeding);
278
297 int ComputeSeparatrices(vtkPolyData* criticalPoints, vtkPolyData* separatrices,
298 vtkPolyData* surfaces, vtkDataSet* dataset, vtkPoints* interestPoints, int integrationStepUnit,
299 double dist, double stepSize, int maxNumSteps, bool computeSurfaces, bool useIterativeSeeding);
300
317 int ComputeSurface(int numberOfSeparatingSurfaces, bool isBackward, double normal[3],
318 double zeroPos[3], vtkPolyData* streamSurfaces, vtkDataSet* dataset, int integrationStepUnit,
319 double dist, double stepSize, int maxNumSteps, bool useIterativeSeeding);
320
325 enum CriticalType2D
326 {
327 DEGENERATE_2D = -1,
328 SINK_2D = 0,
329 SADDLE_2D = 1,
330 SOURCE_2D = 2,
331 CENTER_2D = 3
332 };
333
338 enum CriticalTypeDetailed2D
339 {
340 // DEGENERATE2D = -1,
341 ATTRACTING_NODE_2D = 0,
342 ATTRACTING_FOCUS_2D = 1,
343 NODE_SADDLE_2D = 2,
344 REPELLING_NODE_2D = 3,
345 REPELLING_FOCUS_2D = 4,
346 CENTER_DETAILED_2D = 5
347 };
348
353 enum CriticalType3D
354 {
355 DEGENERATE_3D = -1,
356 SINK_3D = 0,
357 SADDLE_1_3D = 1,
358 SADDLE_2_3D = 2,
359 SOURCE_3D = 3,
360 CENTER_3D = 4
361 };
362
367 enum CriticalTypeDetailed3D
368 {
369 ATTRACTING_NODE_3D = 0,
370 ATTRACTING_FOCUS_3D = 1,
371 NODE_SADDLE_1_3D = 2,
372 FOCUS_SADDLE_1_3D = 3,
373 NODE_SADDLE_2_3D = 4,
374 FOCUS_SADDLE_2_3D = 5,
375 REPELLING_NODE_3D = 6,
376 REPELLING_FOCUS_3D = 7,
377 CENTER_DETAILED_3D = 8
378 };
379
387 static int Classify2D(int countComplex, int countPos, int countNeg);
388
397 static int ClassifyDetailed2D(int countComplex, int countPos, int countNeg);
398
407 static int Classify3D(int countComplex, int countPos, int countNeg);
408
418 static int ClassifyDetailed3D(int countComplex, int countPos, int countNeg);
419
423 int MaxNumSteps = 100;
424
428 double IntegrationStepSize = 1;
429
433 double SeparatrixDistance = 1;
434
438 bool UseIterativeSeeding = false;
439
443 bool ComputeSurfaces = false;
444
448 const char* NameOfVectorArray;
449
454 bool ExcludeBoundary = false;
455
459 int Dimension = 2;
460
468 int IntegrationStepUnit = vtkStreamTracer::CELL_LENGTH_UNIT;
469
476 bool UseBoundarySwitchPoints = false;
477
483
492 double VectorAngleThreshold = 1;
493
501 double OffsetAwayFromBoundary = 1e-3;
502
506 double EpsilonCriticalPoint = 1e-10;
507
508 vtkNew<vtkStreamSurface> StreamSurface;
509};
510VTK_ABI_NAMESPACE_END
511#endif
abstract class to specify dataset behavior
Definition vtkDataSet.h:166
A general filter for gradient estimation.
topologically and geometrically regular array of data
a simple class to control print indentation
Definition vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Allocate and hold a VTK object.
Definition vtkNew.h:167
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.
Advect a stream surface in a vector field.
@ INTERPOLATOR_WITH_DATASET_POINT_LOCATOR
dataset represents arbitrary combinations of all possible cell types
Extract the topological skeleton as output datasets.
void SetInterpolatorTypeToCellLocator()
Set the velocity field interpolator type to the one involving a cell locator.
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
static vtkVectorFieldTopology * New()
void SetInterpolatorTypeToDataSetPointLocator()
Set the velocity field interpolator type to the one involving a dataset point locator.
~vtkVectorFieldTopology() override
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void SetInterpolatorType(int interpType)
Set the type of the velocity field interpolator to determine whether INTERPOLATOR_WITH_DATASET_POINT_...
int FillOutputPortInformation(int port, vtkInformation *info) override
Fill the output port information objects for this algorithm.