Loading [MathJax]/extensions/tex2jax.js
VTK  9.4.20250324
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
257 int ComputeSeparatricesBoundarySwitchPoints(vtkPolyData* boundarySwitchPoints,
258 vtkPolyData* separatrices, vtkDataObject* field, vtkDataSet* dataset, vtkPoints* interestPoints,
259 int integrationStepUnit, double dist, int maxNumSteps, int& numberOfSeparatingLines);
260
279 int ComputeSeparatricesBoundarySwitchLines(vtkPolyData* boundarySwitchLines,
280 vtkPolyData* separatrices, vtkDataObject* field, vtkDataSet* dataset, int integrationStepUnit,
281 double dist, int maxNumSteps, bool computeSurfaces, bool useIterativeSeeding);
282
304 int ComputeSeparatrices(vtkPolyData* criticalPoints, vtkPolyData* separatrices,
305 vtkPolyData* surfaces, vtkDataObject* field, vtkDataSet* dataset, vtkPoints* interestPoints,
306 int integrationStepUnit, double dist, double stepSize, int maxNumSteps, bool computeSurfaces,
307 bool useIterativeSeeding, int& numberOfSeparatingLines, int& numberOfSeparatingSurfaces);
308
326 int ComputeSurface(int numberOfSeparatingSurfaces, bool isBackward, double normal[3],
327 double zeroPos[3], vtkPolyData* streamSurfaces, vtkDataObject* field, vtkDataSet* dataset,
328 int integrationStepUnit, double dist, double stepSize, int maxNumSteps,
329 bool useIterativeSeeding);
330
335 enum CriticalType2D
336 {
337 DEGENERATE_2D = -1,
338 SINK_2D = 0,
339 SADDLE_2D = 1,
340 SOURCE_2D = 2,
341 CENTER_2D = 3
342 };
343
348 enum CriticalTypeDetailed2D
349 {
350 // DEGENERATE2D = -1,
351 ATTRACTING_NODE_2D = 0,
352 ATTRACTING_FOCUS_2D = 1,
353 NODE_SADDLE_2D = 2,
354 REPELLING_NODE_2D = 3,
355 REPELLING_FOCUS_2D = 4,
356 CENTER_DETAILED_2D = 5
357 };
358
363 enum CriticalType3D
364 {
365 DEGENERATE_3D = -1,
366 SINK_3D = 0,
367 SADDLE_1_3D = 1,
368 SADDLE_2_3D = 2,
369 SOURCE_3D = 3,
370 CENTER_3D = 4
371 };
372
377 enum CriticalTypeDetailed3D
378 {
379 ATTRACTING_NODE_3D = 0,
380 ATTRACTING_FOCUS_3D = 1,
381 NODE_SADDLE_1_3D = 2,
382 FOCUS_SADDLE_1_3D = 3,
383 NODE_SADDLE_2_3D = 4,
384 FOCUS_SADDLE_2_3D = 5,
385 REPELLING_NODE_3D = 6,
386 REPELLING_FOCUS_3D = 7,
387 CENTER_DETAILED_3D = 8
388 };
389
397 static int Classify2D(int countComplex, int countPos, int countNeg);
398
407 static int ClassifyDetailed2D(int countComplex, int countPos, int countNeg);
408
417 static int Classify3D(int countComplex, int countPos, int countNeg);
418
428 static int ClassifyDetailed3D(int countComplex, int countPos, int countNeg);
429
438 static void CopyBoundarySwitchLinesArray(vtkDataSet* source, vtkPolyData* target);
439
443 int MaxNumSteps = 100;
444
448 double IntegrationStepSize = 1;
449
453 double SeparatrixDistance = 1;
454
458 bool UseIterativeSeeding = false;
459
463 bool ComputeSurfaces = false;
464
468 const char* NameOfVectorArray;
469
474 bool ExcludeBoundary = false;
475
479 int Dimension = 2;
480
488 int IntegrationStepUnit = vtkStreamTracer::CELL_LENGTH_UNIT;
489
496 bool UseBoundarySwitchPoints = false;
497
503
512 double VectorAngleThreshold = 1;
513
521 double OffsetAwayFromBoundary = 1e-3;
522
526 double EpsilonCriticalPoint = 1e-10;
527
528 vtkNew<vtkStreamSurface> StreamSurface;
529};
530VTK_ABI_NAMESPACE_END
531#endif
general representation of visualization data
abstract class to specify dataset behavior
Definition vtkDataSet.h:165
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 *)