VTK
|
Streamline generator. More...
#include <vtkStreamTracer.h>
Classes | |
struct | IntervalInformation |
Public Types | |
enum | Units { LENGTH_UNIT = 1, CELL_LENGTH_UNIT = 2 } |
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_LENGTH = 4, OUT_OF_STEPS = 5, STAGNATION = 6 } |
enum | { FORWARD, BACKWARD, BOTH } |
enum | { INTERPOLATOR_WITH_DATASET_POINT_LOCATOR, INTERPOLATOR_WITH_CELL_LOCATOR } |
typedef vtkPolyDataAlgorithm | Superclass |
Static Public Member Functions | |
static int | IsTypeOf (const char *type) |
static vtkStreamTracer * | SafeDownCast (vtkObjectBase *o) |
static vtkStreamTracer * | New () |
Protected Member Functions | |
virtual vtkObjectBase * | NewInstanceInternal () const |
vtkStreamTracer () | |
~vtkStreamTracer () | |
virtual vtkExecutive * | CreateDefaultExecutive () |
void | AddInput (vtkDataObject *) |
virtual int | RequestData (vtkInformation *, vtkInformationVector **, vtkInformationVector *) |
virtual int | FillInputPortInformation (int, vtkInformation *) |
void | CalculateVorticity (vtkGenericCell *cell, double pcoords[3], vtkDoubleArray *cellVectors, double vorticity[3]) |
void | Integrate (vtkPointData *inputData, vtkPolyData *output, vtkDataArray *seedSource, vtkIdList *seedIds, vtkIntArray *integrationDirections, double lastPoint[3], vtkAbstractInterpolatedVelocityField *func, int maxCellSize, int vecType, const char *vecFieldName, double &propagation, vtkIdType &numSteps) |
void | SimpleIntegrate (double seed[3], double lastPoint[3], double stepSize, vtkAbstractInterpolatedVelocityField *func) |
int | CheckInputs (vtkAbstractInterpolatedVelocityField *&func, int *maxCellSize) |
void | GenerateNormals (vtkPolyData *output, double *firstNormal, const char *vecName) |
void | ConvertIntervals (double &step, double &minStep, double &maxStep, int direction, double cellLength) |
int | SetupOutput (vtkInformation *inInfo, vtkInformation *outInfo) |
void | InitializeSeeds (vtkDataArray *&seeds, vtkIdList *&seedIds, vtkIntArray *&integrationDirections, vtkDataSet *source) |
Static Protected Member Functions | |
static double | ConvertToLength (double interval, int unit, double cellLength) |
static double | ConvertToLength (IntervalInformation &interval, double cellLength) |
Protected Attributes | |
bool | GenerateNormalsInIntegrate |
double | StartPosition [3] |
double | TerminalSpeed |
double | LastUsedStepSize |
double | MaximumPropagation |
double | MinimumIntegrationStep |
double | MaximumIntegrationStep |
double | InitialIntegrationStep |
int | IntegrationStepUnit |
int | IntegrationDirection |
vtkInitialValueProblemSolver * | Integrator |
double | MaximumError |
vtkIdType | MaximumNumberOfSteps |
bool | ComputeVorticity |
double | RotationScale |
vtkAbstractInterpolatedVelocityField * | InterpolatorPrototype |
vtkCompositeDataSet * | InputData |
bool | HasMatchingPointAttributes |
Static Protected Attributes | |
static const double | EPSILON |
Friends | |
class | PStreamTracerUtils |
Streamline generator.
vtkStreamTracer is a filter that integrates a vector field to generate streamlines. The integration is performed using a specified integrator, by default Runge-Kutta2.
vtkStreamTracer produces polylines as the output, with each cell (i.e., polyline) representing a streamline. The attribute values associated with each streamline are stored in the cell data, whereas those associated with streamline-points are stored in the point data.
vtkStreamTracer supports forward (the default), backward, and combined (i.e., BOTH) integration. The length of a streamline is governed by specifying a maximum value either in physical arc length or in (local) cell length. Otherwise, the integration terminates upon exiting the flow field domain, or if the particle speed is reduced to a value less than a specified terminal speed, or when a maximum number of steps is completed. The specific reason for the termination is stored in a cell array named ReasonForTermination.
Note that normalized vectors are adopted in streamline integration, which achieves high numerical accuracy/smoothness of flow lines that is particularly guaranteed for Runge-Kutta45 with adaptive step size and error control). In support of this feature, the underlying step size is ALWAYS in arc length unit (LENGTH_UNIT) while the 'real' time interval (virtual for steady flows) that a particle actually takes to trave in a single step is obtained by dividing the arc length by the LOCAL speed. The overall elapsed time (i.e., the life span) of the particle is the sum of those individual step-wise time intervals.
The quality of streamline integration can be controlled by setting the initial integration step (InitialIntegrationStep), particularly for Runge-Kutta2 and Runge-Kutta4 (with a fixed step size), and in the case of Runge-Kutta45 (with an adaptive step size and error control) the minimum integration step, the maximum integration step, and the maximum error. These steps are in either LENGTH_UNIT or CELL_LENGTH_UNIT while the error is in physical arc length. For the former two integrators, there is a trade-off between integration speed and streamline quality.
The integration time, vorticity, rotation and angular velocity are stored in point data arrays named "IntegrationTime", "Vorticity", "Rotation" and "AngularVelocity", respectively (vorticity, rotation and angular velocity are computed only when ComputeVorticity is on). All point data attributes in the source dataset are interpolated on the new streamline points.
vtkStreamTracer supports integration through any type of dataset. Thus if the dataset contains 2D cells like polygons or triangles, the integration is constrained to lie on the surface defined by 2D cells.
The starting point, or the so-called 'seed', of a streamline may be set 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, traces will be generated from each point in the source that is inside the dataset.
Definition at line 102 of file vtkStreamTracer.h.
typedef vtkPolyDataAlgorithm vtkStreamTracer::Superclass |
Definition at line 105 of file vtkStreamTracer.h.
Enumerator | |
---|---|
LENGTH_UNIT | |
CELL_LENGTH_UNIT |
Definition at line 144 of file vtkStreamTracer.h.
Enumerator | |
---|---|
RUNGE_KUTTA2 | |
RUNGE_KUTTA4 | |
RUNGE_KUTTA45 | |
NONE | |
UNKNOWN |
Definition at line 150 of file vtkStreamTracer.h.
Enumerator | |
---|---|
OUT_OF_DOMAIN | |
NOT_INITIALIZED | |
UNEXPECTED_VALUE | |
OUT_OF_LENGTH | |
OUT_OF_STEPS | |
STAGNATION |
Definition at line 159 of file vtkStreamTracer.h.
anonymous enum |
Enumerator | |
---|---|
FORWARD | |
BACKWARD | |
BOTH |
Definition at line 258 of file vtkStreamTracer.h.
anonymous enum |
Enumerator | |
---|---|
INTERPOLATOR_WITH_DATASET_POINT_LOCATOR | |
INTERPOLATOR_WITH_CELL_LOCATOR |
Definition at line 265 of file vtkStreamTracer.h.
|
protected |
|
protected |
|
static |
|
virtual |
Reimplemented in vtkTemporalStreamTracer, vtkPTemporalStreamTracer, and vtkPStreamTracer.
|
static |
|
protectedvirtual |
Reimplemented in vtkTemporalStreamTracer, vtkPTemporalStreamTracer, and vtkPStreamTracer.
vtkStreamTracer* vtkStreamTracer::NewInstance | ( | ) | const |
void vtkStreamTracer::PrintSelf | ( | ostream & | os, |
vtkIndent | indent | ||
) |
|
static |
Construct object to start from position (0,0,0), with forward integration, terminal speed 1.0E-12, vorticity computation on, integration step size 0.5 (in cell length unit), maximum number of steps 2000, using Runge-Kutta2, and maximum propagation 1.0 (in arc length unit).
|
virtual |
Specify the starting point (seed) of a streamline in the global coordinate system. Search must be performed to find the initial cell from which to start integration.
|
virtual |
Specify the starting point (seed) of a streamline in the global coordinate system. Search must be performed to find the initial cell from which to start integration.
|
virtual |
Specify the starting point (seed) of a streamline in the global coordinate system. Search must be performed to find the initial cell from which to start integration.
|
virtual |
Specify the starting point (seed) of a streamline in the global coordinate system. Search must be performed to find the initial cell from which to start integration.
|
virtual |
Specify the starting point (seed) of a streamline in the global coordinate system. Search must be performed to find the initial cell from which to start integration.
void vtkStreamTracer::SetSourceData | ( | vtkDataSet * | source | ) |
Specify the source object used to generate starting points (seeds). Old style. Do not use.
vtkDataSet* vtkStreamTracer::GetSource | ( | ) |
Specify the source object used to generate starting points (seeds). Old style. Do not use.
void vtkStreamTracer::SetSourceConnection | ( | vtkAlgorithmOutput * | algOutput | ) |
Specify the source object used to generate starting points (seeds). New style.
void vtkStreamTracer::SetIntegrator | ( | vtkInitialValueProblemSolver * | ) |
Set/get the integrator type to be used for streamline generation. The object passed is not actually used but is cloned with NewInstance in the process of integration (prototype pattern). The default is Runge-Kutta2. The integrator can also be changed using SetIntegratorType. The recognized solvers are: RUNGE_KUTTA2 = 0 RUNGE_KUTTA4 = 1 RUNGE_KUTTA45 = 2
|
virtual |
Set/get the integrator type to be used for streamline generation. The object passed is not actually used but is cloned with NewInstance in the process of integration (prototype pattern). The default is Runge-Kutta2. 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 for streamline generation. The object passed is not actually used but is cloned with NewInstance in the process of integration (prototype pattern). The default is Runge-Kutta2. 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 for streamline generation. The object passed is not actually used but is cloned with NewInstance in the process of integration (prototype pattern). The default is Runge-Kutta2. The integrator can also be changed using SetIntegratorType. The recognized solvers are: RUNGE_KUTTA2 = 0 RUNGE_KUTTA4 = 1 RUNGE_KUTTA45 = 2
|
inline |
Set/get the integrator type to be used for streamline generation. The object passed is not actually used but is cloned with NewInstance in the process of integration (prototype pattern). The default is Runge-Kutta2. The integrator can also be changed using SetIntegratorType. The recognized solvers are: RUNGE_KUTTA2 = 0 RUNGE_KUTTA4 = 1 RUNGE_KUTTA45 = 2
Definition at line 181 of file vtkStreamTracer.h.
|
inline |
Set/get the integrator type to be used for streamline generation. The object passed is not actually used but is cloned with NewInstance in the process of integration (prototype pattern). The default is Runge-Kutta2. The integrator can also be changed using SetIntegratorType. The recognized solvers are: RUNGE_KUTTA2 = 0 RUNGE_KUTTA4 = 1 RUNGE_KUTTA45 = 2
Definition at line 183 of file vtkStreamTracer.h.
|
inline |
Set/get the integrator type to be used for streamline generation. The object passed is not actually used but is cloned with NewInstance in the process of integration (prototype pattern). The default is Runge-Kutta2. The integrator can also be changed using SetIntegratorType. The recognized solvers are: RUNGE_KUTTA2 = 0 RUNGE_KUTTA4 = 1 RUNGE_KUTTA45 = 2
Definition at line 185 of file vtkStreamTracer.h.
void vtkStreamTracer::SetInterpolatorTypeToDataSetPointLocator | ( | ) |
Set the velocity field interpolator type to the one involving a dataset point locator.
void vtkStreamTracer::SetInterpolatorTypeToCellLocator | ( | ) |
Set the velocity field interpolator type to the one involving a cell locator.
|
virtual |
Specify the maximum length of a streamline expressed in LENGTH_UNIT.
|
virtual |
Specify the maximum length of a streamline expressed in LENGTH_UNIT.
void vtkStreamTracer::SetIntegrationStepUnit | ( | int | unit | ) |
Specify a uniform integration step unit for MinimumIntegrationStep, InitialIntegrationStep, and MaximumIntegrationStep. NOTE: The valid unit is now limited to only LENGTH_UNIT (1) and CELL_LENGTH_UNIT (2), EXCLUDING the previously-supported TIME_UNIT.
|
inline |
Specify a uniform integration step unit for MinimumIntegrationStep, InitialIntegrationStep, and MaximumIntegrationStep. NOTE: The valid unit is now limited to only LENGTH_UNIT (1) and CELL_LENGTH_UNIT (2), EXCLUDING the previously-supported TIME_UNIT.
Definition at line 209 of file vtkStreamTracer.h.
|
virtual |
Specify the Initial step size used for line integration, expressed in: LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 (either the starting size for an adaptive integrator, e.g., RK45, or the constant / fixed size for non-adaptive ones, i.e., RK2 and RK4)
|
virtual |
Specify the Initial step size used for line integration, expressed in: LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 (either the starting size for an adaptive integrator, e.g., RK45, or the constant / fixed size for non-adaptive ones, i.e., RK2 and RK4)
|
virtual |
Specify the Minimum step size used for line integration, expressed in: LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 (Only valid for an adaptive integrator, e.g., RK45)
|
virtual |
Specify the Minimum step size used for line integration, expressed in: LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 (Only valid for an adaptive integrator, e.g., RK45)
|
virtual |
Specify the Maximum step size used for line integration, expressed in: LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 (Only valid for an adaptive integrator, e.g., RK45)
|
virtual |
Specify the Maximum step size used for line integration, expressed in: LENGTH_UNIT = 1 CELL_LENGTH_UNIT = 2 (Only valid for an adaptive integrator, e.g., RK45)
|
virtual |
Specify the maximum error tolerated throughout streamline integration.
|
virtual |
Specify the maximum error tolerated throughout streamline integration.
|
virtual |
Specify the maximum number of steps for integrating a streamline.
|
virtual |
Specify the maximum number of steps for integrating a streamline.
|
virtual |
Specify the terminal speed value, below which integration is terminated.
|
virtual |
Specify the terminal speed value, below which integration is terminated.
|
virtual |
Specify whether the streamline is integrated in the upstream or downstream direction.
|
virtual |
Specify whether the streamline is integrated in the upstream or downstream direction.
|
inline |
Specify whether the streamline is integrated in the upstream or downstream direction.
Definition at line 277 of file vtkStreamTracer.h.
|
inline |
Specify whether the streamline is integrated in the upstream or downstream direction.
Definition at line 279 of file vtkStreamTracer.h.
|
inline |
Specify whether the streamline is integrated in the upstream or downstream direction.
Definition at line 281 of file vtkStreamTracer.h.
|
virtual |
Turn on/off vorticity computation at streamline points (necessary for generating proper stream-ribbons using the vtkRibbonFilter.
|
virtual |
Turn on/off vorticity computation at streamline points (necessary for generating proper stream-ribbons using the vtkRibbonFilter.
|
virtual |
This can be used to scale the rate with which the streamribbons twist. The default is 1.
|
virtual |
This can be used to scale the rate with which the streamribbons twist. The default is 1.
void vtkStreamTracer::SetInterpolatorPrototype | ( | vtkAbstractInterpolatedVelocityField * | ivf | ) |
The object used to interpolate the velocity field during integration is of the same class as this prototype.
void vtkStreamTracer::SetInterpolatorType | ( | int | interpType | ) |
Set the type of the velocity field interpolator to determine whether vtkInterpolatedVelocityField (INTERPOLATOR_WITH_DATASET_POINT_LOCATOR) or vtkCellLocatorInterpolatedVelocityField (INTERPOLATOR_WITH_CELL_LOCATOR) is employed for locating cells during streamline integration. The latter (adopting vtkAbstractCellLocator sub-classes such as vtkCellLocator and vtkModifiedBSPTree) is more robust then the former (through vtkDataSet / vtkPointSet::FindCell() coupled with vtkPointLocator).
|
protectedvirtual |
|
inlineprotected |
Definition at line 323 of file vtkStreamTracer.h.
|
protectedvirtual |
Reimplemented in vtkTemporalStreamTracer, vtkPTemporalStreamTracer, and vtkPStreamTracer.
|
protectedvirtual |
Reimplemented in vtkTemporalStreamTracer.
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
staticprotected |
|
staticprotected |
|
protected |
|
protected |
|
friend |
Definition at line 404 of file vtkStreamTracer.h.
|
protected |
Definition at line 351 of file vtkStreamTracer.h.
|
protected |
Definition at line 354 of file vtkStreamTracer.h.
|
staticprotected |
Definition at line 356 of file vtkStreamTracer.h.
|
protected |
Definition at line 357 of file vtkStreamTracer.h.
|
protected |
Definition at line 359 of file vtkStreamTracer.h.
|
protected |
Definition at line 368 of file vtkStreamTracer.h.
|
protected |
Definition at line 369 of file vtkStreamTracer.h.
|
protected |
Definition at line 370 of file vtkStreamTracer.h.
|
protected |
Definition at line 371 of file vtkStreamTracer.h.
|
protected |
Definition at line 387 of file vtkStreamTracer.h.
|
protected |
Definition at line 388 of file vtkStreamTracer.h.
|
protected |
Definition at line 391 of file vtkStreamTracer.h.
|
protected |
Definition at line 393 of file vtkStreamTracer.h.
|
protected |
Definition at line 394 of file vtkStreamTracer.h.
|
protected |
Definition at line 396 of file vtkStreamTracer.h.
|
protected |
Definition at line 397 of file vtkStreamTracer.h.
|
protected |
Definition at line 399 of file vtkStreamTracer.h.
|
protected |
Definition at line 401 of file vtkStreamTracer.h.
|
protected |
Definition at line 402 of file vtkStreamTracer.h.