Main Page | Class Hierarchy | Alphabetical List | Class List | File List | Class Members | File Members | Related Pages

vtkStreamTracer Class Reference

#include <vtkStreamTracer.h>

Inheritance diagram for vtkStreamTracer:

Inheritance graph
[legend]
Collaboration diagram for vtkStreamTracer:

Collaboration graph
[legend]
List of all members.

Detailed Description

Streamline generator.

Date
2003/01/09 19:21:05
Revision
1.13

vtkStreamTracer is a filter that integrates a vector field to generate streamlines. The integration is performed using the provided integrator. The default is second order Runge-Kutta.

vtkStreamTracer generate polylines as output. Each cell (polyline) corresponds to one streamline. The values associated with each streamline are stored in the cell data whereas the values associated with points are stored in point data.

Note that vtkStreamTracer can integrate both forward and backward. The length of the streamline is controlled by specifying either a maximum value in the units of length, cell length or elapsed time (the elapsed time is the time each particle would have traveled if flow were steady). Otherwise, the integration terminates after exiting the dataset or if the particle speed is reduced to a value less than the terminal speed or when a maximum number of steps is reached. The reason for the termination is stored in a cell array named ReasonForTermination.

The quality of integration can be controlled by setting integration step (InitialIntegrationStep) and in the case of adaptive solvers the maximum error, the minimum integration step and the maximum integration step. All of these can have units of length, cell length or elapsed time.

The integration time, vorticity, rotation and angular velocity are stored in point arrays named "IntegrationTime", "Vorticity", "Rotation" and "AngularVelocity" respectively (vorticity, rotation and angular velocity are computed only when ComputeVorticity is on). All point attributes in the source data set are interpolated on the new streamline points.

vtkStreamTracer integrates through any type of dataset. As a result, if the dataset contains 2D cells such as polygons or triangles, the integration is constrained to lie on the surface defined by the 2D cells.

The starting point of traces may be defined in two different ways. Starting from global x-y-z "position" allows you to start a single trace at a specified x-y-z coordinate. If you specify a source object, a trace will be generated for each point in the source that is inside the dataset.

See also:
vtkRibbonFilter vtkRuledSurfaceFilter vtkInitialValueProblemSolver vtkRungeKutta2 vtkRungeKutta4 vtkRungeKutta45
Tests:
vtkStreamTracer (Tests)

Definition at line 86 of file vtkStreamTracer.h.

Public Types

typedef vtkDataSetToPolyDataFilter Superclass
enum  Units { TIME_UNIT, LENGTH_UNIT, CELL_LENGTH_UNIT }
enum  Solvers {
  RUNGE_KUTTA2, RUNGE_KUTTA4, RUNGE_KUTTA45, NONE,
  UNKNOWN
}
enum  ReasonForTermination {
  OUT_OF_DOMAIN = vtkInitialValueProblemSolver::OUT_OF_DOMAIN, NOT_INITIALIZED = vtkInitialValueProblemSolver::NOT_INITIALIZED, UNEXPECTED_VALUE = vtkInitialValueProblemSolver::UNEXPECTED_VALUE, OUT_OF_TIME = 4,
  OUT_OF_STEPS = 5, STAGNATION = 6
}
enum  { FORWARD, BACKWARD, BOTH }

Public Member Functions

