VTK  9.5.20251101
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
144 void SetInterpolatorType(int interpType);
145
150
155
156protected:
159
160 int FillInputPortInformation(int port, vtkInformation* info) override;
161 int FillOutputPortInformation(int port, vtkInformation* info) override;
163
164private:
166 void operator=(const vtkVectorFieldTopology&) = delete;
167
171 int Validate();
172
180 int ImageDataPrepare(vtkDataSet* dataSetInput, vtkUnstructuredGrid* tridataset);
181
189 int UnstructuredGridPrepare(vtkDataSet* dataSetInput, vtkUnstructuredGrid* tridataset);
190
196 int RemoveBoundary(vtkSmartPointer<vtkUnstructuredGrid> tridataset);
197
205 int ComputeCriticalPoints2D(
207
215 int ComputeCriticalPoints3D(
217
228 static void InterpolateVector(
229 double x0, double x1, double x, const double v0[3], const double v1[3], double v[3]);
230
237 int ComputeBoundarySwitchPoints(
238 vtkPolyData* boundarySwitchPoints, vtkUnstructuredGrid* tridataset);
239
256 int ComputeSeparatricesBoundarySwitchPoints(vtkPolyData* boundarySwitchPoints,
257 vtkPolyData* separatrices, vtkDataObject* field, vtkDataSet* dataset, vtkPoints* interestPoints,
258 int integrationStepUnit, double dist, int maxNumSteps, int& numberOfSeparatingLines);
259
277 int ComputeSeparatricesBoundarySwitchLines(vtkPolyData* boundarySwitchLines,
278 vtkPolyData* separatrices, vtkDataObject* field, vtkDataSet* dataset, int integrationStepUnit,
279 double dist, int maxNumSteps, bool computeSurfaces, bool useIterativeSeeding);
280
302 int ComputeSeparatrices(vtkPolyData* criticalPoints, vtkPolyData* separatrices,
303 vtkPolyData* surfaces, vtkDataObject* field, vtkDataSet* dataset, vtkPoints* interestPoints,
304 int integrationStepUnit, double dist, double stepSize, int maxNumSteps, bool computeSurfaces,
305 bool useIterativeSeeding, int& numberOfSeparatingLines, int& numberOfSeparatingSurfaces);
306
324 int ComputeSurface(int numberOfSeparatingSurfaces, bool isBackward, double normal[3],
325 double zeroPos[3], vtkPolyData* streamSurfaces, vtkDataObject* field, vtkDataSet* dataset,
326 int integrationStepUnit, double dist, double stepSize, int maxNumSteps,
327 bool useIterativeSeeding);
328
333 enum CriticalType2D
334 {
335 DEGENERATE_2D = -1,
336 SINK_2D = 0,
337 SADDLE_2D = 1,
338 SOURCE_2D = 2,
339 CENTER_2D = 3
340 };
341
346 enum CriticalTypeDetailed2D
347 {
348 // DEGENERATE2D = -1,
349 ATTRACTING_NODE_2D = 0,
350 ATTRACTING_FOCUS_2D = 1,
351 NODE_SADDLE_2D = 2,
352 REPELLING_NODE_2D = 3,
353 REPELLING_FOCUS_2D = 4,
354 CENTER_DETAILED_2D = 5
355 };
356
361 enum CriticalType3D
362 {
363 DEGENERATE_3D = -1,
364 SINK_3D = 0,
365 SADDLE_1_3D = 1,
366 SADDLE_2_3D = 2,
367 SOURCE_3D = 3,
368 CENTER_3D = 4
369 };
370
375 enum CriticalTypeDetailed3D
376 {
377 ATTRACTING_NODE_3D = 0,
378 ATTRACTING_FOCUS_3D = 1,
379 NODE_SADDLE_1_3D = 2,
380 FOCUS_SADDLE_1_3D = 3,
381 NODE_SADDLE_2_3D = 4,
382 FOCUS_SADDLE_2_3D = 5,
383 REPELLING_NODE_3D = 6,
384 REPELLING_FOCUS_3D = 7,
385 CENTER_DETAILED_3D = 8
386 };
387
395 static int Classify2D(int countComplex, int countPos, int countNeg);
396
405 static int ClassifyDetailed2D(int countComplex, int countPos, int countNeg);
406
415 static int Classify3D(int countComplex, int countPos, int countNeg);
416
426 static int ClassifyDetailed3D(int countComplex, int countPos, int countNeg);
427
436 static void CopyBoundarySwitchLinesArray(vtkDataSet* source, vtkPolyData* target);
437
441 int MaxNumSteps = 100;
442
446 double IntegrationStepSize = 1;
447
451 double SeparatrixDistance = 1;
452
456 bool UseIterativeSeeding = false;
457
461 bool ComputeSurfaces = false;
462
466 const char* NameOfVectorArray;
467
472 bool ExcludeBoundary = false;
473
477 int Dimension = 2;
478
486 int IntegrationStepUnit = vtkStreamTracer::CELL_LENGTH_UNIT;
487
494 bool UseBoundarySwitchPoints = false;
495
501
510 double VectorAngleThreshold = 1;
511
519 double OffsetAwayFromBoundary = 1e-3;
520
524 double EpsilonCriticalPoint = 1e-10;
525
526 vtkNew<vtkStreamSurface> StreamSurface;
527};
528VTK_ABI_NAMESPACE_END
529#endif
general representation of visualization data
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 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 *)