VTK
dox/Filters/FlowPaths/vtkStreamTracer.h
Go to the documentation of this file.
00001 /*=========================================================================
00002 
00003   Program:   Visualization Toolkit
00004   Module:    vtkStreamTracer.h
00005 
00006   Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
00007   All rights reserved.
00008   See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
00009 
00010      This software is distributed WITHOUT ANY WARRANTY; without even
00011      the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
00012      PURPOSE.  See the above copyright notice for more information.
00013 
00014 =========================================================================*/
00085 #ifndef __vtkStreamTracer_h
00086 #define __vtkStreamTracer_h
00087 
00088 #include "vtkFiltersFlowPathsModule.h" // For export macro
00089 #include "vtkPolyDataAlgorithm.h"
00090 
00091 #include "vtkInitialValueProblemSolver.h" // Needed for constants
00092 
00093 class vtkCompositeDataSet;
00094 class vtkDataArray;
00095 class vtkDoubleArray;
00096 class vtkExecutive;
00097 class vtkGenericCell;
00098 class vtkIdList;
00099 class vtkIntArray;
00100 class vtkAbstractInterpolatedVelocityField;
00101 
00102 class VTKFILTERSFLOWPATHS_EXPORT vtkStreamTracer : public vtkPolyDataAlgorithm
00103 {
00104 public:
00105   vtkTypeMacro(vtkStreamTracer,vtkPolyDataAlgorithm);
00106   void PrintSelf(ostream& os, vtkIndent indent);
00107 
00113   static vtkStreamTracer *New();
00114 
00116 
00119   vtkSetVector3Macro(StartPosition, double);
00120   vtkGetVector3Macro(StartPosition, double);
00122 
00124 
00126   void SetSourceData(vtkDataSet *source);
00127   vtkDataSet *GetSource();
00129 
00132   void SetSourceConnection(vtkAlgorithmOutput* algOutput);
00133 
00134 //BTX
00135   // The previously-supported TIME_UNIT is excluded in this current
00136   // enumeration definition because the underlying step size is ALWAYS in
00137   // arc length unit (LENGTH_UNIT) while the 'real' time interval (virtual
00138   // for steady flows) that a particle actually takes to trave in a single
00139   // step is obtained by dividing the arc length by the LOCAL speed. The
00140   // overall elapsed time (i.e., the life span) of the particle is the sum
00141   // of those individual step-wise time intervals. The arc-length-to-time
00142   // conversion only occurs for vorticity computation and for generating a
00143   // point data array named 'IntegrationTime'.
00144   enum Units
00145   {
00146     LENGTH_UNIT = 1,
00147     CELL_LENGTH_UNIT = 2
00148   };
00149 
00150   enum Solvers
00151   {
00152     RUNGE_KUTTA2,
00153     RUNGE_KUTTA4,
00154     RUNGE_KUTTA45,
00155     NONE,
00156     UNKNOWN
00157   };
00158 
00159   enum ReasonForTermination
00160   {
00161     OUT_OF_DOMAIN = vtkInitialValueProblemSolver::OUT_OF_DOMAIN,
00162     NOT_INITIALIZED = vtkInitialValueProblemSolver::NOT_INITIALIZED ,
00163     UNEXPECTED_VALUE = vtkInitialValueProblemSolver::UNEXPECTED_VALUE,
00164     OUT_OF_LENGTH = 4,
00165     OUT_OF_STEPS = 5,
00166     STAGNATION = 6
00167   };
00168 //ETX
00169 
00171 
00177   void SetIntegrator(vtkInitialValueProblemSolver *);
00178   vtkGetObjectMacro ( Integrator, vtkInitialValueProblemSolver );
00179   void SetIntegratorType(int type);
00180   int GetIntegratorType();
00181   void SetIntegratorTypeToRungeKutta2()
00182     {this->SetIntegratorType(RUNGE_KUTTA2);};
00183   void SetIntegratorTypeToRungeKutta4()
00184     {this->SetIntegratorType(RUNGE_KUTTA4);};
00185   void SetIntegratorTypeToRungeKutta45()
00186     {this->SetIntegratorType(RUNGE_KUTTA45);};
00188 
00191   void SetInterpolatorTypeToDataSetPointLocator();
00192 
00195   void SetInterpolatorTypeToCellLocator();
00196 
00198 
00199   vtkSetMacro(MaximumPropagation, double);
00200   vtkGetMacro(MaximumPropagation, double);
00202 
00204 
00208   void SetIntegrationStepUnit( int unit );
00209   int  GetIntegrationStepUnit() { return this->IntegrationStepUnit; }
00211 
00213 
00217   vtkSetMacro(InitialIntegrationStep, double);
00218   vtkGetMacro(InitialIntegrationStep, double);
00220 
00222 
00225   vtkSetMacro(MinimumIntegrationStep, double);
00226   vtkGetMacro(MinimumIntegrationStep, double);
00228 
00230 
00233   vtkSetMacro(MaximumIntegrationStep, double);
00234   vtkGetMacro(MaximumIntegrationStep, double);
00236 
00238 
00240   vtkSetMacro(MaximumError, double);
00241   vtkGetMacro(MaximumError, double);
00243 
00245 
00246   vtkSetMacro(MaximumNumberOfSteps, vtkIdType);
00247   vtkGetMacro(MaximumNumberOfSteps, vtkIdType);
00249 
00251 
00253   vtkSetMacro(TerminalSpeed, double);
00254   vtkGetMacro(TerminalSpeed, double);
00256 
00257 //BTX
00258   enum
00259   {
00260     FORWARD,
00261     BACKWARD,
00262     BOTH
00263   };
00264 
00265   enum
00266   {
00267     INTERPOLATOR_WITH_DATASET_POINT_LOCATOR,
00268     INTERPOLATOR_WITH_CELL_LOCATOR
00269   };
00270 //ETX
00271 
00273 
00275   vtkSetClampMacro(IntegrationDirection, int, FORWARD, BOTH);
00276   vtkGetMacro(IntegrationDirection, int);
00277   void SetIntegrationDirectionToForward()
00278     {this->SetIntegrationDirection(FORWARD);};
00279   void SetIntegrationDirectionToBackward()
00280     {this->SetIntegrationDirection(BACKWARD);};
00281   void SetIntegrationDirectionToBoth()
00282     {this->SetIntegrationDirection(BOTH);};
00284 
00286 
00288   vtkSetMacro(ComputeVorticity, bool);
00289   vtkGetMacro(ComputeVorticity, bool);
00291 
00293 
00295   vtkSetMacro(RotationScale, double);
00296   vtkGetMacro(RotationScale, double);
00298 
00301   void SetInterpolatorPrototype( vtkAbstractInterpolatedVelocityField * ivf );
00302 
00312   void SetInterpolatorType( int interpType );
00313 
00314 protected:
00315 
00316   vtkStreamTracer();
00317   ~vtkStreamTracer();
00318 
00319   // Create a default executive.
00320   virtual vtkExecutive* CreateDefaultExecutive();
00321 
00322   // hide the superclass' AddInput() from the user and the compiler
00323   void AddInput(vtkDataObject *)
00324     { vtkErrorMacro( << "AddInput() must be called with a vtkDataSet not a vtkDataObject."); };
00325 
00326   virtual int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *);
00327   virtual int FillInputPortInformation(int, vtkInformation *);
00328 
00329   void CalculateVorticity( vtkGenericCell* cell, double pcoords[3],
00330                            vtkDoubleArray* cellVectors, double vorticity[3] );
00331   void Integrate(vtkPointData *inputData,
00332                  vtkPolyData* output,
00333                  vtkDataArray* seedSource,
00334                  vtkIdList* seedIds,
00335                  vtkIntArray* integrationDirections,
00336                  double lastPoint[3],
00337                  vtkAbstractInterpolatedVelocityField* func,
00338                  int maxCellSize,
00339                  int vecType,
00340                  const char *vecFieldName,
00341                  double& propagation,
00342                  vtkIdType& numSteps);
00343   void SimpleIntegrate(double seed[3],
00344                        double lastPoint[3],
00345                        double stepSize,
00346                        vtkAbstractInterpolatedVelocityField* func);
00347   int CheckInputs(vtkAbstractInterpolatedVelocityField*& func,
00348                   int* maxCellSize);
00349   void GenerateNormals(vtkPolyData* output, double* firstNormal, const char *vecName);
00350 
00351   bool GenerateNormalsInIntegrate;
00352 
00353   // starting from global x-y-z position
00354   double StartPosition[3];
00355 
00356   static const double EPSILON;
00357   double TerminalSpeed;
00358 
00359   double LastUsedStepSize;
00360 
00361 //BTX
00362   struct IntervalInformation
00363   {
00364     double Interval;
00365     int Unit;
00366   };
00367 
00368   double MaximumPropagation;
00369   double MinimumIntegrationStep;
00370   double MaximumIntegrationStep;
00371   double InitialIntegrationStep;
00372 
00373   void ConvertIntervals( double& step, double& minStep, double& maxStep,
00374                         int direction, double cellLength );
00375   static double ConvertToLength( double interval, int unit, double cellLength );
00376   static double ConvertToLength( IntervalInformation& interval, double cellLength );
00377 
00378 //ETX
00379 
00380   int SetupOutput(vtkInformation* inInfo,
00381                   vtkInformation* outInfo);
00382   void InitializeSeeds(vtkDataArray*& seeds,
00383                        vtkIdList*& seedIds,
00384                        vtkIntArray*& integrationDirections,
00385                        vtkDataSet *source);
00386 
00387   int IntegrationStepUnit;
00388   int IntegrationDirection;
00389 
00390   // Prototype showing the integrator type to be set by the user.
00391   vtkInitialValueProblemSolver* Integrator;
00392 
00393   double MaximumError;
00394   vtkIdType MaximumNumberOfSteps;
00395 
00396   bool ComputeVorticity;
00397   double RotationScale;
00398 
00399   vtkAbstractInterpolatedVelocityField * InterpolatorPrototype;
00400 
00401   vtkCompositeDataSet* InputData;
00402   bool HasMatchingPointAttributes; //does the point data in the multiblocks have the same attributes?
00403 
00404   friend class PStreamTracerUtils;
00405 
00406 private:
00407   vtkStreamTracer(const vtkStreamTracer&);  // Not implemented.
00408   void operator=(const vtkStreamTracer&);  // Not implemented.
00409 };
00410 
00411 
00412 #endif