82#ifndef vtkStreamTracer_h
83#define vtkStreamTracer_h
85#include "vtkFiltersFlowPathsModule.h"
91VTK_ABI_NAMESPACE_BEGIN
107VTK_ABI_NAMESPACE_BEGIN
156 vtkSetVector3Macro(StartPosition,
double);
157 vtkGetVector3Macro(StartPosition,
double);
210 FIXED_REASONS_FOR_TERMINATION_COUNT
256 vtkSetMacro(MaximumPropagation,
double);
257 vtkGetMacro(MaximumPropagation,
double);
277 vtkSetMacro(InitialIntegrationStep,
double);
278 vtkGetMacro(InitialIntegrationStep,
double);
288 vtkSetMacro(MinimumIntegrationStep,
double);
289 vtkGetMacro(MinimumIntegrationStep,
double);
299 vtkSetMacro(MaximumIntegrationStep,
double);
300 vtkGetMacro(MaximumIntegrationStep,
double);
307 vtkSetMacro(MaximumError,
double);
308 vtkGetMacro(MaximumError,
double);
329 vtkSetMacro(TerminalSpeed,
double);
330 vtkGetMacro(TerminalSpeed,
double);
337 vtkGetMacro(SurfaceStreamlines,
bool);
338 vtkSetMacro(SurfaceStreamlines,
bool);
339 vtkBooleanMacro(SurfaceStreamlines,
bool);
352 INTERPOLATOR_WITH_CELL_LOCATOR
363 vtkSetClampMacro(IntegrationDirection,
int, FORWARD, BOTH);
364 vtkGetMacro(IntegrationDirection,
int);
376 vtkSetMacro(ComputeVorticity,
bool);
377 vtkGetMacro(ComputeVorticity,
bool);
385 vtkSetMacro(RotationScale,
double);
386 vtkGetMacro(RotationScale,
double);
416 vtkGetMacro(ForceSerialExecution,
bool);
417 vtkSetMacro(ForceSerialExecution,
bool);
418 vtkBooleanMacro(ForceSerialExecution,
bool);
441 double& step,
double& minStep,
double& maxStep,
int direction,
double cellLength);
464 vtkSetMacro(UseLocalSeedSource,
bool);
465 vtkGetMacro(UseLocalSeedSource,
bool);
466 vtkBooleanMacro(UseLocalSeedSource,
bool);
479 vtkErrorMacro(<<
"AddInput() must be called with a vtkDataSet not a vtkDataObject.");
488 const char* vecFieldName,
double& propagation,
vtkIdType& numSteps,
double& integrationTime,
489 std::vector<CustomTerminationCallbackType>& customTerminationCallback,
490 std::vector<void*>& customTerminationClientData, std::vector<int>& customReasonForTermination);
499 double StartPosition[3];
555 friend class PStreamTracerUtils;
An abstract class for obtaining the interpolated velocity values at a point.
Proxy object to connect input/output ports.
abstract superclass for composite (multi-block or AMR) datasets
abstract superclass for arrays of numeric data
general representation of visualization data
helps manage arrays from multiple vtkDataSetAttributes.
represent and manipulate attribute data in a dataset
abstract class to specify dataset behavior
dynamic, self-adjusting array of double
Superclass for all pipeline executives in VTK.
provides thread-safe access to cells
list of point or cell ids
a simple class to control print indentation
Integrate a set of ordinary differential equations (initial value problem) in time.
dynamic, self-adjusting array of int
represent and manipulate point attribute data
represent and manipulate 3D points
Superclass for algorithms that produce only polydata as output.
concrete dataset represents vertices, lines, polygons, and triangle strips
void SetIntegratorTypeToRungeKutta45()
Set/get the integrator type to be used for streamline generation.
int FillInputPortInformation(int, vtkInformation *) override
Fill the input port information objects for this algorithm.
int SetupOutput(vtkInformation *inInfo, vtkInformation *outInfo)
std::vector< void * > CustomTerminationClientData
vtkDataSetAttributesFieldList InputPD
void SetSourceData(vtkDataSet *source)
Specify the source object used to generate starting points (seeds).
vtkDataSet * GetSource()
Specify the source object used to generate starting points (seeds).
double InitialIntegrationStep
vtkAbstractInterpolatedVelocityField * InterpolatorPrototype
double MaximumPropagation
void SetInterpolatorTypeToCellLocator()
Set the velocity field interpolator type to one that uses a cell locator to perform spatial searching...
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods to obtain type information and print object state.
void CalculateVorticity(vtkGenericCell *cell, double pcoords[3], vtkDoubleArray *cellVectors, double vorticity[3])
Helper methods to generate normals on streamlines.
double MinimumIntegrationStep
void SetIntegratorTypeToRungeKutta4()
Set/get the integrator type to be used for streamline generation.
@ INTERPOLATOR_WITH_DATASET_POINT_LOCATOR
void SetIntegrator(vtkInitialValueProblemSolver *)
Set/get the integrator type to be used for streamline generation.
void SetSourceConnection(vtkAlgorithmOutput *algOutput)
Specify the source object used to generate starting points (seeds).
std::vector< int > CustomReasonForTermination
int CheckInputs(vtkAbstractInterpolatedVelocityField *&func, int *maxCellSize)
void ConvertIntervals(double &step, double &minStep, double &maxStep, int direction, double cellLength)
The following methods should not be called by the user.
bool ForceSerialExecution
void GenerateNormals(vtkPolyData *output, double *firstNormal, const char *vecName)
Helper methods to generate normals on streamlines.
static const double EPSILON
vtkIdType MaximumNumberOfSteps
void SetIntegrationDirectionToForward()
Specify whether the streamline is integrated in the upstream or downstream direction,...
std::vector< CustomTerminationCallbackType > CustomTerminationCallback
static vtkStreamTracer * New()
Construct the object to start from position (0,0,0), with forward integration, terminal speed 1....
bool HasMatchingPointAttributes
vtkCompositeDataSet * InputData
void SetInterpolatorType(int interpType)
Set the type of the velocity field interpolator to determine whether INTERPOLATOR_WITH_DATASET_POINT_...
double MaximumIntegrationStep
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
bool GenerateNormalsInIntegrate
vtkExecutive * CreateDefaultExecutive() override
Create a default executive.
void Integrate(vtkPointData *inputData, vtkPolyData *output, vtkDataArray *seedSource, vtkIdList *seedIds, vtkIntArray *integrationDirections, vtkAbstractInterpolatedVelocityField *func, int maxCellSize, int vecType, const char *vecFieldName, double &propagation, vtkIdType &numSteps, double &integrationTime, std::vector< CustomTerminationCallbackType > &customTerminationCallback, std::vector< void * > &customTerminationClientData, std::vector< int > &customReasonForTermination)
void SetIntegrationDirectionToBackward()
Specify whether the streamline is integrated in the upstream or downstream direction,...
void SetInterpolatorTypeToDataSetPointLocator()
Set the velocity field interpolator type to one that uses a point locator to perform local spatial se...
int GetIntegratorType()
Set/get the integrator type to be used for streamline generation.
void AddCustomTerminationCallback(CustomTerminationCallbackType callback, void *clientdata, int reasonForTermination)
Adds a custom termination callback.
void InitializeSeeds(vtkDataArray *&seeds, vtkIdList *&seedIds, vtkIntArray *&integrationDirections, vtkDataSet *source)
void SetIntegratorTypeToRungeKutta2()
Set/get the integrator type to be used for streamline generation.
void SetIntegrationDirectionToBoth()
Specify whether the streamline is integrated in the upstream or downstream direction,...
double SimpleIntegrate(double seed[3], double lastPoint[3], double stepSize, vtkAbstractInterpolatedVelocityField *func)
~vtkStreamTracer() override
void AddInput(vtkDataObject *)
vtkInitialValueProblemSolver * Integrator
void SetInterpolatorPrototype(vtkAbstractInterpolatedVelocityField *ivf)
The object used to interpolate the velocity field during integration is of the same class as this pro...
void SetIntegrationStepUnit(int unit)
Specify a uniform integration step unit for MinimumIntegrationStep, InitialIntegrationStep,...
void SetIntegratorType(int type)
Set/get the integrator type to be used for streamline generation.
int GetIntegrationStepUnit()
boost::graph_traits< vtkGraph * >::vertex_descriptor source(boost::graph_traits< vtkGraph * >::edge_descriptor e, vtkGraph *)
bool(* CustomTerminationCallbackType)(void *clientdata, vtkPoints *points, vtkDataArray *velocity, int integrationDirection)
Used to specify custom conditions which are evaluated to determine whether a streamline should be ter...