virtual const char * GetClassName ()
virtual int IsA (const char *type)
void PrintSelf (ostream &os, vtkIndent indent)
void AddInput (vtkDataSet *in)
virtual void SetStartPosition (float, float, float)
virtual void SetStartPosition (float[3])
virtual float * GetStartPosition ()
virtual void GetStartPosition (float &, float &, float &)
virtual void GetStartPosition (float[3])
void SetSource (vtkDataSet *source)
vtkDataSetGetSource ()
void SetIntegrator (vtkInitialValueProblemSolver *)
virtual vtkInitialValueProblemSolverGetIntegrator ()
void SetIntegratorType (int type)
int GetIntegratorType ()
void SetIntegratorTypeToRungeKutta2 ()
void SetIntegratorTypeToRungeKutta4 ()
void SetIntegratorTypeToRungeKutta45 ()
void SetMaximumPropagation (int unit, float max)
void SetMaximumPropagation (float max)
void SetMaximumPropagationUnit (int unit)
int GetMaximumPropagationUnit ()
float GetMaximumPropagation ()
void SetMaximumPropagationUnitToTimeUnit ()
void SetMaximumPropagationUnitToLengthUnit ()
void SetMaximumPropagationUnitToCellLengthUnit ()
void SetMinimumIntegrationStep (int unit, float step)
void SetMinimumIntegrationStepUnit (int unit)
void SetMinimumIntegrationStep (float step)
int GetMinimumIntegrationStepUnit ()
float GetMinimumIntegrationStep ()
void SetMinimumIntegrationStepUnitToTimeUnit ()
void SetMinimumIntegrationStepUnitToLengthUnit ()
void SetMinimumIntegrationStepUnitToCellLengthUnit ()
void SetMaximumIntegrationStep (int unit, float step)
void SetMaximumIntegrationStepUnit (int unit)
void SetMaximumIntegrationStep (float step)
int GetMaximumIntegrationStepUnit ()
float GetMaximumIntegrationStep ()
void SetMaximumIntegrationStepUnitToTimeUnit ()
void SetMaximumIntegrationStepUnitToLengthUnit ()
void SetMaximumIntegrationStepUnitToCellLengthUnit ()
void SetInitialIntegrationStep (int unit, float step)
void SetInitialIntegrationStepUnit (int unit)
void SetInitialIntegrationStep (float step)
int GetInitialIntegrationStepUnit ()
float GetInitialIntegrationStep ()
void SetInitialIntegrationStepUnitToTimeUnit ()
void SetInitialIntegrationStepUnitToLengthUnit ()
void SetInitialIntegrationStepUnitToCellLengthUnit ()
virtual void SetMaximumError (float)
virtual float GetMaximumError ()
virtual void SetMaximumNumberOfSteps (vtkIdType)
virtual vtkIdType GetMaximumNumberOfSteps ()
virtual void SetTerminalSpeed (float)
virtual float GetTerminalSpeed ()
virtual void SetIntegrationDirection (int)
virtual int GetIntegrationDirection ()
void SetIntegrationDirectionToForward ()
void SetIntegrationDirectionToBackward ()
void SetIntegrationDirectionToBoth ()
virtual void SetComputeVorticity (int)
virtual int GetComputeVorticity ()
virtual void ComputeVorticityOn ()
virtual void ComputeVorticityOff ()
virtual void SetRotationScale (float)
virtual float GetRotationScale ()

Static Public Member Functions

int IsTypeOf (const char *type)
vtkStreamTracerSafeDownCast (vtkObject *o)
vtkStreamTracerNew ()

Protected Member Functions

 vtkStreamTracer ()
 ~vtkStreamTracer ()
void AddInput (vtkDataObject *)
void Execute ()
void CalculateVorticity (vtkGenericCell *cell, float pcoords[3], vtkFloatArray *cellVectors, float vorticity[3])
void Integrate (vtkPolyData *output, vtkDataArray *seedSource, vtkIdList *seedIds, vtkIntArray *integrationDirections, float lastPoint[3])
int CheckInputs (vtkInterpolatedVelocityField *&func, int *maxCellSize)
virtual void SetInputVectorsSelection (const char *)
void SetIntervalInformation (int unit, float interval, IntervalInformation &currentValues)
void SetIntervalInformation (int unit, IntervalInformation &currentValues)
void ConvertIntervals (float &step, float &minStep, float &maxStep, int direction, float cellLength, float speed)
void InitializeSeeds (vtkDataArray *&seeds, vtkIdList *&seedIds, vtkIntArray *&integrationDirections)

Static Protected Member Functions

float ConvertToTime (IntervalInformation &interval, float cellLength, float speed)
float ConvertToLength (IntervalInformation &interval, float cellLength, float speed)
float ConvertToCellLength (IntervalInformation &interval, float cellLength, float speed)
float ConvertToUnit (IntervalInformation &interval, int unit, float cellLength, float speed)

Protected Attributes

char * InputVectorsSelection
float StartPosition [3]
float TerminalSpeed
IntervalInformation MaximumPropagation
IntervalInformation MinimumIntegrationStep
IntervalInformation MaximumIntegrationStep
IntervalInformation InitialIntegrationStep
int IntegrationDirection
vtkInitialValueProblemSolverIntegrator
float MaximumError
vtkIdType MaximumNumberOfSteps
int ComputeVorticity
float RotationScale

