VTK  9.1.0
vtkStreamTracer.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkStreamTracer.h
5 
6  Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7  All rights reserved.
8  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9 
10  This software is distributed WITHOUT ANY WARRANTY; without even
11  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12  PURPOSE. See the above copyright notice for more information.
13 
14 =========================================================================*/
195 #ifndef vtkStreamTracer_h
196 #define vtkStreamTracer_h
197 
198 #include "vtkFiltersFlowPathsModule.h" // For export macro
199 #include "vtkPolyDataAlgorithm.h"
200 
201 #include "vtkInitialValueProblemSolver.h" // Needed for constants
202 
204 class vtkCompositeDataSet;
205 class vtkDataArray;
207 class vtkDoubleArray;
208 class vtkExecutive;
209 class vtkGenericCell;
210 class vtkIdList;
211 class vtkIntArray;
212 class vtkPoints;
213 
214 #include <vector> // for std::vector
215 
216 class VTKFILTERSFLOWPATHS_EXPORT vtkStreamTracer : public vtkPolyDataAlgorithm
217 {
218 public:
220  void PrintSelf(ostream& os, vtkIndent indent) override;
221 
229  static vtkStreamTracer* New();
230 
232 
237  vtkSetVector3Macro(StartPosition, double);
238  vtkGetVector3Macro(StartPosition, double);
240 
242 
251 
257 
258  // The previously-supported TIME_UNIT is excluded in this current
259  // enumeration definition because the underlying step size is ALWAYS in
260  // arc length unit (LENGTH_UNIT) while the 'real' time interval (virtual
261  // for steady flows) that a particle actually takes to trave in a single
262  // step is obtained by dividing the arc length by the LOCAL speed. The
263  // overall elapsed time (i.e., the life span) of the particle is the sum
264  // of those individual step-wise time intervals. The arc-length-to-time
265  // conversion only occurs for vorticity computation and for generating a
266  // point data array named 'IntegrationTime'.
267  enum Units
268  {
269  LENGTH_UNIT = 1,
270  CELL_LENGTH_UNIT = 2
271  };
272 
273  enum Solvers
274  {
279  UNKNOWN
280  };
281 
283  {
287  OUT_OF_LENGTH = 4,
288  OUT_OF_STEPS = 5,
289  STAGNATION = 6,
290  FIXED_REASONS_FOR_TERMINATION_COUNT
291  };
292 
294 
305  vtkGetObjectMacro(Integrator, vtkInitialValueProblemSolver);
308  void SetIntegratorTypeToRungeKutta2() { this->SetIntegratorType(RUNGE_KUTTA2); }
309  void SetIntegratorTypeToRungeKutta4() { this->SetIntegratorType(RUNGE_KUTTA4); }
310  void SetIntegratorTypeToRungeKutta45() { this->SetIntegratorType(RUNGE_KUTTA45); }
312 
318 
324 
326 
329  vtkSetMacro(MaximumPropagation, double);
330  vtkGetMacro(MaximumPropagation, double);
332 
339  void SetIntegrationStepUnit(int unit);
340  int GetIntegrationStepUnit() { return this->IntegrationStepUnit; }
341 
343 
350  vtkSetMacro(InitialIntegrationStep, double);
351  vtkGetMacro(InitialIntegrationStep, double);
353 
355 
361  vtkSetMacro(MinimumIntegrationStep, double);
362  vtkGetMacro(MinimumIntegrationStep, double);
364 
366 
372  vtkSetMacro(MaximumIntegrationStep, double);
373  vtkGetMacro(MaximumIntegrationStep, double);
375 
377 
380  vtkSetMacro(MaximumError, double);
381  vtkGetMacro(MaximumError, double);
383 
385 
388  vtkSetMacro(MaximumNumberOfSteps, vtkIdType);
389  vtkGetMacro(MaximumNumberOfSteps, vtkIdType);
391 
393 
396  vtkSetMacro(TerminalSpeed, double);
397  vtkGetMacro(TerminalSpeed, double);
399 
401 
404  vtkGetMacro(SurfaceStreamlines, bool);
405  vtkSetMacro(SurfaceStreamlines, bool);
406  vtkBooleanMacro(SurfaceStreamlines, bool);
408 
409  enum
410  {
413  BOTH
414  };
415 
416  enum
417  {
419  INTERPOLATOR_WITH_CELL_LOCATOR
420  };
421 
423 
427  vtkSetClampMacro(IntegrationDirection, int, FORWARD, BOTH);
428  vtkGetMacro(IntegrationDirection, int);
429  void SetIntegrationDirectionToForward() { this->SetIntegrationDirection(FORWARD); }
430  void SetIntegrationDirectionToBackward() { this->SetIntegrationDirection(BACKWARD); }
431  void SetIntegrationDirectionToBoth() { this->SetIntegrationDirection(BOTH); }
433 
435 
440  vtkSetMacro(ComputeVorticity, bool);
441  vtkGetMacro(ComputeVorticity, bool);
443 
445 
449  vtkSetMacro(RotationScale, double);
450  vtkGetMacro(RotationScale, double);
452 
454 
462  vtkSetMacro(UseLocalSeedSource, bool);
463  vtkGetMacro(UseLocalSeedSource, bool);
464  vtkBooleanMacro(UseLocalSeedSource, bool);
466 
472 
482  void SetInterpolatorType(int interpType);
483 
493  typedef bool (*CustomTerminationCallbackType)(
494  void* clientdata, vtkPoints* points, vtkDataArray* velocity, int integrationDirection);
504  CustomTerminationCallbackType callback, void* clientdata, int reasonForTermination);
505 
506 protected:
508  ~vtkStreamTracer() override;
509 
510  // Create a default executive.
512 
513  // hide the superclass' AddInput() from the user and the compiler
515  {
516  vtkErrorMacro(<< "AddInput() must be called with a vtkDataSet not a vtkDataObject.");
517  }
518 
521 
523  vtkGenericCell* cell, double pcoords[3], vtkDoubleArray* cellVectors, double vorticity[3]);
524  void Integrate(vtkPointData* inputData, vtkPolyData* output, vtkDataArray* seedSource,
525  vtkIdList* seedIds, vtkIntArray* integrationDirections, double lastPoint[3],
526  vtkAbstractInterpolatedVelocityField* func, int maxCellSize, int vecType,
527  const char* vecFieldName, double& propagation, vtkIdType& numSteps, double& integrationTime);
528  double SimpleIntegrate(double seed[3], double lastPoint[3], double stepSize,
530  int CheckInputs(vtkAbstractInterpolatedVelocityField*& func, int* maxCellSize);
531  void GenerateNormals(vtkPolyData* output, double* firstNormal, const char* vecName);
532 
534 
535  // starting from global x-y-z position
536  double StartPosition[3];
537 
538  static const double EPSILON;
540 
542 
544  {
545  double Interval;
546  int Unit;
547  };
548 
553 
555  double& step, double& minStep, double& maxStep, int direction, double cellLength);
556  static double ConvertToLength(double interval, int unit, double cellLength);
557  static double ConvertToLength(IntervalInformation& interval, double cellLength);
558 
559  int SetupOutput(vtkInformation* inInfo, vtkInformation* outInfo);
560  void InitializeSeeds(vtkDataArray*& seeds, vtkIdList*& seedIds,
561  vtkIntArray*& integrationDirections, vtkDataSet* source);
562 
565 
566  // Prototype showing the integrator type to be set by the user.
568 
569  double MaximumError;
571 
574 
575  // Compute streamlines only on surface.
577 
578  // Only relevant for the parallel version of this filter (see vtkPStreamTracer)
579  bool UseLocalSeedSource = true;
580 
582 
584  bool
585  HasMatchingPointAttributes; // does the point data in the multiblocks have the same attributes?
586  std::vector<CustomTerminationCallbackType> CustomTerminationCallback;
587  std::vector<void*> CustomTerminationClientData;
588  std::vector<int> CustomReasonForTermination;
589 
590  friend class PStreamTracerUtils;
591 
592 private:
593  vtkStreamTracer(const vtkStreamTracer&) = delete;
594  void operator=(const vtkStreamTracer&) = delete;
595 };
596 
597 #endif
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
Definition: vtkDataArray.h:159
general representation of visualization data
represent and manipulate attribute data in a dataset
abstract class to specify dataset behavior
Definition: vtkDataSet.h:166
dynamic, self-adjusting array of double
Superclass for all pipeline executives in VTK.
Definition: vtkExecutive.h:76
provides thread-safe access to cells
list of point or cell ids
Definition: vtkIdList.h:140
a simple class to control print indentation
Definition: vtkIndent.h:113
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Integrate a set of ordinary differential equations (initial value problem) in time.
dynamic, self-adjusting array of int
Definition: vtkIntArray.h:149
represent and manipulate point attribute data
Definition: vtkPointData.h:142
represent and manipulate 3D points
Definition: vtkPoints.h:143
Superclass for algorithms that produce only polydata as output.
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:195
Streamline generator.
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
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, double &integrationTime)
@ INTERPOLATOR_WITH_DATASET_POINT_LOCATOR
void SetSourceData(vtkDataSet *source)
Specify the source object used to generate starting points (seeds).
double InitialIntegrationStep
vtkAbstractInterpolatedVelocityField * InterpolatorPrototype
void SetInterpolatorTypeToCellLocator()
Set the velocity field interpolator type to the one involving a cell locator.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void CalculateVorticity(vtkGenericCell *cell, double pcoords[3], vtkDoubleArray *cellVectors, double vorticity[3])
double MinimumIntegrationStep
void SetIntegratorTypeToRungeKutta4()
Set/get the integrator type to be used for streamline generation.
vtkDataSet * GetSource()
Specify the source object used to generate starting points (seeds).
static double ConvertToLength(double interval, int unit, double cellLength)
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
static double ConvertToLength(IntervalInformation &interval, double cellLength)
int CheckInputs(vtkAbstractInterpolatedVelocityField *&func, int *maxCellSize)
void ConvertIntervals(double &step, double &minStep, double &maxStep, int direction, double cellLength)
void GenerateNormals(vtkPolyData *output, double *firstNormal, const char *vecName)
static const double EPSILON
vtkIdType MaximumNumberOfSteps
void SetIntegrationDirectionToForward()
Specify whether the streamline is integrated in the upstream or downstream direction.
std::vector< CustomTerminationCallbackType > CustomTerminationCallback
bool HasMatchingPointAttributes
vtkCompositeDataSet * InputData
vtkExecutive * CreateDefaultExecutive() override
Create a default executive.
void SetInterpolatorType(int interpType)
Set the type of the velocity field interpolator to determine whether vtkInterpolatedVelocityField (IN...
double MaximumIntegrationStep
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
bool GenerateNormalsInIntegrate
void SetIntegrationDirectionToBackward()
Specify whether the streamline is integrated in the upstream or downstream direction.
void SetInterpolatorTypeToDataSetPointLocator()
Set the velocity field interpolator type to the one involving a dataset point locator.
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.
static vtkStreamTracer * New()
Construct object to start from position (0,0,0), with forward integration, terminal speed 1....
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()
@ points
Definition: vtkX3D.h:452
@ direction
Definition: vtkX3D.h:266
@ type
Definition: vtkX3D.h:522
boost::graph_traits< vtkGraph * >::vertex_descriptor source(boost::graph_traits< vtkGraph * >::edge_descriptor e, vtkGraph *)
int vtkIdType
Definition: vtkType.h:332