VTK
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
vtkTemporalPathLineFilter Class Reference

Generate a Polydata Pointset from any Dataset. More...

#include <vtkTemporalPathLineFilter.h>

Inherits vtkPolyDataAlgorithm.

Public Member Functions

void Flush ()
 
void SetSelectionConnection (vtkAlgorithmOutput *algOutput)
 
void SetSelectionData (vtkDataSet *input)
 
virtual void SetMaskPoints (int)
 
virtual int GetMaskPoints ()
 
virtual void SetMaxTrackLength (unsigned int)
 
virtual unsigned int GetMaxTrackLength ()
 
virtual void SetIdChannelArray (const char *)
 
virtual char * GetIdChannelArray ()
 
virtual void SetMaxStepDistance (double, double, double)
 
virtual void SetMaxStepDistance (double[3])
 
virtual double * GetMaxStepDistance ()
 
virtual void GetMaxStepDistance (double &, double &, double &)
 
virtual void GetMaxStepDistance (double[3])
 
virtual void SetKeepDeadTrails (int)
 
virtual int GetKeepDeadTrails ()
 

Protected Member Functions

 vtkTemporalPathLineFilter ()
 
 ~vtkTemporalPathLineFilter ()
 
virtual int FillInputPortInformation (int port, vtkInformation *info)
 
virtual int FillOutputPortInformation (int port, vtkInformation *info)
 
TrailPointer GetTrail (vtkIdType i)
 
void IncrementTrail (TrailPointer trail, vtkDataSet *input, vtkIdType i)
 
virtual int RequestInformation (vtkInformation *, vtkInformationVector **, vtkInformationVector *)
 