Static Protected Attributes

const float EPSILON


Member Typedef Documentation

typedef vtkDataSetToPolyDataFilter vtkStreamTracer::Superclass
 

Reimplemented from vtkDataSetToPolyDataFilter.

Reimplemented in vtkPStreamTracer.

Definition at line 89 of file vtkStreamTracer.h.


Member Enumeration Documentation

enum vtkStreamTracer::Units
 

Enumeration values:
TIME_UNIT 
LENGTH_UNIT 
CELL_LENGTH_UNIT 

Definition at line 113 of file vtkStreamTracer.h.

enum vtkStreamTracer::Solvers
 

Enumeration values:
RUNGE_KUTTA2 
RUNGE_KUTTA4 
RUNGE_KUTTA45 
NONE 
UNKNOWN 

Definition at line 120 of file vtkStreamTracer.h.

enum vtkStreamTracer::ReasonForTermination
 

Enumeration values:
OUT_OF_DOMAIN 
NOT_INITIALIZED 
UNEXPECTED_VALUE 
OUT_OF_TIME 
OUT_OF_STEPS 
STAGNATION 

Definition at line 129 of file vtkStreamTracer.h.

anonymous enum
 

Enumeration values:
FORWARD 
BACKWARD 
BOTH 

Definition at line 248 of file vtkStreamTracer.h.


Constructor & Destructor Documentation

vtkStreamTracer::vtkStreamTracer  )  [protected]
 

vtkStreamTracer::~vtkStreamTracer  )  [protected]
 


Member Function Documentation

virtual const char* vtkStreamTracer::GetClassName  )  [virtual]
 

Reimplemented from vtkDataSetToPolyDataFilter.

Reimplemented in vtkPStreamTracer.

int vtkStreamTracer::IsTypeOf const char *  type  )  [static]
 

Return 1 if this class type is the same type of (or a subclass of) the named class. Returns 0 otherwise. This method works in combination with vtkTypeRevisionMacro found in vtkSetGet.h.

Reimplemented from vtkDataSetToPolyDataFilter.

Reimplemented in vtkPStreamTracer.

virtual int vtkStreamTracer::IsA const char *  type  )  [virtual]
 

Return 1 if this class is the same type of (or a subclass of) the named class. Returns 0 otherwise. This method works in combination with vtkTypeRevisionMacro found in vtkSetGet.h.

Reimplemented from vtkDataSetToPolyDataFilter.

Reimplemented in vtkPStreamTracer.

vtkStreamTracer* vtkStreamTracer::SafeDownCast vtkObject o  )  [static]
 

Reimplemented from vtkDataSetToPolyDataFilter.

Reimplemented in vtkPStreamTracer.

void vtkStreamTracer::PrintSelf ostream &  os,
vtkIndent  indent
[virtual]
 

Methods invoked by print to print information about the object including superclasses. Typically not called by the user (use Print() instead) but used in the hierarchical print process to combine the output of several classes.

Reimplemented from vtkDataSetToPolyDataFilter.

Reimplemented in vtkPStreamTracer.

vtkStreamTracer* vtkStreamTracer::New  )  [static]
 

Construct object to start from position (0,0,0), integrate forward, terminal speed 1.0E-12, vorticity computation on, integration step length 0.5 (unit cell length), maximum number of steps 2000, using 2nd order Runge Kutta and maximum propagation 1.0 (unit length).

Reimplemented from vtkObject.

Reimplemented in vtkPStreamTracer.

virtual void vtkStreamTracer::SetStartPosition float  ,
float  ,
float 
[virtual]
 

Specify the start of the streamline in the global coordinate system. Search must be performed to find initial cell to start integration from.

virtual void vtkStreamTracer::SetStartPosition float  [3]  )  [virtual]
 

Specify the start of the streamline in the global coordinate system. Search must be performed to find initial cell to start integration from.

virtual float* vtkStreamTracer::GetStartPosition  )  [virtual]
 

Specify the start of the streamline in the global coordinate system. Search must be performed to find initial cell to start integration from.

virtual void vtkStreamTracer::GetStartPosition float &  ,
float &  ,
float & 
[virtual]
 

Specify the start of the streamline in the global coordinate system. Search must be performed to find initial cell to start integration from.

virtual void vtkStreamTracer::GetStartPosition float  [3]  )  [virtual]
 

