VTK
|
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 double SimulationTime; 00074 // These are useful to track for debugging etc 00075 int ErrorCode; 00076 float age; 00077 // these are needed across time steps to compute vorticity 00078 float rotation; 00079 float angularVel; 00080 float time; 00081 float speed; 00082 00083 vtkIdType PointId; //once the partice is added, PointId is valid 00084 } ParticleInformation; 00085 00086 typedef std::vector<ParticleInformation> ParticleVector; 00087 typedef ParticleVector::iterator ParticleIterator; 00088 typedef std::list<ParticleInformation> ParticleDataList; 00089 typedef ParticleDataList::iterator ParticleListIterator; 00090 }; 00091 //ETX 00092 00093 class VTKFILTERSFLOWPATHS_EXPORT vtkParticleTracerBase : public vtkPolyDataAlgorithm 00094 { 00095 public: 00096 enum Solvers 00097 { 00098 RUNGE_KUTTA2, 00099 RUNGE_KUTTA4, 00100 RUNGE_KUTTA45, 00101 NONE, 00102 UNKNOWN 00103 }; 00104 00105 vtkTypeMacro(vtkParticleTracerBase,vtkPolyDataAlgorithm) 00106 void PrintSelf(ostream& os, vtkIndent indent); 00107 void PrintParticleHistories(); 00108 00110 00112 vtkGetMacro(ComputeVorticity, bool); 00113 void SetComputeVorticity(bool); 00115 00117 00119 vtkGetMacro(TerminalSpeed, double); 00120 void SetTerminalSpeed(double); 00122 00124 00126 vtkGetMacro(RotationScale, double); 00127 void SetRotationScale(double); 00129 00131 00133 vtkSetMacro(IgnorePipelineTime, int); 00134 vtkGetMacro(IgnorePipelineTime, int); 00135 vtkBooleanMacro(IgnorePipelineTime, int); 00137 00139 00146 vtkGetMacro(ForceReinjectionEveryNSteps,int); 00147 void SetForceReinjectionEveryNSteps(int); 00149 00151 00155 void SetTerminationTime(double t); 00156 vtkGetMacro(TerminationTime,double); 00158 00159 void SetIntegrator(vtkInitialValueProblemSolver *); 00160 vtkGetObjectMacro ( Integrator, vtkInitialValueProblemSolver ); 00161 00162 void SetIntegratorType(int type); 00163 int GetIntegratorType(); 00164 00166 00170 vtkGetMacro(StartTime, double); 00171 void SetStartTime(double t); 00173 00175 00181 vtkSetMacro(StaticSeeds,int); 00182 vtkGetMacro(StaticSeeds,int); 00184 00186 00192 vtkSetMacro(StaticMesh,int); 00193 vtkGetMacro(StaticMesh,int); 00195 00197 00201 virtual void SetParticleWriter(vtkAbstractParticleWriter *pw); 00202 vtkGetObjectMacro(ParticleWriter, vtkAbstractParticleWriter); 00204 00206 00208 vtkSetStringMacro(ParticleFileName); 00209 vtkGetStringMacro(ParticleFileName); 00211 00213 00215 vtkSetMacro(EnableParticleWriting,int); 00216 vtkGetMacro(EnableParticleWriting,int); 00217 vtkBooleanMacro(EnableParticleWriting,int); 00219 00221 00223 vtkSetMacro(DisableResetCache,int); 00224 vtkGetMacro(DisableResetCache,int); 00225 vtkBooleanMacro(DisableResetCache,int); 00227 00229 00230 void AddSourceConnection(vtkAlgorithmOutput* input); 00231 void RemoveAllSources(); 00233 00234 protected: 00235 vtkSmartPointer<vtkPolyData> Output; //managed by child classes 00236 vtkSmartPointer<vtkPointData> ProtoPD; 00237 vtkIdType UniqueIdCounter;// global Id counter used to give particles a stamp 00238 vtkParticleTracerBaseNamespace::ParticleDataList ParticleHistories; 00239 vtkSmartPointer<vtkPointData> ParticlePointData; //the current particle point data consistent 00240 //with particle history 00241 //Everything related to time 00242 int IgnorePipelineTime; //whether to use the pipeline time for termination 00243 int DisableResetCache; //whether to enable ResetCache() method 00244 00245 00246 vtkParticleTracerBase(); 00247 virtual ~vtkParticleTracerBase(); 00248 00249 // 00250 // Make sure the pipeline knows what type we expect as input 00251 // 00252 virtual int FillInputPortInformation(int port, vtkInformation* info); 00253 00254 // 00255 // The usual suspects 00256 // 00257 virtual int ProcessRequest(vtkInformation* request, 00258 vtkInformationVector** inputVector, 00259 vtkInformationVector* outputVector); 00260 00261 // 00262 // Store any information we need in the output and fetch what we can 00263 // from the input 00264 // 00265 virtual int RequestInformation(vtkInformation* request, 00266 vtkInformationVector** inputVector, 00267 vtkInformationVector* outputVector); 00268 00269 // 00270 // Compute input time steps given the output step 00271 // 00272 virtual int RequestUpdateExtent(vtkInformation* request, 00273 vtkInformationVector** inputVector, 00274 vtkInformationVector* outputVector); 00275 00276 // 00277 // what the pipeline calls for each time step 00278 // 00279 virtual int RequestData(vtkInformation* request, 00280 vtkInformationVector** inputVector, 00281 vtkInformationVector* outputVector); 00282 00283 // 00284 // these routines are internally called to actually generate the output 00285 // 00286 virtual int ProcessInput(vtkInformationVector** inputVector); 00287 00288 // This is the main part of the algorithm: 00289 // * move all the particles one step 00290 // * Reinject particles (by adding them to this->ParticleHistories) 00291 // either at the beginning or at the end of each step (modulo this->ForceReinjectionEveryNSteps) 00292 // * Output a polydata representing the moved particles 00293 // Note that if the starting and the ending time coincide, the polydata is still valid. 00294 virtual vtkPolyData* Execute(vtkInformationVector** inputVector); 00295 00296 // the RequestData will call these methods in turn 00297 virtual void Initialize(){} //the first iteration 00298 virtual int OutputParticles(vtkPolyData* poly)=0; //every iteration 00299 virtual void Finalize(){} //the last iteration 00300 00301 // 00302 // Initialization of input (vector-field) geometry 00303 // 00304 int InitializeInterpolator(); 00305 int UpdateDataCache(vtkDataObject *td); 00306 00308 00310 void TestParticles( 00311 vtkParticleTracerBaseNamespace::ParticleVector &candidates, 00312 vtkParticleTracerBaseNamespace::ParticleVector &passed, 00313 int &count); 00315 00316 void TestParticles( 00317 vtkParticleTracerBaseNamespace::ParticleVector &candidates, std::vector<int> &passed); 00318 00320 00324 virtual void AssignSeedsToProcessors(double time, 00325 vtkDataSet *source, int sourceID, int ptId, 00326 vtkParticleTracerBaseNamespace::ParticleVector &localSeedPoints, 00327 int &localAssignedCount); 00329 00331 00333 virtual void AssignUniqueIds( 00334 vtkParticleTracerBaseNamespace::ParticleVector &localSeedPoints); 00336 00338 00340 void UpdateParticleList( 00341 vtkParticleTracerBaseNamespace::ParticleVector &candidates); 00343 00346 virtual void UpdateParticleListFromOtherProcesses(){} 00347 00349 00350 void IntegrateParticle( 00351 vtkParticleTracerBaseNamespace::ParticleListIterator &it, 00352 double currenttime, double terminationtime, 00353 vtkInitialValueProblemSolver* integrator); 00355 00356 // if the particle is added to send list, then returns value is 1, 00357 // if it is kept on this process after a retry return value is 0 00358 virtual bool SendParticleToAnotherProcess( 00359 vtkParticleTracerBaseNamespace::ParticleInformation &, 00360 vtkParticleTracerBaseNamespace::ParticleInformation &, vtkPointData*) 00361 { 00362 return true; 00363 } 00364 00366 00370 bool ComputeDomainExitLocation( 00371 double pos[4], double p2[4], double intersection[4], 00372 vtkGenericCell *cell); 00374 00375 // 00376 // Scalar arrays that are generated as each particle is updated 00377 // 00378 void CreateProtoPD(vtkDataObject* input); 00379 00380 vtkFloatArray* GetParticleAge(vtkPointData*); 00381 vtkIntArray* GetParticleIds(vtkPointData*); 00382 vtkCharArray* GetParticleSourceIds(vtkPointData*); 00383 vtkIntArray* GetInjectedPointIds(vtkPointData*); 00384 vtkIntArray* GetInjectedStepIds(vtkPointData*); 00385 vtkIntArray* GetErrorCodeArr(vtkPointData*); 00386 vtkFloatArray* GetParticleVorticity(vtkPointData*); 00387 vtkFloatArray* GetParticleRotation(vtkPointData*); 00388 vtkFloatArray* GetParticleAngularVel(vtkPointData*); 00389 00390 // utility function we use to test if a point is inside any of our local datasets 00391 bool InsideBounds(double point[]); 00392 00393 void CalculateVorticity( vtkGenericCell* cell, double pcoords[3], 00394 vtkDoubleArray* cellVectors, double vorticity[3] ); 00395 00396 //------------------------------------------------------ 00397 00398 00399 double GetCacheDataTime(int i); 00400 double GetCacheDataTime(); 00401 00402 virtual void ResetCache(); 00403 void AddParticle(vtkParticleTracerBaseNamespace::ParticleInformation &info, double* velocity); 00404 00406 00409 virtual bool IsPointDataValid(vtkDataObject* input); 00410 bool IsPointDataValid(vtkCompositeDataSet* input, std::vector<std::string>& arrayNames); 00411 void GetPointDataArrayNames(vtkDataSet* input, std::vector<std::string>& names); 00413 00414 vtkGetMacro(ReinjectionCounter, int); 00415 vtkGetMacro(CurrentTimeValue, double); 00416 00419 virtual void InitializeExtraPointDataArrays(vtkPointData* vtkNotUsed(outputPD)) {} 00420 00421 virtual void AppendToExtraPointDataArrays(vtkParticleTracerBaseNamespace::ParticleInformation &) {} 00422 private: 00424 void SetInterpolatorPrototype(vtkAbstractInterpolatedVelocityField*) {} 00425 00427 00434 bool RetryWithPush( 00435 vtkParticleTracerBaseNamespace::ParticleInformation &info, double* point1,double delT, int subSteps); 00437 00438 bool SetTerminationTimeNoModify(double t); 00439 00440 //Parameters of tracing 00441 vtkInitialValueProblemSolver* Integrator; 00442 double IntegrationStep; 00443 double MaximumError; 00444 bool ComputeVorticity; 00445 double RotationScale; 00446 double TerminalSpeed; 00447 00448 // A counter to keep track of how many times we reinjected 00449 int ReinjectionCounter; 00450 00451 // Important for Caching of Cells/Ids/Weights etc 00452 int AllFixedGeometry; 00453 int StaticMesh; 00454 int StaticSeeds; 00455 00456 std::vector<double> InputTimeValues; 00457 double StartTime; 00458 double TerminationTime; 00459 double CurrentTimeValue; 00460 00461 int StartTimeStep; //InputTimeValues[StartTimeStep] <= StartTime <= InputTimeValues[StartTimeStep+1] 00462 int CurrentTimeStep; 00463 int TerminationTimeStep; //computed from start time 00464 bool FirstIteration; 00465 00466 //Innjection parameters 00467 int ForceReinjectionEveryNSteps; 00468 vtkTimeStamp ParticleInjectionTime; 00469 bool HasCache; 00470 00471 // Particle writing to disk 00472 vtkAbstractParticleWriter *ParticleWriter; 00473 char *ParticleFileName; 00474 int EnableParticleWriting; 00475 00476 00477 // The main lists which are held during operation- between time step updates 00478 vtkParticleTracerBaseNamespace::ParticleVector LocalSeeds; 00479 00480 00481 // The velocity interpolator 00482 vtkSmartPointer<vtkTemporalInterpolatedVelocityField> Interpolator; 00483 vtkAbstractInterpolatedVelocityField * InterpolatorPrototype; 00484 00485 // Data for time step CurrentTimeStep-1 and CurrentTimeStep 00486 vtkSmartPointer<vtkMultiBlockDataSet> CachedData[2]; 00487 00488 // Cache bounds info for each dataset we will use repeatedly 00489 typedef struct { 00490 double b[6]; 00491 } bounds; 00492 std::vector<bounds> CachedBounds[2]; 00493 00494 // temporary variables used by Exeucte(), for convenience only 00495 00496 vtkSmartPointer<vtkPoints> OutputCoordinates; 00497 vtkSmartPointer<vtkFloatArray> ParticleAge; 00498 vtkSmartPointer<vtkIntArray> ParticleIds; 00499 vtkSmartPointer<vtkCharArray> ParticleSourceIds; 00500 vtkSmartPointer<vtkIntArray> InjectedPointIds; 00501 vtkSmartPointer<vtkIntArray> InjectedStepIds; 00502 vtkSmartPointer<vtkIntArray> ErrorCode; 00503 vtkSmartPointer<vtkFloatArray> ParticleVorticity; 00504 vtkSmartPointer<vtkFloatArray> ParticleRotation; 00505 vtkSmartPointer<vtkFloatArray> ParticleAngularVel; 00506 vtkSmartPointer<vtkDoubleArray> CellVectors; 00507 vtkSmartPointer<vtkPointData> OutputPointData; 00508 vtkSmartPointer<vtkDataSet> DataReferenceT[2]; 00509 vtkSmartPointer<vtkCellArray> ParticleCells; 00510 00511 vtkParticleTracerBase(const vtkParticleTracerBase&); // Not implemented. 00512 void operator=(const vtkParticleTracerBase&); // Not implemented. 00513 vtkTimeStamp ExecuteTime; 00514 00515 unsigned int NumberOfParticles(); 00516 00517 friend class ParticlePathFilterInternal; 00518 friend class StreaklineFilterInternal; 00519 00520 static const double Epsilon; 00521 00522 }; 00523 00524 #endif