VTK  9.3.20240328
vtkLagrangianParticleTracker.h
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2 // SPDX-License-Identifier: BSD-3-Clause
79 #ifndef vtkLagrangianParticleTracker_h
80 #define vtkLagrangianParticleTracker_h
81 
82 #include "vtkBoundingBox.h" // For cached bounds
83 #include "vtkDataObjectAlgorithm.h"
84 #include "vtkFiltersFlowPathsModule.h" // For export macro
85 #include "vtkSmartPointer.h" // For smart pointer
86 
87 #include <atomic> // for atomic
88 #include <mutex> // for mutexes
89 #include <queue> // for particle queue
90 
91 VTK_ABI_NAMESPACE_BEGIN
92 class vtkBoundingBox;
93 class vtkCellArray;
94 class vtkDataSet;
95 class vtkDoubleArray;
96 class vtkIdList;
97 class vtkInformation;
103 class vtkPointData;
104 class vtkPoints;
105 class vtkPolyData;
106 class vtkPolyLine;
107 struct IntegratingFunctor;
109 
110 class VTKFILTERSFLOWPATHS_EXPORT vtkLagrangianParticleTracker : public vtkDataObjectAlgorithm
111 {
112 public:
114  void PrintSelf(ostream& os, vtkIndent indent) override;
116 
118  {
119  STEP_CUR_CELL_LENGTH = 1,
120  STEP_CUR_CELL_VEL_DIR = 3,
121  STEP_CUR_CELL_DIV_THEO = 5
122  } CellLengthComputation;
123 
125 
132 
134 
141 
143 
148  vtkSetMacro(GeneratePolyVertexInteractionOutput, bool);
149  vtkGetMacro(GeneratePolyVertexInteractionOutput, bool);
151 
153 
166  vtkSetMacro(CellLengthComputationMode, int);
167  vtkGetMacro(CellLengthComputationMode, int);
169 
171 
174  vtkSetMacro(StepFactor, double);
175  vtkGetMacro(StepFactor, double);
177 
179 
182  vtkSetMacro(StepFactorMin, double);
183  vtkGetMacro(StepFactorMin, double);
185 
187 
190  vtkSetMacro(StepFactorMax, double);
191  vtkGetMacro(StepFactorMax, double);
193 
195 
198  vtkSetMacro(MaximumNumberOfSteps, int);
199  vtkGetMacro(MaximumNumberOfSteps, int);
201 
203 
207  vtkSetMacro(MaximumIntegrationTime, double);
208  vtkGetMacro(MaximumIntegrationTime, double);
210 
212 
218  vtkSetMacro(AdaptiveStepReintegration, bool);
219  vtkGetMacro(AdaptiveStepReintegration, bool);
220  vtkBooleanMacro(AdaptiveStepReintegration, bool);
222 
224 
228  vtkSetMacro(GenerateParticlePathsOutput, bool);
229  vtkGetMacro(GenerateParticlePathsOutput, bool);
230  vtkBooleanMacro(GenerateParticlePathsOutput, bool);
232 
234 
243 
248 
250 
259 
264 
269 
274 
279 
283  int RequestData(vtkInformation* request, vtkInformationVector** inputVector,
284  vtkInformationVector* outputVector) override;
285 
290  vtkMTimeType GetMTime() override;
291 
297 
298 protected:
301 
302  virtual bool InitializeFlow(vtkDataObject* flow, vtkBoundingBox* bounds);
303  virtual bool InitializeParticles(const vtkBoundingBox* bounds, vtkDataSet* seeds,
304  std::queue<vtkLagrangianParticle*>& particles, vtkPointData* seedData);
305  virtual void GenerateParticles(const vtkBoundingBox* bounds, vtkDataSet* seeds,
306  vtkDataArray* initialVelocities, vtkDataArray* initialIntegrationTimes, vtkPointData* seedData,
307  int nVar, std::queue<vtkLagrangianParticle*>& particles);
308  virtual bool UpdateSurfaceCacheIfNeeded(vtkDataObject*& surfaces);
309  virtual void InitializeSurface(vtkDataObject*& surfaces);
310 
314  virtual bool InitializePathsOutput(
315  vtkPointData* seedData, vtkIdType numberOfSeeds, vtkPolyData*& particlePathsOutput);
316 
321  vtkPointData* seedData, vtkDataObject* surfaces, vtkDataObject*& interractionOutput);
322 
323  virtual bool FinalizeOutputs(vtkPolyData* particlePathsOutput, vtkDataObject* interactionOutput);
324 
325  static void InsertPolyVertexCell(vtkPolyData* polydata);
326  static void InsertVertexCells(vtkPolyData* polydata);
327 
328  virtual void GetParticleFeed(std::queue<vtkLagrangianParticle*>& particleQueue);
329 
334  std::queue<vtkLagrangianParticle*>&, vtkPolyData* particlePathsOutput,
335  vtkPolyLine* particlePath, vtkDataObject* interactionOutput);
336 
340  void InsertPathOutputPoint(vtkLagrangianParticle* particle, vtkPolyData* particlePathsOutput,
341  vtkIdList* particlePathPointId, bool prev = false);
342 
347  unsigned int interactedSurfaceFlatIndex, vtkDataObject* interactionOutput);
348 
354 
358  bool ComputeNextStep(vtkInitialValueProblemSolver* integrator, double* xprev, double* xnext,
359  double t, double& delT, double& delTActual, double minStep, double maxStep, double cellLength,
360  int& integrationRes, vtkLagrangianParticle* particle);
361 
366  virtual void DeleteParticle(vtkLagrangianParticle* particle);
367 
370 
372  double StepFactor;
378  bool GenerateParticlePathsOutput = true;
380  std::atomic<vtkIdType> ParticleCounter;
381  std::atomic<vtkIdType> IntegratedParticleCounter;
384 
385  // internal parameters use for step computation
388 
389  // Cache related parameters
393  bool FlowCacheInvalid = true;
396  bool SurfaceCacheInvalid = true;
397 
398  std::mutex ProgressMutex;
399  friend struct IntegratingFunctor;
400 
402 
403 private:
405  void operator=(const vtkLagrangianParticleTracker&) = delete;
406 };
407 
408 VTK_ABI_NAMESPACE_END
409 #endif
Proxy object to connect input/output ports.
Fast, simple class for representing and operating on 3D bounds.
object to represent cell connectivity
Definition: vtkCellArray.h:285
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:154
Superclass for algorithms that produce only data object as output.
general representation of visualization data
abstract class to specify dataset behavior
Definition: vtkDataSet.h:165
dynamic, self-adjusting array of double
list of point or cell ids
Definition: vtkIdList.h:132
a simple class to control print indentation
Definition: vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Integrate a set of ordinary differential equations (initial value problem) in time.
vtkFunctionSet abstract implementation to be used in the vtkLagrangianParticleTracker integrator.
Filter to inject and track particles in a flow.
vtkDataObject * GetSurface()
Specify the source object used to compute surface interaction with Note that this method does not con...
void SetSourceConnection(vtkAlgorithmOutput *algOutput)
Specify the source object used to generate particle initial position (seeds).
virtual bool UpdateSurfaceCacheIfNeeded(vtkDataObject *&surfaces)
virtual vtkIdType GetNewParticleId()
Get an unique id for a particle This method is thread safe.
void InsertPathOutputPoint(vtkLagrangianParticle *particle, vtkPolyData *particlePathsOutput, vtkIdList *particlePathPointId, bool prev=false)
This method is thread safe.
int FillInputPortInformation(int port, vtkInformation *info) override
Declare input port type.
virtual void GetParticleFeed(std::queue< vtkLagrangianParticle * > &particleQueue)
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
Process inputs to integrate particle and generate output.
vtkSmartPointer< vtkInitialValueProblemSolver > Integrator
bool ComputeNextStep(vtkInitialValueProblemSolver *integrator, double *xprev, double *xnext, double t, double &delT, double &delTActual, double minStep, double maxStep, double cellLength, int &integrationRes, vtkLagrangianParticle *particle)
This method is thread safe.
virtual void InitializeSurface(vtkDataObject *&surfaces)
virtual void DeleteParticle(vtkLagrangianParticle *particle)
This method is thread safe Call the ParticleAboutToBeDeleted model method and delete the particle.
virtual void GenerateParticles(const vtkBoundingBox *bounds, vtkDataSet *seeds, vtkDataArray *initialVelocities, vtkDataArray *initialIntegrationTimes, vtkPointData *seedData, int nVar, std::queue< vtkLagrangianParticle * > &particles)
static vtkLagrangianParticleTracker * New()
virtual bool InitializeParticles(const vtkBoundingBox *bounds, vtkDataSet *seeds, std::queue< vtkLagrangianParticle * > &particles, vtkPointData *seedData)
vtkInitialValueProblemSolver * GetIntegrator()
Set/Get the integrator.
static void InsertPolyVertexCell(vtkPolyData *polydata)
double ComputeCellLength(vtkLagrangianParticle *particle)
Computes the cell length for the next step using the method set by CellLengthComputationMode.
vtkLagrangianBasicIntegrationModel * GetIntegrationModel()
Set/Get the integration model.
void SetIntegrator(vtkInitialValueProblemSolver *integrator)
Set/Get the integrator.
void SetSurfaceData(vtkDataObject *source)
Specify the source object used to compute surface interaction with Note that this method does not con...
static void InsertVertexCells(vtkPolyData *polydata)
vtkMTimeType GetMTime() override
Get the tracker modified time taking into account the integration model and the integrator.
std::atomic< vtkIdType > IntegratedParticleCounter
void SetSurfaceConnection(vtkAlgorithmOutput *algOutput)
Specify the object used to compute surface interaction with.
vtkLagrangianThreadedData * SerialThreadedData
void SetIntegrationModel(vtkLagrangianBasicIntegrationModel *integrationModel)
Set/Get the integration model.
virtual int Integrate(vtkInitialValueProblemSolver *integrator, vtkLagrangianParticle *, std::queue< vtkLagrangianParticle * > &, vtkPolyData *particlePathsOutput, vtkPolyLine *particlePath, vtkDataObject *interactionOutput)
This method is thread safe.
virtual bool InitializeFlow(vtkDataObject *flow, vtkBoundingBox *bounds)
int FillOutputPortInformation(int port, vtkInformation *info) override
Declare output port type.
int RequestDataObject(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
Create outputs objects.
virtual bool InitializePathsOutput(vtkPointData *seedData, vtkIdType numberOfSeeds, vtkPolyData *&particlePathsOutput)
This method is thread safe.
virtual bool InitializeInteractionOutput(vtkPointData *seedData, vtkDataObject *surfaces, vtkDataObject *&interractionOutput)
This method is thread safe.
vtkDataObject * GetSource()
Specify the source object used to generate particle initial position (seeds).
void SetSourceData(vtkDataObject *source)
Specify the source object used to generate particle initial position (seeds).
void InsertInteractionOutputPoint(vtkLagrangianParticle *particle, unsigned int interactedSurfaceFlatIndex, vtkDataObject *interactionOutput)
This method is thread safe.
virtual bool FinalizeOutputs(vtkPolyData *particlePathsOutput, vtkDataObject *interactionOutput)
~vtkLagrangianParticleTracker() override
vtkSmartPointer< vtkLagrangianBasicIntegrationModel > IntegrationModel
Basis class for Lagrangian particles.
Composite dataset that organizes datasets into blocks.
composite dataset to encapsulates pieces of dataset.
represent and manipulate point attribute data
Definition: vtkPointData.h:139
represent and manipulate 3D points
Definition: vtkPoints.h:138
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:180
cell represents a set of 1D lines
Definition: vtkPolyLine.h:139
@ info
Definition: vtkX3D.h:376
@ port
Definition: vtkX3D.h:447
struct to hold a user data
boost::graph_traits< vtkGraph * >::vertex_descriptor source(boost::graph_traits< vtkGraph * >::edge_descriptor e, vtkGraph *)
int vtkIdType
Definition: vtkType.h:315
vtkTypeUInt32 vtkMTimeType
Definition: vtkType.h:270