Specify the start of the streamline in the global coordinate system. Search must be performed to find initial cell to start integration from.

void vtkStreamTracer::SetSource vtkDataSet source  ) 
 

Specify the source object used to generate starting points.

vtkDataSet* vtkStreamTracer::GetSource  ) 
 

Specify the source object used to generate starting points.

void vtkStreamTracer::SetIntegrator vtkInitialValueProblemSolver  ) 
 

Set/get the integrator type to be used in the stream line calculation. The object passed is not actually used but is cloned with NewInstance in the process of integration (prototype pattern). The default is 2nd order Runge Kutta. The integrator can also be changed using SetIntegratorType. The recognized solvers are: RUNGE_KUTTA2 = 0 RUNGE_KUTTA4 = 1 RUNGE_KUTTA45 = 2

virtual vtkInitialValueProblemSolver* vtkStreamTracer::GetIntegrator  )  [virtual]
 

Set/get the integrator type to be used in the stream line calculation. The object passed is not actually used but is cloned with NewInstance in the process of integration (prototype pattern). The default is 2nd order Runge Kutta. The integrator can also be changed using SetIntegratorType. The recognized solvers are: RUNGE_KUTTA2 = 0 RUNGE_KUTTA4 = 1 RUNGE_KUTTA45 = 2

void vtkStreamTracer::SetIntegratorType int  type  ) 
 

Set/get the integrator type to be used in the stream line calculation. The object passed is not actually used but is cloned with NewInstance in the process of integration (prototype pattern). The default is 2nd order Runge Kutta. The integrator can also be changed using SetIntegratorType. The recognized solvers are: RUNGE_KUTTA2 = 0 RUNGE_KUTTA4 = 1 RUNGE_KUTTA45 = 2

int vtkStreamTracer::GetIntegratorType  ) 
 

Set/get the integrator type to be used in the stream line calculation. The object passed is not actually used but is cloned with NewInstance in the process of integration (prototype pattern). The default is 2nd order Runge Kutta. The integrator can also be changed using SetIntegratorType. The recognized solvers are: RUNGE_KUTTA2 = 0 RUNGE_KUTTA4 = 1 RUNGE_KUTTA45 = 2

void vtkStreamTracer::SetIntegratorTypeToRungeKutta2  )  [inline]
 

Set/get the integrator type to be used in the stream line calculation. The object passed is not actually used but is cloned with NewInstance in the process of integration (prototype pattern). The default is 2nd order Runge Kutta. The integrator can also be changed using SetIntegratorType. The recognized solvers are: RUNGE_KUTTA2 = 0 RUNGE_KUTTA4 = 1 RUNGE_KUTTA45 = 2

Definition at line 151 of file vtkStreamTracer.h.

void vtkStreamTracer::SetIntegratorTypeToRungeKutta4  )  [inline]
 

Set/get the integrator type to be used in the stream line calculation. The object passed is not actually used but is cloned with NewInstance in the process of integration (prototype pattern). The default is 2nd order Runge Kutta. The integrator can also be changed using SetIntegratorType. The recognized solvers are: RUNGE_KUTTA2 = 0 RUNGE_KUTTA4 = 1 RUNGE_KUTTA45 = 2

Definition at line 153 of file vtkStreamTracer.h.

void vtkStreamTracer::SetIntegratorTypeToRungeKutta45  )  [inline]
 

Set/get the integrator type to be used in the stream line calculation. The object passed is not actually used but is cloned with NewInstance in the process of integration (prototype pattern). The default is 2nd order Runge Kutta. The integrator can also be changed using SetIntegratorType. The recognized solvers are: RUNGE_KUTTA2 = 0 RUNGE_KUTTA4 = 1 RUNGE_KUTTA45 = 2

Definition at line 155 of file vtkStreamTracer.h.

void vtkStreamTracer::SetMaximumPropagation int  unit,
float  max
 

Specify the maximum length of the streamlines expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2

void vtkStreamTracer::SetMaximumPropagation float  max  ) 
 

Specify the maximum length of the streamlines expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2

void vtkStreamTracer::SetMaximumPropagationUnit int  unit  ) 
 

Specify the maximum length of the streamlines expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2

int vtkStreamTracer::GetMaximumPropagationUnit  ) 
 

Specify the maximum length of the streamlines expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2