virtual int RequestData (vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
 

Protected Attributes

int NumberOfTimeSteps
 
int MaskPoints
 
unsigned int MaxTrackLength
 
unsigned int LastTrackLength
 
int FirstTime
 
char * IdChannelArray
 
double MaxStepDistance [3]
 
double LatestTime
 
int KeepDeadTrails
 
int UsingSelection
 
vtkSmartPointer< vtkCellArray > PolyLines
 
vtkSmartPointer< vtkCellArray > Vertices
 
vtkSmartPointer< vtkPoints > LineCoordinates
 
vtkSmartPointer< vtkPoints > VertexCoordinates
 
vtkSmartPointer< vtkFloatArray > TrailId
 
vtkSmartPointer
< vtkTemporalPathLineFilterInternals > 
Internals
 
std::set< vtkIdType > SelectionIds
 
typedef vtkPolyDataAlgorithm Superclass
 
static vtkTemporalPathLineFilterNew ()
 
static int IsTypeOf (const char *type)
 
static vtkTemporalPathLineFilterSafeDownCast (vtkObjectBase *o)
 
virtual int IsA (const char *type)
 
vtkTemporalPathLineFilterNewInstance () const
 
void PrintSelf (ostream &os, vtkIndent indent)
 
virtual vtkObjectBase * NewInstanceInternal () const
 

Detailed Description

Generate a Polydata Pointset from any Dataset.

vtkTemporalPathLineFilter takes any dataset as input, it extracts the point locations of all cells over time to build up a polyline trail. The point number (index) is used as the 'key' if the points are randomly changing their respective order in the points list, then you should specify a scalar that represents the unique ID. This is intended to handle the output of a filter such as the TemporalStreamTracer.

See Also
vtkTemporalStreamTracer
Thanks:
John Bidiscombe of CSCS - Swiss National Supercomputing Centre for creating and contributing this class.

Definition at line 54 of file vtkTemporalPathLineFilter.h.

Member Typedef Documentation

typedef vtkPolyDataAlgorithm vtkTemporalPathLineFilter::Superclass

Standard Type-Macro

Definition at line 59 of file vtkTemporalPathLineFilter.h.

Constructor & Destructor Documentation

vtkTemporalPathLineFilter::vtkTemporalPathLineFilter ( )
protected
vtkTemporalPathLineFilter::~vtkTemporalPathLineFilter ( )
protected

Member Function Documentation

static vtkTemporalPathLineFilter* vtkTemporalPathLineFilter::New ( )
static

Standard Type-Macro

static int vtkTemporalPathLineFilter::IsTypeOf ( const char *  type)
static

Standard Type-Macro

virtual int vtkTemporalPathLineFilter::IsA ( const char *  type)
virtual

Standard Type-Macro

static vtkTemporalPathLineFilter* vtkTemporalPathLineFilter::SafeDownCast ( vtkObjectBase *  o)
static

Standard Type-Macro

virtual vtkObjectBase* vtkTemporalPathLineFilter::NewInstanceInternal ( ) const
protectedvirtual

Standard Type-Macro

vtkTemporalPathLineFilter* vtkTemporalPathLineFilter::NewInstance ( ) const

Standard Type-Macro

void vtkTemporalPathLineFilter::PrintSelf ( ostream &  os,
vtkIndent  indent 
)

Standard Type-Macro

virtual void vtkTemporalPathLineFilter::SetMaskPoints ( int  )
virtual

Set the number of particles to track as a ratio of the input example: setting MaskPoints to 10 will track every 10th point

virtual int vtkTemporalPathLineFilter::GetMaskPoints ( )
virtual

Set the number of particles to track as a ratio of the input example: setting MaskPoints to 10 will track every 10th point

virtual void vtkTemporalPathLineFilter::SetMaxTrackLength ( unsigned  int)
virtual

If the Particles being traced animate for a long time, the trails or traces will become long and stringy. Setting the MaxTraceTimeLength will limit how much of the trace is displayed. Tracks longer then the Max will disappear and the trace will apppear like a snake of fixed length which progresses as the particle moves

virtual unsigned int vtkTemporalPathLineFilter::GetMaxTrackLength ( )
virtual

If the Particles being traced animate for a long time, the trails or traces will become long and stringy. Setting the MaxTraceTimeLength will limit how much of the trace is displayed. Tracks longer then the Max will disappear and the trace will apppear like a snake of fixed length which progresses as the particle moves

virtual void vtkTemporalPathLineFilter::SetIdChannelArray ( const char *  )
virtual

Specify the name of a scalar array which will be used to fetch the index of each point. This is necessary only if the particles change position (Id order) on each time step. The Id can be used to identify particles at each step and hence track them properly. If this array is NULL, the global point ids are used. If an Id array cannot otherwise be found, the point index is used as the ID.

virtual char* vtkTemporalPathLineFilter::GetIdChannelArray ( )
virtual

Specify the name of a scalar array which will be used to fetch the index of each point. This is necessary only if the particles change position (Id order) on each time step. The Id can be used to identify particles at each step and hence track them properly. If this array is NULL, the global point ids are used. If an Id array cannot otherwise be found, the point index is used as the ID.

virtual void vtkTemporalPathLineFilter::SetMaxStepDistance ( double  ,
double  ,
double   
)
virtual

If a particle disappears from one end of a simulation and reappears on the other side, the track left will be unrepresentative. Set a MaxStepDistance{x,y,z} which acts as a threshold above which if a step occurs larger than the value (for the dimension), the track will be dropped and restarted after the step. (ie the part before the wrap around will be dropped and the newer part kept).

virtual void vtkTemporalPathLineFilter::SetMaxStepDistance ( double  [3])
virtual

If a particle disappears from one end of a simulation and reappears on the other side, the track left will be unrepresentative. Set a MaxStepDistance{x,y,z} which acts as a threshold above which if a step occurs larger than the value (for the dimension), the track will be dropped and restarted after the step. (ie the part before the wrap around will be dropped and the newer part kept).

virtual double* vtkTemporalPathLineFilter::GetMaxStepDistance ( )
virtual

If a particle disappears from one end of a simulation and reappears on the other side, the track left will be unrepresentative. Set a MaxStepDistance{x,y,z} which acts as a threshold above which if a step occurs larger than the value (for the dimension), the track will be dropped and restarted after the step. (ie the part before the wrap around will be dropped and the newer part kept).

virtual void vtkTemporalPathLineFilter::GetMaxStepDistance ( double &  ,
double &  ,
double &   
)
virtual

If a particle disappears from one end of a simulation and reappears on the other side, the track left will be unrepresentative. Set a MaxStepDistance{x,y,z} which acts as a threshold above which if a step occurs larger than the value (for the dimension), the track will be dropped and restarted after the step. (ie the part before the wrap around will be dropped and the newer part kept).

virtual void vtkTemporalPathLineFilter::GetMaxStepDistance ( double  [3])
virtual

If a particle disappears from one end of a simulation and reappears on the other side, the track left will be unrepresentative. Set a MaxStepDistance{x,y,z} which acts as a threshold above which if a step occurs larger than the value (for the dimension), the track will be dropped and restarted after the step. (ie the part before the wrap around will be dropped and the newer part kept).

virtual void vtkTemporalPathLineFilter::SetKeepDeadTrails ( int  )
virtual

When a particle 'disappears', the trail belonging to it is removed from the list. When this flag is enabled, dead trails will persist until the next time the list is cleared. Use carefully as it may cause excessive memory consumption if left on by mistake.

virtual int vtkTemporalPathLineFilter::GetKeepDeadTrails ( )
virtual

When a particle 'disappears', the trail belonging to it is removed from the list. When this flag is enabled, dead trails will persist until the next time the list is cleared. Use carefully as it may cause excessive memory consumption if left on by mistake.

void vtkTemporalPathLineFilter::Flush ( )

Flush will wipe any existing data so that traces can be restarted from whatever time step is next supplied.

void vtkTemporalPathLineFilter::SetSelectionConnection ( vtkAlgorithmOutput *  algOutput)

Set a second input which is a selection. Particles with the same Id in the selection as the primary input will be chosen for pathlines Note that you must have the same IdChannelArray in the selection as the input

void vtkTemporalPathLineFilter::SetSelectionData ( vtkDataSet *  input)

Set a second input which is a selection. Particles with the same Id in the selection as the primary input will be chosen for pathlines Note that you must have the same IdChannelArray in the selection as the input

virtual int vtkTemporalPathLineFilter::FillInputPortInformation ( int  port,
vtkInformation *  info 
)
protectedvirtual
virtual int vtkTemporalPathLineFilter::FillOutputPortInformation ( int  port,
vtkInformation *  info 
)
protectedvirtual
virtual int vtkTemporalPathLineFilter::RequestInformation ( vtkInformation *  ,
vtkInformationVector **  ,
vtkInformationVector *   
)
protectedvirtual

The necessary parts of the standard pipeline update mechanism

virtual int vtkTemporalPathLineFilter::RequestData ( vtkInformation *  request,
vtkInformationVector **  inputVector,
vtkInformationVector *  outputVector 
)
protectedvirtual

The necessary parts of the standard pipeline update mechanism

TrailPointer vtkTemporalPathLineFilter::GetTrail ( vtkIdType  i)
protected
void vtkTemporalPathLineFilter::IncrementTrail ( TrailPointer  trail,
vtkDataSet *  input,
vtkIdType  i 
)
protected

Member Data Documentation

int vtkTemporalPathLineFilter::NumberOfTimeSteps
protected

Definition at line 155 of file vtkTemporalPathLineFilter.h.

int vtkTemporalPathLineFilter::MaskPoints
protected

Definition at line 156 of file vtkTemporalPathLineFilter.h.

unsigned int vtkTemporalPathLineFilter::MaxTrackLength
protected

Definition at line 157 of file vtkTemporalPathLineFilter.h.

unsigned int vtkTemporalPathLineFilter::LastTrackLength
protected

Definition at line 158 of file vtkTemporalPathLineFilter.h.

int vtkTemporalPathLineFilter::FirstTime
protected

Definition at line 159 of file vtkTemporalPathLineFilter.h.

char* vtkTemporalPathLineFilter::IdChannelArray
protected

Definition at line 160 of file vtkTemporalPathLineFilter.h.

double vtkTemporalPathLineFilter::MaxStepDistance[3]
protected

Definition at line 161 of file vtkTemporalPathLineFilter.h.

double vtkTemporalPathLineFilter::LatestTime
protected

Definition at line 162 of file vtkTemporalPathLineFilter.h.

int vtkTemporalPathLineFilter::KeepDeadTrails
protected

Definition at line 163 of file vtkTemporalPathLineFilter.h.

int vtkTemporalPathLineFilter::UsingSelection
protected

Definition at line 164 of file vtkTemporalPathLineFilter.h.

vtkSmartPointer<vtkCellArray> vtkTemporalPathLineFilter::PolyLines
protected

Definition at line 167 of file vtkTemporalPathLineFilter.h.

vtkSmartPointer<vtkCellArray> vtkTemporalPathLineFilter::Vertices
protected

Definition at line 168 of file vtkTemporalPathLineFilter.h.

vtkSmartPointer<vtkPoints> vtkTemporalPathLineFilter::LineCoordinates
protected

Definition at line 169 of file vtkTemporalPathLineFilter.h.

vtkSmartPointer<vtkPoints> vtkTemporalPathLineFilter::VertexCoordinates
protected

Definition at line 170 of file vtkTemporalPathLineFilter.h.

vtkSmartPointer<vtkFloatArray> vtkTemporalPathLineFilter::TrailId
protected

Definition at line 171 of file vtkTemporalPathLineFilter.h.

vtkSmartPointer<vtkTemporalPathLineFilterInternals> vtkTemporalPathLineFilter::Internals
protected

Definition at line 172 of file vtkTemporalPathLineFilter.h.

std::set<vtkIdType> vtkTemporalPathLineFilter::SelectionIds
protected

Definition at line 173 of file vtkTemporalPathLineFilter.h.


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