VTK  9.6.20260516
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
17
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
144 VTK_DEPRECATED_IN_9_7_0("Use SetCellLocator(vtkAbstractCellLocator*) instead")
145 void SetInterpolatorType(int interpType);
146
150 VTK_DEPRECATED_IN_9_7_0("Use SetCellLocatorToStaticCellLocator() instead")
152
156 VTK_DEPRECATED_IN_9_7_0("Use SetCellLocatorToJumpAndWalkCellLocator() instead")
158
160
168 vtkGetObjectMacro(CellLocator, vtkAbstractCellLocator);
170
171protected:
174
175 int FillInputPortInformation(int port, vtkInformation* info) override;
176 int FillOutputPortInformation(int port, vtkInformation* info) override;
178
179private:
181 void operator=(const vtkVectorFieldTopology&) = delete;
182
186 int Validate();
187
195 int ImageDataPrepare(vtkDataSet* dataSetInput, vtkUnstructuredGrid* tridataset);
196
204 int UnstructuredGridPrepare(vtkDataSet* dataSetInput, vtkUnstructuredGrid* tridataset);
205
211 int RemoveBoundary(vtkSmartPointer<vtkUnstructuredGrid> tridataset);
212
220 int ComputeCriticalPoints2D(
222
230 int ComputeCriticalPoints3D(
232
243 static void InterpolateVector(
244 double x0, double x1, double x, const double v0[3], const double v1[3], double v[3]);
245
252 int ComputeBoundarySwitchPoints(
253 vtkPolyData* boundarySwitchPoints, vtkUnstructuredGrid* tridataset);
254
271 int ComputeSeparatricesBoundarySwitchPoints(vtkPolyData* boundarySwitchPoints,
272 vtkPolyData* separatrices, vtkDataObject* field, vtkDataSet* dataset, vtkPoints* interestPoints,
273 int integrationStepUnit, double dist, int maxNumSteps, int& numberOfSeparatingLines);
274
292 int ComputeSeparatricesBoundarySwitchLines(vtkPolyData* boundarySwitchLines,
293 vtkPolyData* separatrices, vtkDataObject* field, vtkDataSet* dataset, int integrationStepUnit,
294 double dist, int maxNumSteps, bool computeSurfaces, bool useIterativeSeeding);
295
317 int ComputeSeparatrices(vtkPolyData* criticalPoints, vtkPolyData* separatrices,
318 vtkPolyData* surfaces, vtkDataObject* field, vtkDataSet* dataset, vtkPoints* interestPoints,
319 int integrationStepUnit, double dist, double stepSize, int maxNumSteps, bool computeSurfaces,
320 bool useIterativeSeeding, int& numberOfSeparatingLines, int& numberOfSeparatingSurfaces);
321
339 int ComputeSurface(int numberOfSeparatingSurfaces, bool isBackward, double normal[3],
340 double zeroPos[3], vtkPolyData* streamSurfaces, vtkDataObject* field, vtkDataSet* dataset,
341 int integrationStepUnit, double dist, double stepSize, int maxNumSteps,
342 bool useIterativeSeeding);
343
348 enum CriticalType2D
349 {
350 DEGENERATE_2D = -1,
351 SINK_2D = 0,
352 SADDLE_2D = 1,
353 SOURCE_2D = 2,
354 CENTER_2D = 3
355 };
356
361 enum CriticalTypeDetailed2D
362 {
363 // DEGENERATE2D = -1,
364 ATTRACTING_NODE_2D = 0,
365 ATTRACTING_FOCUS_2D = 1,
366 NODE_SADDLE_2D = 2,
367 REPELLING_NODE_2D = 3,
368 REPELLING_FOCUS_2D = 4,
369 CENTER_DETAILED_2D = 5
370 };
371
376 enum CriticalType3D
377 {
378 DEGENERATE_3D = -1,
379 SINK_3D = 0,
380 SADDLE_1_3D = 1,
381 SADDLE_2_3D = 2,
382 SOURCE_3D = 3,
383 CENTER_3D = 4
384 };
385
390 enum CriticalTypeDetailed3D
391 {
392 ATTRACTING_NODE_3D = 0,
393 ATTRACTING_FOCUS_3D = 1,
394 NODE_SADDLE_1_3D = 2,
395 FOCUS_SADDLE_1_3D = 3,
396 NODE_SADDLE_2_3D = 4,
397 FOCUS_SADDLE_2_3D = 5,
398 REPELLING_NODE_3D = 6,
399 REPELLING_FOCUS_3D = 7,
400 CENTER_DETAILED_3D = 8
401 };
402
410 static int Classify2D(int countComplex, int countPos, int countNeg);
411
420 static int ClassifyDetailed2D(int countComplex, int countPos, int countNeg);
421
430 static int Classify3D(int countComplex, int countPos, int countNeg);
431
441 static int ClassifyDetailed3D(int countComplex, int countPos, int countNeg);
442
451 static void CopyBoundarySwitchLinesArray(vtkDataSet* source, vtkPolyData* target);
452
456 int MaxNumSteps = 100;
457
461 double IntegrationStepSize = 1;
462
466 double SeparatrixDistance = 1;
467
471 bool UseIterativeSeeding = false;
472
476 bool ComputeSurfaces = false;
477
481 const char* NameOfVectorArray;
482
487 bool ExcludeBoundary = false;
488
492 int Dimension = 2;
493
501 int IntegrationStepUnit = vtkStreamTracer::CELL_LENGTH_UNIT;
502
509 bool UseBoundarySwitchPoints = false;
510
511 vtkAbstractCellLocator* CellLocator = nullptr;
512
521 double VectorAngleThreshold = 1;
522
530 double OffsetAwayFromBoundary = 1e-3;
531
535 double EpsilonCriticalPoint = 1e-10;
536
537 vtkNew<vtkStreamSurface> StreamSurface;
538};
539VTK_ABI_NAMESPACE_END
540#endif
an abstract base class for locators which find cells
general representation of visualization data
abstract class to specify dataset behavior
Definition vtkDataSet.h:166
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.
represent and manipulate 3D points
Definition vtkPoints.h:140
concrete dataset represents vertices, lines, polygons, and triangle strips
Hold a reference to a vtkObjectBase instance.
Advect a stream surface in a vector field.
dataset represents arbitrary combinations of all possible cell types
virtual void SetCellLocator(vtkAbstractCellLocator *)
Set / get the cell locator used to perform the FindCell() operation for vtkPointSet.
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.
virtual void SetCellLocatorToStaticCellLocator()
Set / get the cell locator used to perform the FindCell() operation for vtkPointSet.
virtual void SetCellLocatorToJumpAndWalkCellLocator()
Set / get the cell locator used to perform the FindCell() operation for vtkPointSet.
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 vtkInterpolatedVelocityField (IN...
int FillOutputPortInformation(int port, vtkInformation *info) override
Fill the output port information objects for this algorithm.
boost::graph_traits< vtkGraph * >::vertex_descriptor source(boost::graph_traits< vtkGraph * >::edge_descriptor e, vtkGraph *)
boost::graph_traits< vtkGraph * >::vertex_descriptor target(boost::graph_traits< vtkGraph * >::edge_descriptor e, vtkGraph *)
#define VTK_DEPRECATED_IN_9_7_0(reason)
#define vtkGradientFilter