float vtkStreamTracer::GetMaximumPropagation  ) 
 

Specify the maximum length of the streamlines expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2

void vtkStreamTracer::SetMaximumPropagationUnitToTimeUnit  )  [inline]
 

Specify the maximum length of the streamlines expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2

Definition at line 167 of file vtkStreamTracer.h.

void vtkStreamTracer::SetMaximumPropagationUnitToLengthUnit  )  [inline]
 

Specify the maximum length of the streamlines expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2

Definition at line 169 of file vtkStreamTracer.h.

void vtkStreamTracer::SetMaximumPropagationUnitToCellLengthUnit  )  [inline]
 

Specify the maximum length of the streamlines expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2

Definition at line 171 of file vtkStreamTracer.h.

void vtkStreamTracer::SetMinimumIntegrationStep int  unit,
float  step
 

Specify the minimum step used in the integration expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 Only valid when using adaptive integrators.

void vtkStreamTracer::SetMinimumIntegrationStepUnit int  unit  ) 
 

Specify the minimum step used in the integration expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 Only valid when using adaptive integrators.

void vtkStreamTracer::SetMinimumIntegrationStep float  step  ) 
 

Specify the minimum step used in the integration expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 Only valid when using adaptive integrators.

int vtkStreamTracer::GetMinimumIntegrationStepUnit  ) 
 

Specify the minimum step used in the integration expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 Only valid when using adaptive integrators.

float vtkStreamTracer::GetMinimumIntegrationStep  ) 
 

Specify the minimum step used in the integration expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 Only valid when using adaptive integrators.

void vtkStreamTracer::SetMinimumIntegrationStepUnitToTimeUnit  )  [inline]
 

Specify the minimum step used in the integration expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 Only valid when using adaptive integrators.

Definition at line 184 of file vtkStreamTracer.h.

void vtkStreamTracer::SetMinimumIntegrationStepUnitToLengthUnit  )  [inline]
 

Specify the minimum step used in the integration expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 Only valid when using adaptive integrators.

Definition at line 186 of file vtkStreamTracer.h.

void vtkStreamTracer::SetMinimumIntegrationStepUnitToCellLengthUnit  )  [inline]
 

Specify the minimum step used in the integration expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 Only valid when using adaptive integrators.

Definition at line 188 of file vtkStreamTracer.h.

void vtkStreamTracer::SetMaximumIntegrationStep int  unit,
float  step
 

Specify the maximum step used in the integration expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 Only valid when using adaptive integrators.

void vtkStreamTracer::SetMaximumIntegrationStepUnit int  unit  ) 
 

Specify the maximum step used in the integration expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 Only valid when using adaptive integrators.

void vtkStreamTracer::SetMaximumIntegrationStep float  step  ) 
 

Specify the maximum step used in the integration expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 Only valid when using adaptive integrators.

int vtkStreamTracer::GetMaximumIntegrationStepUnit  ) 
 

Specify the maximum step used in the integration expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 Only valid when using adaptive integrators.

float vtkStreamTracer::GetMaximumIntegrationStep  ) 
 

Specify the maximum step used in the integration expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 Only valid when using adaptive integrators.

void vtkStreamTracer::SetMaximumIntegrationStepUnitToTimeUnit  )  [inline]
 

Specify the maximum step used in the integration expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 Only valid when using adaptive integrators.

Definition at line 201 of file vtkStreamTracer.h.

void vtkStreamTracer::SetMaximumIntegrationStepUnitToLengthUnit  )  [inline]
 

Specify the maximum step used in the integration expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 Only valid when using adaptive integrators.

Definition at line 203 of file vtkStreamTracer.h.

void vtkStreamTracer::SetMaximumIntegrationStepUnitToCellLengthUnit  )  [inline]
 

Specify the maximum step used in the integration expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 Only valid when using adaptive integrators.

Definition at line 205 of file vtkStreamTracer.h.

void vtkStreamTracer::SetInitialIntegrationStep int  unit,
float  step
 

Specify the initial step used in the integration expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 If the integrator is not adaptive, this is the actual step used.

void vtkStreamTracer::SetInitialIntegrationStepUnit int  unit  ) 
 

Specify the initial step used in the integration expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 If the integrator is not adaptive, this is the actual step used.

