VTK
dox/Filters/FlowPaths/vtkParticleTracerBase.h
Go to the documentation of this file.
00001 /*=========================================================================
00002 
00003   Program:   Visualization Toolkit
00004   Module:    vtkParticleTracerBase.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 =========================================================================*/
00027 #ifndef __vtkParticleTracerBase_h
00028 #define __vtkParticleTracerBase_h
00029 
00030 #include "vtkFiltersFlowPathsModule.h" // For export macro
00031 #include "vtkSmartPointer.h" // For protected ivars.
00032 #include "vtkPolyDataAlgorithm.h"
00033 //BTX
00034 #include <vector> // STL Header
00035 #include <list>   // STL Header
00036 //ETX
00037 
00038 class vtkAbstractInterpolatedVelocityField;
00039 class vtkAbstractParticleWriter;
00040 class vtkCellArray;
00041 class vtkCharArray;
00042 class vtkCompositeDataSet;
00043 class vtkDataArray;
00044 class vtkDataSet;
00045 class vtkDoubleArray;
00046 class vtkFloatArray;
00047 class vtkGenericCell;
00048 class vtkInitialValueProblemSolver;
00049 class vtkIntArray;
00050 class vtkMultiBlockDataSet;
00051 class vtkMultiProcessController;
00052 class vtkPointData;
00053 class vtkPoints;
00054 class vtkPolyData;
00055 class vtkTemporalInterpolatedVelocityField;
00056 
00057 //BTX
00058 namespace vtkParticleTracerBaseNamespace
00059 {
00060   typedef struct { double x[4]; } Position;
00061   typedef struct {
00062     // These are used during iteration
00063     Position      CurrentPosition;
00064     int           CachedDataSetId[2];
00065     vtkIdType     CachedCellId[2];
00066     int           LocationState;
00067     // These are computed scalars we might display
00068     int           SourceID;
00069     int           TimeStepAge;
00070     int           InjectedPointId;
00071     int           InjectedStepId;
00072     int           UniqueParticleId;
00073     // These are useful to track for debugging etc
00074     int           ErrorCode;
00075     float         age;
00076     // these are needed across time steps to compute vorticity
00077     float         rotation;
00078     float         angularVel;
00079     float         time;
00080     float         speed;
00081 
00082     vtkIdType     PointId; //once the partice is added, PointId is valid
00083   } ParticleInformation;
00084 
00085   typedef std::vector<ParticleInformation>  ParticleVector;
00086   typedef ParticleVector::iterator             ParticleIterator;
00087   typedef std::list<ParticleInformation>    ParticleDataList;
00088   typedef ParticleDataList::iterator           ParticleListIterator;
00089 };
00090 //ETX
00091 
00092 class VTKFILTERSFLOWPATHS_EXPORT vtkParticleTracerBase : public vtkPolyDataAlgorithm
00093 {
00094 public:
00095   enum Solvers
00096   {
00097     RUNGE_KUTTA2,
00098     RUNGE_KUTTA4,
00099     RUNGE_KUTTA45,
00100     NONE,
00101     UNKNOWN
00102   };
00103 
00104   vtkTypeMacro(vtkParticleTracerBase,vtkPolyDataAlgorithm)
00105   void PrintSelf(ostream& os, vtkIndent indent);
00106   void PrintParticleHistories();
00107 
00109 
00111   vtkGetMacro(ComputeVorticity, bool);
00112   void SetComputeVorticity(bool);
00114 
00116 
00118   vtkGetMacro(TerminalSpeed, double);
00119   void SetTerminalSpeed(double);
00121 
00123 
00125   vtkGetMacro(RotationScale, double);
00126   void SetRotationScale(double);
00128 
00130 
00132   vtkSetMacro(IgnorePipelineTime, int);
00133   vtkGetMacro(IgnorePipelineTime, int);
00134   vtkBooleanMacro(IgnorePipelineTime, int);
00136 
00138 
00145   vtkGetMacro(ForceReinjectionEveryNSteps,int);
00146   void SetForceReinjectionEveryNSteps(int);
00148 
00150 
00154   void SetTerminationTime(double t);
00155   vtkGetMacro(TerminationTime,double);
00157 
00158   void SetIntegrator(vtkInitialValueProblemSolver *);
00159   vtkGetObjectMacro ( Integrator, vtkInitialValueProblemSolver );
00160 
00161   void SetIntegratorType(int type);
00162   int GetIntegratorType();
00163 
00165 
00169   vtkGetMacro(StartTime, double);
00170   void SetStartTime(double t);
00172 
00174 
00180   vtkSetMacro(StaticSeeds,int);
00181   vtkGetMacro(StaticSeeds,int);
00183 
00185 
00191   vtkSetMacro(StaticMesh,int);
00192   vtkGetMacro(StaticMesh,int);
00194 
00196 
00200   virtual void SetParticleWriter(vtkAbstractParticleWriter *pw);
00201   vtkGetObjectMacro(ParticleWriter, vtkAbstractParticleWriter);
00203 
00205 
00207   vtkSetStringMacro(ParticleFileName);
00208   vtkGetStringMacro(ParticleFileName);
00210 
00212 
00214   vtkSetMacro(EnableParticleWriting,int);
00215   vtkGetMacro(EnableParticleWriting,int);
00216   vtkBooleanMacro(EnableParticleWriting,int);
00218 
00220 
00222   vtkSetMacro(DisableResetCache,int);
00223   vtkGetMacro(DisableResetCache,int);
00224   vtkBooleanMacro(DisableResetCache,int);
00226 
00228 
00229   void AddSourceConnection(vtkAlgorithmOutput* input);
00230   void RemoveAllSources();
00232 
00233  protected:
00234   vtkSmartPointer<vtkPolyData> Output; //managed by child classes
00235   vtkSmartPointer<vtkPointData> ProtoPD;
00236   vtkIdType UniqueIdCounter;// global Id counter used to give particles a stamp
00237   vtkParticleTracerBaseNamespace::ParticleDataList  ParticleHistories;
00238   vtkSmartPointer<vtkPointData>     ParticlePointData; //the current particle point data consistent
00239                                                        //with particle history
00240   int           ReinjectionCounter;
00241 
00242   //Everything related to time
00243   int IgnorePipelineTime; //whether to use the pipeline time for termination
00244   int DisableResetCache; //whether to enable ResetCache() method
00245 
00246 
00247   vtkParticleTracerBase();
00248   virtual ~vtkParticleTracerBase();
00249 
00250   //
00251   // Make sure the pipeline knows what type we expect as input
00252   //
00253   virtual int FillInputPortInformation(int port, vtkInformation* info);
00254 
00255   //
00256   // The usual suspects
00257   //
00258   virtual int ProcessRequest(vtkInformation* request,
00259                              vtkInformationVector** inputVector,
00260                              vtkInformationVector* outputVector);
00261 
00262   //
00263   // Store any information we need in the output and fetch what we can
00264   // from the input
00265   //
00266   virtual int RequestInformation(vtkInformation* request,
00267                                  vtkInformationVector** inputVector,
00268                                  vtkInformationVector* outputVector);
00269 
00270   //
00271   // Compute input time steps given the output step
00272   //
00273   virtual int RequestUpdateExtent(vtkInformation* request,
00274                                   vtkInformationVector** inputVector,
00275                                   vtkInformationVector* outputVector);
00276 
00277   //
00278   // what the pipeline calls for each time step
00279   //
00280   virtual int RequestData(vtkInformation* request,
00281                           vtkInformationVector** inputVector,
00282                           vtkInformationVector* outputVector);
00283 
00284   //
00285   // these routines are internally called to actually generate the output
00286   //
00287   virtual int ProcessInput(vtkInformationVector** inputVector);
00288 
00289   // This is the main part of the algorithm:
00290   //  * move all the particles one step
00291   //  * Reinject particles (by adding them to this->ParticleHistories)
00292   //    either at the beginning or at the end of each step (modulo this->ForceReinjectionEveryNSteps)
00293   //  * Output a polydata representing the moved particles
00294   // Note that if the starting and the ending time coincide, the polydata is still valid.
00295   virtual vtkPolyData* Execute(vtkInformationVector** inputVector);
00296 
00297   // the RequestData will call these methods in turn
00298   virtual void Initialize(){} //the first iteration
00299   virtual int OutputParticles(vtkPolyData* poly)=0; //every iteration
00300   virtual void Finalize(){} //the last iteration
00301 
00302   //
00303   // Initialization of input (vector-field) geometry
00304   //
00305   int InitializeInterpolator();
00306   int UpdateDataCache(vtkDataObject *td);
00307 
00309 
00311   void TestParticles(
00312     vtkParticleTracerBaseNamespace::ParticleVector &candidates,
00313     vtkParticleTracerBaseNamespace::ParticleVector &passed,
00314     int &count);
00316 
00317   void TestParticles(
00318     vtkParticleTracerBaseNamespace::ParticleVector &candidates, std::vector<int> &passed);
00319 
00321 
00325   virtual void AssignSeedsToProcessors(double time,
00326     vtkDataSet *source, int sourceID, int ptId,
00327     vtkParticleTracerBaseNamespace::ParticleVector &localSeedPoints,
00328     int &localAssignedCount);
00330 
00332 
00334   virtual void AssignUniqueIds(
00335     vtkParticleTracerBaseNamespace::ParticleVector &localSeedPoints);
00337 
00339 
00341   void UpdateParticleList(
00342     vtkParticleTracerBaseNamespace::ParticleVector &candidates);
00344 
00347   virtual void UpdateParticleListFromOtherProcesses(){}
00348 
00350 
00351   void IntegrateParticle(
00352     vtkParticleTracerBaseNamespace::ParticleListIterator &it,
00353     double currenttime, double terminationtime,
00354     vtkInitialValueProblemSolver* integrator);
00356 
00357   // if the particle is added to send list, then returns value is 1,
00358   // if it is kept on this process after a retry return value is 0
00359   virtual bool SendParticleToAnotherProcess(
00360     vtkParticleTracerBaseNamespace::ParticleInformation &,
00361     vtkParticleTracerBaseNamespace::ParticleInformation &, vtkPointData*)
00362   {
00363     return true;
00364   }
00365 
00367 
00371   bool ComputeDomainExitLocation(
00372     double pos[4], double p2[4], double intersection[4],
00373     vtkGenericCell *cell);
00375 
00376   //
00377   // Scalar arrays that are generated as each particle is updated
00378   //
00379   void CreateProtoPD(vtkDataObject* input);
00380 
00381   vtkFloatArray*    GetParticleAge(vtkPointData*);
00382   vtkIntArray*      GetParticleIds(vtkPointData*);
00383   vtkCharArray*     GetParticleSourceIds(vtkPointData*);
00384   vtkIntArray*      GetInjectedPointIds(vtkPointData*);
00385   vtkIntArray*      GetInjectedStepIds(vtkPointData*);
00386   vtkIntArray*      GetErrorCodeArr(vtkPointData*);
00387   vtkFloatArray*    GetParticleVorticity(vtkPointData*);
00388   vtkFloatArray*    GetParticleRotation(vtkPointData*);
00389   vtkFloatArray*    GetParticleAngularVel(vtkPointData*);
00390 
00391 
00392   // utility function we use to test if a point is inside any of our local datasets
00393   bool InsideBounds(double point[]);
00394 
00395 
00396 
00397   void CalculateVorticity( vtkGenericCell* cell, double pcoords[3],
00398                            vtkDoubleArray* cellVectors, double vorticity[3] );
00399 
00400   //------------------------------------------------------
00401 
00402 
00403   double GetCacheDataTime(int i);
00404   double GetCacheDataTime();
00405 
00406   virtual void ResetCache();
00407   void AddParticle(vtkParticleTracerBaseNamespace::ParticleInformation &info, double* velocity);
00408 
00410 
00413   virtual bool IsPointDataValid(vtkDataObject* input);
00414   bool IsPointDataValid(vtkCompositeDataSet* input, std::vector<std::string>& arrayNames);
00415   void GetPointDataArrayNames(vtkDataSet* input, std::vector<std::string>& names);
00417 
00418 private:
00420   void SetInterpolatorPrototype(vtkAbstractInterpolatedVelocityField*) {}
00421 
00423 
00430   bool RetryWithPush(
00431     vtkParticleTracerBaseNamespace::ParticleInformation &info, double* point1,double delT, int subSteps);
00433 
00434   //Parameters of tracing
00435   vtkInitialValueProblemSolver* Integrator;
00436   double IntegrationStep;
00437   double MaximumError;
00438   bool ComputeVorticity;
00439   double RotationScale;
00440   double TerminalSpeed;
00441 
00442   // Important for Caching of Cells/Ids/Weights etc
00443   int           AllFixedGeometry;
00444   int           StaticMesh;
00445   int           StaticSeeds;
00446 
00447   std::vector<double>  InputTimeValues;
00448   double StartTime;
00449   double TerminationTime;
00450   double CurrentTime;
00451 
00452   int  StartTimeStep; //InputTimeValues[StartTimeStep] <= StartTime <= InputTimeValues[StartTimeStep+1]
00453   int  CurrentTimeStep;
00454   int  TerminationTimeStep; //computed from start time
00455   bool FirstIteration;
00456 
00457   //Innjection parameters
00458   int           ForceReinjectionEveryNSteps;
00459   vtkTimeStamp  ParticleInjectionTime;
00460   bool          HasCache;
00461 
00462   // Particle writing to disk
00463   vtkAbstractParticleWriter *ParticleWriter;
00464   char                      *ParticleFileName;
00465   int                        EnableParticleWriting;
00466 
00467 
00468   // The main lists which are held during operation- between time step updates
00469   vtkParticleTracerBaseNamespace::ParticleVector    LocalSeeds;
00470 
00471 
00472   // The velocity interpolator
00473   vtkSmartPointer<vtkTemporalInterpolatedVelocityField>  Interpolator;
00474   vtkAbstractInterpolatedVelocityField * InterpolatorPrototype;
00475 
00476   // Data for time step CurrentTimeStep-1 and CurrentTimeStep
00477   vtkSmartPointer<vtkMultiBlockDataSet> CachedData[2];
00478 
00479   // Cache bounds info for each dataset we will use repeatedly
00480   typedef struct {
00481     double b[6];
00482   } bounds;
00483   std::vector<bounds> CachedBounds[2];
00484 
00485   // temporary variables used by Exeucte(), for convenience only
00486 
00487   vtkSmartPointer<vtkPoints> OutputCoordinates;
00488   vtkSmartPointer<vtkFloatArray>    ParticleAge;
00489   vtkSmartPointer<vtkIntArray>      ParticleIds;
00490   vtkSmartPointer<vtkCharArray>     ParticleSourceIds;
00491   vtkSmartPointer<vtkIntArray>      InjectedPointIds;
00492   vtkSmartPointer<vtkIntArray>      InjectedStepIds;
00493   vtkSmartPointer<vtkIntArray>      ErrorCode;
00494   vtkSmartPointer<vtkFloatArray>    ParticleVorticity;
00495   vtkSmartPointer<vtkFloatArray>    ParticleRotation;
00496   vtkSmartPointer<vtkFloatArray>    ParticleAngularVel;
00497   vtkSmartPointer<vtkDoubleArray>   CellVectors;
00498   vtkSmartPointer<vtkPointData>     OutputPointData;
00499   vtkSmartPointer<vtkDataSet>       DataReferenceT[2];
00500   vtkSmartPointer<vtkCellArray>     ParticleCells;
00501 
00502   vtkParticleTracerBase(const vtkParticleTracerBase&);  // Not implemented.
00503   void operator=(const vtkParticleTracerBase&);  // Not implemented.
00504   vtkTimeStamp ExecuteTime;
00505 
00506   unsigned int NumberOfParticles();
00507 
00508   friend class ParticlePathFilterInternal;
00509   friend class StreaklineFilterInternal;
00510 
00511   static const double Epsilon;
00512 
00513 };
00514 
00515 #endif