VTK  9.0.1
vtkLagrangianParticle.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkLagrangianParticle.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 =========================================================================*/
29 #ifndef vtkLagrangianParticle_h
30 #define vtkLagrangianParticle_h
31 
32 #include "vtkFiltersFlowPathsModule.h" // For export macro
33 #include "vtkNew.h" // For vtkNew
34 #include "vtkSystemIncludes.h" // For PrintSelf signature and vtkType
35 
36 #include <vector>
37 
40 class vtkDataSet;
41 class vtkGenericCell;
42 class vtkIdList;
43 class vtkPointData;
45 
46 class VTKFILTERSFLOWPATHS_EXPORT vtkLagrangianParticle
47 {
48 public:
65  typedef enum ParticleTermination
66  {
67  PARTICLE_TERMINATION_NOT_TERMINATED = 0,
73  PARTICLE_TERMINATION_OUT_OF_TIME
74  } ParticleTermination;
75 
87  typedef enum SurfaceInteraction
88  {
89  SURFACE_INTERACTION_NO_INTERACTION = 0,
94  SURFACE_INTERACTION_OTHER
95  } SurfaceInteraction;
96 
106  vtkLagrangianParticle(int numberOfVariables, vtkIdType seedId, vtkIdType particleId,
107  vtkIdType seedArrayTupleIndex, double integrationTime, vtkPointData* seedData, int weightsSize,
108  int numberOfTrackedUserData);
109 
114  static vtkLagrangianParticle* NewInstance(int numberOfVariables, vtkIdType seedId,
115  vtkIdType particleId, vtkIdType seedArrayTupleIndex, double integrationTime,
116  vtkPointData* seedData, int weightsSize, int numberOfTrackedUserData,
117  vtkIdType numberOfSteps = 0, double previousIntegrationTime = 0);
118 
125  vtkLagrangianParticle* NewParticle(vtkIdType particleId);
126 
130  vtkLagrangianParticle* CloneParticle();
131 
135  virtual ~vtkLagrangianParticle();
136 
138 
142  inline double* GetPrevEquationVariables() { return this->PrevEquationVariables.data(); }
144 
146 
158  inline double* GetEquationVariables() { return this->EquationVariables.data(); }
160 
162 
167  inline double* GetNextEquationVariables() { return this->NextEquationVariables.data(); }
169 
171 
176  inline double* GetPrevPosition() { return this->PrevEquationVariables.data(); }
178 
180 
185  inline double* GetPosition() { return this->EquationVariables.data(); }
187 
189 
194  inline double* GetNextPosition() { return this->NextEquationVariables.data(); }
196 
198 
203  inline double* GetPrevVelocity() { return this->PrevVelocity; }
205 
207 
212  inline double* GetVelocity() { return this->Velocity; }
214 
216 
221  inline double* GetNextVelocity() { return this->NextVelocity; }
223 
225 
230  inline double* GetPrevUserVariables() { return this->PrevUserVariables; }
232 
234 
239  inline double* GetUserVariables() { return this->UserVariables; }
241 
243 
248  inline double* GetNextUserVariables() { return this->NextUserVariables; }
250 
252 
256  inline std::vector<double>& GetPrevTrackedUserData() { return this->PrevTrackedUserData; }
258 
260 
270  inline std::vector<double>& GetTrackedUserData() { return this->TrackedUserData; }
272 
274 
278  inline std::vector<double>& GetNextTrackedUserData() { return this->NextTrackedUserData; }
280 
282 
289  inline vtkLagrangianThreadedData* GetThreadedData() { return this->ThreadedData; }
290  inline void SetThreadedData(vtkLagrangianThreadedData* threadedData)
291  {
292  this->ThreadedData = threadedData;
293  }
295 
302  virtual void MoveToNextPosition();
303 
307  virtual vtkIdType GetId();
308 
310 
314  virtual void SetParentId(vtkIdType parentId);
315  virtual vtkIdType GetParentId();
317 
322  virtual vtkIdType GetSeedId();
323 
327  virtual int GetNumberOfVariables();
328 
332  virtual int GetNumberOfUserVariables();
333 
337  virtual vtkPointData* GetSeedData();
338 
343  virtual vtkIdType GetSeedArrayTupleIndex() const;
344 
349  double* GetLastWeights();
350 
354  vtkIdType GetLastCellId();
355 
359  double* GetLastCellPosition();
360 
364  vtkDataSet* GetLastDataSet();
365 
369  vtkAbstractCellLocator* GetLastLocator();
370 
374  vtkIdType GetLastSurfaceCellId();
375 
379  vtkDataSet* GetLastSurfaceDataSet();
380 
384  void SetLastCell(vtkAbstractCellLocator* locator, vtkDataSet* dataset, vtkIdType cellId,
385  double lastCellPosition[3]);
386 
390  void SetLastSurfaceCell(vtkDataSet* dataset, vtkIdType cellId);
391 
395  virtual vtkIdType GetNumberOfSteps();
396 
398 
403  virtual void SetTermination(int termination);
404  virtual int GetTermination();
406 
408 
413  virtual void SetInteraction(int interaction);
414  virtual int GetInteraction();
416 
418 
421  virtual void SetUserFlag(int flag);
422  virtual int GetUserFlag();
424 
426 
431  virtual void SetPInsertPreviousPosition(bool val);
432  virtual bool GetPInsertPreviousPosition();
434 
436 
441  virtual void SetPManualShift(bool val);
442  virtual bool GetPManualShift();
444 
448  virtual double& GetStepTimeRef();
449 
453  virtual double GetIntegrationTime();
454 
458  virtual double GetPrevIntegrationTime();
459 
468  virtual void SetIntegrationTime(double time);
469 
473  double GetPositionVectorMagnitude();
474 
478  virtual void PrintSelf(ostream& os, vtkIndent indent);
479 
480 protected:
482  vtkLagrangianParticle() = delete;
483  void operator=(const vtkLagrangianParticle&) = delete;
484 
485  std::vector<double> PrevEquationVariables;
486  double* PrevVelocity;
488 
489  std::vector<double> EquationVariables;
490  double* Velocity;
491  double* UserVariables;
492 
493  std::vector<double> NextEquationVariables;
494  double* NextVelocity;
496 
497  std::vector<double> PrevTrackedUserData;
498  std::vector<double> TrackedUserData;
499  std::vector<double> NextTrackedUserData;
500 
501  vtkLagrangianThreadedData* ThreadedData = nullptr;
502 
509 
513  double LastCellPosition[3];
515  std::vector<double> LastWeights;
516 
517  double StepTime;
522  int UserFlag;
526 
527  // Parallel related flags
530 };
531 
532 #endif
533 // VTK-HeaderTest-Exclude: vtkLagrangianParticle.h
vtkLagrangianParticle::GetNextEquationVariables
double * GetNextEquationVariables()
Get a pointer to the particle variables array at its next position.
Definition: vtkLagrangianParticle.h:167
vtkLagrangianParticle::PrevVelocity
double * PrevVelocity
Definition: vtkLagrangianParticle.h:486
vtkLagrangianParticle::GetTrackedUserData
std::vector< double > & GetTrackedUserData()
Get a reference to TrackedUserData.
Definition: vtkLagrangianParticle.h:270
vtkLagrangianParticle::PrevIntegrationTime
double PrevIntegrationTime
Definition: vtkLagrangianParticle.h:519
vtkLagrangianParticle::SeedId
vtkIdType SeedId
Definition: vtkLagrangianParticle.h:505
vtkLagrangianParticle::PARTICLE_TERMINATION_OUT_OF_STEPS
@ PARTICLE_TERMINATION_OUT_OF_STEPS
Definition: vtkLagrangianParticle.h:72
vtkLagrangianParticle::LastDataSet
vtkDataSet * LastDataSet
Definition: vtkLagrangianParticle.h:511
vtkPointData
represent and manipulate point attribute data
Definition: vtkPointData.h:31
vtkLagrangianParticle::GetNextUserVariables
double * GetNextUserVariables()
Get a pointer to the next user variables.
Definition: vtkLagrangianParticle.h:248
vtkIdType
int vtkIdType
Definition: vtkType.h:338
vtkLagrangianParticle::NumberOfVariables
int NumberOfVariables
Definition: vtkLagrangianParticle.h:525
vtkLagrangianParticle::GetPrevPosition
double * GetPrevPosition()
Get a pointer to the previous particle position.
Definition: vtkLagrangianParticle.h:176
vtkLagrangianParticle::PARTICLE_TERMINATION_FLIGHT_TERMINATED
@ PARTICLE_TERMINATION_FLIGHT_TERMINATED
Definition: vtkLagrangianParticle.h:69
vtkLagrangianParticle::Id
vtkIdType Id
Definition: vtkLagrangianParticle.h:503
vtkLagrangianParticle::GetEquationVariables
double * GetEquationVariables()
Get a pointer to the particle variables array.
Definition: vtkLagrangianParticle.h:158
vtkLagrangianParticle::LastSurfaceCellId
vtkIdType LastSurfaceCellId
Definition: vtkLagrangianParticle.h:524
vtkLagrangianParticle::PARTICLE_TERMINATION_OUT_OF_DOMAIN
@ PARTICLE_TERMINATION_OUT_OF_DOMAIN
Definition: vtkLagrangianParticle.h:71
vtkLagrangianParticle::PrevTrackedUserData
std::vector< double > PrevTrackedUserData
Definition: vtkLagrangianParticle.h:497
vtkLagrangianParticle::EquationVariables
std::vector< double > EquationVariables
Definition: vtkLagrangianParticle.h:489
vtkLagrangianParticle::GetPrevTrackedUserData
std::vector< double > & GetPrevTrackedUserData()
Get a reference to PrevTrackedUserData See GetTrackedUserData for an explanation on how to use it.
Definition: vtkLagrangianParticle.h:256
vtkLagrangianParticle::GetPrevEquationVariables
double * GetPrevEquationVariables()
Get a pointer to Particle variables at its previous position See GetEquationVariables for content des...
Definition: vtkLagrangianParticle.h:142
vtkLagrangianParticle::UserVariables
double * UserVariables
Definition: vtkLagrangianParticle.h:491
vtkX3D::time
@ time
Definition: vtkX3D.h:503
vtkLagrangianParticle::NextEquationVariables
std::vector< double > NextEquationVariables
Definition: vtkLagrangianParticle.h:493
vtkLagrangianParticle::GetNextPosition
double * GetNextPosition()
Get a pointer to the next particle position.
Definition: vtkLagrangianParticle.h:194
vtkLagrangianParticle::GetNextVelocity
double * GetNextVelocity()
Get a pointer to the next particle velocity.
Definition: vtkLagrangianParticle.h:221
vtkLagrangianParticle::Velocity
double * Velocity
Definition: vtkLagrangianParticle.h:490
vtkLagrangianParticle::SurfaceInteraction
SurfaceInteraction
An enum to inform about a surface interaction SURFACE_INTERACTION_NO_INTERACTION = 0,...
Definition: vtkLagrangianParticle.h:87
vtkLagrangianParticle::GetVelocity
double * GetVelocity()
Get a pointer to the particle velocity.
Definition: vtkLagrangianParticle.h:212
vtkLagrangianParticle::GetPosition
double * GetPosition()
Get a pointer to the particle position.
Definition: vtkLagrangianParticle.h:185
vtkLagrangianParticle::GetNextTrackedUserData
std::vector< double > & GetNextTrackedUserData()
Get a reference to NextTrackedUserData See GetTrackedUserData for an explanation on how to use it.
Definition: vtkLagrangianParticle.h:278
vtkLagrangianParticle::PrevUserVariables
double * PrevUserVariables
Definition: vtkLagrangianParticle.h:487
vtkLagrangianParticle::ParentId
vtkIdType ParentId
Definition: vtkLagrangianParticle.h:504
vtkLagrangianParticle::TrackedUserData
std::vector< double > TrackedUserData
Definition: vtkLagrangianParticle.h:498
vtkLagrangianParticle::LastWeights
std::vector< double > LastWeights
Definition: vtkLagrangianParticle.h:515
vtkLagrangianParticle::Interaction
int Interaction
Definition: vtkLagrangianParticle.h:521
vtkLagrangianParticle::SeedData
vtkPointData * SeedData
Definition: vtkLagrangianParticle.h:508
vtkLagrangianParticle::UserFlag
int UserFlag
Definition: vtkLagrangianParticle.h:522
vtkLagrangianParticle::PARTICLE_TERMINATION_SURF_TERMINATED
@ PARTICLE_TERMINATION_SURF_TERMINATED
Definition: vtkLagrangianParticle.h:68
vtkIndent
a simple class to control print indentation
Definition: vtkIndent.h:33
vtkLagrangianParticle::SURFACE_INTERACTION_BREAK
@ SURFACE_INTERACTION_BREAK
Definition: vtkLagrangianParticle.h:91
vtkLagrangianThreadedData
struct to hold a user data
Definition: vtkLagrangianThreadedData.h:34
vtkLagrangianParticle::PARTICLE_TERMINATION_SURF_BREAK
@ PARTICLE_TERMINATION_SURF_BREAK
Definition: vtkLagrangianParticle.h:70
vtkLagrangianParticle::NextUserVariables
double * NextUserVariables
Definition: vtkLagrangianParticle.h:495
vtkIdList
list of point or cell ids
Definition: vtkIdList.h:30
vtkLagrangianParticle::StepTime
double StepTime
Definition: vtkLagrangianParticle.h:517
vtkLagrangianParticle::PrevEquationVariables
std::vector< double > PrevEquationVariables
Definition: vtkLagrangianParticle.h:485
vtkLagrangianParticle::NextVelocity
double * NextVelocity
Definition: vtkLagrangianParticle.h:494
vtkLagrangianParticle::NumberOfSteps
vtkIdType NumberOfSteps
Definition: vtkLagrangianParticle.h:506
vtkLagrangianParticle::LastCellId
vtkIdType LastCellId
Definition: vtkLagrangianParticle.h:512
vtkLagrangianParticle::PInsertPreviousPosition
bool PInsertPreviousPosition
Definition: vtkLagrangianParticle.h:528
vtkLagrangianParticle::GetUserVariables
double * GetUserVariables()
Get a pointer to the user variables.
Definition: vtkLagrangianParticle.h:239
vtkLagrangianParticle::Termination
int Termination
Definition: vtkLagrangianParticle.h:520
vtkAbstractCellLocator
an abstract base class for locators which find cells
Definition: vtkAbstractCellLocator.h:48
vtkLagrangianParticle::SetThreadedData
void SetThreadedData(vtkLagrangianThreadedData *threadedData)
Definition: vtkLagrangianParticle.h:290
vtkLagrangianParticle::SURFACE_INTERACTION_TERMINATED
@ SURFACE_INTERACTION_TERMINATED
Definition: vtkLagrangianParticle.h:90
vtkLagrangianParticle
Basis class for Lagrangian particles.
Definition: vtkLagrangianParticle.h:46
vtkLagrangianParticle::GetPrevUserVariables
double * GetPrevUserVariables()
Get a pointer to the previous user variables.
Definition: vtkLagrangianParticle.h:230
vtkLagrangianParticle::LastSurfaceDataSet
vtkDataSet * LastSurfaceDataSet
Definition: vtkLagrangianParticle.h:523
vtkDataSet
abstract class to specify dataset behavior
Definition: vtkDataSet.h:56
vtkNew.h
vtkLagrangianParticle::LastLocator
vtkAbstractCellLocator * LastLocator
Definition: vtkLagrangianParticle.h:510
vtkLagrangianParticle::ParticleTermination
ParticleTermination
An enum to inform about a reason for termination PARTICLE_TERMINATION_NOT_TERMINATED = 0,...
Definition: vtkLagrangianParticle.h:65
vtkLagrangianParticle::IntegrationTime
double IntegrationTime
Definition: vtkLagrangianParticle.h:518
vtkGenericCell
provides thread-safe access to cells
Definition: vtkGenericCell.h:36
vtkLagrangianParticle::GetPrevVelocity
double * GetPrevVelocity()
Get a pointer to the previous particle velocity.
Definition: vtkLagrangianParticle.h:203
vtkLagrangianParticle::SURFACE_INTERACTION_PASS
@ SURFACE_INTERACTION_PASS
Definition: vtkLagrangianParticle.h:93
vtkLagrangianParticle::WeightsSize
int WeightsSize
Definition: vtkLagrangianParticle.h:514
vtkLagrangianParticle::GetThreadedData
vtkLagrangianThreadedData * GetThreadedData()
Get/Set a pointer to a vtkLagrangianThreadedData that is considered to be local to the thread.
Definition: vtkLagrangianParticle.h:289
vtkLagrangianParticle::NextTrackedUserData
std::vector< double > NextTrackedUserData
Definition: vtkLagrangianParticle.h:499
vtkLagrangianParticle::PManualShift
bool PManualShift
Definition: vtkLagrangianParticle.h:529
vtkLagrangianParticle::SURFACE_INTERACTION_BOUNCE
@ SURFACE_INTERACTION_BOUNCE
Definition: vtkLagrangianParticle.h:92
vtkLagrangianParticle::SeedArrayTupleIndex
vtkIdType SeedArrayTupleIndex
Definition: vtkLagrangianParticle.h:507
vtkSystemIncludes.h
vtkBilinearQuadIntersection
Class to perform non planar quad intersection.
Definition: vtkBilinearQuadIntersection.h:30