void vtkStreamTracer::SetInitialIntegrationStep float  step  ) 
 

Specify the initial step used in the integration expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 If the integrator is not adaptive, this is the actual step used.

int vtkStreamTracer::GetInitialIntegrationStepUnit  ) 
 

Specify the initial step used in the integration expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 If the integrator is not adaptive, this is the actual step used.

float vtkStreamTracer::GetInitialIntegrationStep  ) 
 

Specify the initial step used in the integration expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 If the integrator is not adaptive, this is the actual step used.

void vtkStreamTracer::SetInitialIntegrationStepUnitToTimeUnit  )  [inline]
 

Specify the initial step used in the integration expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 If the integrator is not adaptive, this is the actual step used.

Definition at line 218 of file vtkStreamTracer.h.

void vtkStreamTracer::SetInitialIntegrationStepUnitToLengthUnit  )  [inline]
 

Specify the initial step used in the integration expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 If the integrator is not adaptive, this is the actual step used.

Definition at line 220 of file vtkStreamTracer.h.

void vtkStreamTracer::SetInitialIntegrationStepUnitToCellLengthUnit  )  [inline]
 

Specify the initial step used in the integration expressed in one of the: TIME_UNIT = 0 LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 If the integrator is not adaptive, this is the actual step used.

Definition at line 222 of file vtkStreamTracer.h.

virtual void vtkStreamTracer::SetMaximumError float   )  [virtual]
 

Specify the maximum error in the integration. This value is passed to the integrator. Therefore, it's meaning depends on the integrator used.

virtual float vtkStreamTracer::GetMaximumError  )  [virtual]
 

Specify the maximum error in the integration. This value is passed to the integrator. Therefore, it's meaning depends on the integrator used.

virtual void vtkStreamTracer::SetMaximumNumberOfSteps vtkIdType   )  [virtual]
 

Specify the maximum number of steps used in the integration.

virtual vtkIdType vtkStreamTracer::GetMaximumNumberOfSteps  )  [virtual]
 

Specify the maximum number of steps used in the integration.

virtual void vtkStreamTracer::SetTerminalSpeed float   )  [virtual]
 

If at any point, the speed is below this value, the integration is terminated.

virtual float vtkStreamTracer::GetTerminalSpeed  )  [virtual]
 

If at any point, the speed is below this value, the integration is terminated.

virtual void vtkStreamTracer::SetIntegrationDirection int   )  [virtual]
 

Specify whether the streamtrace will be generated in the upstream or downstream direction.

virtual int vtkStreamTracer::GetIntegrationDirection  )  [virtual]
 

Specify whether the streamtrace will be generated in the upstream or downstream direction.

void vtkStreamTracer::SetIntegrationDirectionToForward  )  [inline]
 

Specify whether the streamtrace will be generated in the upstream or downstream direction.

Definition at line 261 of file vtkStreamTracer.h.

void vtkStreamTracer::SetIntegrationDirectionToBackward  )  [inline]
 

Specify whether the streamtrace will be generated in the upstream or downstream direction.

Definition at line 263 of file vtkStreamTracer.h.

void vtkStreamTracer::SetIntegrationDirectionToBoth  )  [inline]
 

Specify whether the streamtrace will be generated in the upstream or downstream direction.

Definition at line 265 of file vtkStreamTracer.h.

virtual void vtkStreamTracer::SetComputeVorticity int   )  [virtual]
 

Turn on/off calculation of vorticity at streamline points (necessary for generating proper streamribbons using the vtkRibbonFilter.

virtual int vtkStreamTracer::GetComputeVorticity  )  [virtual]
 

Turn on/off calculation of vorticity at streamline points (necessary for generating proper streamribbons using the vtkRibbonFilter.

virtual void vtkStreamTracer::ComputeVorticityOn  )  [virtual]
 

Turn on/off calculation of vorticity at streamline points (necessary for generating proper streamribbons using the vtkRibbonFilter.

virtual void vtkStreamTracer::ComputeVorticityOff  )  [virtual]
 

Turn on/off calculation of vorticity at streamline points (necessary for generating proper streamribbons using the vtkRibbonFilter.

virtual void vtkStreamTracer::SetRotationScale float   )  [virtual]
 

This can be used to scale the rate with which the streamribbons twist. The default is 1.

virtual float vtkStreamTracer::GetRotationScale  )  [virtual]
 

This can be used to scale the rate with which the streamribbons twist. The default is 1.

void vtkStreamTracer::AddInput vtkDataSet in  ) 
 

Add a dataset to the list inputs

void vtkStreamTracer::AddInput vtkDataObject  )  [inline, protected, virtual]
 

Reimplemented from vtkProcessObject.

Definition at line 293 of file vtkStreamTracer.h.

void vtkStreamTracer::Execute  )  [protected, virtual]
 

This method is the old style execute method

Reimplemented from vtkSource.

Reimplemented in vtkPStreamTracer.

void vtkStreamTracer::CalculateVorticity vtkGenericCell cell,
float  pcoords[3],
vtkFloatArray cellVectors,
float  vorticity[3]
[protected]
 

void vtkStreamTracer::Integrate vtkPolyData output,
vtkDataArray seedSource,
vtkIdList seedIds,
vtkIntArray integrationDirections,
float  lastPoint[3]
[protected]
 

int vtkStreamTracer::CheckInputs vtkInterpolatedVelocityField *&  func,
int *  maxCellSize
[protected]
 

virtual void vtkStreamTracer::SetInputVectorsSelection const char *   )  [protected, virtual]
 

Referenced by vtkPStreamTracer::SelectInputVectors().

void vtkStreamTracer::SetIntervalInformation int  unit,
float  interval,
IntervalInformation currentValues
[protected]
 

void vtkStreamTracer::SetIntervalInformation int  unit,
IntervalInformation currentValues
[protected]
 

float vtkStreamTracer::ConvertToTime IntervalInformation interval,
float  cellLength,
float  speed
[static, protected]
 

float vtkStreamTracer::ConvertToLength IntervalInformation interval,
float  cellLength,
float  speed
[static, protected]
 

float vtkStreamTracer::ConvertToCellLength IntervalInformation interval,
float  cellLength,
float  speed
[static, protected]
 

float vtkStreamTracer::ConvertToUnit IntervalInformation interval,
int  unit,
float  cellLength,
float  speed
[static, protected]
 

void vtkStreamTracer::ConvertIntervals float &  step,
float &  minStep,
float &  maxStep,
int  direction,
float  cellLength,
float  speed
[protected]
 

void vtkStreamTracer::InitializeSeeds vtkDataArray *&  seeds,
vtkIdList *&  seedIds,
vtkIntArray *&  integrationDirections
[protected]
 


Member Data Documentation

char* vtkStreamTracer::InputVectorsSelection [protected]
 

Definition at line 308 of file vtkStreamTracer.h.

float vtkStreamTracer::StartPosition[3] [protected]
 

Definition at line 312 of file vtkStreamTracer.h.

const float vtkStreamTracer::EPSILON [static, protected]
 

Definition at line 314 of file vtkStreamTracer.h.

float vtkStreamTracer::TerminalSpeed [protected]
 

Definition at line 315 of file vtkStreamTracer.h.

IntervalInformation vtkStreamTracer::MaximumPropagation [protected]
 

Definition at line 324 of file vtkStreamTracer.h.

IntervalInformation vtkStreamTracer::MinimumIntegrationStep [protected]
 

Definition at line 325 of file vtkStreamTracer.h.

IntervalInformation vtkStreamTracer::MaximumIntegrationStep [protected]
 

Definition at line 326 of file vtkStreamTracer.h.

IntervalInformation vtkStreamTracer::InitialIntegrationStep [protected]
 

Definition at line 327 of file vtkStreamTracer.h.

int vtkStreamTracer::IntegrationDirection [protected]
 

Definition at line 348 of file vtkStreamTracer.h.

vtkInitialValueProblemSolver* vtkStreamTracer::Integrator [protected]
 

Definition at line 351 of file vtkStreamTracer.h.

float vtkStreamTracer::MaximumError [protected]
 

Definition at line 353 of file vtkStreamTracer.h.

vtkIdType vtkStreamTracer::MaximumNumberOfSteps [protected]
 

Definition at line 354 of file vtkStreamTracer.h.

int vtkStreamTracer::ComputeVorticity [protected]
 

Definition at line 356 of file vtkStreamTracer.h.

float vtkStreamTracer::RotationScale [protected]
 

Definition at line 357 of file vtkStreamTracer.h.


The documentation for this class was generated from the following file: